fmm_tree.txx 23 KB

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