pvfmm.hpp 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220
  1. /**
  2. * \file pvfmm.hpp
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 1-2-2014
  5. * \brief This file contains wrapper functions for PvFMM.
  6. */
  7. #include <mpi.h>
  8. #include <vector>
  9. #include <cmath>
  10. #include <pvfmm_common.hpp>
  11. #include <cheb_node.hpp>
  12. #include <mpi_node.hpp>
  13. #include <fmm_tree.hpp>
  14. #include <fmm_node.hpp>
  15. #include <fmm_cheb.hpp>
  16. #include <fmm_pts.hpp>
  17. #include <vector.hpp>
  18. #include <parUtils.h>
  19. #ifndef _PVFMM_HPP_
  20. #define _PVFMM_HPP_
  21. namespace pvfmm{
  22. typedef FMM_Node<Cheb_Node<double> > ChebFMM_Node;
  23. typedef FMM_Cheb<ChebFMM_Node> ChebFMM;
  24. typedef FMM_Tree<ChebFMM> ChebFMM_Tree;
  25. typedef ChebFMM_Node::NodeData ChebFMM_Data;
  26. typedef void (*ChebFn)(double* , int , double*);
  27. ChebFMM_Tree* ChebFMM_CreateTree(int cheb_deg, int data_dim, ChebFn fn_ptr, std::vector<double>& trg_coord, MPI_Comm& comm,
  28. double tol=1e-6, int max_pts=100, BoundaryType bndry=FreeSpace, int init_depth=0){
  29. int np, myrank;
  30. MPI_Comm_size(comm, &np);
  31. MPI_Comm_rank(comm, &myrank);
  32. ChebFMM_Data tree_data;
  33. tree_data.cheb_deg=cheb_deg;
  34. tree_data.data_dof=data_dim;
  35. tree_data.input_fn=fn_ptr;
  36. tree_data.tol=tol;
  37. bool adap=true;
  38. tree_data.dim=COORD_DIM;
  39. tree_data.max_depth=MAX_DEPTH;
  40. tree_data.max_pts=max_pts;
  41. { // Set points for initial tree.
  42. std::vector<double> coord;
  43. size_t N=pow(8.0,init_depth);
  44. N=(N<np?np:N)*max_pts;
  45. size_t NN=ceil(pow((double)N,1.0/3.0));
  46. size_t N_total=NN*NN*NN;
  47. size_t start= myrank *N_total/np;
  48. size_t end =(myrank+1)*N_total/np;
  49. for(size_t i=start;i<end;i++){
  50. coord.push_back(((double)((i/ 1 )%NN)+0.5)/NN);
  51. coord.push_back(((double)((i/ NN )%NN)+0.5)/NN);
  52. coord.push_back(((double)((i/(NN*NN))%NN)+0.5)/NN);
  53. }
  54. tree_data.pt_coord=coord;
  55. }
  56. // Set target points.
  57. tree_data.trg_coord=trg_coord;
  58. ChebFMM_Tree* tree=new ChebFMM_Tree(comm);
  59. tree->Initialize(&tree_data);
  60. tree->InitFMM_Tree(adap,bndry);
  61. return tree;
  62. }
  63. void ChebFMM_Evaluate(ChebFMM_Tree* tree, std::vector<double>& trg_val, size_t loc_size=0){
  64. tree->RunFMM();
  65. Vector<double> trg_value;
  66. Vector<size_t> trg_scatter;
  67. {// Collect data from each node to trg_value and trg_scatter.
  68. std::vector<double> trg_value_;
  69. std::vector<size_t> trg_scatter_;
  70. std::vector<ChebFMM_Node*>& nodes=tree->GetNodeList();
  71. for(size_t i=0;i<nodes.size();i++){
  72. if(nodes[i]->IsLeaf() && !nodes[i]->IsGhost()){
  73. Vector<double>& trg_value=nodes[i]->trg_value;
  74. Vector<size_t>& trg_scatter=nodes[i]->trg_scatter;
  75. for(size_t j=0;j<trg_value.Dim();j++) trg_value_.push_back(trg_value[j]);
  76. for(size_t j=0;j<trg_scatter.Dim();j++) trg_scatter_.push_back(trg_scatter[j]);
  77. }
  78. }
  79. trg_value=trg_value_;
  80. trg_scatter=trg_scatter_;
  81. }
  82. par::ScatterReverse(trg_value,trg_scatter,*tree->Comm(),loc_size);
  83. trg_val.assign(&trg_value[0],&trg_value[0]+trg_value.Dim());;
  84. }
  85. typedef FMM_Node<MPI_Node<double> > PtFMM_Node;
  86. typedef FMM_Pts<PtFMM_Node> PtFMM;
  87. typedef FMM_Tree<PtFMM> PtFMM_Tree;
  88. typedef PtFMM_Node::NodeData PtFMM_Data;
  89. PtFMM_Tree* PtFMM_CreateTree(std::vector<double>& src_coord, std::vector<double>& src_value,
  90. std::vector<double>& surf_coord, std::vector<double>& surf_value,
  91. std::vector<double>& trg_coord, MPI_Comm& comm, int max_pts=100,
  92. BoundaryType bndry=FreeSpace, int init_depth=0){
  93. int np, myrank;
  94. MPI_Comm_size(comm, &np);
  95. MPI_Comm_rank(comm, &myrank);
  96. PtFMM_Data tree_data;
  97. bool adap=true;
  98. tree_data.dim=COORD_DIM;
  99. tree_data.max_depth=MAX_DEPTH;
  100. tree_data.max_pts=max_pts;
  101. // Set source points.
  102. tree_data. src_coord= src_coord;
  103. tree_data. src_value= src_value;
  104. tree_data.surf_coord=surf_coord;
  105. tree_data.surf_value=surf_value;
  106. // Set target points.
  107. tree_data.trg_coord=trg_coord;
  108. tree_data. pt_coord=trg_coord;
  109. PtFMM_Tree* tree=new PtFMM_Tree(comm);
  110. tree->Initialize(&tree_data);
  111. tree->InitFMM_Tree(adap,bndry);
  112. return tree;
  113. }
  114. PtFMM_Tree* PtFMM_CreateTree(std::vector<double>& src_coord, std::vector<double>& src_value,
  115. std::vector<double>& trg_coord, MPI_Comm& comm, int max_pts=100,
  116. BoundaryType bndry=FreeSpace, int init_depth=0){
  117. std::vector<double> surf_coord;
  118. std::vector<double> surf_value;
  119. return PtFMM_CreateTree(src_coord, src_value, surf_coord,surf_value, trg_coord, comm, max_pts, bndry, init_depth);
  120. }
  121. void PtFMM_Evaluate(PtFMM_Tree* tree, std::vector<double>& trg_val, size_t loc_size=0, std::vector<double>* src_val=NULL, std::vector<double>* surf_val=NULL){
  122. if(src_val){
  123. std::vector<size_t> src_scatter_;
  124. std::vector<PtFMM_Node*>& nodes=tree->GetNodeList();
  125. for(size_t i=0;i<nodes.size();i++){
  126. if(nodes[i]->IsLeaf() && !nodes[i]->IsGhost()){
  127. Vector<size_t>& src_scatter=nodes[i]->src_scatter;
  128. for(size_t j=0;j<src_scatter.Dim();j++) src_scatter_.push_back(src_scatter[j]);
  129. }
  130. }
  131. Vector<double> src_value=*src_val;
  132. Vector<size_t> src_scatter=src_scatter_;
  133. par::ScatterForward(src_value,src_scatter,*tree->Comm());
  134. size_t indx=0;
  135. for(size_t i=0;i<nodes.size();i++){
  136. if(nodes[i]->IsLeaf() && !nodes[i]->IsGhost()){
  137. Vector<double>& src_value_=nodes[i]->src_value;
  138. for(size_t j=0;j<src_value_.Dim();j++){
  139. src_value_[j]=src_value[indx];
  140. indx++;
  141. }
  142. }
  143. }
  144. }
  145. if(surf_val){
  146. std::vector<size_t> surf_scatter_;
  147. std::vector<PtFMM_Node*>& nodes=tree->GetNodeList();
  148. for(size_t i=0;i<nodes.size();i++){
  149. if(nodes[i]->IsLeaf() && !nodes[i]->IsGhost()){
  150. Vector<size_t>& surf_scatter=nodes[i]->surf_scatter;
  151. for(size_t j=0;j<surf_scatter.Dim();j++) surf_scatter_.push_back(surf_scatter[j]);
  152. }
  153. }
  154. Vector<double> surf_value=*surf_val;
  155. Vector<size_t> surf_scatter=surf_scatter_;
  156. par::ScatterForward(surf_value,surf_scatter,*tree->Comm());
  157. size_t indx=0;
  158. for(size_t i=0;i<nodes.size();i++){
  159. if(nodes[i]->IsLeaf() && !nodes[i]->IsGhost()){
  160. Vector<double>& surf_value_=nodes[i]->surf_value;
  161. for(size_t j=0;j<surf_value_.Dim();j++){
  162. surf_value_[j]=surf_value[indx];
  163. indx++;
  164. }
  165. }
  166. }
  167. }
  168. tree->RunFMM();
  169. Vector<double> trg_value;
  170. Vector<size_t> trg_scatter;
  171. {
  172. std::vector<double> trg_value_;
  173. std::vector<size_t> trg_scatter_;
  174. std::vector<PtFMM_Node*>& nodes=tree->GetNodeList();
  175. for(size_t i=0;i<nodes.size();i++){
  176. if(nodes[i]->IsLeaf() && !nodes[i]->IsGhost()){
  177. Vector<double>& trg_value=nodes[i]->trg_value;
  178. Vector<size_t>& trg_scatter=nodes[i]->trg_scatter;
  179. for(size_t j=0;j<trg_value.Dim();j++) trg_value_.push_back(trg_value[j]);
  180. for(size_t j=0;j<trg_scatter.Dim();j++) trg_scatter_.push_back(trg_scatter[j]);
  181. }
  182. }
  183. trg_value=trg_value_;
  184. trg_scatter=trg_scatter_;
  185. }
  186. par::ScatterReverse(trg_value,trg_scatter,*tree->Comm(),loc_size);
  187. trg_val.assign(&trg_value[0],&trg_value[0]+trg_value.Dim());;
  188. }
  189. }//end namespace
  190. #endif //_PVFMM_HPP_