PoissonConvolve3D Class Reference

#include <PoissonConvolve3D.hh>

List of all members.


Detailed Description

class PoissonConvolve3D

A solver for Poisson's equation in 3-D Cartesian coordinates {x,y,z}, using convolution of the Green's function (whch has infinite boundary conditions). The convolution is sped up by using an FFT on a grid twice as large (in each dimension) as the solution (physical) grid.

This class uses the Poisson solver on a rectangular 3-d grid; inside the grid, it calculates E by applying the analytic derivative of the interpolation function to the values of phi at the surrounding 8 grid points. Outside the grid, by default it uses the analytic value of E for a single point charge, with its charge equal to the total charge in the grid, and its position at the mean position of the charge in the grid. More complicated approximations, based on multiple point charges, can be used (not discussed here).

Units are determined by the code using this class. Normally the charge values will not be integers, but rather integers times the positron charge divided by epsilon_0 (the permittivity of free space), the grid sizes will have units of length; the result is that phi() and getE() will then have the appropriate units.

This class solves: del^2 phi = -4*pi*rho On a specified grid in {x,y,z}. rho(x,y,z) is the charge density. The -4*pi and the conversion from charge to charge density is handled internally in putCharge(); they are NOT handled in setRhs() (which is intended primarily for testing).

This class includes an approximation for use outside the grid. The default is simply to use the total charge put into the grid (via putCharge()), and the weighted mean position of that charge. Better approximations can be computed externally as a vector of point charges and used via setApproximation().

Note: array[] contains CHARGE between calls to zeroRhs() and solve(); it contains PHI between calls to solve() and zeroRhs(); status reflects its contents, ans assert()-s are used to ensure accesses are valid..

Accuracy for a 65x65x65 grid of unit size: ... MORE TO COME ...

Public Types

enum  PutMode { LOWEST, NEAREST, PRORATED }
 enum PutMode determines how putCharge() places charge into the grid. PRORATED is the default. LOWEST places charge into the grid point lower than (x,y,z). NEAREST places charge into the nearest grid point to (x,y,z). PRORATED places charge into the 8 nearest grid points, pro-rated by proximity to (x,y,z). For 1M particles into a 129^3 grid, NEAREST and LOWER were less than 10% faster than PRORATED, giving slightly more ragged results. More...

Public Member Functions

 PoissonConvolve3D (int _nx, int _ny, int _nz, float _xMin, float _xMax, float _yMin, float _yMax, float _zMin, float _zMax)
 Constructor. nx,ny,nz are the # grid points in the 3 dimensions -- NOTE: these are the size of the PHYSICAL grid, the actual grid is doubled in each dimension. They can be any integers >= 9, but powers of 2 have a performance advantage. Note that 512x512x512 allocates too much memory for a 32-bit program on an Intel processor. xMin,xMax,yMin,yMax,zMin,zMax are the grid limits in 3 dimensions. That implies array[index(0,iy,iz)] is the boundary at x=xMin; array[index(nx-1,iy,iz)] is the boundary at x=xMax; etc. 1.0e-4 is appropriate.
virtual ~PoissonConvolve3D ()
 Destructor.
void updatePosition (float _xMin, float _xMax, float _yMin, float _yMax, float _zMin, float _zMax)
 updatePosition() updates the position of the boundaries. Usually followed by a call to setApproximation(), with the approximation charge positions updated similarly, or by a series of calls to putCharge() and then defaultApproximation().
float phi (float x, float y, float z)
 phi() works for all {x,y,z}, using the solution inside the grid and the approximation outside the grid. Note that phi() is not necessarily continuous or differentiable, but it should be approximately continuous. Do not differentiate this function, use getE() instead (it is MUCH more intelligent about the grid and its bounadries). Intended primarily for testing.
float phiGrid (int ix, int iy, int iz) const
 phiGrid() returns the value of the solution to Poisson's equation (i.e. the electrostatic potential) on grid points. ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 Intended primarily for testing.
void setPhi (int ix, int iy, int iz, float v)
 setPhi() will set a value in array[]. Intended primarily for testing.
void addPhi (int ix, int iy, int iz, float v)
 addPhi() will add v to a value in array[]. Intended primarily for testing.
float rhs (int ix, int iy, int iz) const
 rhs() returns the charge at a grid point of the right-hand-side of Poisson's equation. NOTE: ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1
void setRhs (int ix, int iy, int iz, float v)
 setRhs() will set an entry in rhs[]. Normally putCharge() is better. ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 This directly sets the rhs array entry, so you probably need to multiply by -4*pi/(dx*dy*dz) to get the proper R.H.S. for Poisson's equation; putCharge() does that internally. Intended primarily for testing.
bool putCharge (float x, float y, float z, float c)
 putCharge() puts a charge into the rhs grid according to the current putMode. Does nothing if outside the grid. Multiplies by -1.0/(dx*dy*dz) to convert the charge to the R.H.S. of Poisson's equation, which is minus the charge density. NOTE: dividing by epsilon0 is up to the user's choice of units. Accumulates the total charge and mean position of charge actually placed into the grid (for use by defaultApproximation()). Returns true if the charge was placed within the grid, false if not.
bool putCharge (Charge &c)
void setPutMode (PutMode m)
 setPutMode() sets the mode of putCharge(). See PutMode above.
