utils.txx 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421
  1. /**
  2. * \file mat_utils.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 1-1-2013
  5. */
  6. template <class FMM_Mat_t>
  7. void CheckFMMOutput(pvfmm::FMM_Tree<FMM_Mat_t>* mytree, pvfmm::Kernel<typename FMM_Mat_t::Real_t>* mykernel){
  8. if(mykernel==NULL) return;
  9. // Find out number of OMP thereads.
  10. int np=omp_get_max_threads();
  11. // Find out my identity in the default communicator
  12. int myrank, p;
  13. MPI_Comm c1=MPI_COMM_WORLD;
  14. MPI_Comm_rank(c1, &myrank);
  15. MPI_Comm_size(c1,&p);
  16. typedef typename FMM_Mat_t::FMMData FMM_Data_t;
  17. typedef typename FMM_Mat_t::FMMNode_t FMMNode_t;
  18. typedef typename FMM_Mat_t::Real_t Real_t;
  19. // Read source data.
  20. std::vector<Real_t> src_coord;
  21. std::vector<Real_t> src_value;
  22. FMMNode_t* n=static_cast<FMMNode_t*>(mytree->PreorderFirst());
  23. while(n!=NULL){
  24. if(n->IsLeaf() && !n->IsGhost()){
  25. pvfmm::Vector<Real_t>& coord_vec=n->src_coord;
  26. pvfmm::Vector<Real_t>& value_vec=n->src_value;
  27. for(size_t i=0;i<coord_vec.Dim();i++) src_coord.push_back(coord_vec[i]);
  28. for(size_t i=0;i<value_vec.Dim();i++) src_value.push_back(value_vec[i]);
  29. }
  30. n=static_cast<FMMNode_t*>(mytree->PreorderNxt(n));
  31. }
  32. long long glb_src_cnt=0, src_cnt=src_coord.size()/3;
  33. MPI_Allreduce(&src_cnt, &glb_src_cnt, 1, MPI_LONG_LONG, MPI_SUM, c1);
  34. long long glb_val_cnt=0, val_cnt=src_value.size();
  35. MPI_Allreduce(&val_cnt, &glb_val_cnt, 1, MPI_LONG_LONG, MPI_SUM, c1);
  36. if(glb_src_cnt==0) return;
  37. int dof=glb_val_cnt/glb_src_cnt/mykernel->ker_dim[0];
  38. int trg_dof=dof*mykernel->ker_dim[1];
  39. // Read target data.
  40. std::vector<Real_t> trg_coord;
  41. std::vector<Real_t> trg_poten_fmm;
  42. long long trg_iter=0;
  43. size_t step_size=1+glb_src_cnt*glb_src_cnt*1e-9/p;
  44. n=static_cast<FMMNode_t*>(mytree->PreorderFirst());
  45. while(n!=NULL){
  46. if(n->IsLeaf() && !n->IsGhost()){
  47. pvfmm::Vector<Real_t>& coord_vec=n->trg_coord;
  48. pvfmm::Vector<Real_t>& poten_vec=n->trg_value;
  49. for(size_t i=0;i<coord_vec.Dim()/3 ;i++){
  50. if(trg_iter%step_size==0){
  51. for(int j=0;j<3 ;j++) trg_coord .push_back(coord_vec[i*3 +j]);
  52. for(int j=0;j<trg_dof ;j++) trg_poten_fmm.push_back(poten_vec[i*trg_dof +j]);
  53. }
  54. trg_iter++;
  55. }
  56. }
  57. n=static_cast<FMMNode_t*>(mytree->PreorderNxt(n));
  58. }
  59. int trg_cnt=trg_coord.size()/3;
  60. int send_cnt=trg_cnt*3;
  61. std::vector<int> recv_cnts(p), recv_disp(p,0);
  62. MPI_Allgather(&send_cnt , 1, MPI_INT,
  63. &recv_cnts[0], 1, MPI_INT, c1);
  64. pvfmm::omp_par::scan(&recv_cnts[0], &recv_disp[0], p);
  65. int glb_trg_cnt=(recv_disp[p-1]+recv_cnts[p-1])/3;
  66. std::vector<Real_t> glb_trg_coord(glb_trg_cnt*3);
  67. MPI_Allgatherv(&trg_coord[0] , send_cnt , pvfmm::par::Mpi_datatype<Real_t>::value(),
  68. &glb_trg_coord[0], &recv_cnts[0], &recv_disp[0], pvfmm::par::Mpi_datatype<Real_t>::value(), c1);
  69. if(glb_trg_cnt==0) return;
  70. //Direct N-Body.
  71. std::vector<Real_t> trg_poten_dir(glb_trg_cnt*trg_dof ,0);
  72. std::vector<Real_t> glb_trg_poten_dir(glb_trg_cnt*trg_dof ,0);
  73. pvfmm::Profile::Tic("N-Body Direct",&c1,false,1);
  74. #pragma omp parallel for
  75. for(int i=0;i<np;i++){
  76. size_t a=(i*glb_trg_cnt)/np;
  77. size_t b=((i+1)*glb_trg_cnt)/np;
  78. mykernel->ker_poten(&src_coord[0], src_cnt, &src_value[0], dof, &glb_trg_coord[a*3], b-a, &trg_poten_dir[a*trg_dof ],NULL);
  79. }
  80. MPI_Allreduce(&trg_poten_dir[0], &glb_trg_poten_dir[0], trg_poten_dir.size(), pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::sum(), c1);
  81. pvfmm::Profile::Toc();
  82. //Compute error.
  83. {
  84. Real_t max_=0;
  85. Real_t max_err=0;
  86. for(size_t i=0;i<trg_poten_fmm.size();i++){
  87. Real_t err=fabs(glb_trg_poten_dir[i+(recv_disp[myrank]/3)*trg_dof]-trg_poten_fmm[i]);
  88. Real_t max=fabs(glb_trg_poten_dir[i+(recv_disp[myrank]/3)*trg_dof]);
  89. if(err>max_err) max_err=err;
  90. if(max>max_) max_=max;
  91. }
  92. Real_t glb_max, glb_max_err;
  93. MPI_Reduce(&max_ , &glb_max , 1, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::max(), 0, c1);
  94. MPI_Reduce(&max_err, &glb_max_err, 1, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::max(), 0, c1);
  95. if(!myrank){
  96. std::cout<<"Maximum Absolute Error ["<<mykernel->ker_name<<"] : "<<std::scientific<<glb_max_err<<'\n';
  97. std::cout<<"Maximum Relative Error ["<<mykernel->ker_name<<"] : "<<std::scientific<<glb_max_err/glb_max<<'\n';
  98. }
  99. }
  100. }
  101. template <class FMMTree_t>
  102. void CheckChebOutput(FMMTree_t* mytree, typename TestFn<typename FMMTree_t::Real_t>::Fn_t fn_poten, int fn_dof, std::string t_name){
  103. typedef typename FMMTree_t::FMM_Node_t FMMNode_t;
  104. typedef typename FMMTree_t::Real_t Real_t;
  105. MPI_Comm c1=*mytree->Comm();
  106. pvfmm::Profile::Tic((std::string("Compute Error ")+t_name).c_str(),&c1,true,1);
  107. int myrank; MPI_Comm_rank(c1, &myrank);
  108. FMMNode_t* r_node=static_cast<FMMNode_t*>(mytree->RootNode());
  109. int dof=r_node->DataDOF()/fn_dof;
  110. std::vector<FMMNode_t*> nodes;
  111. {
  112. FMMNode_t* node=static_cast<FMMNode_t*>(mytree->PreorderFirst());
  113. while(node!=NULL){
  114. if(node->IsLeaf() && !node->IsGhost()) nodes.push_back(node);
  115. node=static_cast<FMMNode_t*>(mytree->PreorderNxt(node));
  116. }
  117. if(nodes.size()==0) return;
  118. }
  119. int cheb_deg=nodes[0]->ChebDeg();
  120. std::vector<Real_t> cheb_nds=pvfmm::cheb_nodes<Real_t>(cheb_deg/3+1, 1);
  121. for(size_t i=0;i<cheb_nds.size();i++) cheb_nds[i]=2.0*cheb_nds[i]-1.0;
  122. std::vector<Real_t> cheb_pts=pvfmm::cheb_nodes<Real_t>(cheb_deg/3+1, COORD_DIM);
  123. int n_pts=cheb_pts.size()/COORD_DIM;
  124. int omp_p=omp_get_max_threads();
  125. std::vector<Real_t> glb_err_avg(dof*fn_dof,0);
  126. { // Determine glb_err_avg.
  127. std::vector<Real_t> err_avg(omp_p*dof*fn_dof,0);
  128. #pragma omp parallel for
  129. for(size_t tid=0;tid<omp_p;tid++){
  130. pvfmm::Vector<Real_t> out; out.SetZero();
  131. pvfmm::Vector<Real_t> fn_out(dof*fn_dof);
  132. for(size_t i=(nodes.size()*tid)/omp_p;i<(nodes.size()*(tid+1))/omp_p;i++){
  133. pvfmm::Vector<Real_t>& cheb_coeff=nodes[i]->ChebData();
  134. cheb_eval(cheb_coeff, cheb_deg, cheb_nds, cheb_nds, cheb_nds, out);
  135. Real_t* c=nodes[i]->Coord();
  136. Real_t s=pow(2.0,-nodes[i]->Depth());
  137. Real_t s3=s*s*s;
  138. for(size_t j=0;j<n_pts;j++){
  139. Real_t coord[3]={c[0]+s*cheb_pts[j*COORD_DIM+0],
  140. c[1]+s*cheb_pts[j*COORD_DIM+1],
  141. c[2]+s*cheb_pts[j*COORD_DIM+2]};
  142. fn_out.SetZero();
  143. for(size_t k=0;k<dof;k++)
  144. fn_poten(coord,1,&fn_out[k*fn_dof]);
  145. for(size_t k=0;k<dof*fn_dof;k++){
  146. Real_t err=out[n_pts*k+j]-fn_out[k];
  147. err_avg[tid*dof*fn_dof+k]+=err*s3;
  148. }
  149. }
  150. }
  151. for(size_t k=0;k<dof*fn_dof;k++)
  152. err_avg[tid*dof*fn_dof+k]/=n_pts;
  153. }
  154. for(size_t tid=1;tid<omp_p;tid++)
  155. for(size_t k=0;k<dof*fn_dof;k++)
  156. err_avg[k]+=err_avg[tid*dof*fn_dof+k];
  157. MPI_Allreduce(&err_avg[0], &glb_err_avg[0], dof*fn_dof, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::sum(), c1);
  158. }
  159. { // Write error to file.
  160. int nn_x=1;
  161. int nn_y=201;
  162. int nn_z=201;
  163. std::vector<Real_t> x(nn_x,1.0/2.0);
  164. std::vector<Real_t> y(nn_y);
  165. std::vector<Real_t> z(nn_z);
  166. //for(int i=0;i<nn_x;i++)
  167. // x[i]=((Real_t)i)/(nn_x-1);
  168. for(int i=0;i<nn_y;i++)
  169. y[i]=((Real_t)i)/(nn_y-1);
  170. for(int i=0;i<nn_z;i++)
  171. z[i]=((Real_t)i)/(nn_z-1);
  172. FMMNode_t* r_node=static_cast<FMMNode_t*>(mytree->RootNode());
  173. int fn_dof=r_node->DataDOF()/dof;
  174. Real_t* fn_out=new Real_t[dof*fn_dof];
  175. pvfmm::Matrix<Real_t> M_out (nn_z*fn_dof*dof,nn_y*nn_x,NULL,true); M_out .SetZero();
  176. pvfmm::Matrix<Real_t> M_out_err(nn_z*fn_dof*dof,nn_y*nn_x,NULL,true); M_out_err.SetZero();
  177. r_node->ReadVal(x,y,z,M_out[0],false);
  178. for(int l0=0;l0<dof;l0++)
  179. for(int l1=0;l1<fn_dof;l1++)
  180. for(int i=0;i<nn_x;i++)
  181. for(int j=0;j<nn_y;j++)
  182. for(int k=0;k<nn_z;k++){
  183. int l=l0*fn_dof+l1;
  184. Real_t ch_coord[3]={x[i],y[j],z[k]};
  185. Real_t cheb_val=M_out[k+l*nn_z][i+j*nn_x];
  186. if(fabs(cheb_val)>1.0e-20){
  187. fn_poten(ch_coord,1,fn_out);
  188. Real_t err=(cheb_val-fn_out[l])-glb_err_avg[l];
  189. M_out_err[k+l*nn_z][i+j*nn_x]=err;
  190. }
  191. }
  192. delete[] fn_out;
  193. pvfmm::Matrix<Real_t> M_global (nn_z*fn_dof*dof,nn_y*nn_x,NULL,true);
  194. pvfmm::Matrix<Real_t> M_global_err(nn_z*fn_dof*dof,nn_y*nn_x,NULL,true);
  195. MPI_Reduce(&M_out [0][0], &M_global [0][0], nn_x*nn_y*nn_z*fn_dof*dof, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::sum(), 0, c1);
  196. MPI_Reduce(&M_out_err[0][0], &M_global_err[0][0], nn_x*nn_y*nn_z*fn_dof*dof, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::sum(), 0, c1);
  197. //std::string fname;
  198. //fname=std::string("result/" )+t_name+std::string(".mat");
  199. //if(!myrank) M_global .Write(fname.c_str());
  200. //fname=std::string("result/err_")+t_name+std::string(".mat");
  201. //if(!myrank) M_global_err.Write(fname.c_str());
  202. }
  203. std::vector<Real_t> max (omp_p,0), l2 (omp_p,0);
  204. std::vector<Real_t> max_err(omp_p,0), l2_err(omp_p,0);
  205. #pragma omp parallel for
  206. for(size_t tid=0;tid<omp_p;tid++){
  207. pvfmm::Vector<Real_t> out; out.SetZero();
  208. pvfmm::Vector<Real_t> fn_out(dof*fn_dof);
  209. for(size_t i=(nodes.size()*tid)/omp_p;i<(nodes.size()*(tid+1))/omp_p;i++){
  210. pvfmm::Vector<Real_t>& cheb_coeff=nodes[i]->ChebData();
  211. cheb_eval(cheb_coeff, cheb_deg, cheb_nds, cheb_nds, cheb_nds, out);
  212. Real_t* c=nodes[i]->Coord();
  213. Real_t s=pow(2.0,-nodes[i]->Depth());
  214. Real_t s3=s*s*s;
  215. for(size_t j=0;j<n_pts;j++){
  216. Real_t coord[3]={c[0]+s*cheb_pts[j*COORD_DIM+0],
  217. c[1]+s*cheb_pts[j*COORD_DIM+1],
  218. c[2]+s*cheb_pts[j*COORD_DIM+2]};
  219. fn_out.SetZero();
  220. for(size_t k=0;k<dof;k++)
  221. fn_poten(coord,1,&fn_out[k*fn_dof]);
  222. for(size_t k=0;k<dof*fn_dof;k++){
  223. Real_t err=out[n_pts*k+j]-fn_out[k]-glb_err_avg[k];
  224. if(fabs(fn_out[k])>max [tid]) max [tid]=fabs(fn_out[k]);
  225. if(fabs(err )>max_err[tid]) max_err[tid]=fabs(err );
  226. l2[tid]+=fn_out[k]*fn_out[k]*s3;
  227. l2_err[tid]+=err*err*s3;
  228. }
  229. }
  230. }
  231. l2 [tid]/=n_pts;
  232. l2_err[tid]/=n_pts;
  233. }
  234. for(size_t tid=1;tid<omp_p;tid++){
  235. if(max [tid]>max [0]) max [0]=max [tid];
  236. if(max_err[tid]>max_err[0]) max_err[0]=max_err[tid];
  237. l2 [0]+=l2 [tid];
  238. l2_err[0]+=l2_err[tid];
  239. }
  240. Real_t global_l2, global_l2_err;
  241. Real_t global_max, global_max_err;
  242. MPI_Reduce(&l2 [0], &global_l2 , 1, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::sum(), 0, c1);
  243. MPI_Reduce(&l2_err [0], &global_l2_err , 1, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::sum(), 0, c1);
  244. MPI_Reduce(&max [0], &global_max , 1, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::max(), 0, c1);
  245. MPI_Reduce(&max_err[0], &global_max_err, 1, pvfmm::par::Mpi_datatype<Real_t>::value(), pvfmm::par::Mpi_datatype<Real_t>::max(), 0, c1);
  246. if(!myrank){
  247. std::cout<<"Absolute L2 Error ["<<t_name<<"] : "<<std::scientific<<sqrt(global_l2_err)<<'\n';
  248. std::cout<<"Relative L2 Error ["<<t_name<<"] : "<<std::scientific<<sqrt(global_l2_err/global_l2)<<'\n';
  249. std::cout<<"Maximum Absolute Error ["<<t_name<<"]: "<<std::scientific<<global_max_err<<'\n';
  250. std::cout<<"Maximum Relative Error ["<<t_name<<"]: "<<std::scientific<<global_max_err/global_max<<'\n';
  251. }
  252. pvfmm::Profile::Toc();
  253. }
  254. template <class Real_t>
  255. std::vector<Real_t> point_distrib(DistribType dist_type, size_t N, MPI_Comm comm){
  256. int np, myrank;
  257. MPI_Comm_size(comm, &np);
  258. MPI_Comm_rank(comm, &myrank);
  259. std::vector<Real_t> coord;
  260. srand48(myrank+1);
  261. switch(dist_type){
  262. case UnifGrid:
  263. {
  264. size_t NN=(size_t)round(pow((double)N,1.0/3.0));
  265. size_t N_total=NN*NN*NN;
  266. size_t start= myrank *N_total/np;
  267. size_t end =(myrank+1)*N_total/np;
  268. for(size_t i=start;i<end;i++){
  269. coord.push_back(((Real_t)((i/ 1 )%NN)+0.5)/NN);
  270. coord.push_back(((Real_t)((i/ NN )%NN)+0.5)/NN);
  271. coord.push_back(((Real_t)((i/(NN*NN))%NN)+0.5)/NN);
  272. }
  273. }
  274. break;
  275. case RandUnif:
  276. {
  277. size_t N_total=N;
  278. size_t start= myrank *N_total/np;
  279. size_t end =(myrank+1)*N_total/np;
  280. size_t N_local=end-start;
  281. coord.resize(N_local*3);
  282. for(size_t i=0;i<N_local*3;i++)
  283. coord[i]=((Real_t)drand48());
  284. }
  285. break;
  286. case RandGaus:
  287. {
  288. Real_t e=2.7182818284590452;
  289. Real_t log_e=log(e);
  290. size_t N_total=N;
  291. size_t start= myrank *N_total/np;
  292. size_t end =(myrank+1)*N_total/np;
  293. for(size_t i=start*3;i<end*3;i++){
  294. Real_t y=-1;
  295. while(y<=0 || y>=1){
  296. Real_t r1=sqrt(-2*log(drand48())/log_e)*cos(2*M_PI*drand48());
  297. Real_t r2=pow(0.5,i*10/N_total);
  298. y=0.5+r1*r2;
  299. }
  300. coord.push_back(y);
  301. }
  302. }
  303. break;
  304. case RandElps:
  305. {
  306. size_t N_total=N;
  307. size_t start= myrank *N_total/np;
  308. size_t end =(myrank+1)*N_total/np;
  309. size_t N_local=end-start;
  310. coord.resize(N_local*3);
  311. const Real_t r=0.25;
  312. const Real_t center[3]={0.5,0.5,0.5};
  313. for(size_t i=0;i<N_local;i++){
  314. Real_t* y=&coord[i*3];
  315. Real_t phi=2*M_PI*drand48();
  316. Real_t theta=M_PI*drand48();
  317. y[0]=center[0]+0.25*r*sin(theta)*cos(phi);
  318. y[1]=center[1]+0.25*r*sin(theta)*sin(phi);
  319. y[2]=center[2]+r*cos(theta);
  320. }
  321. }
  322. break;
  323. case RandSphr:
  324. {
  325. size_t N_total=N;
  326. size_t start= myrank *N_total/np;
  327. size_t end =(myrank+1)*N_total/np;
  328. size_t N_local=end-start;
  329. coord.resize(N_local*3);
  330. const Real_t center[3]={0.5,0.5,0.5};
  331. for(size_t i=0;i<N_local;i++){
  332. Real_t* y=&coord[i*3];
  333. Real_t r=1;
  334. while(r>0.5 || r==0){
  335. y[0]=drand48(); y[1]=drand48(); y[2]=drand48();
  336. r=sqrt((y[0]-center[0])*(y[0]-center[0])
  337. +(y[1]-center[1])*(y[1]-center[1])
  338. +(y[2]-center[2])*(y[2]-center[2]));
  339. y[0]=center[0]+0.1*(y[0]-center[0])/r;
  340. y[1]=center[1]+0.1*(y[1]-center[1])/r;
  341. y[2]=center[2]+0.1*(y[2]-center[2])/r;
  342. }
  343. }
  344. }
  345. break;
  346. default:
  347. break;
  348. }
  349. return coord;
  350. }
  351. void commandline_option_start(int argc, char** argv, const char* help_text){
  352. char help[]="--help";
  353. for(int i=0;i<argc;i++){
  354. if(!strcmp(argv[i],help)){
  355. if(help_text!=NULL) std::cout<<help_text<<'\n';
  356. std::cout<<"Usage:\n";
  357. std::cout<<" "<<argv[0]<<" [options]\n\n";
  358. }
  359. }
  360. }
  361. const char* commandline_option(int argc, char** argv, const char* opt, const char* def_val, bool required, const char* err_msg){
  362. char help[]="--help";
  363. for(int i=0;i<argc;i++){
  364. if(!strcmp(argv[i],help)){
  365. std::cout<<" "<<err_msg<<'\n';
  366. return def_val;
  367. }
  368. }
  369. for(int i=0;i<argc;i++){
  370. if(!strcmp(argv[i],opt)){
  371. return argv[(i+1)%argc];
  372. }
  373. }
  374. if(required){
  375. std::cout<<"Missing: required option\n"<<" "<<err_msg<<"\n\n";
  376. std::cout<<"To see usage options\n"<<" "<<argv[0]<<" --help\n\n";
  377. exit(0);
  378. }
  379. return def_val;
  380. }
  381. void commandline_option_end(int argc, char** argv){
  382. char help[]="--help";
  383. for(int i=0;i<argc;i++){
  384. if(!strcmp(argv[i],help)){
  385. std::cout<<"\n";
  386. exit(0);
  387. }
  388. }
  389. }