OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
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 #ifdef OPAL_DKS
36 #include "DKSOPAL.h"
37 #endif
38 
39 class PartBunch;
40 
42 
44 public:
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 #ifdef OPAL_DKS
53  DKSOPAL dksbase;
54 #endif
55 
56  // given a charge-density field rho and a set of mesh spacings hr,
57  // compute the scalar potential with image charges at -z
58  void computePotential(Field_t &rho, Vector_t hr, double zshift);
59 
60  // given a charge-density field rho and a set of mesh spacings hr,
61  // compute the scalar potential in open space
62  void computePotential(Field_t &rho, Vector_t hr);
63 
64  void computePotentialDKS(Field_t &rho);
65  // compute the green's function for a Poisson problem and put it in in grntm_m
66  // uses grnIField_m to eliminate excess calculation in greenFunction()
67  // given mesh information in nr and hr
68  void greensFunction();
69 
72 
74  //compute the integrated Green function as described in <A HREF="http://prst-ab.aps.org/abstract/PRSTAB/v9/i4/e044204">Three-dimensional quasistatic model for high brightness beam dynamics simulation</A> by Qiang et al.
76 
78  void shiftedIntGreensFunction(double zshift);
79 
80  double getXRangeMin(unsigned short level) {return 1.0;}
81  double getXRangeMax(unsigned short level) {return 1.0;}
82  double getYRangeMin(unsigned short level) {return 1.0;}
83  double getYRangeMax(unsigned short level) {return 1.0;}
84  double getZRangeMin(unsigned short level) {return 1.0;}
85  double getZRangeMax(unsigned short level) {return 1.0;}
86  void test(PartBunchBase<double, 3> *bunch) { }
87 
88  Inform &print(Inform &os) const;
89 private:
90  void initializeFields();
92  void mirrorRhoField(Field_t & ggrn2);// FIELDASSIGNOPTIMIZATION;
93 
94  // rho2_m is the charge-density field with mesh doubled in each dimension
96 
97  // real field with layout of complex field: domain3_m
99 
100  // rho2tr_m is the Fourier transformed charge-density field
101  // domain3_m and mesh3_ are used
104 
105  // grntr_m is the Fourier transformed Green's function
106  // domain3_m and mesh3_ are used
108 
109 #ifdef OPAL_DKS
110  //pointer for Fourier transformed Green's function on GPU
111  void * grntr_m_ptr;
112  void * rho2_m_ptr;
113  void * tmpgreen_ptr;
114  void *rho2real_m_ptr;
115  void *rho2tr_m_ptr;
116 
117  //stream id for calculating greens function
118  int streamGreens;
119  int streamFFT;
120 
121 #endif
122 
123  // Fields used to eliminate excess calculation in greensFunction()
124  // mesh2_m and layout2_m are used
126 
127  // the FFT object
128  std::unique_ptr<FFT_t> fft_m;
129 
130  // mesh and layout objects for rho_m
133 
134  // mesh and layout objects for rho2_m
135  std::unique_ptr<Mesh_t> mesh2_m;
136  std::unique_ptr<FieldLayout_t> layout2_m;
137 
138  //
139  std::unique_ptr<Mesh_t> mesh3_m;
140  std::unique_ptr<FieldLayout_t> layout3_m;
141 
142  // mesh and layout for integrated greens function
143  std::unique_ptr<Mesh_t> mesh4_m;
144  std::unique_ptr<FieldLayout_t> layout4_m;
145 
146  // tmp
148 
149  // domains for the various fields
150  NDIndex<3> domain_m; // original domain, gridsize
151  // mesh and gridsize defined outside of FFT class, given as
152  // parameter to the constructor (mesh and layout object).
153  NDIndex<3> domain2_m; // doubled gridsize (2*Nx,2*Ny,2*Nz)
154  NDIndex<3> domain3_m; // field for the complex values of the RC transformation
156  // (2*Nx,Ny,2*Nz)
158 
159  // mesh spacing and size values
162 
166 
167  bool bcz_m;
170  /*
171  IpplTimings::TimerRef IntGreensFunctionTimer1_m;
172  IpplTimings::TimerRef IntGreensFunctionTimer2_m;
173  IpplTimings::TimerRef IntGreensFunctionTimer3_m;
174  IpplTimings::TimerRef IntGreensMirrorTimer1_m;
175 
176  IpplTimings::TimerRef ShIntGreensFunctionTimer1_m;
177  IpplTimings::TimerRef ShIntGreensFunctionTimer2_m;
178  IpplTimings::TimerRef ShIntGreensFunctionTimer3_m;
179  IpplTimings::TimerRef ShIntGreensFunctionTimer4_m;
180  IpplTimings::TimerRef IntGreensMirrorTimer2_m;
181 
182  IpplTimings::TimerRef GreensFunctionTimer1_m;
183  IpplTimings::TimerRef GreensFunctionTimer4_m;
184  */
186 };
187 
189  return fs.print(os);
190 }
191 
192 
193 
194 #endif
195 
196 /***************************************************************************
197  * $RCSfile: FFTPoissonSolver.hh,v $ $Author: adelmann $
198  * $Revision: 1.1.1.1 $ $Date: 2001/08/08 11:21:48 $
199  ***************************************************************************/
std::ostream & operator<<(std::ostream &os, const Attribute &attr)
Definition: Attribute.cpp:167
std::unique_ptr< Mesh_t > mesh2_m
double getXRangeMin(unsigned short level)
void shiftedIntGreensFunction(double zshift)
compute the shifted integrated Green function as described in Three-dimensional quasistatic model for...
NDIndex< 3 > domain2_m
std::unique_ptr< FieldLayout_t > layout4_m
std::unique_ptr< Mesh_t > mesh3_m
double getZRangeMax(unsigned short level)
FRONT * fs
Definition: hypervolume.cpp:59
Inform & print(Inform &os) const
BConds< Vector_t, 3, Mesh_t, Center_t > vbc_m
IpplTimings::TimerRef GreensFunctionTimer_m
double getXRangeMax(unsigned short level)
IpplTimings::TimerRef ComputePotential_m
NDIndex< 3 > domainFFTConstruct_m
Vektor< int, 3 > nr_m
Definition: BCond.h:34
void integratedGreensFunction()
compute the integrated Green function as described in Three-dimensional quasistatic model for high br...
double getZRangeMin(unsigned short level)
std::unique_ptr< Mesh_t > mesh4_m
void integratedGreensFunctionDKS()
Uses DKS to offload the computation of Greens function on the GPU.
FFTPoissonSolver(PartBunch &bunch, std::string greensFuntion)
double getYRangeMax(unsigned short level)
#define FIELDASSIGNOPTIMIZATION
BConds< double, 3, Mesh_t, Center_t > bc_m
for defining the boundary conditions
double getYRangeMin(unsigned short level)
IField_t grnIField_m[3]
NDIndex< 3 > domain3_m
std::unique_ptr< FFT_t > fft_m
Particle Bunch.
Definition: PartBunch.h:30
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176
FieldLayout_t * layout_m
NDIndex< 3 > domain4_m
void computePotentialDKS(Field_t &rho)
void computePotential(Field_t &rho, Vector_t hr, double zshift)
Definition: Inform.h:41
std::unique_ptr< FieldLayout_t > layout3_m
std::unique_ptr< FieldLayout_t > layout2_m
NDIndex< 3 > domain_m
void test(PartBunchBase< double, 3 > *bunch)