fmm_cheb.txx 41 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084
  1. /**
  2. * \file fmm_cheb.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 3-07-2011
  5. * \brief This file contains the implementation of the FMM_Cheb class.
  6. */
  7. #include <omp.h>
  8. #include <sstream>
  9. #include <iostream>
  10. #include <cstdlib>
  11. #include <cmath>
  12. #ifdef PVFMM_HAVE_SYS_STAT_H
  13. #include <sys/stat.h>
  14. #endif
  15. #include <dtypes.h>
  16. #include <parUtils.h>
  17. #include <cheb_utils.hpp>
  18. #include <mem_mgr.hpp>
  19. #include <profile.hpp>
  20. namespace pvfmm{
  21. template <class FMMNode>
  22. FMM_Cheb<FMMNode>::~FMM_Cheb() {
  23. if(this->mat!=NULL){
  24. int rank;
  25. MPI_Comm_rank(this->comm,&rank);
  26. if(!rank){
  27. FILE* f=fopen(this->mat_fname.c_str(),"r");
  28. if(f==NULL) { //File does not exists.
  29. { // Delete easy to compute matrices.
  30. Mat_Type type=W_Type;
  31. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  32. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  33. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  34. M.Resize(0,0);
  35. }
  36. type=V_Type;
  37. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  38. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  39. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  40. M.Resize(0,0);
  41. }
  42. type=V1_Type;
  43. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  44. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  45. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  46. M.Resize(0,0);
  47. }
  48. type=U2U_Type;
  49. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  50. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  51. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  52. M.Resize(0,0);
  53. }
  54. type=D2D_Type;
  55. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  56. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  57. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  58. M.Resize(0,0);
  59. }
  60. type=D2T_Type;
  61. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  62. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  63. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  64. M.Resize(0,0);
  65. }
  66. }
  67. this->mat->Save2File(this->mat_fname.c_str());
  68. }else fclose(f);
  69. }
  70. }
  71. }
  72. template <class FMMNode>
  73. void FMM_Cheb<FMMNode>::Initialize(int mult_order, int cheb_deg_, const MPI_Comm& comm_, const Kernel<Real_t>* kernel_){
  74. Profile::Tic("InitFMM_Cheb",&comm_,true);{
  75. int rank;
  76. MPI_Comm_rank(comm_,&rank);
  77. int dim=3; //Only supporting 3D
  78. cheb_deg=cheb_deg_;
  79. if(this->mat_fname.size()==0){
  80. std::stringstream st;
  81. st<<PVFMM_PRECOMP_DATA_PATH;
  82. if(!st.str().size()){ // look in PVFMM_DIR
  83. char* pvfmm_dir = getenv ("PVFMM_DIR");
  84. if(pvfmm_dir) st<<pvfmm_dir;
  85. }
  86. #ifndef STAT_MACROS_BROKEN
  87. if(st.str().size()){ // check if the path is a directory
  88. struct stat stat_buff;
  89. if(stat(st.str().c_str(), &stat_buff) || !S_ISDIR(stat_buff.st_mode)){
  90. std::cout<<"error: path not found: "<<st.str()<<'\n';
  91. exit(0);
  92. }
  93. }
  94. #endif
  95. if(st.str().size()) st<<'/';
  96. st<<"Precomp_"<<kernel_->ker_name.c_str()<<"_q"<<cheb_deg<<"_m"<<mult_order;
  97. if(sizeof(Real_t)==8) st<<"";
  98. else if(sizeof(Real_t)==4) st<<"_f";
  99. else st<<"_t"<<sizeof(Real_t);
  100. st<<".data";
  101. this->mat_fname=st.str();
  102. }
  103. if(!rank){
  104. FILE* f=fopen(this->mat_fname.c_str(),"r");
  105. if(f==NULL) { //File does not exists.
  106. std::cout<<"Could not find precomputed data file for "<<kernel_->ker_name.c_str()<<" kernel with q="<<cheb_deg<<" and m="<<mult_order<<".\n";
  107. std::cout<<"This data will be computed and stored for future use at:\n"<<this->mat_fname<<'\n';
  108. std::cout<<"This may take a while...\n";
  109. }else fclose(f);
  110. }
  111. //this->mat->LoadFile(this->mat_fname.c_str(), this->comm);
  112. FMM_Pts<FMMNode>::Initialize(mult_order, comm_, kernel_);
  113. this->mat->RelativeTrgCoord()=cheb_nodes<Real_t>(ChebDeg(),dim);
  114. Profile::Tic("PrecompD2T",&this->comm,false,4);
  115. this->PrecompAll(D2T_Type);
  116. Profile::Toc();
  117. //Volume solver.
  118. Profile::Tic("PrecompS2M",&this->comm,false,4);
  119. this->PrecompAll(S2U_Type);
  120. Profile::Toc();
  121. Profile::Tic("PrecompX",&this->comm,false,4);
  122. this->PrecompAll(X_Type);
  123. Profile::Toc();
  124. Profile::Tic("PrecompW",&this->comm,false,4);
  125. this->PrecompAll(W_Type);
  126. Profile::Toc();
  127. Profile::Tic("PrecompU0",&this->comm,false,4);
  128. this->PrecompAll(U0_Type);
  129. Profile::Toc();
  130. Profile::Tic("PrecompU1",&this->comm,false,4);
  131. this->PrecompAll(U1_Type);
  132. Profile::Toc();
  133. Profile::Tic("PrecompU2",&this->comm,false,4);
  134. this->PrecompAll(U2_Type);
  135. Profile::Toc();
  136. Profile::Tic("Save2File",&this->comm,false,4);
  137. if(!rank){
  138. FILE* f=fopen(this->mat_fname.c_str(),"r");
  139. if(f==NULL) { //File does not exists.
  140. { // Delete easy to compute matrices.
  141. Mat_Type type=W_Type;
  142. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  143. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  144. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  145. M.Resize(0,0);
  146. }
  147. type=V_Type;
  148. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  149. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  150. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  151. M.Resize(0,0);
  152. }
  153. type=V1_Type;
  154. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  155. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  156. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  157. M.Resize(0,0);
  158. }
  159. type=U2U_Type;
  160. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  161. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  162. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  163. M.Resize(0,0);
  164. }
  165. type=D2D_Type;
  166. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  167. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  168. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  169. M.Resize(0,0);
  170. }
  171. type=D2T_Type;
  172. for(int l=-BC_LEVELS;l<MAX_DEPTH;l++)
  173. for(size_t indx=0;indx<this->interac_list.ListCount(type);indx++){
  174. Matrix<Real_t>& M=this->mat->Mat(l, type, indx);
  175. M.Resize(0,0);
  176. }
  177. }
  178. this->mat->Save2File(this->mat_fname.c_str());
  179. }else fclose(f);
  180. }
  181. Profile::Toc();
  182. Profile::Tic("Recompute",&this->comm,false,4);
  183. { // Recompute matrices.
  184. this->PrecompAll(W_Type);
  185. this->PrecompAll(V_Type);
  186. this->PrecompAll(V1_Type);
  187. this->PrecompAll(U2U_Type);
  188. this->PrecompAll(D2D_Type);
  189. this->PrecompAll(D2T_Type);
  190. }
  191. Profile::Toc();
  192. }Profile::Toc();
  193. }
  194. template <class Real_t>
  195. Permutation<Real_t> cheb_perm(size_t q, size_t p_indx, const Permutation<Real_t>& ker_perm, const Vector<Real_t>* scal_exp=NULL){
  196. int dim=3; //Only supporting 3D
  197. int dof=ker_perm.Dim();
  198. int coeff_cnt=((q+1)*(q+2)*(q+3))/6;
  199. int n3=(int)pow((Real_t)(q+1),dim);
  200. Permutation<Real_t> P0(n3*dof);
  201. if(scal_exp && p_indx==Scaling){ // Set level-by-level scaling
  202. assert(dof==scal_exp->Dim());
  203. Vector<Real_t> scal(scal_exp->Dim());
  204. for(size_t i=0;i<scal.Dim();i++){
  205. scal[i]=pow(2.0,(*scal_exp)[i]);
  206. }
  207. for(int j=0;j<dof;j++){
  208. for(int i=0;i<n3;i++){
  209. P0.scal[j*n3+i]*=scal[j];
  210. }
  211. }
  212. }
  213. { // Set P0.scal
  214. for(int j=0;j<dof;j++){
  215. for(int i=0;i<n3;i++){
  216. P0.scal[j*n3+i]*=ker_perm.scal[j];
  217. }
  218. }
  219. }
  220. if(p_indx==ReflecX || p_indx==ReflecY || p_indx==ReflecZ){
  221. for(int j=0;j<dof;j++)
  222. for(int i=0;i<n3;i++){
  223. size_t x[3]={i%(q+1), (i/(q+1))%(q+1), i/(q+1)/(q+1)};
  224. P0.scal[i+n3*j]*=(x[p_indx-ReflecX]%2?-1.0:1.0);
  225. }
  226. }
  227. { // Set P0.perm
  228. int indx[3]={0,1,2};
  229. if(p_indx==SwapXY) {indx[0]=1; indx[1]=0; indx[2]=2;}
  230. if(p_indx==SwapXZ) {indx[0]=2; indx[1]=1; indx[2]=0;}
  231. for(int j=0;j<dof;j++)
  232. for(int i=0;i<n3;i++){
  233. size_t x[3]={i%(q+1), (i/(q+1))%(q+1), i/(q+1)/(q+1)};
  234. P0.perm[i+n3*j]=x[indx[0]]+(x[indx[1]]+x[indx[2]]*(q+1))*(q+1)
  235. +n3*ker_perm.perm[j];
  236. }
  237. }
  238. std::vector<size_t> coeff_map(n3*dof,0);
  239. {
  240. int indx=0;
  241. for(int j=0;j<dof;j++)
  242. for(int i=0;i<n3;i++){
  243. size_t x[3]={i%(q+1), (i/(q+1))%(q+1), i/(q+1)/(q+1)};
  244. if(x[0]+x[1]+x[2]<=q){
  245. coeff_map[i+n3*j]=indx;
  246. indx++;
  247. }
  248. }
  249. }
  250. Permutation<Real_t> P=Permutation<Real_t>(coeff_cnt*dof);
  251. {
  252. int indx=0;
  253. for(int j=0;j<dof;j++)
  254. for(int i=0;i<n3;i++){
  255. size_t x[3]={i%(q+1), (i/(q+1))%(q+1), i/(q+1)/(q+1)};
  256. if(x[0]+x[1]+x[2]<=q){
  257. P.perm[indx]=coeff_map[P0.perm[i+n3*j]];
  258. P.scal[indx]= P0.scal[i+n3*j] ;
  259. indx++;
  260. }
  261. }
  262. }
  263. return P;
  264. }
  265. template <class FMMNode>
  266. Permutation<typename FMMNode::Real_t>& FMM_Cheb<FMMNode>::PrecompPerm(Mat_Type type, Perm_Type perm_indx){
  267. //int dim=3; //Only supporting 3D
  268. Real_t eps=1e-10;
  269. //Check if the matrix already exists.
  270. Permutation<Real_t>& P_ = FMM_Pts<FMMNode>::PrecompPerm(type, perm_indx);
  271. if(P_.Dim()!=0) return P_;
  272. size_t q=cheb_deg;
  273. size_t m=this->MultipoleOrder();
  274. size_t p_indx=perm_indx % C_Perm;
  275. //Compute the matrix.
  276. Permutation<Real_t> P;
  277. switch (type){
  278. case S2U_Type:
  279. {
  280. Vector<Real_t> scal_exp;
  281. Permutation<Real_t> ker_perm;
  282. if(perm_indx<C_Perm){
  283. ker_perm=this->kernel->k_s2m->perm_vec[0 +p_indx];
  284. scal_exp=this->kernel->k_s2m->src_scal;
  285. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]-=COORD_DIM;
  286. P=cheb_perm(q, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  287. }else{
  288. ker_perm=this->kernel->k_m2m->perm_vec[0 +p_indx];
  289. scal_exp=this->kernel->k_m2m->src_scal;
  290. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]=-scal_exp[i];
  291. // Check that target perm for the two kernels agree.
  292. Permutation<Real_t> ker_perm0=this->kernel->k_s2m->perm_vec[C_Perm+p_indx];
  293. Permutation<Real_t> ker_perm1=this->kernel->k_m2m->perm_vec[C_Perm+p_indx];
  294. assert(ker_perm0.Dim()==ker_perm1.Dim());
  295. if(ker_perm0.Dim()>0 && fabs(ker_perm0.scal[0]-ker_perm1.scal[0])>eps){
  296. for(size_t i=0;i<ker_perm0.Dim();i++){
  297. ker_perm0.scal[i]*=-1.0;
  298. }
  299. for(size_t i=0;i<ker_perm.Dim();i++){
  300. ker_perm.scal[i]*=-1.0;
  301. }
  302. }
  303. for(size_t i=0;i<ker_perm0.Dim();i++){
  304. assert( (ker_perm0.perm[i]-ker_perm1.perm[i])== 0);
  305. assert(fabs(ker_perm0.scal[i]-ker_perm1.scal[i])<eps);
  306. }
  307. Real_t s=0;
  308. // Check that target scaling for the two kernels agree.
  309. const Vector<Real_t>& scal_exp0=this->kernel->k_s2m->trg_scal;
  310. const Vector<Real_t>& scal_exp1=this->kernel->k_m2m->trg_scal;
  311. assert(scal_exp0.Dim()>0 && scal_exp0.Dim()==scal_exp1.Dim());
  312. if(scal_exp0.Dim()>0){
  313. s=(scal_exp0[0]-scal_exp1[0]);
  314. for(size_t i=1;i<scal_exp0.Dim();i++){
  315. assert(fabs(s-(scal_exp0[i]-scal_exp1[i]))<eps);
  316. // In general this is not necessary, but to allow this, we must
  317. // also change src_scal accordingly.
  318. }
  319. }
  320. // Apply the difference in scaling of the two kernels.
  321. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]+=s;
  322. P=equiv_surf_perm(m, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  323. }
  324. break;
  325. }
  326. case D2T_Type:
  327. {
  328. Vector<Real_t> scal_exp;
  329. Permutation<Real_t> ker_perm;
  330. if(perm_indx<C_Perm){ // Source permutation
  331. ker_perm=this->kernel->k_l2l->perm_vec[C_Perm+p_indx];
  332. scal_exp=this->kernel->k_l2l->trg_scal;
  333. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]=-scal_exp[i];
  334. // Check that source perm for the two kernels agree.
  335. Permutation<Real_t> ker_perm0=this->kernel->k_l2t->perm_vec[0 +p_indx];
  336. Permutation<Real_t> ker_perm1=this->kernel->k_l2l->perm_vec[0 +p_indx];
  337. assert(ker_perm0.Dim()==ker_perm1.Dim());
  338. if(ker_perm0.Dim()>0 && fabs(ker_perm0.scal[0]-ker_perm1.scal[0])>eps){
  339. for(size_t i=0;i<ker_perm0.Dim();i++){
  340. ker_perm0.scal[i]*=-1.0;
  341. }
  342. for(size_t i=0;i<ker_perm.Dim();i++){
  343. ker_perm.scal[i]*=-1.0;
  344. }
  345. }
  346. for(size_t i=0;i<ker_perm0.Dim();i++){
  347. assert( (ker_perm0.perm[i]-ker_perm1.perm[i])== 0);
  348. assert(fabs(ker_perm0.scal[i]-ker_perm1.scal[i])<eps);
  349. }
  350. Real_t s=0;
  351. // Check that source scaling for the two kernels agree.
  352. const Vector<Real_t>& scal_exp0=this->kernel->k_l2t->src_scal;
  353. const Vector<Real_t>& scal_exp1=this->kernel->k_l2l->src_scal;
  354. assert(scal_exp0.Dim()>0 && scal_exp0.Dim()==scal_exp1.Dim());
  355. if(scal_exp0.Dim()>0){
  356. s=(scal_exp0[0]-scal_exp1[0]);
  357. for(size_t i=1;i<scal_exp0.Dim();i++){
  358. assert(fabs(s-(scal_exp0[i]-scal_exp1[i]))<eps);
  359. // In general this is not necessary, but to allow this, we must
  360. // also change trg_scal accordingly.
  361. }
  362. }
  363. // Apply the difference in scaling of the two kernels.
  364. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]+=s;
  365. P=equiv_surf_perm(m, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  366. }else{ // Target permutation
  367. ker_perm=this->kernel->k_l2t->perm_vec[C_Perm+p_indx];
  368. scal_exp=this->kernel->k_l2t->trg_scal;
  369. P=cheb_perm(q, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  370. }
  371. break;
  372. }
  373. case U0_Type:
  374. {
  375. Vector<Real_t> scal_exp;
  376. Permutation<Real_t> ker_perm;
  377. if(perm_indx<C_Perm){
  378. ker_perm=this->kernel->k_s2t->perm_vec[0 +p_indx];
  379. scal_exp=this->kernel->k_s2t->src_scal;
  380. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]-=COORD_DIM;
  381. }else{
  382. ker_perm=this->kernel->k_s2t->perm_vec[C_Perm+p_indx];
  383. scal_exp=this->kernel->k_s2t->trg_scal;
  384. }
  385. P=cheb_perm(q, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  386. break;
  387. }
  388. case U1_Type:
  389. {
  390. P=PrecompPerm(U0_Type, perm_indx);
  391. break;
  392. }
  393. case U2_Type:
  394. {
  395. P=PrecompPerm(U0_Type, perm_indx);
  396. break;
  397. }
  398. case W_Type:
  399. {
  400. Vector<Real_t> scal_exp;
  401. Permutation<Real_t> ker_perm;
  402. if(perm_indx<C_Perm){ // Source permutation
  403. ker_perm=this->kernel->k_m2t->perm_vec[0 +p_indx];
  404. scal_exp=this->kernel->k_m2t->src_scal;
  405. P=equiv_surf_perm(m, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  406. }else{ // Target permutation
  407. ker_perm=this->kernel->k_m2t->perm_vec[C_Perm+p_indx];
  408. scal_exp=this->kernel->k_m2t->trg_scal;
  409. P=cheb_perm(q, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  410. }
  411. break;
  412. }
  413. case X_Type:
  414. {
  415. Vector<Real_t> scal_exp;
  416. Permutation<Real_t> ker_perm;
  417. if(perm_indx<C_Perm){
  418. ker_perm=this->kernel->k_s2l->perm_vec[0 +p_indx];
  419. scal_exp=this->kernel->k_s2l->src_scal;
  420. for(size_t i=0;i<scal_exp.Dim();i++) scal_exp[i]-=COORD_DIM;
  421. P=cheb_perm(q, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  422. }else{
  423. ker_perm=this->kernel->k_s2l->perm_vec[C_Perm+p_indx];
  424. scal_exp=this->kernel->k_s2l->trg_scal;
  425. P=equiv_surf_perm(m, p_indx, ker_perm, (this->ScaleInvar()?&scal_exp:NULL));
  426. }
  427. break;
  428. }
  429. default:
  430. return FMM_Pts<FMMNode>::PrecompPerm(type, perm_indx);
  431. break;
  432. }
  433. //Save the matrix for future use.
  434. #pragma omp critical (PRECOMP_MATRIX_PTS)
  435. if(P_.Dim()==0){ P_=P;}
  436. return P_;
  437. }
  438. template <class FMMNode>
  439. Matrix<typename FMMNode::Real_t>& FMM_Cheb<FMMNode>::Precomp(int level, Mat_Type type, size_t mat_indx){
  440. if(this->ScaleInvar()) level=0;
  441. //Check if the matrix already exists.
  442. //Compute matrix from symmetry class (if possible).
  443. Matrix<Real_t>& M_ = this->mat->Mat(level, type, mat_indx);
  444. if(M_.Dim(0)!=0 && M_.Dim(1)!=0) return M_;
  445. else{ //Compute matrix from symmetry class (if possible).
  446. size_t class_indx = this->interac_list.InteracClass(type, mat_indx);
  447. if(class_indx!=mat_indx){
  448. Matrix<Real_t>& M0 = this->Precomp(level, type, class_indx);
  449. Permutation<Real_t>& Pr = this->interac_list.Perm_R(level, type, mat_indx);
  450. Permutation<Real_t>& Pc = this->interac_list.Perm_C(level, type, mat_indx);
  451. if(Pr.Dim()>0 && Pc.Dim()>0 && M0.Dim(0)>0 && M0.Dim(1)>0) return M_;
  452. }
  453. }
  454. int myrank, np;
  455. MPI_Comm_rank(this->comm, &myrank);
  456. MPI_Comm_size(this->comm,&np);
  457. size_t progress=0, class_count=0;
  458. { // Determine precomputation progress.
  459. size_t mat_cnt=this->interac_list.ListCount((Mat_Type)type);
  460. for(size_t i=0; i<mat_cnt; i++){
  461. size_t indx=this->interac_list.InteracClass((Mat_Type)type,i);
  462. if(indx==i){
  463. class_count++;
  464. if(i<mat_indx) progress++;
  465. }
  466. }
  467. }
  468. //Compute the matrix.
  469. Matrix<Real_t> M;
  470. int n_src=((cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3))/6;
  471. switch (type){
  472. case S2U_Type:
  473. {
  474. if(this->MultipoleOrder()==0) break;
  475. Real_t r=pow(0.5,level);
  476. Real_t c[3]={0,0,0};
  477. // Coord of upward check surface
  478. std::vector<Real_t> uc_coord=u_check_surf(this->MultipoleOrder(),c,level);
  479. size_t n_uc=uc_coord.size()/3;
  480. // Evaluate potential at check surface.
  481. Matrix<Real_t> M_s2c(n_src*this->kernel->k_s2m->ker_dim[0],n_uc*this->kernel->k_s2m->ker_dim[1]); //source 2 check
  482. Matrix<Real_t> M_s2c_local(n_src*this->kernel->k_s2m->ker_dim[0],n_uc*this->kernel->k_s2m->ker_dim[1]);
  483. {
  484. M_s2c.SetZero();
  485. M_s2c_local.SetZero();
  486. size_t cnt_done=0;
  487. #pragma omp parallel for schedule(dynamic)
  488. for(size_t i=myrank;i<n_uc;i+=np){
  489. std::vector<Real_t> M_=cheb_integ(cheb_deg, &uc_coord[i*3], r, *this->kernel->k_s2m);
  490. #ifdef __VERBOSE__
  491. #pragma omp critical
  492. if(!myrank){
  493. cnt_done++;
  494. std::cout<<"\r Progress: "<<(100*progress*n_uc+100*cnt_done*np)/(class_count*n_uc)<<"% "<<std::flush;
  495. }
  496. #endif
  497. for(int k=0; k<this->kernel->k_s2m->ker_dim[1]; k++)
  498. for(size_t j=0; j<(size_t)M_s2c.Dim(0); j++)
  499. M_s2c_local[j][i*this->kernel->k_s2m->ker_dim[1]+k] = M_[j+k*M_s2c.Dim(0)];
  500. }
  501. if(!myrank) std::cout<<"\r \r"<<std::flush;
  502. MPI_Allreduce(M_s2c_local[0], M_s2c[0], M_s2c.Dim(0)*M_s2c.Dim(1), par::Mpi_datatype<Real_t>::value(), par::Mpi_datatype<Real_t>::sum(), this->comm);
  503. }
  504. Matrix<Real_t>& M_c2e0 = this->Precomp(level, UC2UE0_Type, 0);
  505. Matrix<Real_t>& M_c2e1 = this->Precomp(level, UC2UE1_Type, 0);
  506. M=(M_s2c*M_c2e0)*M_c2e1;
  507. break;
  508. }
  509. case D2T_Type:
  510. {
  511. if(this->MultipoleOrder()==0) break;
  512. Matrix<Real_t>& M_s2t=FMM_Pts<FMMNode>::Precomp(level, type, mat_indx);
  513. int n_trg=M_s2t.Dim(1)/this->kernel->k_l2t->ker_dim[1];
  514. // Compute Chebyshev approx from target potential.
  515. M.Resize(M_s2t.Dim(0), n_src*this->kernel->k_l2t->ker_dim [1]);
  516. #pragma omp parallel for schedule(dynamic)
  517. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++){
  518. Matrix<Real_t> M_trg(n_trg,this->kernel->k_l2t->ker_dim[1],M_s2t[j],false);
  519. M_trg=M_trg.Transpose();
  520. cheb_approx<Real_t,Real_t>(M_s2t[j],cheb_deg,this->kernel->k_l2t->ker_dim[1],M[j]);
  521. }
  522. #pragma omp critical (PRECOMP_MATRIX_PTS)
  523. {
  524. M_s2t.Resize(0,0);
  525. }
  526. break;
  527. }
  528. case U0_Type:
  529. {
  530. // Coord of target points
  531. Real_t s=pow(0.5,level);
  532. int* coord=this->interac_list.RelativeCoord(type,mat_indx);
  533. Real_t coord_diff[3]={(coord[0]-1)*s*0.5,(coord[1]-1)*s*0.5,(coord[2]-1)*s*0.5};
  534. std::vector<Real_t>& rel_trg_coord = this->mat->RelativeTrgCoord();
  535. size_t n_trg = rel_trg_coord.size()/3;
  536. std::vector<Real_t> trg_coord(n_trg*3);
  537. for(size_t j=0;j<n_trg;j++){
  538. trg_coord[j*3+0]=rel_trg_coord[j*3+0]*s-coord_diff[0];
  539. trg_coord[j*3+1]=rel_trg_coord[j*3+1]*s-coord_diff[1];
  540. trg_coord[j*3+2]=rel_trg_coord[j*3+2]*s-coord_diff[2];
  541. }
  542. // Evaluate potential at target points.
  543. Matrix<Real_t> M_s2t(n_src*this->kernel->k_s2t->ker_dim [0], n_trg*this->kernel->k_s2t->ker_dim [1]);
  544. Matrix<Real_t> M_s2t_local(n_src*this->kernel->k_s2t->ker_dim [0], n_trg*this->kernel->k_s2t->ker_dim [1]);
  545. {
  546. M_s2t.SetZero();
  547. M_s2t_local.SetZero();
  548. size_t cnt_done=0;
  549. #pragma omp parallel for schedule(dynamic)
  550. for(size_t i=myrank;i<n_trg;i+=np){
  551. std::vector<Real_t> s2t=cheb_integ(cheb_deg, &trg_coord[i*3], (Real_t)(s*2.0), *this->kernel->k_s2t);
  552. #ifdef __VERBOSE__
  553. #pragma omp critical
  554. if(!myrank){
  555. cnt_done++;
  556. std::cout<<"\r Progress: "<<(100*progress*n_trg+100*cnt_done*np)/(class_count*n_trg)<<"% "<<std::flush;
  557. }
  558. #endif
  559. for(int k=0; k<this->kernel->k_s2t->ker_dim[1]; k++)
  560. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++)
  561. M_s2t_local[j][i*this->kernel->k_s2t->ker_dim[1]+k] = s2t[j+k*M_s2t.Dim(0)];
  562. }
  563. if(!myrank) std::cout<<"\r \r"<<std::flush;
  564. MPI_Allreduce(M_s2t_local[0], M_s2t[0], M_s2t.Dim(0)*M_s2t.Dim(1), par::Mpi_datatype<Real_t>::value(), par::Mpi_datatype<Real_t>::sum(), this->comm);
  565. }
  566. // Compute Chebyshev approx from target potential.
  567. M.Resize(M_s2t.Dim(0), n_src*this->kernel->k_s2t->ker_dim [1]);
  568. #pragma omp parallel for schedule(dynamic)
  569. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++){
  570. Matrix<Real_t> M_trg(n_trg,this->kernel->k_s2t->ker_dim[1],M_s2t[j],false);
  571. M_trg=M_trg.Transpose();
  572. cheb_approx<Real_t,Real_t>(M_s2t[j],cheb_deg,this->kernel->k_s2t->ker_dim[1],M[j]);
  573. }
  574. break;
  575. }
  576. case U1_Type:
  577. {
  578. // Coord of target points
  579. Real_t s=pow(0.5,level);
  580. int* coord=this->interac_list.RelativeCoord(type,mat_indx);
  581. Real_t coord_diff[3]={coord[0]*s,coord[1]*s,coord[2]*s};
  582. std::vector<Real_t>& rel_trg_coord = this->mat->RelativeTrgCoord();
  583. size_t n_trg = rel_trg_coord.size()/3;
  584. std::vector<Real_t> trg_coord(n_trg*3);
  585. for(size_t j=0;j<n_trg;j++){
  586. trg_coord[j*3+0]=rel_trg_coord[j*3+0]*s-coord_diff[0];
  587. trg_coord[j*3+1]=rel_trg_coord[j*3+1]*s-coord_diff[1];
  588. trg_coord[j*3+2]=rel_trg_coord[j*3+2]*s-coord_diff[2];
  589. }
  590. // Evaluate potential at target points.
  591. Matrix<Real_t> M_s2t(n_src*this->kernel->k_s2t->ker_dim [0], n_trg*this->kernel->k_s2t->ker_dim [1]);
  592. Matrix<Real_t> M_s2t_local(n_src*this->kernel->k_s2t->ker_dim [0], n_trg*this->kernel->k_s2t->ker_dim [1]);
  593. {
  594. M_s2t.SetZero();
  595. M_s2t_local.SetZero();
  596. size_t cnt_done=0;
  597. #pragma omp parallel for schedule(dynamic)
  598. for(size_t i=myrank;i<n_trg;i+=np){
  599. std::vector<Real_t> s2t=cheb_integ(cheb_deg, &trg_coord[i*3], s, *this->kernel->k_s2t);
  600. #ifdef __VERBOSE__
  601. #pragma omp critical
  602. if(!myrank){
  603. cnt_done++;
  604. std::cout<<"\r Progress: "<<(100*progress*n_trg+100*cnt_done*np)/(class_count*n_trg)<<"% "<<std::flush;
  605. }
  606. #endif
  607. for(int k=0; k<this->kernel->k_s2t->ker_dim[1]; k++)
  608. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++)
  609. M_s2t_local[j][i*this->kernel->k_s2t->ker_dim[1]+k] = s2t[j+k*M_s2t.Dim(0)];
  610. }
  611. if(!myrank) std::cout<<"\r \r"<<std::flush;
  612. MPI_Allreduce(M_s2t_local[0], M_s2t[0], M_s2t.Dim(0)*M_s2t.Dim(1), par::Mpi_datatype<Real_t>::value(), par::Mpi_datatype<Real_t>::sum(), this->comm);
  613. }
  614. // Compute Chebyshev approx from target potential.
  615. M.Resize(M_s2t.Dim(0), n_src*this->kernel->k_s2t->ker_dim [1]);
  616. #pragma omp parallel for schedule(dynamic)
  617. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++){
  618. Matrix<Real_t> M_trg(n_trg,this->kernel->k_s2t->ker_dim[1],M_s2t[j],false);
  619. M_trg=M_trg.Transpose();
  620. cheb_approx<Real_t,Real_t>(M_s2t[j],cheb_deg,this->kernel->k_s2t->ker_dim[1],M[j]);
  621. }
  622. break;
  623. }
  624. case U2_Type:
  625. {
  626. // Coord of target points
  627. Real_t s=pow(0.5,level);
  628. int* coord=this->interac_list.RelativeCoord(type,mat_indx);
  629. Real_t coord_diff[3]={(coord[0]+1)*s*0.25,(coord[1]+1)*s*0.25,(coord[2]+1)*s*0.25};
  630. std::vector<Real_t>& rel_trg_coord = this->mat->RelativeTrgCoord();
  631. size_t n_trg = rel_trg_coord.size()/3;
  632. std::vector<Real_t> trg_coord(n_trg*3);
  633. for(size_t j=0;j<n_trg;j++){
  634. trg_coord[j*3+0]=rel_trg_coord[j*3+0]*s-coord_diff[0];
  635. trg_coord[j*3+1]=rel_trg_coord[j*3+1]*s-coord_diff[1];
  636. trg_coord[j*3+2]=rel_trg_coord[j*3+2]*s-coord_diff[2];
  637. }
  638. // Evaluate potential at target points.
  639. Matrix<Real_t> M_s2t(n_src*this->kernel->k_s2t->ker_dim [0], n_trg*this->kernel->k_s2t->ker_dim [1]);
  640. Matrix<Real_t> M_s2t_local(n_src*this->kernel->k_s2t->ker_dim [0], n_trg*this->kernel->k_s2t->ker_dim [1]);
  641. {
  642. M_s2t.SetZero();
  643. M_s2t_local.SetZero();
  644. size_t cnt_done=0;
  645. #pragma omp parallel for schedule(dynamic)
  646. for(size_t i=myrank;i<n_trg;i+=np){
  647. std::vector<Real_t> s2t=cheb_integ(cheb_deg, &trg_coord[i*3], (Real_t)(s*0.5), *this->kernel->k_s2t);
  648. #ifdef __VERBOSE__
  649. #pragma omp critical
  650. if(!myrank){
  651. cnt_done++;
  652. std::cout<<"\r Progress: "<<(100*progress*n_trg+100*cnt_done*np)/(class_count*n_trg)<<"% "<<std::flush;
  653. }
  654. #endif
  655. for(int k=0; k<this->kernel->k_s2t->ker_dim[1]; k++)
  656. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++)
  657. M_s2t_local[j][i*this->kernel->k_s2t->ker_dim[1]+k] = s2t[j+k*M_s2t.Dim(0)];
  658. }
  659. if(!myrank) std::cout<<"\r \r"<<std::flush;
  660. MPI_Allreduce(M_s2t_local[0], M_s2t[0], M_s2t.Dim(0)*M_s2t.Dim(1), par::Mpi_datatype<Real_t>::value(), par::Mpi_datatype<Real_t>::sum(), this->comm);
  661. }
  662. // Compute Chebyshev approx from target potential.
  663. M.Resize(M_s2t.Dim(0), n_src*this->kernel->k_s2t->ker_dim [1]);
  664. #pragma omp parallel for schedule(dynamic)
  665. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++){
  666. Matrix<Real_t> M_trg(n_trg,this->kernel->k_s2t->ker_dim[1],M_s2t[j],false);
  667. M_trg=M_trg.Transpose();
  668. cheb_approx<Real_t,Real_t>(M_s2t[j],cheb_deg,this->kernel->k_s2t->ker_dim[1],M[j]);
  669. }
  670. break;
  671. }
  672. case W_Type:
  673. {
  674. if(this->MultipoleOrder()==0) break;
  675. Matrix<Real_t>& M_s2t=FMM_Pts<FMMNode>::Precomp(level, type, mat_indx);
  676. int n_trg=M_s2t.Dim(1)/this->kernel->k_m2t->ker_dim[1];
  677. // Compute Chebyshev approx from target potential.
  678. M.Resize(M_s2t.Dim(0), n_src*this->kernel->k_m2t->ker_dim [1]);
  679. #pragma omp parallel for schedule(dynamic)
  680. for(size_t j=0; j<(size_t)M_s2t.Dim(0); j++){
  681. Matrix<Real_t> M_trg(n_trg,this->kernel->k_m2t->ker_dim[1],M_s2t[j],false);
  682. M_trg=M_trg.Transpose();
  683. cheb_approx<Real_t,Real_t>(M_s2t[j],cheb_deg,this->kernel->k_m2t->ker_dim[1],M[j]);
  684. }
  685. #pragma omp critical (PRECOMP_MATRIX_PTS)
  686. {
  687. M_s2t.Resize(0,0);
  688. }
  689. break;
  690. }
  691. case X_Type:
  692. {
  693. if(this->MultipoleOrder()==0) break;
  694. // Coord of target points
  695. Real_t s=pow(0.5,level-1);
  696. int* coord=this->interac_list.RelativeCoord(type,mat_indx);
  697. Real_t c[3]={-(coord[0]-1)*s*0.25,-(coord[1]-1)*s*0.25,-(coord[2]-1)*s*0.25};
  698. std::vector<Real_t> trg_coord=d_check_surf(this->MultipoleOrder(),c,level);
  699. size_t n_trg=trg_coord.size()/3;
  700. // Evaluate potential at target points.
  701. Matrix<Real_t> M_xs2c(n_src*this->kernel->k_s2l->ker_dim[0], n_trg*this->kernel->k_s2l->ker_dim[1]);
  702. Matrix<Real_t> M_xs2c_local(n_src*this->kernel->k_s2l->ker_dim[0], n_trg*this->kernel->k_s2l->ker_dim[1]);
  703. {
  704. M_xs2c.SetZero();
  705. M_xs2c_local.SetZero();
  706. size_t cnt_done=0;
  707. #pragma omp parallel for schedule(dynamic)
  708. for(size_t i=myrank;i<n_trg;i+=np){
  709. std::vector<Real_t> M_=cheb_integ(cheb_deg, &trg_coord[i*3], s, *this->kernel->k_s2l);
  710. #ifdef __VERBOSE__
  711. #pragma omp critical
  712. if(!myrank){
  713. cnt_done++;
  714. std::cout<<"\r Progress: "<<(100*progress*n_trg+100*cnt_done*np)/(class_count*n_trg)<<"% "<<std::flush;
  715. }
  716. #endif
  717. for(int k=0; k<this->kernel->k_s2l->ker_dim[1]; k++)
  718. for(size_t j=0; j<(size_t)M_xs2c.Dim(0); j++)
  719. M_xs2c_local[j][i*this->kernel->k_s2l->ker_dim[1]+k] = M_[j+k*M_xs2c.Dim(0)];
  720. }
  721. if(!myrank) std::cout<<"\r \r"<<std::flush;
  722. MPI_Allreduce(M_xs2c_local[0], M_xs2c[0], M_xs2c.Dim(0)*M_xs2c.Dim(1), par::Mpi_datatype<Real_t>::value(), par::Mpi_datatype<Real_t>::sum(), this->comm);
  723. }
  724. M=M_xs2c;
  725. break;
  726. }
  727. default:
  728. {
  729. return FMM_Pts<FMMNode>::Precomp(level, type, mat_indx);
  730. break;
  731. }
  732. }
  733. //Save the matrix for future use.
  734. #pragma omp critical (PRECOMP_MATRIX_PTS)
  735. if(M_.Dim(0)==0 && M_.Dim(1)==0){ M_=M;}
  736. return M_;
  737. }
  738. template <class FMMNode>
  739. void FMM_Cheb<FMMNode>::CollectNodeData(FMMTree_t* tree, std::vector<FMMNode*>& node, std::vector<Matrix<Real_t> >& buff, std::vector<Vector<FMMNode_t*> >& n_list, std::vector<std::vector<Vector<Real_t>* > > vec_list){
  740. if(vec_list.size()<6) vec_list.resize(6);
  741. size_t n_coeff=(cheb_deg+1)*(cheb_deg+2)*(cheb_deg+3)/6;
  742. if(node.size()==0) return;
  743. {// 4. cheb_in
  744. int indx=4;
  745. size_t vec_sz=this->kernel->ker_dim[0]*n_coeff;
  746. for(size_t i=0;i<node.size();i++){
  747. if(node[i]->IsLeaf()){
  748. Vector<Real_t>& data_vec=node[i]->ChebData();
  749. if(data_vec.Dim()!=vec_sz) data_vec.ReInit(vec_sz);
  750. vec_list[indx].push_back(&data_vec);
  751. }
  752. }
  753. }
  754. {// 5. cheb_out
  755. int indx=5;
  756. size_t vec_sz=this->kernel->ker_dim[1]*n_coeff;
  757. for(size_t i=0;i<node.size();i++){
  758. if(node[i]->IsLeaf() && !node[i]->IsGhost()){
  759. Vector<Real_t>& data_vec=((FMMData*)node[i]->FMMData())->cheb_out;
  760. if(data_vec.Dim()!=vec_sz) data_vec.ReInit(vec_sz);
  761. vec_list[indx].push_back(&data_vec);
  762. }
  763. }
  764. }
  765. FMM_Pts<FMMNode_t>::CollectNodeData(tree, node, buff, n_list, vec_list);
  766. #pragma omp parallel for
  767. for(size_t i=0;i<node.size();i++){
  768. node[i]->pt_cnt[0]+=2*n_coeff;
  769. node[i]->pt_cnt[1]+=2*n_coeff;
  770. }
  771. }
  772. template <class FMMNode>
  773. void FMM_Cheb<FMMNode>::Source2UpSetup(SetupData<Real_t>& setup_data, FMMTree_t* tree, std::vector<Matrix<Real_t> >& buff, std::vector<Vector<FMMNode_t*> >& n_list, int level, bool device){
  774. if(!this->MultipoleOrder()) return;
  775. FMM_Pts<FMMNode>::Source2UpSetup(setup_data, tree, buff, n_list, level, device);
  776. { // Set setup_data
  777. setup_data.level=level;
  778. setup_data.kernel=this->kernel->k_s2m;
  779. setup_data.interac_type.resize(1);
  780. setup_data.interac_type[0]=S2U_Type;
  781. setup_data. input_data=&buff[4];
  782. setup_data.output_data=&buff[0];
  783. Vector<FMMNode_t*>& nodes_in =n_list[4];
  784. Vector<FMMNode_t*>& nodes_out=n_list[0];
  785. setup_data.nodes_in .clear();
  786. setup_data.nodes_out.clear();
  787. for(size_t i=0;i<nodes_in .Dim();i++) if(nodes_in [i]->Depth()==level || level==-1) setup_data.nodes_in .push_back(nodes_in [i]);
  788. for(size_t i=0;i<nodes_out.Dim();i++) if(nodes_out[i]->Depth()==level || level==-1) setup_data.nodes_out.push_back(nodes_out[i]);
  789. }
  790. std::vector<void*>& nodes_in =setup_data.nodes_in ;
  791. std::vector<void*>& nodes_out=setup_data.nodes_out;
  792. std::vector<Vector<Real_t>*>& input_vector=setup_data. input_vector; input_vector.clear();
  793. std::vector<Vector<Real_t>*>& output_vector=setup_data.output_vector; output_vector.clear();
  794. for(size_t i=0;i<nodes_in .size();i++) input_vector.push_back(&( ((FMMNode*)nodes_in [i]) )->ChebData() );
  795. for(size_t i=0;i<nodes_out.size();i++) output_vector.push_back(&((FMMData*)((FMMNode*)nodes_out[i])->FMMData())->upward_equiv);
  796. this->SetupInterac(setup_data,device);
  797. }
  798. template <class FMMNode>
  799. void FMM_Cheb<FMMNode>::Source2Up (SetupData<Real_t>& setup_data, bool device){
  800. if(!this->MultipoleOrder()) return;
  801. //Add Source2Up contribution.
  802. this->EvalList(setup_data, device);
  803. }
  804. template <class FMMNode>
  805. void FMM_Cheb<FMMNode>::X_ListSetup(SetupData<Real_t>& setup_data, FMMTree_t* tree, std::vector<Matrix<Real_t> >& buff, std::vector<Vector<FMMNode_t*> >& n_list, int level, bool device){
  806. if(!this->MultipoleOrder()) return;
  807. FMM_Pts<FMMNode>::X_ListSetup(setup_data, tree, buff, n_list, level, device);
  808. { // Set setup_data
  809. setup_data.level=level;
  810. setup_data.kernel=this->kernel->k_s2l;
  811. setup_data.interac_type.resize(1);
  812. setup_data.interac_type[0]=X_Type;
  813. setup_data. input_data=&buff[4];
  814. setup_data.output_data=&buff[1];
  815. Vector<FMMNode_t*>& nodes_in =n_list[4];
  816. Vector<FMMNode_t*>& nodes_out=n_list[1];
  817. setup_data.nodes_in .clear();
  818. setup_data.nodes_out.clear();
  819. for(size_t i=0;i<nodes_in .Dim();i++) if(nodes_in [i]->Depth()==level-1 || level==-1) setup_data.nodes_in .push_back(nodes_in [i]);
  820. for(size_t i=0;i<nodes_out.Dim();i++) if(nodes_out[i]->Depth()==level || level==-1) setup_data.nodes_out.push_back(nodes_out[i]);
  821. }
  822. std::vector<void*>& nodes_in =setup_data.nodes_in ;
  823. std::vector<void*>& nodes_out=setup_data.nodes_out;
  824. std::vector<Vector<Real_t>*>& input_vector=setup_data. input_vector; input_vector.clear();
  825. std::vector<Vector<Real_t>*>& output_vector=setup_data.output_vector; output_vector.clear();
  826. for(size_t i=0;i<nodes_in .size();i++) input_vector.push_back(&( ((FMMNode*)nodes_in [i]) )->ChebData() );
  827. for(size_t i=0;i<nodes_out.size();i++) output_vector.push_back(&((FMMData*)((FMMNode*)nodes_out[i])->FMMData())->dnward_equiv);
  828. this->SetupInterac(setup_data,device);
  829. { // Resize device buffer
  830. size_t n=setup_data.output_data->Dim(0)*setup_data.output_data->Dim(1)*sizeof(Real_t);
  831. if(this->dev_buffer.Dim()<n) this->dev_buffer.ReInit(n);
  832. }
  833. }
  834. template <class FMMNode>
  835. void FMM_Cheb<FMMNode>::X_List (SetupData<Real_t>& setup_data, bool device){
  836. if(!this->MultipoleOrder()) return;
  837. //Add X_List contribution.
  838. FMM_Pts<FMMNode>::X_List(setup_data, device);
  839. this->EvalList(setup_data, device);
  840. }
  841. template <class FMMNode>
  842. void FMM_Cheb<FMMNode>::W_ListSetup(SetupData<Real_t>& setup_data, FMMTree_t* tree, std::vector<Matrix<Real_t> >& buff, std::vector<Vector<FMMNode_t*> >& n_list, int level, bool device){
  843. if(!this->MultipoleOrder()) return;
  844. { // Set setup_data
  845. setup_data.level=level;
  846. setup_data.kernel=this->kernel->k_m2t;
  847. setup_data.interac_type.resize(1);
  848. setup_data.interac_type[0]=W_Type;
  849. setup_data. input_data=&buff[0];
  850. setup_data.output_data=&buff[5];
  851. Vector<FMMNode_t*>& nodes_in =n_list[0];
  852. Vector<FMMNode_t*>& nodes_out=n_list[5];
  853. setup_data.nodes_in .clear();
  854. setup_data.nodes_out.clear();
  855. for(size_t i=0;i<nodes_in .Dim();i++) if(nodes_in [i]->Depth()==level+1 || level==-1) setup_data.nodes_in .push_back(nodes_in [i]);
  856. for(size_t i=0;i<nodes_out.Dim();i++) if(nodes_out[i]->Depth()==level || level==-1) setup_data.nodes_out.push_back(nodes_out[i]);
  857. }
  858. std::vector<void*>& nodes_in =setup_data.nodes_in ;
  859. std::vector<void*>& nodes_out=setup_data.nodes_out;
  860. std::vector<Vector<Real_t>*>& input_vector=setup_data. input_vector; input_vector.clear();
  861. std::vector<Vector<Real_t>*>& output_vector=setup_data.output_vector; output_vector.clear();
  862. for(size_t i=0;i<nodes_in .size();i++) input_vector.push_back(&((FMMData*)((FMMNode*)nodes_in [i])->FMMData())->upward_equiv);
  863. for(size_t i=0;i<nodes_out.size();i++) output_vector.push_back(&((FMMData*)((FMMNode*)nodes_out[i])->FMMData())->cheb_out );
  864. this->SetupInterac(setup_data,device);
  865. { // Resize device buffer
  866. size_t n=setup_data.output_data->Dim(0)*setup_data.output_data->Dim(1)*sizeof(Real_t);
  867. if(this->dev_buffer.Dim()<n) this->dev_buffer.ReInit(n);
  868. }
  869. }
  870. template <class FMMNode>
  871. void FMM_Cheb<FMMNode>::W_List (SetupData<Real_t>& setup_data, bool device){
  872. if(!this->MultipoleOrder()) return;
  873. //Add W_List contribution.
  874. this->EvalList(setup_data, device);
  875. }
  876. template <class FMMNode>
  877. void FMM_Cheb<FMMNode>::U_ListSetup(SetupData<Real_t>& setup_data, FMMTree_t* tree, std::vector<Matrix<Real_t> >& buff, std::vector<Vector<FMMNode_t*> >& n_list, int level, bool device){
  878. FMM_Pts<FMMNode>::U_ListSetup(setup_data, tree, buff, n_list, level, device);
  879. { // Set setup_data
  880. setup_data.level=level;
  881. setup_data.kernel=this->kernel->k_s2t;
  882. setup_data.interac_type.resize(3);
  883. setup_data.interac_type[0]=U0_Type;
  884. setup_data.interac_type[1]=U1_Type;
  885. setup_data.interac_type[2]=U2_Type;
  886. setup_data. input_data=&buff[4];
  887. setup_data.output_data=&buff[5];
  888. Vector<FMMNode_t*>& nodes_in =n_list[4];
  889. Vector<FMMNode_t*>& nodes_out=n_list[5];
  890. setup_data.nodes_in .clear();
  891. setup_data.nodes_out.clear();
  892. for(size_t i=0;i<nodes_in .Dim();i++) if((level-1<=nodes_in [i]->Depth() && nodes_in [i]->Depth()<=level+1) || level==-1) setup_data.nodes_in .push_back(nodes_in [i]);
  893. for(size_t i=0;i<nodes_out.Dim();i++) if(( nodes_out[i]->Depth()==level ) || level==-1) setup_data.nodes_out.push_back(nodes_out[i]);
  894. }
  895. std::vector<void*>& nodes_in =setup_data.nodes_in ;
  896. std::vector<void*>& nodes_out=setup_data.nodes_out;
  897. std::vector<Vector<Real_t>*>& input_vector=setup_data. input_vector; input_vector.clear();
  898. std::vector<Vector<Real_t>*>& output_vector=setup_data.output_vector; output_vector.clear();
  899. for(size_t i=0;i<nodes_in .size();i++) input_vector.push_back(&( ((FMMNode*)nodes_in [i]) )->ChebData());
  900. for(size_t i=0;i<nodes_out.size();i++) output_vector.push_back(&((FMMData*)((FMMNode*)nodes_out[i])->FMMData())->cheb_out );
  901. this->SetupInterac(setup_data,device);
  902. { // Resize device buffer
  903. size_t n=setup_data.output_data->Dim(0)*setup_data.output_data->Dim(1)*sizeof(Real_t);
  904. if(this->dev_buffer.Dim()<n) this->dev_buffer.ReInit(n);
  905. }
  906. }
  907. template <class FMMNode>
  908. void FMM_Cheb<FMMNode>::U_List (SetupData<Real_t>& setup_data, bool device){
  909. //Add U_List contribution.
  910. FMM_Pts<FMMNode>::U_List(setup_data, device);
  911. this->EvalList(setup_data, device);
  912. }
  913. template <class FMMNode>
  914. void FMM_Cheb<FMMNode>::Down2TargetSetup(SetupData<Real_t>& setup_data, FMMTree_t* tree, std::vector<Matrix<Real_t> >& buff, std::vector<Vector<FMMNode_t*> >& n_list, int level, bool device){
  915. if(!this->MultipoleOrder()) return;
  916. { // Set setup_data
  917. setup_data.level=level;
  918. setup_data.kernel=this->kernel->k_l2t;
  919. setup_data.interac_type.resize(1);
  920. setup_data.interac_type[0]=D2T_Type;
  921. setup_data. input_data=&buff[1];
  922. setup_data.output_data=&buff[5];
  923. Vector<FMMNode_t*>& nodes_in =n_list[1];
  924. Vector<FMMNode_t*>& nodes_out=n_list[5];
  925. setup_data.nodes_in .clear();
  926. setup_data.nodes_out.clear();
  927. for(size_t i=0;i<nodes_in .Dim();i++) if(nodes_in [i]->Depth()==level || level==-1) setup_data.nodes_in .push_back(nodes_in [i]);
  928. for(size_t i=0;i<nodes_out.Dim();i++) if(nodes_out[i]->Depth()==level || level==-1) setup_data.nodes_out.push_back(nodes_out[i]);
  929. }
  930. std::vector<void*>& nodes_in =setup_data.nodes_in ;
  931. std::vector<void*>& nodes_out=setup_data.nodes_out;
  932. std::vector<Vector<Real_t>*>& input_vector=setup_data. input_vector; input_vector.clear();
  933. std::vector<Vector<Real_t>*>& output_vector=setup_data.output_vector; output_vector.clear();
  934. for(size_t i=0;i<nodes_in .size();i++) input_vector.push_back(&((FMMData*)((FMMNode*)nodes_in [i])->FMMData())->dnward_equiv);
  935. for(size_t i=0;i<nodes_out.size();i++) output_vector.push_back(&((FMMData*)((FMMNode*)nodes_out[i])->FMMData())->cheb_out );
  936. this->SetupInterac(setup_data,device);
  937. }
  938. template <class FMMNode>
  939. void FMM_Cheb<FMMNode>::Down2Target (SetupData<Real_t>& setup_data, bool device){
  940. if(!this->MultipoleOrder()) return;
  941. //Add Down2Target contribution.
  942. this->EvalList(setup_data, device);
  943. }
  944. template <class FMMNode>
  945. void FMM_Cheb<FMMNode>::PostProcessing(std::vector<FMMNode_t*>& nodes){
  946. size_t n=nodes.size();
  947. #pragma omp parallel
  948. {
  949. int omp_p=omp_get_num_threads();
  950. int pid = omp_get_thread_num();
  951. size_t a=(pid*n)/omp_p;
  952. size_t b=((pid+1)*n)/omp_p;
  953. //For each node, compute interaction from point sources.
  954. for(size_t i=a;i<b;i++){
  955. Vector<Real_t>& trg_coord=nodes[i]->trg_coord;
  956. Vector<Real_t>& trg_value=nodes[i]->trg_value;
  957. Vector<Real_t>& cheb_out =((FMMData*)nodes[i]->FMMData())->cheb_out;
  958. //Initialize target potential.
  959. size_t trg_cnt=trg_coord.Dim()/COORD_DIM;
  960. //trg_value.assign(trg_cnt*dof*this->kernel->ker_dim[1],0);
  961. //Sample local expansion at target points.
  962. if(trg_cnt>0 && cheb_out.Dim()>0){
  963. Real_t* c=nodes[i]->Coord();
  964. Real_t scale=pow(2.0,nodes[i]->Depth()+1);
  965. std::vector<Real_t> rel_coord(COORD_DIM*trg_cnt);
  966. for(size_t j=0;j<trg_cnt;j++) for(int k=0;k<COORD_DIM;k++)
  967. rel_coord[j*COORD_DIM+k]=(trg_coord[j*COORD_DIM+k]-c[k])*scale-1.0;
  968. cheb_eval(cheb_out, cheb_deg, rel_coord, trg_value);
  969. }
  970. }
  971. }
  972. }
  973. }//end namespace