void zeroRhs ()
 zeroRhs() will zero rhs.
G4ThreeVector getE (float x, float y, float z)
 getE() returns the value of the E field at a given point. Uses the solution inside the grid and the approximation outside. It handles grid cells and boundaries intelligently.
G4ThreeVector approxE (double x, double y, double z)
 approxE() evaluates the approximation E field. {x,y,z} can be anywhere, but normally are outside the grid or on its boundary.
float approxPhi (double x, double y, double z)
 approxPhi() evaluates the approximation phi. {x,y,z} can be anywhere, but normally are outside the grid or on its boundary.
G4ThreeVector getPosition (int ix, int iy, int iz)
 getPosition() returns the position corresponding to grid point indexes. Intended primarily for testing.
void setApproximation (std::vector< Charge > v)
 setApproximation() sets the approximation to use outside the grid, and on its boundary. v is a vector of point charges. Updates the boundary values of phi.
void defaultApproximation ()
 defaultApproximation() calls setApproximation() with the total charge and mean position in rhs (from putCharge()). Updates the boundary values of phi, and the interior if interior=true.
void computeGreensFunction ()
 computeGreensFunction() does just that. Must be called whenever the physical size of the grid is changed.
void solve ()
 solve() will solve Poisson's equation
void printRhs ()
 printRhs() prints the rhs. Intended primarily for testing.
void printPhi ()
 printPhi() prints phi. Intended primarily for testing.
int getNx () const
 accessors.
int getNy () const
int getNz () const

Protected Types

enum  Status { UNINITIALIZED, CHARGE, PHI }

Protected Member Functions

int index (int ix, int iy, int iz) const

Protected Attributes

int nx
int ny
int nz
float xMin
float xMax
float yMin
float yMax
float zMin
float zMax
Status status
float dx
float dy
float dz
float * array
fftwf_complex * transf
fftwf_complex * greensFunction
std::vector< Chargeapprox
float sumX
float sumY
float sumZ
float totalCharge
float vol
PutMode putMode
int nc
fftwf_plan gfPlan
fftwf_plan plan1
fftwf_plan plan2


Member Enumeration Documentation

enum PutMode determines how putCharge() places charge into the grid. PRORATED is the default. LOWEST places charge into the grid point lower than (x,y,z). NEAREST places charge into the nearest grid point to (x,y,z). PRORATED places charge into the 8 nearest grid points, pro-rated by proximity to (x,y,z). For 1M particles into a 129^3 grid, NEAREST and LOWER were less than 10% faster than PRORATED, giving slightly more ragged results.

Enumerator:
LOWEST 
NEAREST 
PRORATED 
00120 {LOWEST, NEAREST, PRORATED};

enum PoissonConvolve3D::Status [protected]

Enumerator:
UNINITIALIZED 
CHARGE 
PHI 


Constructor & Destructor Documentation

PoissonConvolve3D::PoissonConvolve3D ( int  _nx,
int  _ny,
int  _nz,
float  _xMin,
float  _xMax,
float  _yMin,
float  _yMax,
float  _zMin,
float  _zMax 
)

Constructor. nx,ny,nz are the # grid points in the 3 dimensions -- NOTE: these are the size of the PHYSICAL grid, the actual grid is doubled in each dimension. They can be any integers >= 9, but powers of 2 have a performance advantage. Note that 512x512x512 allocates too much memory for a 32-bit program on an Intel processor. xMin,xMax,yMin,yMax,zMin,zMax are the grid limits in 3 dimensions. That implies array[index(0,iy,iz)] is the boundary at x=xMin; array[index(nx-1,iy,iz)] is the boundary at x=xMax; etc. 1.0e-4 is appropriate.

References array, gfPlan, greensFunction, nc, nx, ny, nz, plan1, plan2, PRORATED, putMode, status, sumX, sumY, sumZ, totalCharge, transf, UNINITIALIZED, and updatePosition().

00015                                                           : approx() {
00016         nx=_nx; ny=_ny; nz=_nz;
00017         nc = 4*nx*ny*(nz+1); // double size, except last dimension
00018         status = UNINITIALIZED;
00019         array = (float*)fftwf_malloc(sizeof(float)*8*nx*ny*nz);
00020         transf = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*nc);
00021         greensFunction = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*nc);
00022         sumX = 0.0; sumY = 0.0; sumZ = 0.0; totalCharge = 0.0;
00023         putMode = PRORATED;
00024 
00025         fftwf_set_timelimit(10.0); // otherwise can take up to 226 sec (!)
00026 
00027         gfPlan = fftwf_plan_dft_r2c_3d(2*nx,2*ny,2*nz,array,greensFunction,
00028                                                         FFTW_ESTIMATE);
00029         plan1 = fftwf_plan_dft_r2c_3d(2*nx,2*ny,2*nz,array,transf,
00030                                                         FFTW_ESTIMATE);
00031         plan2 = fftwf_plan_dft_c2r_3d(2*nx,2*ny,2*nz,transf,array,
00032                                                         FFTW_ESTIMATE);
00033 
00034         updatePosition(_xMin,_xMax,_yMin,_yMax,_zMin,_zMax);
00035 }

virtual PoissonConvolve3D::~PoissonConvolve3D (  )  [inline, virtual]

Destructor.

References approx, array, gfPlan, greensFunction, plan1, plan2, and transf.

00158                                      {
00159                 if(array) fftwf_free(array);
00160                 if(transf) fftwf_free(transf);
00161                 if(greensFunction) fftwf_free(greensFunction);
00162                 fftwf_destroy_plan(gfPlan);
00163                 fftwf_destroy_plan(plan1);
00164                 fftwf_destroy_plan(plan2);
00165                 approx.clear();
00166         }


Member Function Documentation

int PoissonConvolve3D::index ( int  ix,
int  iy,
int  iz 
) const [inline, protected]

References ny, and nz.

Referenced by addPhi(), computeGreensFunction(), getE(), phi(), phiGrid(), printPhi(), printRhs(), putCharge(), rhs(), setPhi(), setRhs(), and zeroRhs().

00138                                                 {
00139                 assert(ix>=0&&ix<2*nx && iy>=0&&iy<2*ny && iz>=0&&iz<2*nz);
00140                 return iz+2*nz*(iy+2*ny*ix);  // C order
00141         }

void PoissonConvolve3D::updatePosition ( float  _xMin,
float  _xMax,
float  _yMin,
float  _yMax,
float  _zMin,
float  _zMax 
) [inline]

updatePosition() updates the position of the boundaries. Usually followed by a call to setApproximation(), with the approximation charge positions updated similarly, or by a series of calls to putCharge() and then defaultApproximation().

References computeGreensFunction(), dx, dy, dz, nx, ny, nz, vol, xMax, xMin, yMax, yMin, zMax, and zMin.

Referenced by Bunch::computeField(), Bunch::init(), and PoissonConvolve3D().

00173                                                                             {
00174                 xMin=_xMin; xMax=_xMax;
00175                 yMin=_yMin; yMax=_yMax;
00176                 zMin=_zMin; zMax=_zMax;
00177                 dx = (xMax-xMin)/(nx-1);
00178                 dy = (yMax-yMin)/(ny-1);
00179                 dz = (zMax-zMin)/(nz-1);
00180                 vol = dx*dy*dz;
00181                 computeGreensFunction();
00182         }

float PoissonConvolve3D::phi ( float  x,
float  y,
float  z 
)

phi() works for all {x,y,z}, using the solution inside the grid and the approximation outside the grid. Note that phi() is not necessarily continuous or differentiable, but it should be approximately continuous. Do not differentiate this function, use getE() instead (it is MUCH more intelligent about the grid and its bounadries). Intended primarily for testing.

References approx, array, dx, dy, dz, index(), M_PI, nx, ny, nz, PHI, status, xMin, yMin, and zMin.

Referenced by approxPhi().

00038 {
00039         assert(status==PHI);
00040         int ix = (int)floor((x-xMin)/dx);
00041         int iy = (int)floor((y-yMin)/dy);
00042         int iz = (int)floor((z-zMin)/dz);
00043         if(ix < 0 || iy < 0 || iz < 0 || ix+1 >= nx || iy+1 >= ny || 
00044                                                                 iz+1 >= nz) {
00045                 // outside grid, use the approximation
00046                 float v=0.0;
00047                 G4ThreeVector pos(x,y,z);
00048                 for(unsigned i=0; i<approx.size(); ++i) {
00049                         G4ThreeVector pt = approx[i].getPosition();
00050                         pt -= pos;
00051                         double r = pt.mag();
00052                         v += approx[i].c / r / (4.0*M_PI);
00053                 }
00054                 return v;
00055         }
00056 
00057         // within grid, interpolate linearly among 8 surrounding grid points
00058         float fx = (x - (xMin+ix*dx))/dx;
00059         assert(-0.01 <= fx && fx <= 1.01);
00060         float fy = (y - (yMin+iy*dy))/dy;
00061         assert(-0.01 <= fy && fy <= 1.01);
00062         float fz = (z - (zMin+iz*dz))/dz;
00063         assert(-0.01 <= fz && fz <= 1.01);
00064         float gx=1.0-fx, gy=1.0-fy, gz=1.0-fz;
00065 
00066         float v =
00067                 gx*gy*gz*array[index(ix,iy,iz)] +
00068                 fx*gy*gz*array[index(ix+1,iy,iz)] +
00069                 gx*fy*gz*array[index(ix,iy+1,iz)] +
00070                 fx*fy*gz*array[index(ix+1,iy+1,iz)] +
00071                 gx*gy*fz*array[index(ix,iy,iz+1)] +
00072                 fx*gy*fz*array[index(ix+1,iy,iz+1)] +
00073                 gx*fy*fz*array[index(ix,iy+1,iz+1)] +
00074                 fx*fy*fz*array[index(ix+1,iy+1,iz+1)];
00075 
00076         return v;
00077 }

float PoissonConvolve3D::phiGrid ( int  ix,
int  iy,
int  iz 
) const [inline]

phiGrid() returns the value of the solution to Poisson's equation (i.e. the electrostatic potential) on grid points. ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 Intended primarily for testing.

References array, index(), PHI, and status.

00198                 { assert(status==PHI);  return array[index(ix,iy,iz)]; }

void PoissonConvolve3D::setPhi ( int  ix,
int  iy,
int  iz,
float  v 
) [inline]

setPhi() will set a value in array[]. Intended primarily for testing.

References array, index(), PHI, and status.

