OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
Public Member Functions | Private Member Functions | Private Attributes | List of all members
FFTPoissonSolver Class Reference

#include <FFTPoissonSolver.h>

Inheritance diagram for FFTPoissonSolver:
Inheritance graph
[legend]
Collaboration diagram for FFTPoissonSolver:
Collaboration graph
[legend]

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 computePotentialDKS (Field_t &rho)
 
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 integratedGreensFunctionDKS ()
 Uses DKS to offload the computation of Greens function on the GPU. 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 level)
 
double getXRangeMax (unsigned short level)
 
double getYRangeMin (unsigned short level)
 
double getYRangeMax (unsigned short level)
 
double getZRangeMin (unsigned short level)
 
double getZRangeMax (unsigned short level)
 
void test (PartBunchBase< double, 3 > *bunch)
 
Informprint (Inform &os) const
 
- Public Member Functions inherited from PoissonSolver
virtual void solve (AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true)
 
virtual void hasToRegrid ()
 
virtual ~PoissonSolver ()
 

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_tfft_m
 
Mesh_tmesh_m
 
FieldLayout_tlayout_m
 
std::unique_ptr< Mesh_tmesh2_m
 
std::unique_ptr< FieldLayout_tlayout2_m
 
std::unique_ptr< Mesh_tmesh3_m
 
std::unique_ptr< FieldLayout_tlayout3_m
 
std::unique_ptr< Mesh_tmesh4_m
 
std::unique_ptr< FieldLayout_tlayout4_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
 

Detailed Description

Definition at line 43 of file FFTPoissonSolver.h.

Constructor & Destructor Documentation

FFTPoissonSolver::FFTPoissonSolver ( PartBunch bunch,
std::string  greensFuntion 
)

Definition at line 70 of file FFTPoissonSolver.cpp.

References ComputePotential_m, IpplTimings::getTimer(), GreensFunctionTimer_m, initializeFields(), and integratedGreens_m.

Here is the call graph for this function:

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(), GreensFunctionTimer_m, initializeFields(), and integratedGreens_m.

Here is the call graph for this function:

FFTPoissonSolver::~FFTPoissonSolver ( )

Member Function Documentation

void FFTPoissonSolver::computePotential ( Field_t rho,
Vector_t  hr,
double  zshift 
)
virtual

Implements PoissonSolver.

Definition at line 289 of file FFTPoissonSolver.cpp.

References domain_m, fft_m, GreensFunctionTimer_m, grntr_m, hr_m, imgrho2tr_m, K, nr_m, rho2_m, rho2tr_m, shiftedIntGreensFunction(), IpplTimings::startTimer(), and IpplTimings::stopTimer().

Here is the call graph for this function:

void FFTPoissonSolver::computePotential ( Field_t rho,
Vector_t  hr 
)
virtual
void FFTPoissonSolver::computePotentialDKS ( Field_t rho)
double FFTPoissonSolver::getXRangeMax ( unsigned short  level)
inlinevirtual

Implements PoissonSolver.

Definition at line 81 of file FFTPoissonSolver.h.

double FFTPoissonSolver::getXRangeMin ( unsigned short  level)
inlinevirtual

Implements PoissonSolver.

Definition at line 80 of file FFTPoissonSolver.h.

double FFTPoissonSolver::getYRangeMax ( unsigned short  level)
inlinevirtual

Implements PoissonSolver.

Definition at line 83 of file FFTPoissonSolver.h.

double FFTPoissonSolver::getYRangeMin ( unsigned short  level)
inlinevirtual

Implements PoissonSolver.

Definition at line 82 of file FFTPoissonSolver.h.

double FFTPoissonSolver::getZRangeMax ( unsigned short  level)
inlinevirtual

Implements PoissonSolver.

Definition at line 85 of file FFTPoissonSolver.h.

double FFTPoissonSolver::getZRangeMin ( unsigned short  level)
inlinevirtual

Implements PoissonSolver.

Definition at line 84 of file FFTPoissonSolver.h.

void FFTPoissonSolver::greensFunction ( )

Definition at line 426 of file FFTPoissonSolver.cpp.

References fft_m, grnIField_m, grntr_m, hr_m, and rho2_m.

Referenced by computePotential().

void FFTPoissonSolver::initializeFields ( )
private
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 459 of file FFTPoissonSolver.cpp.

References atan(), fft_m, grntr_m, hr_m, K, layout4_m, log(), mirrorRhoField(), nr_m, rho2_m, sqrt(), and tmpgreen_m.

Referenced by computePotential().

Here is the call graph for this function:

void FFTPoissonSolver::integratedGreensFunctionDKS ( )

Uses DKS to offload the computation of Greens function on the GPU.

Definition at line 524 of file FFTPoissonSolver.cpp.

References hr_m, K, layout4_m, and nr_m.

Referenced by computePotentialDKS().

void FFTPoissonSolver::mirrorRhoField ( )
private

Definition at line 641 of file FFTPoissonSolver.cpp.

References bcz_m, K, nr_m, and rho2_m.

Referenced by integratedGreensFunction(), and shiftedIntGreensFunction().

void FFTPoissonSolver::mirrorRhoField ( Field_t ggrn2)
private

Definition at line 666 of file FFTPoissonSolver.cpp.

References bcz_m, K, nr_m, and rho2_m.

