00001 // PoissonConvolve3D.hh 00002 /* 00003 NOTE: fftw must be installed with float libraries. 00004 */ 00005 00006 #ifndef POISSONCONVOLVE3D_HH 00007 #define POISSONCONVOLVE3D_HH 00008 00009 #include <assert.h> 00010 #include <math.h> 00011 #include <vector> 00012 00013 #include "fftw3.h" 00014 00015 #ifdef USE_DOUBLE 00016 #define float double 00017 #define sqrtf sqrt 00018 #define fftwf_complex fftw_complex 00019 #define fftwf_plan fftw_plan 00020 #define fftwf_malloc fftw_malloc 00021 #define fftwf_free fftw_free 00022 #define fftwf_destroy_plan fftw_destroy_plan 00023 #define fftwf_set_timelimit fftw_set_timelimit 00024 #define fftwf_plan_dft_r2c_3d fftw_plan_dft_r2c_3d 00025 #define fftwf_plan_dft_c2r_3d fftw_plan_dft_c2r_3d 00026 #define fftwf_execute fftw_execute 00027 #endif 00028 00029 #ifndef M_PI 00030 #define M_PI 3.14159265358979323846264338327950288 00031 #endif 00032 00033 #ifdef MAIN 00034 ///@ testing - don't require Geant4 headers. 00035 struct G4ThreeVector { 00036 double xx, yy, zz; 00037 G4ThreeVector(double x, double y, double z) { xx=x; yy=y; zz=z; } 00038 double x() { return xx; } 00039 double y() { return yy; } 00040 double z() { return zz; } 00041 double mag() { return sqrt(xx*xx+yy*yy+zz*zz); } 00042 G4ThreeVector& operator += (G4ThreeVector& r) 00043 { xx += r.xx; yy += r.yy; zz += r.zz; return *this; } 00044 G4ThreeVector& operator -= (G4ThreeVector& r) 00045 { xx -= r.xx; yy -= r.yy; zz -= r.zz; return *this; } 00046 G4ThreeVector& operator *= (double v) 00047 { xx *= v; yy*= v; xx*= v; return *this; } 00048 G4ThreeVector& operator /= (double v) 00049 { xx /= v; yy/= v; xx/= v; return *this; } 00050 }; 00051 const double megavolt=1.0, meter=1.0; 00052 #else 00053 #include "G4ThreeVector.hh" 00054 #endif //MAIN 00055 00056 /// Charge represents one point charge for the approximation outside the grid. 00057 struct Charge { 00058 float x, y, z, c; 00059 Charge(float _x, float _y, float _z, float _c) 00060 { x=_x; y=_y; z=_z; c=_c; } 00061 G4ThreeVector getPosition() { return G4ThreeVector(x,y,z); } 00062 float getCharge() { return c; } 00063 }; 00064 00065 /** class PoissonConvolve3D 00066 * 00067 * A solver for Poisson's equation in 3-D Cartesian coordinates {x,y,z}, 00068 * using convolution of the Green's function (whch has infinite boundary 00069 * conditions). The convolution is sped up by using an FFT on a grid 00070 * twice as large (in each dimension) as the solution (physical) grid. 00071 * 00072 * This class uses the Poisson solver on a rectangular 3-d grid; inside 00073 * the grid, it calculates E by applying the analytic derivative of the 00074 * interpolation function to the values of phi at the surrounding 8 00075 * grid points. Outside the grid, by default it uses the analytic value 00076 * of E for a single point charge, with its charge equal to the total 00077 * charge in the grid, and its position at the mean position of the 00078 * charge in the grid. More complicated approximations, based on 00079 * multiple point charges, can be used (not discussed here). 00080 * 00081 * Units are determined by the code using this class. Normally the 00082 * charge values will not be integers, but rather integers times the 00083 * positron charge divided by epsilon_0 (the permittivity of free space), 00084 * the grid sizes will have units of length; the result is that 00085 * phi() and getE() will then have the appropriate units. 00086 * 00087 * This class solves: 00088 * del^2 phi = -4*pi*rho 00089 * On a specified grid in {x,y,z}. rho(x,y,z) is the charge density. 00090 * The -4*pi and the conversion from charge to charge density is handled 00091 * internally in putCharge(); they are NOT handled in setRhs() (which 00092 * is intended primarily for testing). 00093 * 00094 * This class includes an approximation for use outside the grid. The 00095 * default is simply to use the total charge put into the grid (via 00096 * putCharge()), and the weighted mean position of that charge. Better 00097 * approximations can be computed externally as a vector of point 00098 * charges and used via setApproximation(). 00099 * 00100 * Note: array[] contains CHARGE between calls to zeroRhs() and solve(); 00101 * it contains PHI between calls to solve() and zeroRhs(); status 00102 * reflects its contents, ans assert()-s are used to ensure accesses 00103 * are valid.. 00104 * 00105 * 00106 * Accuracy for a 65x65x65 grid of unit size: 00107 * ... MORE TO COME ... 00108 * 00109 **/ 00110 class PoissonConvolve3D { 00111 public: 00112 /// enum PutMode determines how putCharge() places charge into 00113 /// the grid. PRORATED is the default. 00114 /// LOWEST places charge into the grid point lower than (x,y,z). 00115 /// NEAREST places charge into the nearest grid point to (x,y,z). 00116 /// PRORATED places charge into the 8 nearest grid points, pro-rated 00117 /// by proximity to (x,y,z). 00118 /// For 1M particles into a 129^3 grid, NEAREST and LOWER were less 00119 /// than 10% faster than PRORATED, giving slightly more ragged results. 00120 enum PutMode {LOWEST, NEAREST, PRORATED}; 00121 protected: 00122 enum Status {UNINITIALIZED, CHARGE, PHI}; 00123 int nx, ny, nz; 00124 float xMin, xMax, yMin, yMax, zMin, zMax; 00125 Status status; 00126 float dx, dy, dz; // grid spacings 00127 float *array; // phi[index(ix,iy,iz)] or rho[index(ix,iy,iz)] 00128 fftwf_complex *transf; // transform of array (working) 00129 fftwf_complex *greensFunction; 00130 std::vector<Charge> approx; 00131 float sumX, sumY, sumZ, totalCharge; // for defaultApproximation() 00132 float vol; 00133 PutMode putMode; 00134 int nc; 00135 fftwf_plan gfPlan; 00136 fftwf_plan plan1; 00137 fftwf_plan plan2; 00138 int index(int ix, int iy, int iz) const { 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 } 00142 public: 00143 /// Constructor. 00144 /// nx,ny,nz are the # grid points in the 3 dimensions -- NOTE: 00145 /// these are the size of the PHYSICAL grid, the actual grid is 00146 /// doubled in each dimension. They can be any integers >= 9, but 00147 /// powers of 2 have a performance advantage. 00148 /// Note that 512x512x512 allocates too much memory for a 32-bit 00149 /// program on an Intel processor. 00150 /// xMin,xMax,yMin,yMax,zMin,zMax are the grid limits in 3 dimensions. 00151 /// That implies array[index(0,iy,iz)] is the boundary at x=xMin; 00152 /// array[index(nx-1,iy,iz)] is the boundary at x=xMax; etc. 00153 /// 1.0e-4 is appropriate. 00154 PoissonConvolve3D(int _nx, int _ny, int _nz, float _xMin, float _xMax, 00155 float _yMin, float _yMax, float _zMin, float _zMax); 00156 00157 /// Destructor. 00158 virtual ~PoissonConvolve3D() { 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 } 00167 00168 /// updatePosition() updates the position of the boundaries. 00169 /// Usually followed by a call to setApproximation(), with the 00170 /// approximation charge positions updated similarly, or by 00171 /// a series of calls to putCharge() and then defaultApproximation(). 00172 void updatePosition(float _xMin, float _xMax, 00173 float _yMin, float _yMax, float _zMin, float _zMax) { 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 } 00183 00184 /// phi() works for all {x,y,z}, using the solution inside the grid 00185 /// and the approximation outside the grid. 00186 /// Note that phi() is not necessarily continuous or differentiable, 00187 /// but it should be approximately continuous. 00188 /// Do not differentiate this function, use getE() instead (it is 00189 /// MUCH more intelligent about the grid and its bounadries). 00190 /// Intended primarily for testing. 00191 float phi(float x, float y, float z); 00192 00193 /// phiGrid() returns the value of the solution to Poisson's equation 00194 /// (i.e. the electrostatic potential) on grid points. 00195 /// ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 00196 /// Intended primarily for testing. 00197 float phiGrid(int ix, int iy, int iz) const 00198 { assert(status==PHI); return array[index(ix,iy,iz)]; } 00199 00200 /// setPhi() will set a value in array[]. 00201 /// Intended primarily for testing. 00202 void setPhi(int ix, int iy, int iz, float v) 00203 { assert(status==PHI); array[index(ix,iy,iz)] = v; } 00204 00205 /// addPhi() will add v to a value in array[]. 00206 /// Intended primarily for testing. 00207 void addPhi(int ix, int iy, int iz, float v) 00208 { assert(status==PHI); array[index(ix,iy,iz)] += v; } 00209 00210 /// rhs() returns the charge at a grid point of the right-hand-side 00211 /// of Poisson's equation. 00212 /// NOTE: ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 00213 float rhs(int ix, int iy, int iz) const 00214 { assert(status==CHARGE); 00215 return 4.0*M_PI*array[index(ix,iy,iz)]; } 00216 00217 /// setRhs() will set an entry in rhs[]. Normally putCharge() is better. 00218 /// ix,iy,iz MUST be in range: 0..nx-1, 0..ny-1, 0..nz-1 00219 /// This directly sets the rhs array entry, so you probably need 00220 /// to multiply by -4*pi/(dx*dy*dz) to get the proper R.H.S. for 00221 /// Poisson's equation; putCharge() does that internally. 00222 /// Intended primarily for testing. 00223 void setRhs(int ix, int iy, int iz, float v) 00224 { assert(status==CHARGE); array[index(ix,iy,iz)] = v; } 00225 00226 /// putCharge() puts a charge into the rhs grid according to 00227 /// the current putMode. Does nothing if outside the grid. 00228 /// Multiplies by -1.0/(dx*dy*dz) to convert the charge to the R.H.S. 00229 /// of Poisson's equation, which is minus the charge density. 00230 /// NOTE: dividing by epsilon0 is up to the user's choice of units. 00231 /// Accumulates the total charge and mean position of charge actually 00232 /// placed into the grid (for use by defaultApproximation()). 00233 /// Returns true if the charge was placed within the grid, false if not. 00234 bool putCharge(float x, float y, float z, float c); 00235 bool putCharge(Charge &c) { return putCharge(c.x,c.y,c.z,c.c); } 00236 00237 /// setPutMode() sets the mode of putCharge(). See PutMode above. 00238 void setPutMode(PutMode m) { putMode = m; } 00239 00240 /// zeroRhs() will zero rhs. 00241 void zeroRhs(); 00242 00243 /// getE() returns the value of the E field at a given point. 00244 /// Uses the solution inside the grid and the approximation outside. 00245 /// It handles grid cells and boundaries intelligently. 00246 G4ThreeVector getE(float x, float y, float z); 00247 00248 /// approxE() evaluates the approximation E field. 00249 /// {x,y,z} can be anywhere, but normally are outside the grid or on 00250 /// its boundary. 00251 G4ThreeVector approxE(double x, double y, double z); 00252 00253 /// approxPhi() evaluates the approximation phi. 00254 /// {x,y,z} can be anywhere, but normally are outside the grid or on 00255 /// its boundary. 00256 float approxPhi(double x, double y, double z); 00257 00258 /// getPosition() returns the position corresponding to grid point 00259 /// indexes. Intended primarily for testing. 00260 G4ThreeVector getPosition(int ix, int iy, int iz) 00261 { return G4ThreeVector(xMin+ix*dx,yMin+iy*dy,zMin+iz*dz); } 00262 00263 /// setApproximation() sets the approximation to use outside the 00264 /// grid, and on its boundary. v is a vector of point charges. 00265 /// Updates the boundary values of phi. 00266 void setApproximation(std::vector<Charge> v) 00267 { approx = v; } 00268 00269 /// defaultApproximation() calls setApproximation() with the 00270 /// total charge and mean position in rhs (from putCharge()). 00271 /// Updates the boundary values of phi, and the interior if 00272 /// interior=true. 00273 void defaultApproximation(); 00274 00275 /// computeGreensFunction() does just that. Must be called whenever 00276 /// the physical size of the grid is changed. 00277 void computeGreensFunction(); 00278 00279 /// solve() will solve Poisson's equation 00280 void solve(); 00281 00282 /// printRhs() prints the rhs. 00283 /// Intended primarily for testing. 00284 void printRhs(); 00285 00286 /// printPhi() prints phi. 00287 /// Intended primarily for testing. 00288 void printPhi(); 00289 00290 /// accessors. 00291 int getNx() const { return nx; } 00292 int getNy() const { return ny; } 00293 int getNz() const { return nz; } 00294 }; 00295 00296 #endif // POISSONCONVOLVE3D_HH