00203                 { assert(status==PHI);  array[index(ix,iy,iz)] = v; }

void PoissonConvolve3D::addPhi ( int  ix,
int  iy,
int  iz,
float  v 
) [inline]

addPhi() will add v to a value in array[]. Intended primarily for testing.

References array, index(), PHI, and status.

00208                 { assert(status==PHI);  array[index(ix,iy,iz)] += v; }

float PoissonConvolve3D::rhs ( int  ix,
int  iy,
int  iz 
) const [inline]

rhs() returns the charge at a grid point of the right-hand-side of Poisson's equation. NOTE: ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1

References array, CHARGE, index(), M_PI, and status.

00214                 { assert(status==CHARGE);
00215                   return 4.0*M_PI*array[index(ix,iy,iz)]; }

void PoissonConvolve3D::setRhs ( int  ix,
int  iy,
int  iz,
float  v 
) [inline]

setRhs() will set an entry in rhs[]. Normally putCharge() is better. ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 This directly sets the rhs array entry, so you probably need to multiply by -4*pi/(dx*dy*dz) to get the proper R.H.S. for Poisson's equation; putCharge() does that internally. Intended primarily for testing.

References array, CHARGE, index(), and status.

00224                 { assert(status==CHARGE); array[index(ix,iy,iz)] = v; }

bool PoissonConvolve3D::putCharge ( float  x,
float  y,
float  z,
float  c 
)

putCharge() puts a charge into the rhs grid according to the current putMode. Does nothing if outside the grid. Multiplies by -1.0/(dx*dy*dz) to convert the charge to the R.H.S. of Poisson's equation, which is minus the charge density. NOTE: dividing by epsilon0 is up to the user's choice of units. Accumulates the total charge and mean position of charge actually placed into the grid (for use by defaultApproximation()). Returns true if the charge was placed within the grid, false if not.

References array, CHARGE, dx, dy, dz, index(), LOWEST, M_PI, NEAREST, nx, ny, nz, PRORATED, putMode, status, sumX, sumY, sumZ, totalCharge, xMin, yMin, and zMin.

Referenced by Bunch::computeField(), putCharge(), Bunch::testBfield(), and Bunch::testEfield().

00080 {
00081         assert(status==CHARGE);
00082         int ix = (int)floor((x-xMin)/dx);
00083         int iy = (int)floor((y-yMin)/dy);
00084         int iz = (int)floor((z-zMin)/dz);
00085         // (omit boundaries)
00086         if(ix<=0 || iy<=0 || iz<=0 || ix+1>=nx || iy+1>=ny || iz+1>=nz) {
00087                 return false;
00088         }
00089 
00090         sumX += c*x;
00091         sumY += c*y;
00092         sumZ += c*z;
00093         totalCharge += c;
00094 
00095         c /= 4.0*M_PI;
00096 
00097         switch(putMode) {
00098         case LOWEST:
00099                 array[index(ix,iy,iz)] += c;
00100                 break;
00101         case NEAREST:
00102                 if(x > xMin+ix*dx+dx/2.0) ++ix; // could be upper boundary
00103                 if(y > yMin+iy*dy+dy/2.0) ++iy; // but is guaranteed to 
00104                 if(z > zMin+iz*dz+dz/2.0) ++iz; // be in range.
00105                 array[index(ix,iy,iz)] += c;
00106                 break;
00107         case PRORATED:
00108                 { float fx = (x - (xMin+ix*dx))/dx;
00109                   assert(-0.01 <= fx && fx <= 1.01);
00110                   float fy = (y - (yMin+iy*dy))/dy;
00111                   assert(-0.01 <= fy && fy <= 1.01);
00112                   float fz = (z - (zMin+iz*dz))/dz;
00113                   assert(-0.01 <= fz && fz <= 1.01);
00114                   float gx=1.0-fx, gy=1.0-fy, gz=1.0-fz;
00115                   array[index(ix,iy,iz)] += gx*gy*gz*c;
00116                   array[index(ix+1,iy,iz)] += fx*gy*gz*c;
00117                   array[index(ix,iy+1,iz)] += gx*fy*gz*c;
00118                   array[index(ix+1,iy+1,iz)] += fx*fy*gz*c;
00119                   array[index(ix,iy,iz+1)] += gx*gy*fz*c;
00120                   array[index(ix+1,iy,iz+1)] += fx*gy*fz*c;
00121                   array[index(ix,iy+1,iz+1)] += gx*fy*fz*c;
00122                   array[index(ix+1,iy+1,iz+1)] += fx*fy*fz*c;
00123                 }
00124                 break;
00125         }
00126 
00127         return true;
00128 }

bool PoissonConvolve3D::putCharge ( Charge c  )  [inline]

References Charge::c, putCharge(), Charge::x, Charge::y, and Charge::z.

00235 { return putCharge(c.x,c.y,c.z,c.c); }

void PoissonConvolve3D::setPutMode ( PutMode  m  )  [inline]

setPutMode() sets the mode of putCharge(). See PutMode above.

References putMode.

00238 { putMode = m; }

void PoissonConvolve3D::zeroRhs (  ) 

zeroRhs() will zero rhs.

References array, CHARGE, index(), nx, ny, nz, status, sumX, sumY, sumZ, and totalCharge.

Referenced by Bunch::computeField(), Bunch::testBfield(), and Bunch::testEfield().

00131 {
00132         sumX = sumY = sumZ = totalCharge = 0;
00133         for(int i=0; i<2*nx; ++i) {
00134                 for(int j=0; j<2*ny; ++j) {
00135                         for(int k=0; k<2*nz; ++k) {
00136                                 array[index(i,j,k)] = 0.0;
00137                         }
00138                 }
00139         }
00140 
00141         // indicate CHARGE is now in array[]
00142         status = CHARGE;
00143 }

G4ThreeVector PoissonConvolve3D::getE ( float  x,
float  y,
float  z 
)

getE() returns the value of the E field at a given point. Uses the solution inside the grid and the approximation outside. It handles grid cells and boundaries intelligently.

References approxE(), array, dx, dy, dz, index(), nx, ny, nz, PHI, status, xMax, xMin, yMax, yMin, zMax, and zMin.

Referenced by Bunch::addFieldValue(), Bunch::testBfield(), and Bunch::testEfield().

00146 {
00147         assert(status==PHI);
00148         int ix = (int)floor((x-xMin)/dx);
00149         int iy = (int)floor((y-yMin)/dy);
00150         int iz = (int)floor((z-zMin)/dz);
00151         if(ix < 0 || iy < 0 || iz < 0 || ix+1 >= nx || iy+1 >= ny || 
00152                                                                 iz+1 >= nz) {
00153                 return approxE(x,y,z);
00154         }
00155 
00156 #ifdef SIMPLE_LINEAR
00157         // Has the problem that the derivative is constant in each box.
00158         // variables used in mathematica notebook "interpolationOn3DGrid.nb"
00159         // -- it differentiates the 8-point interpolation formula.
00160         double x0=xMin+ix*dx, y0=yMin+iy*dy, z0=zMin+iz*dz;
00161         double x1=x0+dx, y1=y0+dy, z1=z0+dz;
00162         double v000=array[index(ix,iy,iz)];
00163         double v100=array[index(ix+1,iy,iz)];
00164         double v010=array[index(ix,iy+1,iz)];
00165         double v110=array[index(ix+1,iy+1,iz)];
00166         double v001=array[index(ix,iy,iz+1)];
00167         double v101=array[index(ix+1,iy,iz+1)];
00168         double v011=array[index(ix,iy+1,iz+1)];
00169         double v111=array[index(ix+1,iy+1,iz+1)];
00170 
00171         // E = - grad array[]
00172         double Ex=-(1.0/(dx*dy*dz))*(dy*(dz*(-v000+v100)+(v000-v001-v100+v101)
00173                 *(z-z0))+(y-y0)*(dz*(v000-v010-v100+v110)-
00174                 (v000-v001-v010+v011-v100+v101+v110-v111)*(z-z0)));
00175         double Ey=-(1.0/(dx*dy*dz))*(dx*(dz*(-v000+v010)+(v000-v001-v010+v011)
00176                 *(z-z0))+(x-x0)*(dz*(v000-v010-v100+v110)-
00177                 (v000-v001-v010+v011-v100+v101+v110-v111)*(z-z0)));
00178         double Ez=-(1.0/(dx*dy*dz))*(dx*(dy*(-v000+v001)+(v000-v001-v010+v011)
00179                 *(y-y0))+(x-x0)*(dy*(v000-v001-v100+v101)-
00180                 (v000-v001-v010+v011-v100+v101+v110-v111)*(y-y0)));
00181 #else
00182         // linear interpolation between derivatives in adjacent cells.
00183         // Associate the derivative in each grid box with the center of the box;
00184         // interpolate the derivative between the current box and the closer
00185         // adjacent one, or the boundary if this box is an edge.
00186         float fx = (x - (xMin+ix*dx))/dx;
00187         assert(-0.01 <= fx && fx <= 1.01);
00188         float fy = (y - (yMin+iy*dy))/dy;
00189         assert(-0.01 <= fy && fy <= 1.01);
00190         float fz = (z - (zMin+iz*dz))/dz;
00191         assert(-0.01 <= fz && fz <= 1.01);
00192         float gx=1.0-fx, gy=1.0-fy, gz=1.0-fz;
00193 
00194         double v000=array[index(ix,iy,iz)];
00195         double v100=array[index(ix+1,iy,iz)];
00196         double v010=array[index(ix,iy+1,iz)];
00197         double v110=array[index(ix+1,iy+1,iz)];
00198         double v001=array[index(ix,iy,iz+1)];
00199         double v101=array[index(ix+1,iy,iz+1)];
00200         double v011=array[index(ix,iy+1,iz+1)];
00201         double v111=array[index(ix+1,iy+1,iz+1)];
00202 
00203         // partial d/dx
00204         double v0=gy*gz*v000 + fy*gz*v010 + gy*fz*v001 + fy*fz*v011;
00205         double v1=gy*gz*v100 + fy*gz*v110 + gy*fz*v101 + fy*fz*v111;
00206         double Ex=-(v1-v0)/dx; // for this grid box
00207         if(fx > 0.5) {
00208                 if(ix+1 == nx-1) {
00209                         // interpolate between center of box and xMax
00210                         double Ex2=approxE(xMax,y,z).x();
00211                         double f=2.0*fx-1.0;
00212                         Ex = (1.0-f)*Ex + f*Ex2;
00213                 } else {
00214                         // interpolate between this box and next higher one
00215                         double v200=array[index(ix+2,iy,iz)];
00216                         double v210=array[index(ix+2,iy+1,iz)];
00217                         double v201=array[index(ix+2,iy,iz+1)];
00218                         double v211=array[index(ix+2,iy+1,iz+1)];
00219                         double v2=gy*gz*v200 + fy*gz*v210 + gy*fz*v201 + 
00220                                                                 fy*fz*v211;
00221                         double f=fx-0.5;
00222                         Ex = (1.0-f)*Ex - f*(v2-v1)/dx;
00223                 }
00224         } else {
00225                 if(ix == 0) {
00226                         // interpolate between center of box and xMin
00227                         double Ex2=approxE(xMin,y,z).x();
00228                         double f=fx*2.0;
00229                         Ex = (1.0-f)*Ex2 + f*Ex;
00230                 } else {
00231                         // interpolate between this box and next lower one
00232                         double v200=array[index(ix-1,iy,iz)];
00233                         double v210=array[index(ix-1,iy+1,iz)];
00234                         double v201=array[index(ix-1,iy,iz+1)];
00235                         double v211=array[index(ix-1,iy+1,iz+1)];
00236                         double v2=gy*gz*v200 + fy*gz*v210 + gy*fz*v201 + 
00237                                                                 fy*fz*v211;
00238                         double f=fx+0.5;
00239                         Ex = -(1.0-f)*(v0-v2)/dx + f*Ex;
00240                 }
00241         }
00242 
00243         // partial d/dy
00244         v0=gx*gz*v000 + fx*gz*v100 + gx*fz*v001 + fx*fz*v101;
00245         v1=gx*gz*v010 + fx*gz*v110 + gx*fz*v011 + fx*fz*v111;
00246         double Ey=-(v1-v0)/dy; // for this grid box
00247         if(fy > 0.5) {
00248                 if(iy+1 == ny-1) {
00249                         // interpolate between center of box and yMax
00250                         double Ey2=approxE(x,yMax,z).y();
00251                         double f=2.0*fy-1.0;
00252                         Ey = (1.0-f)*Ey + f*Ey2;
00253                 } else {
00254                         // interpolate between this box and next higher one
00255                         double v020=array[index(ix,iy+2,iz)];
00256                         double v120=array[index(ix+1,iy+2,iz)];
00257                         double v021=array[index(ix,iy+2,iz+1)];
00258                         double v121=array[index(ix+1,iy+2,iz+1)];
00259                         double v2=gx*gz*v020 + fx*gz*v120 + gx*fz*v021 + 
00260                                                                 fx*fz*v121;
00261                         double f=fy-0.5;
00262                         Ey = (1.0-f)*Ey - f*(v2-v1)/dy;
00263                 }
00264         } else {
00265                 if(iy == 0) {
00266                         // interpolate between center of box and yMin
00267                         double Ey2=approxE(x,yMin,z).y();
00268                         double f=fy*2.0;
00269                         Ey = (1.0-f)*Ey2 + f*Ey;
00270                 } else {
00271                         // interpolate between this box and next lower one
00272                         double v020=array[index(ix,iy-1,iz)];
00273                         double v120=array[index(ix+1,iy-1,iz)];
00274                         double v021=array[index(ix,iy-1,iz+1)];
00275                         double v121=array[index(ix+1,iy-1,iz+1)];
00276                         double v2=gx*gz*v020 + fx*gz*v120 + gx*fz*v021 + 
00277                                                                 fx*fz*v121;
00278                         double f=fy+0.5;
00279                         Ey = -(1.0-f)*(v0-v2)/dy + f*Ey;
00280                 }
00281         }
00282 
00283         // partial d/dz
00284         v0=gx*gy*v000 + fx*gy*v100 + gx*fy*v010 + fx*fy*v110;
00285         v1=gx*gy*v001 + fx*gy*v101 + gx*fy*v011 + fx*fy*v111;
00286         double Ez=-(v1-v0)/dz; // for this grid box
00287         if(fz > 0.5) {
00288                 if(iz+1 == nz-1) {
00289                         // interpolate between center of box and zMax
00290                         double Ez2=approxE(x,y,zMax).z();
00291                         double f=2.0*fz-1.0;
00292                         Ez = (1.0-f)*Ez + f*Ez2;
00293                 } else {
00294                         // interpolate between this box and next higher one
00295                         double v002=array[index(ix,iy,iz+2)];
00296                         double v102=array[index(ix+1,iy,iz+2)];
00297                         double v012=array[index(ix,iy+1,iz+2)];
00298                         double v112=array[index(ix+1,iy+1,iz+2)];
00299                         double v2=gx*gy*v002 + fx*gy*v102 + gx*fy*v012 + 
00300                                                                 fx*fy*v112;
00301                         double f=fz-0.5;
00302                         Ez = (1.0-f)*Ez - f*(v2-v1)/dz;
00303                 }
00304         } else {
00305                 if(iz == 0) {
00306                         // interpolate between center of box and zMin
00307                         double Ez2=approxE(x,y,zMin).z();
00308                         double f=fz*2.0;
00309                         Ez = (1.0-f)*Ez2 + f*Ez;
00310                 } else {
00311                         // interpolate between this box and next lower one
00312                         double v002=array[index(ix,iy,iz-1)];
00313                         double v102=array[index(ix+1,iy,iz-1)];
00314                         double v012=array[index(ix,iy+1,iz-1)];
00315                         double v112=array[index(ix+1,iy+1,iz-1)];
00316                         double v2=gx*gy*v002 + fx*gy*v102 + gx*fy*v012 + 
00317                                                                 fx*fy*v112;
00318                         double f=fz+0.5;
00319                         Ez = -(1.0-f)*(v0-v2)/dz + f*Ez;
00320                 }
00321         }
00322 #endif
00323 
00324         return G4ThreeVector(Ex,Ey,Ez);
00325 }

