fmm_tree.txx 23 KB

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