Inform & FFTPoissonSolver::print ( Inform os) const

Definition at line 693 of file FFTPoissonSolver.cpp.

References endl(), and hr_m.

Referenced by operator<<().

Here is the call graph for this function:

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 546 of file FFTPoissonSolver.cpp.

References atan(), fft_m, grntr_m, hr_m, K, layout4_m, log(), mesh4_m, mirrorRhoField(), nr_m, rho2_m, sqrt(), and tmpgreen_m.

Referenced by computePotential().

Here is the call graph for this function:

void FFTPoissonSolver::test ( PartBunchBase< double, 3 > *  bunch)
inlinevirtual

Implements PoissonSolver.

Definition at line 86 of file FFTPoissonSolver.h.

Member Data Documentation

BConds<double, 3, Mesh_t, Center_t> FFTPoissonSolver::bc_m
private

for defining the boundary conditions

Definition at line 164 of file FFTPoissonSolver.h.

Referenced by initializeFields().

bool FFTPoissonSolver::bcz_m
private

Definition at line 167 of file FFTPoissonSolver.h.

Referenced by FFTPoissonSolver(), initializeFields(), and mirrorRhoField().

IpplTimings::TimerRef FFTPoissonSolver::ComputePotential_m
private

Definition at line 185 of file FFTPoissonSolver.h.

Referenced by computePotential(), and FFTPoissonSolver().

NDIndex<3> FFTPoissonSolver::domain2_m
private

Definition at line 153 of file FFTPoissonSolver.h.

Referenced by initializeFields().

NDIndex<3> FFTPoissonSolver::domain3_m
private

Definition at line 154 of file FFTPoissonSolver.h.

Referenced by initializeFields().

NDIndex<3> FFTPoissonSolver::domain4_m
private

Definition at line 155 of file FFTPoissonSolver.h.

Referenced by initializeFields().

NDIndex<3> FFTPoissonSolver::domain_m
private

Definition at line 150 of file FFTPoissonSolver.h.

Referenced by computePotential(), and initializeFields().

NDIndex<3> FFTPoissonSolver::domainFFTConstruct_m
private

Definition at line 157 of file FFTPoissonSolver.h.

Referenced by initializeFields().

std::unique_ptr<FFT_t> FFTPoissonSolver::fft_m
private
IpplTimings::TimerRef FFTPoissonSolver::GreensFunctionTimer_m
private

Definition at line 169 of file FFTPoissonSolver.h.

Referenced by computePotential(), computePotentialDKS(), and FFTPoissonSolver().

Field_t FFTPoissonSolver::greentr_m
private

Definition at line 98 of file FFTPoissonSolver.h.

Referenced by initializeFields().

IField_t FFTPoissonSolver::grnIField_m[3]
private

Definition at line 125 of file FFTPoissonSolver.h.

Referenced by greensFunction(), and initializeFields().

CxField_t FFTPoissonSolver::grntr_m
private
Vector_t FFTPoissonSolver::hr_m
private
CxField_t FFTPoissonSolver::imgrho2tr_m
private

Definition at line 103 of file FFTPoissonSolver.h.

Referenced by computePotential(), and initializeFields().

bool FFTPoissonSolver::integratedGreens_m
private

Definition at line 168 of file FFTPoissonSolver.h.

Referenced by computePotential(), and FFTPoissonSolver().

std::unique_ptr<FieldLayout_t> FFTPoissonSolver::layout2_m
private

Definition at line 136 of file FFTPoissonSolver.h.

Referenced by initializeFields().

std::unique_ptr<FieldLayout_t> FFTPoissonSolver::layout3_m
private

Definition at line 140 of file FFTPoissonSolver.h.

Referenced by initializeFields().

std::unique_ptr<FieldLayout_t> FFTPoissonSolver::layout4_m
private
FieldLayout_t* FFTPoissonSolver::layout_m
private

Definition at line 132 of file FFTPoissonSolver.h.

Referenced by initializeFields().

std::unique_ptr<Mesh_t> FFTPoissonSolver::mesh2_m
private

Definition at line 135 of file FFTPoissonSolver.h.

Referenced by initializeFields().

std::unique_ptr<Mesh_t> FFTPoissonSolver::mesh3_m
private

Definition at line 139 of file FFTPoissonSolver.h.

Referenced by initializeFields().

std::unique_ptr<Mesh_t> FFTPoissonSolver::mesh4_m
private

Definition at line 143 of file FFTPoissonSolver.h.

Referenced by initializeFields(), and shiftedIntGreensFunction().

Mesh_t* FFTPoissonSolver::mesh_m
private

Definition at line 131 of file FFTPoissonSolver.h.

Referenced by initializeFields().

Vektor<int, 3> FFTPoissonSolver::nr_m
private
Field_t FFTPoissonSolver::rho2_m
private
CxField_t FFTPoissonSolver::rho2tr_m
private

Definition at line 102 of file FFTPoissonSolver.h.

Referenced by computePotential(), and initializeFields().

Field_t FFTPoissonSolver::tmpgreen_m
private
BConds<Vector_t, 3, Mesh_t, Center_t> FFTPoissonSolver::vbc_m
private

Definition at line 165 of file FFTPoissonSolver.h.

Referenced by initializeFields().


The documentation for this class was generated from the following files: