OPAL (Object Oriented Parallel Accelerator Library)
2021.1.99
OPAL
|
#include <FFTPoissonSolver.h>
Public Types | |
typedef FFT< RCTransform, 3, double > | FFT_t |
Public Member Functions | |
FFTPoissonSolver (PartBunch &bunch, std::string greensFuntion) | |
FFTPoissonSolver (Mesh_t *mesh, FieldLayout_t *fl, std::string greensFunction, std::string bcz) | |
~FFTPoissonSolver () | |
void | computePotential (Field_t &rho, Vector_t hr, double zshift) |
void | computePotential (Field_t &rho, Vector_t hr) |
void | greensFunction () |
void | integratedGreensFunction () |
compute the integrated Green function as described in Three-dimensional quasistatic model for high brightness beam dynamics simulation by Qiang et al. More... | |
void | shiftedIntGreensFunction (double zshift) |
compute the shifted integrated Green function as described in Three-dimensional quasistatic model for high brightness beam dynamics simulation by Qiang et al. More... | |
double | getXRangeMin (unsigned short) |
double | getXRangeMax (unsigned short) |
double | getYRangeMin (unsigned short) |
double | getYRangeMax (unsigned short) |
double | getZRangeMin (unsigned short) |
double | getZRangeMax (unsigned short) |
void | test (PartBunchBase< double, 3 > *) |
Inform & | print (Inform &os) const |
![]() | |
virtual void | solve (AmrScalarFieldContainer_t &, AmrScalarFieldContainer_t &, AmrVectorFieldContainer_t &, unsigned short, unsigned short, bool=true) |
virtual void | hasToRegrid () |
virtual | ~PoissonSolver () |
virtual void | resizeMesh (Vector_t &, Vector_t &, const Vector_t &, const Vector_t &, double) |
Private Member Functions | |
void | initializeFields () |
void | mirrorRhoField () |
void | mirrorRhoField (Field_t &ggrn2) |
Private Attributes | |
Field_t | rho2_m |
Field_t | greentr_m |
CxField_t | rho2tr_m |
CxField_t | imgrho2tr_m |
CxField_t | grntr_m |
IField_t | grnIField_m [3] |
std::unique_ptr< FFT_t > | fft_m |
Mesh_t * | mesh_m |
FieldLayout_t * | layout_m |
std::unique_ptr< Mesh_t > | mesh2_m |
std::unique_ptr< FieldLayout_t > | layout2_m |
std::unique_ptr< Mesh_t > | mesh3_m |
std::unique_ptr< FieldLayout_t > | layout3_m |
std::unique_ptr< Mesh_t > | mesh4_m |
std::unique_ptr< FieldLayout_t > | layout4_m |
Field_t | tmpgreen_m |
NDIndex< 3 > | domain_m |
NDIndex< 3 > | domain2_m |
NDIndex< 3 > | domain3_m |
NDIndex< 3 > | domain4_m |
NDIndex< 3 > | domainFFTConstruct_m |
Vector_t | hr_m |
Vektor< int, 3 > | nr_m |
BConds< double, 3, Mesh_t, Center_t > | bc_m |
for defining the boundary conditions More... | |
BConds< Vector_t, 3, Mesh_t, Center_t > | vbc_m |
bool | bcz_m |
bool | integratedGreens_m |
IpplTimings::TimerRef | GreensFunctionTimer_m |
IpplTimings::TimerRef | ComputePotential_m |
Additional Inherited Members | |
![]() | |
typedef Field< int, 3, Mesh_t, Center_t > | IField_t |
typedef Field< std::complex< double >, 3, Mesh_t, Center_t > | CxField_t |
Definition at line 41 of file FFTPoissonSolver.h.
typedef FFT<RCTransform, 3, double> FFTPoissonSolver::FFT_t |
Definition at line 43 of file FFTPoissonSolver.h.
FFTPoissonSolver::FFTPoissonSolver | ( | PartBunch & | bunch, |
std::string | greensFuntion | ||
) |
Definition at line 70 of file FFTPoissonSolver.cpp.
References ComputePotential_m, IpplTimings::getTimer(), greensFunction(), GreensFunctionTimer_m, initializeFields(), and integratedGreens_m.
FFTPoissonSolver::FFTPoissonSolver | ( | Mesh_t * | mesh, |
FieldLayout_t * | fl, | ||
std::string | greensFunction, | ||
std::string | bcz | ||
) |
Definition at line 50 of file FFTPoissonSolver.cpp.
References bcz_m, ComputePotential_m, IpplTimings::getTimer(), greensFunction(), GreensFunctionTimer_m, initializeFields(), and integratedGreens_m.
FFTPoissonSolver::~FFTPoissonSolver | ( | ) |
Definition at line 90 of file FFTPoissonSolver.cpp.
Implements PoissonSolver.
Definition at line 260 of file FFTPoissonSolver.cpp.
References ComputePotential_m, domain_m, fft_m, greensFunction(), GreensFunctionTimer_m, grntr_m, hr_m, integratedGreens_m, integratedGreensFunction(), rho2_m, rho2tr_m, IpplTimings::startTimer(), and IpplTimings::stopTimer().
Implements PoissonSolver.
Definition at line 212 of file FFTPoissonSolver.cpp.
References domain_m, fft_m, GreensFunctionTimer_m, grntr_m, hr_m, imgrho2tr_m, nr_m, rho2_m, rho2tr_m, shiftedIntGreensFunction(), IpplTimings::startTimer(), and IpplTimings::stopTimer().
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 72 of file FFTPoissonSolver.h.
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 71 of file FFTPoissonSolver.h.
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 74 of file FFTPoissonSolver.h.
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 73 of file FFTPoissonSolver.h.
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 76 of file FFTPoissonSolver.h.
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 75 of file FFTPoissonSolver.h.
void FFTPoissonSolver::greensFunction | ( | ) |
Definition at line 305 of file FFTPoissonSolver.cpp.
References fft_m, grnIField_m, grntr_m, hr_m, and rho2_m.
Referenced by computePotential(), and FFTPoissonSolver().
|
private |
Definition at line 92 of file FFTPoissonSolver.cpp.
References bc_m, bcz_m, domain2_m, domain3_m, domain4_m, domain_m, domainFFTConstruct_m, fft_m, UniformCartesian< Dim, MFLOAT >::get_meshSpacing(), FieldLayout< Dim >::getDomain(), FieldLayout< Dim >::getRequestedDistribution(), greentr_m, grnIField_m, grntr_m, hr_m, imgrho2tr_m, Field< T, Dim, M, C >::initialize(), layout2_m, layout3_m, layout4_m, layout_m, cmp_diff::lt(), mesh2_m, mesh3_m, mesh4_m, mesh_m, nr_m, rho2_m, rho2tr_m, tmpgreen_m, vbc_m, and where().
Referenced by FFTPoissonSolver().
void FFTPoissonSolver::integratedGreensFunction | ( | ) |
compute the integrated Green function as described in Three-dimensional quasistatic model for high brightness beam dynamics simulation by Qiang et al.
If the beam has a longitudinal size >> transverse size the direct Green function at each mesh point is not efficient (needs a lot of mesh points along the transverse size to get a good resolution)
If the charge density function is uniform within each cell the following Green's function can be defined:
\[ \overline{G}(x_i - x_{i'}, y_j - y_{j'}, z_k - z_{k'} cout << I << endl; cout << J << endl; cout << K << endl; cout << IE << endl; cout << JE << endl; cout << KE << endl; ) = \int_{x_{i'} - h_x/2}^{x_{i'} + h_x/2} dx' \int_{y_{j'} - h_y/2}^{y_{j'} + h_y/2} dy' \int_{z_{k'} - h_z/2}^{z_{k'} + h_z/2} dz' G(x_i - x_{i'}, y_j - y_{j'}, z_k - z_{k'}). \]
This integral can be calculated analytically in a closed from:
Definition at line 338 of file FFTPoissonSolver.cpp.
References atan(), fft_m, grntr_m, hr_m, layout4_m, log(), mirrorRhoField(), nr_m, rho2_m, sqrt(), and tmpgreen_m.
Referenced by computePotential().
|
private |
Definition at line 498 of file FFTPoissonSolver.cpp.
References bcz_m, nr_m, and rho2_m.
Referenced by integratedGreensFunction(), and shiftedIntGreensFunction().
|
private |
Definition at line 523 of file FFTPoissonSolver.cpp.
void FFTPoissonSolver::shiftedIntGreensFunction | ( | double | zshift | ) |
compute the shifted integrated Green function as described in Three-dimensional quasistatic model for high brightness beam dynamics simulation by Qiang et al.
(x[0:nr_m[0]-1]^2 + y[0:nr_m[1]-1]^2 + (z_c + z[0:nr_m[2]-1])^2)^{-0.5} (x[nr_m[0]:1]^2 + y[0:nr_m[1]-1]^2 + (z_c + z[0:nr_m[2]-1])^2)^{-0.5} (x[0:nr_m[0]-1]^2 + y[nr_m[1]:1]^2 + (z_c + z[0:nr_m[2]-1])^2)^{-0.5} (x[nr_m[0]:1]^2 + y[nr_m[1]:1]^2 + (z_c + z[0:nr_m[2]-1])^2)^{-0.5}
(x[0:nr_m[0]-1]^2 + y[0:nr_m[1]-1]^2 + (z_c - z[nr_m[2]:1])^2)^{-0.5} (x[nr_m[0]:1]^2 + y[0:nr_m[1]-1]^2 + (z_c - z[nr_m[2]:1])^2)^{-0.5} (x[0:nr_m[0]-1]^2 + y[nr_m[1]:1]^2 + (z_c - z[nr_m[2]:1])^2)^{-0.5} (x[nr_m[0]:1]^2 + y[nr_m[1]:1]^2 + (z_c - z[nr_m[2]:1])^2)^{-0.5}
Definition at line 403 of file FFTPoissonSolver.cpp.
References atan(), fft_m, grntr_m, hr_m, layout4_m, log(), mesh4_m, mirrorRhoField(), nr_m, rho2_m, sqrt(), and tmpgreen_m.
Referenced by computePotential().
|
inlinevirtual |
Implements PoissonSolver.
Definition at line 77 of file FFTPoissonSolver.h.
for defining the boundary conditions
Definition at line 141 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 144 of file FFTPoissonSolver.h.
Referenced by FFTPoissonSolver(), initializeFields(), and mirrorRhoField().
|
private |
Definition at line 162 of file FFTPoissonSolver.h.
Referenced by computePotential(), and FFTPoissonSolver().
|
private |
Definition at line 130 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 131 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 132 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 127 of file FFTPoissonSolver.h.
Referenced by computePotential(), and initializeFields().
|
private |
Definition at line 134 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 105 of file FFTPoissonSolver.h.
Referenced by computePotential(), greensFunction(), initializeFields(), integratedGreensFunction(), and shiftedIntGreensFunction().
|
private |
Definition at line 146 of file FFTPoissonSolver.h.
Referenced by computePotential(), and FFTPoissonSolver().
|
private |
Definition at line 89 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 102 of file FFTPoissonSolver.h.
Referenced by greensFunction(), and initializeFields().
|
private |
Definition at line 98 of file FFTPoissonSolver.h.
Referenced by computePotential(), greensFunction(), initializeFields(), integratedGreensFunction(), and shiftedIntGreensFunction().
|
private |
Definition at line 137 of file FFTPoissonSolver.h.
Referenced by computePotential(), greensFunction(), initializeFields(), integratedGreensFunction(), print(), and shiftedIntGreensFunction().
|
private |
Definition at line 94 of file FFTPoissonSolver.h.
Referenced by computePotential(), and initializeFields().
|
private |
Definition at line 145 of file FFTPoissonSolver.h.
Referenced by computePotential(), and FFTPoissonSolver().
|
private |
Definition at line 113 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 117 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 121 of file FFTPoissonSolver.h.
Referenced by initializeFields(), integratedGreensFunction(), and shiftedIntGreensFunction().
|
private |
Definition at line 109 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 112 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 116 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 120 of file FFTPoissonSolver.h.
Referenced by initializeFields(), and shiftedIntGreensFunction().
|
private |
Definition at line 108 of file FFTPoissonSolver.h.
Referenced by initializeFields().
|
private |
Definition at line 138 of file FFTPoissonSolver.h.
Referenced by computePotential(), initializeFields(), integratedGreensFunction(), mirrorRhoField(), and shiftedIntGreensFunction().
|
private |
Definition at line 86 of file FFTPoissonSolver.h.
Referenced by computePotential(), greensFunction(), initializeFields(), integratedGreensFunction(), mirrorRhoField(), and shiftedIntGreensFunction().
|
private |
Definition at line 93 of file FFTPoissonSolver.h.
Referenced by computePotential(), and initializeFields().
|
private |
Definition at line 124 of file FFTPoissonSolver.h.
Referenced by initializeFields(), integratedGreensFunction(), and shiftedIntGreensFunction().
Definition at line 142 of file FFTPoissonSolver.h.
Referenced by initializeFields().