G4ThreeVector PoissonConvolve3D::approxE ( double  x,
double  y,
double  z 
)

approxE() evaluates the approximation E field. {x,y,z} can be anywhere, but normally are outside the grid or on its boundary.

References approx, dx, dy, dz, E, and M_PI.

Referenced by Bunch::addFieldValue(), and getE().

00328 {
00329         double tol = ((dx<dy?dy:dx)<dz ? (dx<dy?dy:dx) : dz) * 0.05;
00330         G4ThreeVector v(0.0,0.0,0.0);
00331         G4ThreeVector pos(x,y,z);
00332         for(unsigned i=0; i<approx.size(); ++i) {
00333                 G4ThreeVector p = approx[i].getPosition();
00334                 p -= pos; // note - sign
00335                 double r=p.mag();
00336                 if(r < tol) continue; // should never happen
00337                 double denom=4.0*M_PI*r*r*r;
00338                 double c=approx[i].getCharge();
00339                 G4ThreeVector E(c*p.x()/denom,c*p.y()/denom,c*p.z()/denom);
00340                 v -= E; // note - sign above
00341         }
00342         return v;
00343 }

float PoissonConvolve3D::approxPhi ( double  x,
double  y,
double  z 
)

approxPhi() evaluates the approximation phi. {x,y,z} can be anywhere, but normally are outside the grid or on its boundary.

