OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
FFTPoissonSolver.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 /***************************************************************************
3  *
4  *
5  * FFTPoissonSolver.hh
6  *
7  * Open BC in x,y and z.
8  *
9  *
10  *
11  *
12  *
13  *
14  ***************************************************************************/
15 
17 // This class contains methods for solving Poisson's equation for the
18 // space charge portion of the calculation.
20 
21 #ifndef FFT_POISSON_SOLVER_H_
22 #define FFT_POISSON_SOLVER_H_
23 
24 
25 #ifdef dontOPTIMIZE_FIELD_ASSIGNMENT
26 #define FIELDASSIGNOPTIMIZATION __attribute__((optimize(0)))
27 #else
28 #define FIELDASSIGNOPTIMIZATION
29 #endif
30 
31 #include <memory>
33 #include "PoissonSolver.h"
34 
35 #include "FFT/FFT.h"
36 
37 class PartBunch;
38 
40 
42 public:
44 
45  // constructor and destructor
46  FFTPoissonSolver(PartBunch &bunch, std::string greensFuntion);
47 
48  FFTPoissonSolver(Mesh_t *mesh, FieldLayout_t *fl, std::string greensFunction, std::string bcz);
49 
51 
52  // given a charge-density field rho and a set of mesh spacings hr,
53  // compute the scalar potential with image charges at -z
54  void computePotential(Field_t &rho, Vector_t hr, double zshift);
55 
56  // given a charge-density field rho and a set of mesh spacings hr,
57  // compute the scalar potential in open space
58  void computePotential(Field_t &rho, Vector_t hr);
59 
60  // compute the green's function for a Poisson problem and put it in in grntm_m
61  // uses grnIField_m to eliminate excess calculation in greenFunction()
62  // given mesh information in nr and hr
63  void greensFunction();
64 
67 
69  void shiftedIntGreensFunction(double zshift);
70 
71  double getXRangeMin(unsigned short /*level*/) {return 1.0;}
72  double getXRangeMax(unsigned short /*level*/) {return 1.0;}
73  double getYRangeMin(unsigned short /*level*/) {return 1.0;}
74  double getYRangeMax(unsigned short /*level*/) {return 1.0;}
75  double getZRangeMin(unsigned short /*level*/) {return 1.0;}
76  double getZRangeMax(unsigned short /*level*/) {return 1.0;}
77  void test(PartBunchBase<double, 3> */*bunch*/) { }
78 
79  Inform &print(Inform &os) const;
80 private:
81  void initializeFields();
83  void mirrorRhoField(Field_t & ggrn2);// FIELDASSIGNOPTIMIZATION;
84 
85  // rho2_m is the charge-density field with mesh doubled in each dimension
87 
88  // real field with layout of complex field: domain3_m
90 
91  // rho2tr_m is the Fourier transformed charge-density field
92  // domain3_m and mesh3_ are used
95 
96  // grntr_m is the Fourier transformed Green's function
97  // domain3_m and mesh3_ are used
99 
100  // Fields used to eliminate excess calculation in greensFunction()
101  // mesh2_m and layout2_m are used
103 
104  // the FFT object
105  std::unique_ptr<FFT_t> fft_m;
106 
107  // mesh and layout objects for rho_m
110 
111  // mesh and layout objects for rho2_m
112  std::unique_ptr<Mesh_t> mesh2_m;
114 
115  //
116  std::unique_ptr<Mesh_t> mesh3_m;
118 
119  // mesh and layout for integrated greens function
120  std::unique_ptr<Mesh_t> mesh4_m;
122 
123  // tmp
125 
126  // domains for the various fields
127  NDIndex<3> domain_m; // original domain, gridsize
128  // mesh and gridsize defined outside of FFT class, given as
129  // parameter to the constructor (mesh and layout object).
130  NDIndex<3> domain2_m; // doubled gridsize (2*Nx,2*Ny,2*Nz)
131  NDIndex<3> domain3_m; // field for the complex values of the RC transformation
133  // (2*Nx,Ny,2*Nz)
135 
136  // mesh spacing and size values
138  Vektor<int, 3> nr_m;
139 
141  BConds<double, 3, Mesh_t, Center_t> bc_m;
143 
144  bool bcz_m;
147  /*
148  IpplTimings::TimerRef IntGreensFunctionTimer1_m;
149  IpplTimings::TimerRef IntGreensFunctionTimer2_m;
150  IpplTimings::TimerRef IntGreensFunctionTimer3_m;
151  IpplTimings::TimerRef IntGreensMirrorTimer1_m;
152 
153  IpplTimings::TimerRef ShIntGreensFunctionTimer1_m;
154  IpplTimings::TimerRef ShIntGreensFunctionTimer2_m;
155  IpplTimings::TimerRef ShIntGreensFunctionTimer3_m;
156  IpplTimings::TimerRef ShIntGreensFunctionTimer4_m;
157  IpplTimings::TimerRef IntGreensMirrorTimer2_m;
158 
159  IpplTimings::TimerRef GreensFunctionTimer1_m;
160  IpplTimings::TimerRef GreensFunctionTimer4_m;
161  */
163 };
164 
165 inline Inform &operator<<(Inform &os, const FFTPoissonSolver &fs) {
166  return fs.print(os);
167 }
168 
169 
170 
171 #endif
#define FIELDASSIGNOPTIMIZATION
FRONT * fs
Definition: hypervolume.cpp:59
NDIndex< 3 > domainFFTConstruct_m
FieldLayout_t * layout_m
std::unique_ptr< Mesh_t > mesh4_m
void shiftedIntGreensFunction(double zshift)
compute the shifted integrated Green function as described in Three-dimensional quasistatic model for...
std::unique_ptr< FieldLayout_t > layout3_m
void computePotential(Field_t &rho, Vector_t hr, double zshift)
FFT< RCTransform, 3, double > FFT_t
IpplTimings::TimerRef GreensFunctionTimer_m
NDIndex< 3 > domain2_m
double getZRangeMax(unsigned short)
std::unique_ptr< FFT_t > fft_m
NDIndex< 3 > domain_m
FFTPoissonSolver(PartBunch &bunch, std::string greensFuntion)
std::unique_ptr< FieldLayout_t > layout2_m
BConds< Vector_t, 3, Mesh_t, Center_t > vbc_m
IField_t grnIField_m[3]
double getZRangeMin(unsigned short)
NDIndex< 3 > domain3_m
std::unique_ptr< Mesh_t > mesh3_m
double getYRangeMax(unsigned short)
double getXRangeMax(unsigned short)
double getYRangeMin(unsigned short)
Inform & print(Inform &os) const
NDIndex< 3 > domain4_m
std::unique_ptr< Mesh_t > mesh2_m
Vektor< int, 3 > nr_m
std::unique_ptr< FieldLayout_t > layout4_m
BConds< double, 3, Mesh_t, Center_t > bc_m
for defining the boundary conditions
IpplTimings::TimerRef ComputePotential_m
void integratedGreensFunction()
compute the integrated Green function as described in Three-dimensional quasistatic model for high br...
double getXRangeMin(unsigned short)
void test(PartBunchBase< double, 3 > *)
Definition: FFT.h:48
Definition: BCond.h:199
Definition: Centering.h:32
Definition: Inform.h:42