fmm_tree.txx 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676
  1. /**
  2. * \file fmm_tree.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 12-11-2010
  5. * \brief This file contains the implementation of the class FMM_Tree.
  6. */
  7. #include <assert.h>
  8. #include <fmm_node.hpp>
  9. #include <profile.hpp>
  10. namespace pvfmm{
  11. template <class FMM_Mat_t>
  12. void FMM_Tree<FMM_Mat_t>::Initialize(typename FMM_Node_t::NodeData* init_data) {
  13. Profile::Tic("InitTree",this->Comm(),true);{
  14. //Build octree from points.
  15. MPI_Tree<FMM_Node_t>::Initialize(init_data);
  16. Profile::Tic("InitFMMData",this->Comm(),true,5);
  17. { //Initialize FMM data.
  18. std::vector<FMM_Node_t*>& nodes=this->GetNodeList();
  19. #pragma omp parallel for
  20. for(size_t i=0;i<nodes.size();i++){
  21. if(nodes[i]->FMMData()==NULL) nodes[i]->FMMData()=new typename FMM_Mat_t::FMMData;
  22. }
  23. }
  24. Profile::Toc();
  25. }Profile::Toc();
  26. }
  27. template <class FMM_Mat_t>
  28. void FMM_Tree<FMM_Mat_t>::InitFMM_Tree(bool refine, BoundaryType bndry_) {
  29. Profile::Tic("InitFMM_Tree",this->Comm(),true);{
  30. interac_list.Initialize(this->Dim());
  31. bndry=bndry_;
  32. if(refine){
  33. //RefineTree
  34. Profile::Tic("RefineTree",this->Comm(),true,2);
  35. this->RefineTree();
  36. Profile::Toc();
  37. }
  38. //2:1 Balancing
  39. Profile::Tic("2:1Balance",this->Comm(),true,2);
  40. this->Balance21(bndry);
  41. Profile::Toc();
  42. //Redistribute nodes.
  43. Profile::Tic("Redistribute",this->Comm(),true,3);
  44. this->RedistNodes();
  45. Profile::Toc();
  46. }Profile::Toc();
  47. }
  48. template <class FMM_Mat_t>
  49. void FMM_Tree<FMM_Mat_t>::SetupFMM(FMM_Mat_t* fmm_mat_) {
  50. Profile::Tic("SetupFMM",this->Comm(),true);{
  51. bool device=true;
  52. #ifdef __INTEL_OFFLOAD
  53. Profile::Tic("InitLocks",this->Comm(),false,3);
  54. MIC_Lock::init();
  55. Profile::Toc();
  56. #endif
  57. int omp_p=omp_get_max_threads();
  58. fmm_mat=fmm_mat_;
  59. //Construct LET
  60. Profile::Tic("ConstructLET",this->Comm(),false,2);
  61. this->ConstructLET(bndry);
  62. Profile::Toc();
  63. //Set Colleagues (Needed to build U, V, W and X lists.)
  64. Profile::Tic("SetColleagues",this->Comm(),false,3);
  65. this->SetColleagues(bndry);
  66. Profile::Toc();
  67. Profile::Tic("BuildLists",this->Comm(),false,3);
  68. BuildInteracLists();
  69. Profile::Toc();
  70. Profile::Tic("CollectNodeData",this->Comm(),false,3);
  71. //Build node list.
  72. FMM_Node_t* n=dynamic_cast<FMM_Node_t*>(this->PostorderFirst());
  73. std::vector<FMM_Node_t*> all_nodes;
  74. while(n!=NULL){
  75. all_nodes.push_back(n);
  76. n=static_cast<FMM_Node_t*>(this->PostorderNxt(n));
  77. }
  78. //Collect node data into continuous array.
  79. std::vector<Vector<FMM_Node_t*> > node_lists; // TODO: Remove this parameter, not really needed
  80. fmm_mat->CollectNodeData(all_nodes, node_data_buff, node_lists);
  81. Profile::Toc();
  82. setup_data.resize(8*MAX_DEPTH);
  83. precomp_lst.resize(8);
  84. Profile::Tic("UListSetup",this->Comm(),false,3);
  85. for(size_t i=0;i<MAX_DEPTH;i++){
  86. setup_data[i+MAX_DEPTH*0].precomp_data=&precomp_lst[0];
  87. fmm_mat->U_ListSetup(setup_data[i+MAX_DEPTH*0],node_data_buff,node_lists,fmm_mat->Homogen()?(i==0?-1:MAX_DEPTH+1):i, device);
  88. }
  89. Profile::Toc();
  90. Profile::Tic("WListSetup",this->Comm(),false,3);
  91. for(size_t i=0;i<MAX_DEPTH;i++){
  92. setup_data[i+MAX_DEPTH*1].precomp_data=&precomp_lst[1];
  93. fmm_mat->W_ListSetup(setup_data[i+MAX_DEPTH*1],node_data_buff,node_lists,fmm_mat->Homogen()?(i==0?-1:MAX_DEPTH+1):i, device);
  94. }
  95. Profile::Toc();
  96. Profile::Tic("XListSetup",this->Comm(),false,3);
  97. for(size_t i=0;i<MAX_DEPTH;i++){
  98. setup_data[i+MAX_DEPTH*2].precomp_data=&precomp_lst[2];
  99. fmm_mat->X_ListSetup(setup_data[i+MAX_DEPTH*2],node_data_buff,node_lists,fmm_mat->Homogen()?(i==0?-1:MAX_DEPTH+1):i, device);
  100. }
  101. Profile::Toc();
  102. Profile::Tic("VListSetup",this->Comm(),false,3);
  103. for(size_t i=0;i<MAX_DEPTH;i++){
  104. setup_data[i+MAX_DEPTH*3].precomp_data=&precomp_lst[3];
  105. fmm_mat->V_ListSetup(setup_data[i+MAX_DEPTH*3],node_data_buff,node_lists,fmm_mat->Homogen()?(i==0?-1:MAX_DEPTH+1):i, /*device*/ false);
  106. }
  107. Profile::Toc();
  108. Profile::Tic("D2DSetup",this->Comm(),false,3);
  109. for(size_t i=0;i<MAX_DEPTH;i++){
  110. setup_data[i+MAX_DEPTH*4].precomp_data=&precomp_lst[4];
  111. fmm_mat->Down2DownSetup(setup_data[i+MAX_DEPTH*4],node_data_buff,node_lists,i, /*device*/ false);
  112. }
  113. Profile::Toc();
  114. Profile::Tic("D2TSetup",this->Comm(),false,3);
  115. for(size_t i=0;i<MAX_DEPTH;i++){
  116. setup_data[i+MAX_DEPTH*5].precomp_data=&precomp_lst[5];
  117. fmm_mat->Down2TargetSetup(setup_data[i+MAX_DEPTH*5],node_data_buff,node_lists,fmm_mat->Homogen()?(i==0?-1:MAX_DEPTH+1):i, /*device*/ false);
  118. }
  119. Profile::Toc();
  120. Profile::Tic("S2USetup",this->Comm(),false,3);
  121. for(size_t i=0;i<MAX_DEPTH;i++){
  122. setup_data[i+MAX_DEPTH*6].precomp_data=&precomp_lst[6];
  123. fmm_mat->Source2UpSetup(setup_data[i+MAX_DEPTH*6],node_data_buff,node_lists,fmm_mat->Homogen()?(i==0?-1:MAX_DEPTH+1):i, /*device*/ false);
  124. }
  125. Profile::Toc();
  126. Profile::Tic("U2USetup",this->Comm(),false,3);
  127. for(size_t i=0;i<MAX_DEPTH;i++){
  128. setup_data[i+MAX_DEPTH*7].precomp_data=&precomp_lst[7];
  129. fmm_mat->Up2UpSetup(setup_data[i+MAX_DEPTH*7],node_data_buff,node_lists,i, /*device*/ false);
  130. }
  131. Profile::Toc();
  132. #ifdef __INTEL_OFFLOAD
  133. int wait_lock_idx=-1;
  134. wait_lock_idx=MIC_Lock::curr_lock();
  135. #pragma offload target(mic:0)
  136. {MIC_Lock::wait_lock(wait_lock_idx);}
  137. #endif
  138. }Profile::Toc();
  139. }
  140. template <class FMM_Mat_t>
  141. void FMM_Tree<FMM_Mat_t>::ClearFMMData() {
  142. Profile::Tic("ClearFMMData",this->Comm(),true);{
  143. bool device=true;
  144. if(setup_data[0+MAX_DEPTH*1]. input_data!=NULL) setup_data[0+MAX_DEPTH*1]. input_data->SetZero();
  145. if(setup_data[0+MAX_DEPTH*2].output_data!=NULL) setup_data[0+MAX_DEPTH*2].output_data->SetZero();
  146. if(setup_data[0+MAX_DEPTH*0].output_data!=NULL) setup_data[0+MAX_DEPTH*0].output_data->SetZero();
  147. if(device){ // Host2Device
  148. if(setup_data[0+MAX_DEPTH*1]. input_data!=NULL) setup_data[0+MAX_DEPTH*1]. input_data->AllocDevice(true);
  149. if(setup_data[0+MAX_DEPTH*2].output_data!=NULL) setup_data[0+MAX_DEPTH*2].output_data->AllocDevice(true);
  150. if(setup_data[0+MAX_DEPTH*0].output_data!=NULL) setup_data[0+MAX_DEPTH*0].output_data->AllocDevice(true);
  151. #ifdef __INTEL_OFFLOAD
  152. if(!fmm_mat->Homogen()){ // Wait
  153. int wait_lock_idx=-1;
  154. wait_lock_idx=MIC_Lock::curr_lock();
  155. #pragma offload target(mic:0)
  156. {MIC_Lock::wait_lock(wait_lock_idx);}
  157. }
  158. MIC_Lock::init();
  159. #endif
  160. }
  161. }Profile::Toc();
  162. }
  163. template <class FMM_Mat_t>
  164. void FMM_Tree<FMM_Mat_t>::RunFMM() {
  165. Profile::Tic("RunFMM",this->Comm(),true);{
  166. //Upward Pass
  167. Profile::Tic("UpwardPass",this->Comm(),false,2);
  168. UpwardPass();
  169. Profile::Toc();
  170. //Multipole Reduce Broadcast.
  171. Profile::Tic("ReduceBcast",this->Comm(),true,2);
  172. MultipoleReduceBcast();
  173. Profile::Toc();
  174. //Local 2:1 Balancing.
  175. //This can cause load imbalance, always use global 2:1 balance instead.
  176. //Profile::Tic("2:1Balance(local)",this->Comm(),false,3);
  177. //this->Balance21_local(bndry);
  178. //UpwardPass(true);
  179. //Profile::Toc();
  180. //Downward Pass
  181. Profile::Tic("DownwardPass",this->Comm(),true,2);
  182. DownwardPass();
  183. Profile::Toc();
  184. }Profile::Toc();
  185. }
  186. template <class FMM_Mat_t>
  187. void FMM_Tree<FMM_Mat_t>::UpwardPass() {
  188. bool device=true;
  189. //Upward Pass (initialize all leaf nodes)
  190. Profile::Tic("S2U",this->Comm(),false,5);
  191. for(int i=0; i<(fmm_mat->Homogen()?1:MAX_DEPTH); i++){ // Source2Up
  192. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*6],/*device*/ false);
  193. fmm_mat->Source2Up(setup_data[i+MAX_DEPTH*6]);
  194. }
  195. Profile::Toc();
  196. //Upward Pass (level by level)
  197. Profile::Tic("U2U",this->Comm(),false,5);
  198. for(int i=MAX_DEPTH-1; i>=0; i--){ // Up2Up
  199. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*7],/*device*/ false);
  200. fmm_mat->Up2Up(setup_data[i+MAX_DEPTH*7]);
  201. }
  202. Profile::Toc();
  203. }
  204. template <class FMM_Mat_t>
  205. void FMM_Tree<FMM_Mat_t>::BuildInteracLists() {
  206. std::vector<FMM_Node_t*>& n_list=this->GetNodeList();
  207. // Build interaction lists.
  208. int omp_p=omp_get_max_threads();
  209. {
  210. size_t k=n_list.size();
  211. #pragma omp parallel for
  212. for(int j=0;j<omp_p;j++){
  213. size_t a=(k*j)/omp_p;
  214. size_t b=(k*(j+1))/omp_p;
  215. for(size_t i=a;i<b;i++){
  216. FMM_Node_t* n=n_list[i];
  217. n->interac_list.resize(Type_Count);
  218. n->interac_list[S2U_Type]=interac_list.BuildList(n,S2U_Type);
  219. n->interac_list[U2U_Type]=interac_list.BuildList(n,U2U_Type);
  220. n->interac_list[D2D_Type]=interac_list.BuildList(n,D2D_Type);
  221. n->interac_list[D2T_Type]=interac_list.BuildList(n,D2T_Type);
  222. n->interac_list[U0_Type]=interac_list.BuildList(n,U0_Type);
  223. n->interac_list[U1_Type]=interac_list.BuildList(n,U1_Type);
  224. n->interac_list[U2_Type]=interac_list.BuildList(n,U2_Type);
  225. n->interac_list[V_Type]=interac_list.BuildList(n,V_Type);
  226. n->interac_list[V1_Type]=interac_list.BuildList(n,V1_Type);
  227. n->interac_list[W_Type]=interac_list.BuildList(n,W_Type);
  228. n->interac_list[X_Type]=interac_list.BuildList(n,X_Type);
  229. }
  230. }
  231. }
  232. }
  233. template <class FMM_Mat_t>
  234. void FMM_Tree<FMM_Mat_t>::MultipoleReduceBcast() {
  235. int num_p,rank;
  236. MPI_Comm_size(*this->Comm(),&num_p);
  237. MPI_Comm_rank(*this->Comm(),&rank );
  238. if(num_p==1) return;
  239. Profile::Tic("Reduce",this->Comm(),true,3);
  240. std::vector<MortonId> mins=this->GetMins();
  241. size_t bit_mask=1;
  242. size_t max_child=(1UL<<this->Dim());
  243. //Initialize initial send nodes.
  244. std::vector<FMM_Node_t*> send_nodes[2];
  245. //Initialize send_node[0]
  246. FMM_Node_t* tmp_node=static_cast<FMM_Node_t*>(this->RootNode());
  247. assert(!tmp_node->IsGhost());
  248. while(!tmp_node->IsLeaf()){
  249. FMM_Node_t* tmp_node_=NULL;
  250. for(size_t i=0;i<max_child;i++){
  251. tmp_node_=static_cast<FMM_Node_t*>(tmp_node->Child(i));
  252. if(tmp_node_!=NULL) if(!tmp_node_->IsGhost()) break;
  253. }
  254. tmp_node=tmp_node_; assert(tmp_node!=NULL);
  255. }
  256. int n[2];
  257. n[0]=tmp_node->Depth()+1;
  258. send_nodes[0].resize(n[0]);
  259. send_nodes[0][n[0]-1]=tmp_node;
  260. for(int i=n[0]-1;i>0;i--)
  261. send_nodes[0][i-1]=static_cast<FMM_Node_t*>(send_nodes[0][i]->Parent());
  262. //Initialize send_node[1]
  263. tmp_node=static_cast<FMM_Node_t*>(this->RootNode());
  264. while(!tmp_node->IsLeaf()){
  265. FMM_Node_t* tmp_node_=NULL;
  266. for(int i=max_child-1;i>=0;i--){
  267. tmp_node_=static_cast<FMM_Node_t*>(tmp_node->Child(i));
  268. if(tmp_node_!=NULL) if(!tmp_node_->IsGhost()) break;
  269. }
  270. tmp_node=tmp_node_; assert(tmp_node!=NULL);
  271. }
  272. n[1]=tmp_node->Depth()+1;
  273. send_nodes[1].resize(n[1]);
  274. send_nodes[1][n[1]-1]=tmp_node;
  275. for(int i=n[1]-1;i>0;i--)
  276. send_nodes[1][i-1]=static_cast<FMM_Node_t*>(send_nodes[1][i]->Parent());
  277. //Hypercube reduction.
  278. while(bit_mask<(size_t)num_p){
  279. size_t partner=rank^bit_mask; //Partner process id
  280. int merge_indx=((bit_mask & rank)==0?1:0);
  281. bit_mask=bit_mask<<1;
  282. //if(rank >= num_p - (num_p % bit_mask)) break;
  283. //Initialize send data.
  284. size_t s_node_cnt[2]={send_nodes[0].size(),send_nodes[1].size()};
  285. int send_size=2*sizeof(size_t)+(s_node_cnt[0]+s_node_cnt[1])*sizeof(MortonId);
  286. std::vector<PackedData> send_data(s_node_cnt[0]+s_node_cnt[1]);
  287. size_t s_iter=0;
  288. for(int i=0;i<2;i++)
  289. for(size_t j=0;j<s_node_cnt[i];j++){
  290. assert(send_nodes[i][j]!=NULL);
  291. send_data[s_iter]=send_nodes[i][j]->PackMultipole();
  292. send_size+=send_data[s_iter].length+sizeof(size_t);
  293. s_iter++;
  294. }
  295. char* send_buff=new char[send_size];
  296. char* buff_iter=send_buff;
  297. ((size_t*)buff_iter)[0]=s_node_cnt[0];
  298. ((size_t*)buff_iter)[1]=s_node_cnt[1];
  299. buff_iter+=2*sizeof(size_t);
  300. s_iter=0;
  301. for(int i=0;i<2;i++)
  302. for(size_t j=0;j<s_node_cnt[i];j++){
  303. ((MortonId*)buff_iter)[0]=send_nodes[i][j]->GetMortonId();
  304. buff_iter+=sizeof(MortonId);
  305. ((size_t*)buff_iter)[0]=send_data[s_iter].length;
  306. buff_iter+=sizeof(size_t);
  307. mem::memcopy((void*)buff_iter,send_data[s_iter].data,send_data[s_iter].length);
  308. buff_iter+=send_data[s_iter].length;
  309. s_iter++;
  310. }
  311. //Exchange send and recv sizes
  312. int recv_size=0;
  313. MPI_Status status;
  314. char* recv_buff=NULL;
  315. if(partner<(size_t)num_p){
  316. MPI_Sendrecv(&send_size, 1, MPI_INT, partner, 0, &recv_size, 1, MPI_INT, partner, 0, *this->Comm(), &status);
  317. recv_buff=new char[recv_size];
  318. MPI_Sendrecv(send_buff, send_size, MPI_BYTE, partner, 0, recv_buff, recv_size, MPI_BYTE, partner, 0, *this->Comm(), &status);
  319. }
  320. //Need an extra broadcast for incomplete hypercubes.
  321. size_t p0_start=num_p - (num_p % (bit_mask ));
  322. size_t p0_end =num_p - (num_p % (bit_mask>>1));
  323. if(((size_t)rank >= p0_start) && ((size_t)num_p>p0_end) && ((size_t)rank < p0_end) ){
  324. size_t bit_mask0=1;
  325. size_t num_p0=p0_end-p0_start;
  326. while( bit_mask0 < num_p0 ){
  327. if( (bit_mask0<<1) > (num_p - p0_end) ){
  328. size_t partner0=rank^bit_mask0;
  329. if( rank-p0_start < bit_mask0 ){
  330. //Send
  331. MPI_Send(&recv_size, 1, MPI_INT , partner0, 0, *this->Comm());
  332. MPI_Send( recv_buff, recv_size, MPI_BYTE, partner0, 0, *this->Comm());
  333. }else if( rank-p0_start < (bit_mask0<<1) ){
  334. //Receive
  335. if(recv_size>0) delete[] recv_buff;
  336. MPI_Recv(&recv_size, 1, MPI_INT , partner0, 0, *this->Comm(), &status);
  337. recv_buff=new char[recv_size];
  338. MPI_Recv( recv_buff, recv_size, MPI_BYTE, partner0, 0, *this->Comm(), &status);
  339. }
  340. }
  341. bit_mask0=bit_mask0<<1;
  342. }
  343. }
  344. //Construct nodes from received data.
  345. if(recv_size>0){
  346. buff_iter=recv_buff;
  347. size_t r_node_cnt[2]={((size_t*)buff_iter)[0],((size_t*)buff_iter)[1]};
  348. buff_iter+=2*sizeof(size_t);
  349. std::vector<MortonId> r_mid[2];
  350. r_mid[0].resize(r_node_cnt[0]);
  351. r_mid[1].resize(r_node_cnt[1]);
  352. std::vector<FMM_Node_t*> recv_nodes[2];
  353. recv_nodes[0].resize(r_node_cnt[0]);
  354. recv_nodes[1].resize(r_node_cnt[1]);
  355. std::vector<PackedData> recv_data[2];
  356. recv_data[0].resize(r_node_cnt[0]);
  357. recv_data[1].resize(r_node_cnt[1]);
  358. for(int i=0;i<2;i++)
  359. for(size_t j=0;j<r_node_cnt[i];j++){
  360. r_mid[i][j]=((MortonId*)buff_iter)[0];
  361. buff_iter+=sizeof(MortonId);
  362. recv_data[i][j].length=((size_t*)buff_iter)[0];
  363. buff_iter+=sizeof(size_t);
  364. recv_data[i][j].data=(void*)buff_iter;
  365. buff_iter+=recv_data[i][j].length;
  366. }
  367. // Add multipole expansion to existing nodes.
  368. for(size_t i=0;i<r_node_cnt[1-merge_indx];i++){
  369. if(i<send_nodes[merge_indx].size()){
  370. if(r_mid[1-merge_indx][i]==send_nodes[merge_indx][i]->GetMortonId()){
  371. send_nodes[merge_indx][i]->AddMultipole(recv_data[1-merge_indx][i]);
  372. }else break;
  373. }else break;
  374. }
  375. bool new_branch=false;
  376. for(size_t i=0;i<r_node_cnt[merge_indx];i++){
  377. if(i<send_nodes[merge_indx].size() && !new_branch){
  378. if(r_mid[merge_indx][i]==send_nodes[merge_indx][i]->GetMortonId()){
  379. recv_nodes[merge_indx][i]=send_nodes[merge_indx][i];
  380. }else{
  381. new_branch=true;
  382. size_t n_=(i<(size_t)n[merge_indx]?n[merge_indx]:i);
  383. for(size_t j=n_;j<send_nodes[merge_indx].size();j++)
  384. delete send_nodes[merge_indx][j];
  385. if(i<(size_t)n[merge_indx]) n[merge_indx]=i;
  386. }
  387. }
  388. if(i>=send_nodes[merge_indx].size() || new_branch){
  389. recv_nodes[merge_indx][i]=static_cast<FMM_Node_t*>(this->NewNode());
  390. recv_nodes[merge_indx][i]->SetCoord(r_mid[merge_indx][i]);
  391. recv_nodes[merge_indx][i]->InitMultipole(recv_data[merge_indx][i]);
  392. }
  393. }
  394. send_nodes[merge_indx]=recv_nodes[merge_indx];
  395. }
  396. delete[] send_buff;
  397. delete[] recv_buff;
  398. }
  399. for(int i=0;i<2;i++)
  400. for(size_t j=n[i];j<send_nodes[i].size();j++)
  401. delete send_nodes[i][j];
  402. Profile::Toc();
  403. //Now Broadcast nodes to build LET.
  404. Profile::Tic("Broadcast",this->Comm(),true,4);
  405. this->ConstructLET(bndry);
  406. Profile::Toc();
  407. }
  408. template <class FMM_Mat_t>
  409. void FMM_Tree<FMM_Mat_t>::DownwardPass() {
  410. bool device=true;
  411. Profile::Tic("Setup",this->Comm(),true,3);
  412. std::vector<FMM_Node_t*> leaf_nodes;
  413. int max_depth=0;
  414. { // Build leaf node list
  415. int max_depth_loc=0;
  416. std::vector<FMM_Node_t*>& nodes=this->GetNodeList();
  417. for(size_t i=0;i<nodes.size();i++){
  418. FMM_Node_t* n=nodes[i];
  419. if(!n->IsGhost() && n->IsLeaf()) leaf_nodes.push_back(n);
  420. if(n->Depth()>max_depth_loc) max_depth_loc=n->Depth();
  421. }
  422. MPI_Allreduce(&max_depth_loc, &max_depth, 1, MPI_INT, MPI_MAX, *this->Comm());
  423. }
  424. Profile::Toc();
  425. #ifdef __INTEL_OFFLOAD
  426. if(device){ // Host2Device:Src
  427. Profile::Tic("Host2Device:Src",this->Comm(),false,5);
  428. if(setup_data[0+MAX_DEPTH*2]. coord_data!=NULL) setup_data[0+MAX_DEPTH*2]. coord_data->AllocDevice(true);
  429. if(setup_data[0+MAX_DEPTH*2]. input_data!=NULL) setup_data[0+MAX_DEPTH*2]. input_data->AllocDevice(true);
  430. Profile::Toc();
  431. }
  432. #endif
  433. if(bndry==Periodic){ //Add contribution from periodic infinite tiling.
  434. Profile::Tic("BoundaryCondition",this->Comm(),false,5);
  435. fmm_mat->PeriodicBC(dynamic_cast<FMM_Node_t*>(this->RootNode()));
  436. Profile::Toc();
  437. }
  438. for(size_t i=0; i<=(fmm_mat->Homogen()?0:max_depth); i++){ // U,V,W,X-lists
  439. if(!fmm_mat->Homogen()){ // Precomp
  440. std::stringstream level_str;
  441. level_str<<"Level-"<<std::setfill('0')<<std::setw(2)<<i<<"\0";
  442. Profile::Tic(level_str.str().c_str(),this->Comm(),false,5);
  443. Profile::Tic("Precomp",this->Comm(),false,5);
  444. {// Precomp U
  445. Profile::Tic("Precomp-U",this->Comm(),false,10);
  446. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*0],device);
  447. Profile::Toc();
  448. }
  449. {// Precomp W
  450. Profile::Tic("Precomp-W",this->Comm(),false,10);
  451. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*1],device);
  452. Profile::Toc();
  453. }
  454. {// Precomp X
  455. Profile::Tic("Precomp-X",this->Comm(),false,10);
  456. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*2],device);
  457. Profile::Toc();
  458. }
  459. {// Precomp V
  460. Profile::Tic("Precomp-V",this->Comm(),false,10);
  461. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*3], /*device*/ false);
  462. Profile::Toc();
  463. }
  464. Profile::Toc();
  465. }
  466. {// X-List
  467. Profile::Tic("X-List",this->Comm(),false,5);
  468. fmm_mat->X_List(setup_data[i+MAX_DEPTH*2], device);
  469. Profile::Toc();
  470. }
  471. #ifdef __INTEL_OFFLOAD
  472. if(i==0 && device){ // Host2Device:Mult
  473. Profile::Tic("Host2Device:Mult",this->Comm(),false,5);
  474. if(setup_data[0+MAX_DEPTH*1]. input_data!=NULL) setup_data[0+MAX_DEPTH*1]. input_data->AllocDevice(true);
  475. Profile::Toc();
  476. }
  477. if(device) if(i==(fmm_mat->Homogen()?0:max_depth)){ // Device2Host: LocalExp
  478. Profile::Tic("Device2Host:LocExp",this->Comm(),false,5);
  479. if(setup_data[0+MAX_DEPTH*2].output_data!=NULL){
  480. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*2].output_data;
  481. assert(fmm_mat->dev_buffer.Dim()>=output_data.Dim(0)*output_data.Dim(1));
  482. output_data.Device2Host((Real_t*)&fmm_mat->dev_buffer[0]);
  483. }
  484. Profile::Toc();
  485. }
  486. #endif
  487. {// W-List
  488. Profile::Tic("W-List",this->Comm(),false,5);
  489. fmm_mat->W_List(setup_data[i+MAX_DEPTH*1], device);
  490. Profile::Toc();
  491. }
  492. {// U-List
  493. Profile::Tic("U-List",this->Comm(),false,5);
  494. fmm_mat->U_List(setup_data[i+MAX_DEPTH*0], device);
  495. Profile::Toc();
  496. }
  497. {// V-List
  498. Profile::Tic("V-List",this->Comm(),false,5);
  499. fmm_mat->V_List(setup_data[i+MAX_DEPTH*3], /*device*/ false);
  500. Profile::Toc();
  501. }
  502. if(!fmm_mat->Homogen()){ // Wait
  503. #ifdef __INTEL_OFFLOAD
  504. int wait_lock_idx=-1;
  505. if(device) wait_lock_idx=MIC_Lock::curr_lock();
  506. #pragma offload if(device) target(mic:0)
  507. {if(device) MIC_Lock::wait_lock(wait_lock_idx);}
  508. #endif
  509. Profile::Toc();
  510. }
  511. }
  512. #ifdef __INTEL_OFFLOAD
  513. Profile::Tic("D2H_Wait:LocExp",this->Comm(),false,5);
  514. if(device) if(setup_data[0+MAX_DEPTH*2].output_data!=NULL){
  515. Real_t* dev_ptr=(Real_t*)&fmm_mat->dev_buffer[0];
  516. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*2].output_data;
  517. size_t n=output_data.Dim(0)*output_data.Dim(1);
  518. Real_t* host_ptr=output_data[0];
  519. output_data.Device2HostWait();
  520. #pragma omp parallel for
  521. for(size_t i=0;i<n;i++){
  522. host_ptr[i]+=dev_ptr[i];
  523. }
  524. }
  525. Profile::Toc();
  526. Profile::Tic("Device2Host:Trg",this->Comm(),false,5);
  527. if(device) if(setup_data[0+MAX_DEPTH*0].output_data!=NULL){ // Device2Host: Target
  528. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*0].output_data;
  529. assert(fmm_mat->dev_buffer.Dim()>=output_data.Dim(0)*output_data.Dim(1));
  530. output_data.Device2Host((Real_t*)&fmm_mat->dev_buffer[0]);
  531. }
  532. Profile::Toc();
  533. #endif
  534. Profile::Tic("D2D",this->Comm(),false,5);
  535. for(size_t i=0; i<=max_depth; i++){ // Down2Down
  536. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*4],/*device*/ false);
  537. fmm_mat->Down2Down(setup_data[i+MAX_DEPTH*4]);
  538. }
  539. Profile::Toc();
  540. Profile::Tic("D2T",this->Comm(),false,5);
  541. for(int i=0; i<=(fmm_mat->Homogen()?0:max_depth); i++){ // Down2Target
  542. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*5],/*device*/ false);
  543. fmm_mat->Down2Target(setup_data[i+MAX_DEPTH*5]);
  544. }
  545. Profile::Toc();
  546. #ifdef __INTEL_OFFLOAD
  547. Profile::Tic("D2H_Wait:Trg",this->Comm(),false,5);
  548. if(device) if(setup_data[0+MAX_DEPTH*0].output_data!=NULL){
  549. Real_t* dev_ptr=(Real_t*)&fmm_mat->dev_buffer[0];
  550. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*0].output_data;
  551. size_t n=output_data.Dim(0)*output_data.Dim(1);
  552. Real_t* host_ptr=output_data[0];
  553. output_data.Device2HostWait();
  554. #pragma omp parallel for
  555. for(size_t i=0;i<n;i++){
  556. host_ptr[i]+=dev_ptr[i];
  557. }
  558. }
  559. Profile::Toc();
  560. #endif
  561. Profile::Tic("PostProc",this->Comm(),false,5);
  562. fmm_mat->PostProcessing(leaf_nodes);
  563. Profile::Toc();
  564. }
  565. template <class FMM_Mat_t>
  566. void FMM_Tree<FMM_Mat_t>::Copy_FMMOutput() {
  567. std::vector<FMM_Node_t*>& all_nodes=this->GetNodeList();
  568. int omp_p=omp_get_max_threads();
  569. // Copy output to the tree.
  570. {
  571. size_t k=all_nodes.size();
  572. #pragma omp parallel for
  573. for(int j=0;j<omp_p;j++){
  574. size_t a=(k*j)/omp_p;
  575. size_t b=(k*(j+1))/omp_p;
  576. fmm_mat->CopyOutput(&(all_nodes[a]),b-a);
  577. }
  578. }
  579. }
  580. }//end namespace