References approx, dx, dy, dz, M_PI, and phi().

00346 {
00347         double tol = ((dx<dy?dy:dx)<dz ? (dx<dy?dy:dx) : dz) * 0.05;
00348         double phi=0.0;
00349         G4ThreeVector pos(x,y,z);
00350         for(unsigned i=0; i<approx.size(); ++i) {
00351                 G4ThreeVector p = approx[i].getPosition();
00352                 p -= pos; // note - sign
00353                 double r=p.mag();
00354                 if(r < tol) r = tol; // should never happen
00355                 double c=approx[i].getCharge();
00356                 phi += c/r;
00357         }
00358         return phi/(4.0*M_PI);
00359 }

G4ThreeVector PoissonConvolve3D::getPosition ( int  ix,
int  iy,
int  iz 
) [inline]

getPosition() returns the position corresponding to grid point indexes. Intended primarily for testing.

References dx, dy, dz, xMin, yMin, and zMin.

00261                 { return G4ThreeVector(xMin+ix*dx,yMin+iy*dy,zMin+iz*dz); }

void PoissonConvolve3D::setApproximation ( std::vector< Charge v  )  [inline]

setApproximation() sets the approximation to use outside the grid, and on its boundary. v is a vector of point charges. Updates the boundary values of phi.

References approx.

Referenced by Bunch::computeField(), and Bunch::testBfield().

00267                 { approx = v; }

void PoissonConvolve3D::defaultApproximation (  ) 

defaultApproximation() calls setApproximation() with the total charge and mean position in rhs (from putCharge()). Updates the boundary values of phi, and the interior if interior=true.

References approx, sumX, sumY, sumZ, and totalCharge.

Referenced by Bunch::testEfield().

00362 {
00363         Charge c(sumX/totalCharge,sumY/totalCharge,sumZ/totalCharge,
00364                                                                 totalCharge);
00365         approx.clear();
00366         approx.push_back(c);
00367 }

void PoissonConvolve3D::computeGreensFunction (  ) 

computeGreensFunction() does just that. Must be called whenever the physical size of the grid is changed.

References array, dx, dy, dz, gfPlan, index(), nx, ny, nz, status, and UNINITIALIZED.

Referenced by updatePosition().

00370 {
00371         // fill in the physical region
00372         // include the 1/N overall normalization
00373         // (Need one more sample than usual, for the symmetry.)
00374         float norm=0.125f/(float)(nx*ny*nz);
00375         for(int ix=0; ix<=nx; ++ix) {
00376                 float x2=ix*dx*ix*dx;
00377                 for(int iy=0; iy<=ny; ++iy) {
00378                         float y2=iy*dy*iy*dy;
00379                         for(int iz=0; iz<=nz; ++iz) {
00380                                 float z2=iz*dz*iz*dz;
00381                                 float r=sqrtf(x2+y2+z2);
00382                                 if(ix==0 && iy==0 && iz==0)
00383                                         r = (dx+dy+dz)/3.0;
00384                                 array[index(ix,iy,iz)] = norm/r;
00385                         }
00386                 }
00387         }
00388 
00389         // now fill in the un-physical region via symmetry
00390         for(int ix=0; ix<2*nx; ++ix) {
00391                 int jx=ix;
00392                 if(jx > nx) jx=2*nx-jx;
00393                 assert(jx>=0 && jx<=nx);
00394                 for(int iy=0; iy<2*ny; ++iy) {
00395                         int jy=iy;
00396                         if(jy > ny) jy=2*ny-jy;
00397                         assert(jy>=0 && jy<=ny);
00398                         for(int iz=0; iz<2*nz; ++iz) {
00399                                 if(ix<=nx && iy<=ny && iz<=nz) continue;
00400                                 int jz=iz;
00401                                 if(jz > nz) jz=2*nz-jz;
00402                                 assert(jz>=0 && jz<=nz);
00403                                 array[index(ix,iy,iz)] = array[index(jx,jy,jz)];
00404                         }
00405                 }
00406         }
00407 
00408         // compute the FFT of the Green's function
00409         fftwf_execute(gfPlan);
00410 
00411         // indicate we clobbered array[]
00412         status = UNINITIALIZED;
00413 }

