fmm_tree.txx 23 KB

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