fmm_tree.txx 22 KB

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