void PoissonConvolve3D::solve (  ) 

solve() will solve Poisson's equation

References CHARGE, greensFunction, nc, PHI, plan1, plan2, status, and transf.

Referenced by Bunch::computeField(), Bunch::testBfield(), and Bunch::testEfield().

00416 {
00417         assert(status==CHARGE);
00418 
00419         // compute the FFT of array[], which contains the charge,
00420         // into transf[].
00421         fftwf_execute(plan1);
00422 
00423         // multiply by the FFT of the Green's function
00424         for(int i=0; i<nc; ++i) {
00425                 double tmp = transf[i][0]*greensFunction[i][0] -
00426                                         transf[i][1]*greensFunction[i][1];
00427                 transf[i][1] = transf[i][1]*greensFunction[i][0] +
00428                                         transf[i][0]*greensFunction[i][1];
00429                 transf[i][0] = tmp;
00430         }
00431 
00432         // inverse transform transf[] to get phi into array[]
00433         fftwf_execute(plan2); // transf[] is overwritten, but don't care
00434 
00435         // indicate a valid phi in array[]
00436         status = PHI;
00437 }

void PoissonConvolve3D::printRhs (  ) 

printRhs() prints the rhs. Intended primarily for testing.

References array, CHARGE, index(), nx, ny, nz, and status.

00440 {
00441         assert(status==CHARGE);
00442         for(int ix=0; ix<nx; ++ix) {
00443                 for(int iy=0; iy<ny; ++iy) {
00444                         for(int iz=0; iz<nz; ++iz) {
00445                                 int k=index(ix,iy,iz);
00446                                 if(array[k] == 0.0) continue;
00447                                 printf("rhs[%d] = %.3f  ix=%d iy=%d iz=%d\n",
00448                                         k,array[k],ix,iy,iz);
00449                         }
00450                 }
00451         }
00452 }

void PoissonConvolve3D::printPhi (  ) 

printPhi() prints phi. Intended primarily for testing.

References array, index(), nx, ny, nz, PHI, and status.

00455 {
00456         assert(status==PHI);
00457         for(int ix=0; ix<nx; ++ix) {
00458                 for(int iy=0; iy<ny; ++iy) {
00459                         for(int iz=0; iz<nz; ++iz) {
00460                                 int k=index(ix,iy,iz);
00461                                 printf("phi[%d] = %.3f  ix=%d iy=%d iz=%d\n",
00462                                         k,array[k],ix,iy,iz);
00463                         }
00464                 }
00465         }
00466 }

int PoissonConvolve3D::getNx (  )  const [inline]

accessors.

References nx.

00291 { return nx; }

int PoissonConvolve3D::getNy (  )  const [inline]

References ny.

00292 { return ny; }

int PoissonConvolve3D::getNz (  )  const [inline]

References nz.

00293 { return nz; }


Member Data Documentation

int PoissonConvolve3D::nx [protected]

int PoissonConvolve3D::ny [protected]

int PoissonConvolve3D::nz [protected]

float PoissonConvolve3D::xMin [protected]

float PoissonConvolve3D::xMax [protected]

Referenced by getE(), and updatePosition().

float PoissonConvolve3D::yMin [protected]

float PoissonConvolve3D::yMax [protected]

Referenced by getE(), and updatePosition().

float PoissonConvolve3D::zMin [protected]

float PoissonConvolve3D::zMax [protected]

Referenced by getE(), and updatePosition().

float PoissonConvolve3D::dx [protected]

float PoissonConvolve3D::dy [protected]

float PoissonConvolve3D::dz [protected]

float* PoissonConvolve3D::array [protected]

fftwf_complex* PoissonConvolve3D::transf [protected]

fftwf_complex* PoissonConvolve3D::greensFunction [protected]

std::vector<Charge> PoissonConvolve3D::approx [protected]

float PoissonConvolve3D::sumX [protected]

float PoissonConvolve3D::sumY [protected]

float PoissonConvolve3D::sumZ [protected]

float PoissonConvolve3D::totalCharge [protected]

float PoissonConvolve3D::vol [protected]

Referenced by updatePosition().

int PoissonConvolve3D::nc [protected]

Referenced by PoissonConvolve3D(), and solve().

fftwf_plan PoissonConvolve3D::gfPlan [protected]

fftwf_plan PoissonConvolve3D::plan1 [protected]

fftwf_plan PoissonConvolve3D::plan2 [protected]


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