fmm_tree.txx 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687
  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,2);
  43. this->RefineTree();
  44. Profile::Toc();
  45. }
  46. //2:1 Balancing
  47. Profile::Tic("2:1Balance",this->Comm(),true,2);
  48. this->Balance21(bndry);
  49. Profile::Toc();
  50. //Redistribute nodes.
  51. Profile::Tic("Redistribute",this->Comm(),true,3);
  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. //Upward Pass (initialize all leaf nodes)
  200. Profile::Tic("S2U",this->Comm(),false,5);
  201. for(int i=0; i<(fmm_mat->Homogen()?1:MAX_DEPTH); i++){ // Source2Up
  202. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*6],/*device*/ false);
  203. fmm_mat->Source2Up(setup_data[i+MAX_DEPTH*6]);
  204. }
  205. Profile::Toc();
  206. //Upward Pass (level by level)
  207. Profile::Tic("U2U",this->Comm(),false,5);
  208. for(int i=MAX_DEPTH-1; i>=0; i--){ // Up2Up
  209. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*7],/*device*/ false);
  210. fmm_mat->Up2Up(setup_data[i+MAX_DEPTH*7]);
  211. }
  212. Profile::Toc();
  213. }
  214. template <class FMM_Mat_t>
  215. void FMM_Tree<FMM_Mat_t>::BuildInteracLists() {
  216. std::vector<Node_t*>& n_list=this->GetNodeList();
  217. // Build interaction lists.
  218. int omp_p=omp_get_max_threads();
  219. {
  220. size_t k=n_list.size();
  221. #pragma omp parallel for
  222. for(int j=0;j<omp_p;j++){
  223. size_t a=(k*j)/omp_p;
  224. size_t b=(k*(j+1))/omp_p;
  225. for(size_t i=a;i<b;i++){
  226. Node_t* n=n_list[i];
  227. n->interac_list.resize(Type_Count);
  228. n->interac_list[S2U_Type]=interac_list.BuildList(n,S2U_Type);
  229. n->interac_list[U2U_Type]=interac_list.BuildList(n,U2U_Type);
  230. n->interac_list[D2D_Type]=interac_list.BuildList(n,D2D_Type);
  231. n->interac_list[D2T_Type]=interac_list.BuildList(n,D2T_Type);
  232. n->interac_list[U0_Type]=interac_list.BuildList(n,U0_Type);
  233. n->interac_list[U1_Type]=interac_list.BuildList(n,U1_Type);
  234. n->interac_list[U2_Type]=interac_list.BuildList(n,U2_Type);
  235. n->interac_list[V_Type]=interac_list.BuildList(n,V_Type);
  236. n->interac_list[V1_Type]=interac_list.BuildList(n,V1_Type);
  237. n->interac_list[W_Type]=interac_list.BuildList(n,W_Type);
  238. n->interac_list[X_Type]=interac_list.BuildList(n,X_Type);
  239. }
  240. }
  241. }
  242. }
  243. template <class FMM_Mat_t>
  244. void FMM_Tree<FMM_Mat_t>::MultipoleReduceBcast() {
  245. int num_p,rank;
  246. MPI_Comm_size(*this->Comm(),&num_p);
  247. MPI_Comm_rank(*this->Comm(),&rank );
  248. if(num_p==1) return;
  249. Profile::Tic("Reduce",this->Comm(),true,3);
  250. std::vector<MortonId> mins=this->GetMins();
  251. size_t bit_mask=1;
  252. size_t max_child=(1UL<<this->Dim());
  253. //Initialize initial send nodes.
  254. std::vector<Node_t*> send_nodes[2];
  255. //Initialize send_node[0]
  256. Node_t* tmp_node=static_cast<Node_t*>(this->RootNode());
  257. assert(!tmp_node->IsGhost());
  258. while(!tmp_node->IsLeaf()){
  259. Node_t* tmp_node_=NULL;
  260. for(size_t i=0;i<max_child;i++){
  261. tmp_node_=static_cast<Node_t*>(tmp_node->Child(i));
  262. if(tmp_node_!=NULL) if(!tmp_node_->IsGhost()) break;
  263. }
  264. tmp_node=tmp_node_; assert(tmp_node!=NULL);
  265. }
  266. int n[2];
  267. n[0]=tmp_node->Depth()+1;
  268. send_nodes[0].resize(n[0]);
  269. send_nodes[0][n[0]-1]=tmp_node;
  270. for(int i=n[0]-1;i>0;i--)
  271. send_nodes[0][i-1]=static_cast<Node_t*>(send_nodes[0][i]->Parent());
  272. //Initialize send_node[1]
  273. tmp_node=static_cast<Node_t*>(this->RootNode());
  274. while(!tmp_node->IsLeaf()){
  275. Node_t* tmp_node_=NULL;
  276. for(int i=max_child-1;i>=0;i--){
  277. tmp_node_=static_cast<Node_t*>(tmp_node->Child(i));
  278. if(tmp_node_!=NULL) if(!tmp_node_->IsGhost()) break;
  279. }
  280. tmp_node=tmp_node_; assert(tmp_node!=NULL);
  281. }
  282. n[1]=tmp_node->Depth()+1;
  283. send_nodes[1].resize(n[1]);
  284. send_nodes[1][n[1]-1]=tmp_node;
  285. for(int i=n[1]-1;i>0;i--)
  286. send_nodes[1][i-1]=static_cast<Node_t*>(send_nodes[1][i]->Parent());
  287. //Hypercube reduction.
  288. while(bit_mask<(size_t)num_p){
  289. size_t partner=rank^bit_mask; //Partner process id
  290. int merge_indx=((bit_mask & rank)==0?1:0);
  291. bit_mask=bit_mask<<1;
  292. //if(rank >= num_p - (num_p % bit_mask)) break;
  293. //Initialize send data.
  294. size_t s_node_cnt[2]={send_nodes[0].size(),send_nodes[1].size()};
  295. int send_size=2*sizeof(size_t)+(s_node_cnt[0]+s_node_cnt[1])*sizeof(MortonId);
  296. std::vector<PackedData> send_data(s_node_cnt[0]+s_node_cnt[1]);
  297. size_t s_iter=0;
  298. for(int i=0;i<2;i++)
  299. for(size_t j=0;j<s_node_cnt[i];j++){
  300. assert(send_nodes[i][j]!=NULL);
  301. send_data[s_iter]=send_nodes[i][j]->PackMultipole();
  302. send_size+=send_data[s_iter].length+sizeof(size_t);
  303. s_iter++;
  304. }
  305. char* send_buff=mem::aligned_new<char>(send_size);
  306. char* buff_iter=send_buff;
  307. ((size_t*)buff_iter)[0]=s_node_cnt[0];
  308. ((size_t*)buff_iter)[1]=s_node_cnt[1];
  309. buff_iter+=2*sizeof(size_t);
  310. s_iter=0;
  311. for(int i=0;i<2;i++)
  312. for(size_t j=0;j<s_node_cnt[i];j++){
  313. ((MortonId*)buff_iter)[0]=send_nodes[i][j]->GetMortonId();
  314. buff_iter+=sizeof(MortonId);
  315. ((size_t*)buff_iter)[0]=send_data[s_iter].length;
  316. buff_iter+=sizeof(size_t);
  317. mem::memcopy((void*)buff_iter,send_data[s_iter].data,send_data[s_iter].length);
  318. buff_iter+=send_data[s_iter].length;
  319. s_iter++;
  320. }
  321. //Exchange send and recv sizes
  322. int recv_size=0;
  323. MPI_Status status;
  324. char* recv_buff=NULL;
  325. if(partner<(size_t)num_p){
  326. MPI_Sendrecv(&send_size, 1, MPI_INT, partner, 0, &recv_size, 1, MPI_INT, partner, 0, *this->Comm(), &status);
  327. recv_buff=mem::aligned_new<char>(recv_size);
  328. MPI_Sendrecv(send_buff, send_size, MPI_BYTE, partner, 0, recv_buff, recv_size, MPI_BYTE, partner, 0, *this->Comm(), &status);
  329. }
  330. //Need an extra broadcast for incomplete hypercubes.
  331. size_t p0_start=num_p - (num_p % (bit_mask ));
  332. size_t p0_end =num_p - (num_p % (bit_mask>>1));
  333. if(((size_t)rank >= p0_start) && ((size_t)num_p>p0_end) && ((size_t)rank < p0_end) ){
  334. size_t bit_mask0=1;
  335. size_t num_p0=p0_end-p0_start;
  336. while( bit_mask0 < num_p0 ){
  337. if( (bit_mask0<<1) > (num_p - p0_end) ){
  338. size_t partner0=rank^bit_mask0;
  339. if( rank-p0_start < bit_mask0 ){
  340. //Send
  341. MPI_Send(&recv_size, 1, MPI_INT , partner0, 0, *this->Comm());
  342. MPI_Send( recv_buff, recv_size, MPI_BYTE, partner0, 0, *this->Comm());
  343. }else if( rank-p0_start < (bit_mask0<<1) ){
  344. //Receive
  345. if(recv_size>0) mem::aligned_delete<char>(recv_buff);
  346. MPI_Recv(&recv_size, 1, MPI_INT , partner0, 0, *this->Comm(), &status);
  347. recv_buff=mem::aligned_new<char>(recv_size);
  348. MPI_Recv( recv_buff, recv_size, MPI_BYTE, partner0, 0, *this->Comm(), &status);
  349. }
  350. }
  351. bit_mask0=bit_mask0<<1;
  352. }
  353. }
  354. //Construct nodes from received data.
  355. if(recv_size>0){
  356. buff_iter=recv_buff;
  357. size_t r_node_cnt[2]={((size_t*)buff_iter)[0],((size_t*)buff_iter)[1]};
  358. buff_iter+=2*sizeof(size_t);
  359. std::vector<MortonId> r_mid[2];
  360. r_mid[0].resize(r_node_cnt[0]);
  361. r_mid[1].resize(r_node_cnt[1]);
  362. std::vector<Node_t*> recv_nodes[2];
  363. recv_nodes[0].resize(r_node_cnt[0]);
  364. recv_nodes[1].resize(r_node_cnt[1]);
  365. std::vector<PackedData> recv_data[2];
  366. recv_data[0].resize(r_node_cnt[0]);
  367. recv_data[1].resize(r_node_cnt[1]);
  368. for(int i=0;i<2;i++)
  369. for(size_t j=0;j<r_node_cnt[i];j++){
  370. r_mid[i][j]=((MortonId*)buff_iter)[0];
  371. buff_iter+=sizeof(MortonId);
  372. recv_data[i][j].length=((size_t*)buff_iter)[0];
  373. buff_iter+=sizeof(size_t);
  374. recv_data[i][j].data=(void*)buff_iter;
  375. buff_iter+=recv_data[i][j].length;
  376. }
  377. // Add multipole expansion to existing nodes.
  378. for(size_t i=0;i<r_node_cnt[1-merge_indx];i++){
  379. if(i<send_nodes[merge_indx].size()){
  380. if(r_mid[1-merge_indx][i]==send_nodes[merge_indx][i]->GetMortonId()){
  381. send_nodes[merge_indx][i]->AddMultipole(recv_data[1-merge_indx][i]);
  382. }else break;
  383. }else break;
  384. }
  385. bool new_branch=false;
  386. for(size_t i=0;i<r_node_cnt[merge_indx];i++){
  387. if(i<send_nodes[merge_indx].size() && !new_branch){
  388. if(r_mid[merge_indx][i]==send_nodes[merge_indx][i]->GetMortonId()){
  389. recv_nodes[merge_indx][i]=send_nodes[merge_indx][i];
  390. }else{
  391. new_branch=true;
  392. size_t n_=(i<(size_t)n[merge_indx]?n[merge_indx]:i);
  393. for(size_t j=n_;j<send_nodes[merge_indx].size();j++)
  394. delete send_nodes[merge_indx][j];
  395. if(i<(size_t)n[merge_indx]) n[merge_indx]=i;
  396. }
  397. }
  398. if(i>=send_nodes[merge_indx].size() || new_branch){
  399. recv_nodes[merge_indx][i]=static_cast<Node_t*>(this->NewNode());
  400. recv_nodes[merge_indx][i]->SetCoord(r_mid[merge_indx][i]);
  401. recv_nodes[merge_indx][i]->InitMultipole(recv_data[merge_indx][i]);
  402. }
  403. }
  404. send_nodes[merge_indx]=recv_nodes[merge_indx];
  405. }
  406. mem::aligned_delete<char>(send_buff);
  407. mem::aligned_delete<char>(recv_buff);
  408. }
  409. for(int i=0;i<2;i++)
  410. for(size_t j=n[i];j<send_nodes[i].size();j++)
  411. delete send_nodes[i][j];
  412. Profile::Toc();
  413. //Now Broadcast nodes to build LET.
  414. Profile::Tic("Broadcast",this->Comm(),true,4);
  415. this->ConstructLET(bndry);
  416. Profile::Toc();
  417. }
  418. template <class FMM_Mat_t>
  419. void FMM_Tree<FMM_Mat_t>::DownwardPass() {
  420. bool device=true;
  421. Profile::Tic("Setup",this->Comm(),true,3);
  422. std::vector<Node_t*> leaf_nodes;
  423. int max_depth=0;
  424. { // Build leaf node list
  425. int max_depth_loc=0;
  426. std::vector<Node_t*>& nodes=this->GetNodeList();
  427. for(size_t i=0;i<nodes.size();i++){
  428. Node_t* n=nodes[i];
  429. if(!n->IsGhost() && n->IsLeaf()) leaf_nodes.push_back(n);
  430. if(n->Depth()>max_depth_loc) max_depth_loc=n->Depth();
  431. }
  432. MPI_Allreduce(&max_depth_loc, &max_depth, 1, MPI_INT, MPI_MAX, *this->Comm());
  433. }
  434. Profile::Toc();
  435. #ifdef __INTEL_OFFLOAD
  436. if(device){ // Host2Device:Src
  437. Profile::Tic("Host2Device:Src",this->Comm(),false,5);
  438. if(setup_data[0+MAX_DEPTH*2]. coord_data!=NULL) setup_data[0+MAX_DEPTH*2]. coord_data->AllocDevice(true);
  439. if(setup_data[0+MAX_DEPTH*2]. input_data!=NULL) setup_data[0+MAX_DEPTH*2]. input_data->AllocDevice(true);
  440. Profile::Toc();
  441. }
  442. #endif
  443. if(bndry==Periodic){ //Add contribution from periodic infinite tiling.
  444. Profile::Tic("BoundaryCondition",this->Comm(),false,5);
  445. fmm_mat->PeriodicBC(dynamic_cast<Node_t*>(this->RootNode()));
  446. Profile::Toc();
  447. }
  448. for(size_t i=0; i<=(fmm_mat->Homogen()?0:max_depth); i++){ // U,V,W,X-lists
  449. if(!fmm_mat->Homogen()){ // Precomp
  450. std::stringstream level_str;
  451. level_str<<"Level-"<<std::setfill('0')<<std::setw(2)<<i<<"\0";
  452. Profile::Tic(level_str.str().c_str(),this->Comm(),false,5);
  453. Profile::Tic("Precomp",this->Comm(),false,5);
  454. {// Precomp U
  455. Profile::Tic("Precomp-U",this->Comm(),false,10);
  456. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*0],device);
  457. Profile::Toc();
  458. }
  459. {// Precomp W
  460. Profile::Tic("Precomp-W",this->Comm(),false,10);
  461. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*1],device);
  462. Profile::Toc();
  463. }
  464. {// Precomp X
  465. Profile::Tic("Precomp-X",this->Comm(),false,10);
  466. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*2],device);
  467. Profile::Toc();
  468. }
  469. {// Precomp V
  470. Profile::Tic("Precomp-V",this->Comm(),false,10);
  471. fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*3], /*device*/ false);
  472. Profile::Toc();
  473. }
  474. Profile::Toc();
  475. }
  476. {// X-List
  477. Profile::Tic("X-List",this->Comm(),false,5);
  478. fmm_mat->X_List(setup_data[i+MAX_DEPTH*2], device);
  479. Profile::Toc();
  480. }
  481. #ifdef __INTEL_OFFLOAD
  482. if(i==0 && device){ // Host2Device:Mult
  483. Profile::Tic("Host2Device:Mult",this->Comm(),false,5);
  484. if(setup_data[0+MAX_DEPTH*1]. input_data!=NULL) setup_data[0+MAX_DEPTH*1]. input_data->AllocDevice(true);
  485. Profile::Toc();
  486. }
  487. if(device) if(i==(fmm_mat->Homogen()?0:max_depth)){ // Device2Host: LocalExp
  488. Profile::Tic("Device2Host:LocExp",this->Comm(),false,5);
  489. if(setup_data[0+MAX_DEPTH*2].output_data!=NULL){
  490. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*2].output_data;
  491. assert(fmm_mat->dev_buffer.Dim()>=output_data.Dim(0)*output_data.Dim(1));
  492. output_data.Device2Host((Real_t*)&fmm_mat->dev_buffer[0]);
  493. }
  494. Profile::Toc();
  495. }
  496. #endif
  497. {// W-List
  498. Profile::Tic("W-List",this->Comm(),false,5);
  499. fmm_mat->W_List(setup_data[i+MAX_DEPTH*1], device);
  500. Profile::Toc();
  501. }
  502. {// U-List
  503. Profile::Tic("U-List",this->Comm(),false,5);
  504. fmm_mat->U_List(setup_data[i+MAX_DEPTH*0], device);
  505. Profile::Toc();
  506. }
  507. {// V-List
  508. Profile::Tic("V-List",this->Comm(),false,5);
  509. fmm_mat->V_List(setup_data[i+MAX_DEPTH*3], /*device*/ false);
  510. Profile::Toc();
  511. }
  512. if(!fmm_mat->Homogen()){ // Wait
  513. #ifdef __INTEL_OFFLOAD
  514. int wait_lock_idx=-1;
  515. if(device) wait_lock_idx=MIC_Lock::curr_lock();
  516. #pragma offload if(device) target(mic:0)
  517. {if(device) MIC_Lock::wait_lock(wait_lock_idx);}
  518. #endif
  519. Profile::Toc();
  520. }
  521. }
  522. #ifdef __INTEL_OFFLOAD
  523. Profile::Tic("D2H_Wait:LocExp",this->Comm(),false,5);
  524. if(device) if(setup_data[0+MAX_DEPTH*2].output_data!=NULL){
  525. Real_t* dev_ptr=(Real_t*)&fmm_mat->dev_buffer[0];
  526. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*2].output_data;
  527. size_t n=output_data.Dim(0)*output_data.Dim(1);
  528. Real_t* host_ptr=output_data[0];
  529. output_data.Device2HostWait();
  530. #pragma omp parallel for
  531. for(size_t i=0;i<n;i++){
  532. host_ptr[i]+=dev_ptr[i];
  533. }
  534. }
  535. Profile::Toc();
  536. Profile::Tic("Device2Host:Trg",this->Comm(),false,5);
  537. if(device) if(setup_data[0+MAX_DEPTH*0].output_data!=NULL){ // Device2Host: Target
  538. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*0].output_data;
  539. assert(fmm_mat->dev_buffer.Dim()>=output_data.Dim(0)*output_data.Dim(1));
  540. output_data.Device2Host((Real_t*)&fmm_mat->dev_buffer[0]);
  541. }
  542. Profile::Toc();
  543. #endif
  544. Profile::Tic("D2D",this->Comm(),false,5);
  545. for(size_t i=0; i<=max_depth; i++){ // Down2Down
  546. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*4],/*device*/ false);
  547. fmm_mat->Down2Down(setup_data[i+MAX_DEPTH*4]);
  548. }
  549. Profile::Toc();
  550. Profile::Tic("D2T",this->Comm(),false,5);
  551. for(int i=0; i<=(fmm_mat->Homogen()?0:max_depth); i++){ // Down2Target
  552. if(!fmm_mat->Homogen()) fmm_mat->SetupPrecomp(setup_data[i+MAX_DEPTH*5],/*device*/ false);
  553. fmm_mat->Down2Target(setup_data[i+MAX_DEPTH*5]);
  554. }
  555. Profile::Toc();
  556. #ifdef __INTEL_OFFLOAD
  557. Profile::Tic("D2H_Wait:Trg",this->Comm(),false,5);
  558. if(device) if(setup_data[0+MAX_DEPTH*0].output_data!=NULL){
  559. Real_t* dev_ptr=(Real_t*)&fmm_mat->dev_buffer[0];
  560. Matrix<Real_t>& output_data=*setup_data[0+MAX_DEPTH*0].output_data;
  561. size_t n=output_data.Dim(0)*output_data.Dim(1);
  562. Real_t* host_ptr=output_data[0];
  563. output_data.Device2HostWait();
  564. #pragma omp parallel for
  565. for(size_t i=0;i<n;i++){
  566. host_ptr[i]+=dev_ptr[i];
  567. }
  568. }
  569. Profile::Toc();
  570. #endif
  571. Profile::Tic("PostProc",this->Comm(),false,5);
  572. fmm_mat->PostProcessing(leaf_nodes);
  573. Profile::Toc();
  574. }
  575. template <class FMM_Mat_t>
  576. void FMM_Tree<FMM_Mat_t>::Copy_FMMOutput() {
  577. std::vector<Node_t*>& all_nodes=this->GetNodeList();
  578. int omp_p=omp_get_max_threads();
  579. // Copy output to the tree.
  580. {
  581. size_t k=all_nodes.size();
  582. #pragma omp parallel for
  583. for(int j=0;j<omp_p;j++){
  584. size_t a=(k*j)/omp_p;
  585. size_t b=(k*(j+1))/omp_p;
  586. fmm_mat->CopyOutput(&(all_nodes[a]),b-a);
  587. }
  588. }
  589. }
  590. }//end namespace