example2.cpp 3.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102
  1. #include <mpi.h>
  2. #include <omp.h>
  3. #include <iostream>
  4. #include <pvfmm.hpp>
  5. #include <utils.hpp>
  6. //Input function
  7. void fn_input(double* coord, int n, double* out){
  8. double a=-160;
  9. for(int i=0;i<n;i++){
  10. double* c=&coord[i*3];
  11. double r_2=(c[0]-0.5)*(c[0]-0.5)+(c[1]-0.5)*(c[1]-0.5)+(c[2]-0.5)*(c[2]-0.5);
  12. out[i]=(2*a*r_2+3)*2*a*exp(a*r_2);
  13. }
  14. }
  15. //Analytical solution (Expected output)
  16. void fn_output(double* coord, int n, double* out){
  17. double a=-160;
  18. for(int i=0;i<n;i++){
  19. double* c=&coord[i*3];
  20. double r_2=(c[0]-0.5)*(c[0]-0.5)+(c[1]-0.5)*(c[1]-0.5)+(c[2]-0.5)*(c[2]-0.5);
  21. out[i]=-exp(a*r_2);
  22. }
  23. }
  24. void fmm_test(size_t N, int mult_order, int cheb_deg, double tol, MPI_Comm comm){
  25. // Set kernel.
  26. const pvfmm::Kernel<double>& kernel_fn=pvfmm::laplace_potn_d;
  27. // Construct tree.
  28. size_t max_pts=100;
  29. std::vector<double> trg_coord=point_distrib<double>(RandUnif,N,comm);
  30. pvfmm::ChebFMM_Tree* tree=ChebFMM_CreateTree(cheb_deg, kernel_fn.ker_dim[0], fn_input,
  31. trg_coord, comm, tol, max_pts, pvfmm::FreeSpace);
  32. // Load matrices.
  33. pvfmm::ChebFMM matrices;
  34. matrices.Initialize(mult_order, cheb_deg, comm, &kernel_fn);
  35. // FMM Setup
  36. tree->SetupFMM(&matrices);
  37. // Run FMM
  38. std::vector<double> trg_value;
  39. size_t n_trg=trg_coord.size()/COORD_DIM;
  40. pvfmm::ChebFMM_Evaluate(tree, trg_value, n_trg);
  41. // Re-run FMM
  42. tree->ClearFMMData();
  43. pvfmm::ChebFMM_Evaluate(tree, trg_value, n_trg);
  44. {// Check error
  45. std::vector<double> trg_value_(n_trg*kernel_fn.ker_dim[1]);
  46. fn_output(&trg_coord[0],n_trg,&trg_value_[0]);
  47. double max_err=0;
  48. double max_err_glb=0;
  49. for(size_t i=0;i<n_trg;i++){
  50. if(fabs(trg_value_[i]-trg_value[i])>max_err)
  51. max_err=fabs(trg_value_[i]-trg_value[i]);
  52. }
  53. MPI_Reduce(&max_err, &max_err_glb, 1, MPI_DOUBLE, MPI_MAX, 0, comm);
  54. int rank;
  55. MPI_Comm_rank(comm, &rank);
  56. if(!rank) std::cout<<"Maximum Error:"<<max_err_glb<<'\n';
  57. }
  58. // Free memory
  59. delete tree;
  60. }
  61. int main(int argc, char **argv){
  62. MPI_Init(&argc, &argv);
  63. MPI_Comm comm=MPI_COMM_WORLD;
  64. // Read command line options.
  65. commandline_option_start(argc, argv, "\
  66. This example demonstrates solving a volume potential problem,\n\
  67. with Laplace kernel, using the PvFMM library.\n");
  68. commandline_option_start(argc, argv);
  69. omp_set_num_threads( atoi(commandline_option(argc, argv, "-omp", "1", false, "-omp <int> = (1) : Number of OpenMP threads." )));
  70. size_t N=(size_t)strtod(commandline_option(argc, argv, "-N", "1", true, "-N <int> : Number of target points." ),NULL);
  71. int m= strtoul(commandline_option(argc, argv, "-m", "10", false, "-m <int> = (10) : Multipole order (+ve even integer)."),NULL,10);
  72. int q= strtoul(commandline_option(argc, argv, "-q", "14", false, "-q <int> = (14) : Chebyshev order (+ve integer)." ),NULL,10);
  73. double tol= strtod(commandline_option(argc, argv, "-tol", "1e-5", false, "-tol <real> = (1e-5) : Tolerance for adaptive refinement." ),NULL);
  74. commandline_option_end(argc, argv);
  75. pvfmm::Profile::Enable(true);
  76. // Run FMM with above options.
  77. fmm_test(N, m,q, tol, comm);
  78. //Output Profiling results.
  79. pvfmm::Profile::print(&comm);
  80. // Shut down MPI
  81. MPI_Finalize();
  82. return 0;
  83. }