Voro++
container_prd.cc
Go to the documentation of this file.
00001 // Voro++, a 3D cell-based Voronoi library
00002 //
00003 // Author   : Chris H. Rycroft (LBL / UC Berkeley)
00004 // Email    : chr@alum.mit.edu
00005 // Date     : August 30th 2011
00006 
00007 /** \file container_prd.cc
00008  * \brief Function implementations for the container_periodic_base and
00009  * related classes. */
00010 
00011 #include "container_prd.hh"
00012 
00013 namespace voro {
00014 
00015 /** The class constructor sets up the geometry of container, initializing the
00016  * minimum and maximum coordinates in each direction, and setting whether each
00017  * direction is periodic or not. It divides the container into a rectangular
00018  * grid of blocks, and allocates memory for each of these for storing particle
00019  * positions and IDs.
00020  * \param[in] (bx_) The x coordinate of the first unit vector.
00021  * \param[in] (bxy_,by_) The x and y coordinates of the second unit vector.
00022  * \param[in] (bxz_,byz_,bz_) The x, y, and z coordinates of the third unit
00023  *                            vector.
00024  * \param[in] (nx_,ny_,nz_) the number of grid blocks in each of the three
00025  *                       coordinate directions.
00026  * \param[in] init_mem_ the initial memory allocation for each block.
00027  * \param[in] ps_ the number of floating point entries to store for each
00028  *                particle. */
00029 container_periodic_base::container_periodic_base(double bx_,double bxy_,double by_,
00030                 double bxz_,double byz_,double bz_,int nx_,int ny_,int nz_,int init_mem_,int ps_)
00031         : unitcell(bx_,bxy_,by_,bxz_,byz_,bz_), voro_base(nx_,ny_,nz_,bx_/nx_,by_/ny_,bz_/nz_),
00032         ey(int(max_uv_y*ysp+1)), ez(int(max_uv_z*zsp+1)), wy(ny+ey), wz(nz+ez),
00033         oy(ny+2*ey), oz(nz+2*ez), oxyz(nx*oy*oz), id(new int*[oxyz]), p(new double*[oxyz]),
00034         co(new int[oxyz]), mem(new int[oxyz]), img(new char[oxyz]), init_mem(init_mem_), ps(ps_) {
00035         int i,j,k,l;
00036 
00037         // Clear the global arrays
00038         int *pp=co;while(pp<co+oxyz) *(pp++)=0;
00039         pp=mem;while(pp<mem+oxyz) *(pp++)=0;
00040         char *cp=img;while(cp<img+oxyz) *(cp++)=0;
00041 
00042         // Set up memory for the blocks in the primary domain
00043         for(k=ez;k<wz;k++) for(j=ey;j<wy;j++) for(i=0;i<nx;i++) {
00044                 l=i+nx*(j+oy*k);
00045                 mem[l]=init_mem;
00046                 id[l]=new int[init_mem];
00047                 p[l]=new double[ps*init_mem];
00048         }
00049 }
00050 
00051 /** The container destructor frees the dynamically allocated memory. */
00052 container_periodic_base::~container_periodic_base() {
00053         int l;
00054         for(l=0;l<oxyz;l++) if(mem[l]>0) delete [] p[l];
00055         for(l=0;l<oxyz;l++) if(mem[l]>0) delete [] id[l];
00056         delete [] p;
00057         delete [] id;
00058         delete [] mem;
00059         delete [] co;
00060 }
00061 
00062 /** The class constructor sets up the geometry of container.
00063  * \param[in] (bx_) The x coordinate of the first unit vector.
00064  * \param[in] (bxy_,by_) The x and y coordinates of the second unit vector.
00065  * \param[in] (bxz_,byz_,bz_) The x, y, and z coordinates of the third unit
00066  *                            vector.
00067  * \param[in] (nx_,ny_,nz_) the number of grid blocks in each of the three
00068  *                          coordinate directions.
00069  * \param[in] init_mem_ the initial memory allocation for each block. */
00070 container_periodic::container_periodic(double bx_,double bxy_,double by_,double bxz_,double byz_,double bz_,
00071         int nx_,int ny_,int nz_,int init_mem_)
00072         : container_periodic_base(bx_,bxy_,by_,bxz_,byz_,bz_,nx_,ny_,nz_,init_mem_,3),
00073         vc(*this,2*nx_+1,2*ey+1,2*ez+1) {}
00074 
00075 /** The class constructor sets up the geometry of container.
00076  * \param[in] (bx_) The x coordinate of the first unit vector.
00077  * \param[in] (bxy_,by_) The x and y coordinates of the second unit vector.
00078  * \param[in] (bxz_,byz_,bz_) The x, y, and z coordinates of the third unit
00079  *                            vector.
00080  * \param[in] (nx_,ny_,nz_) the number of grid blocks in each of the three
00081  *                          coordinate directions.
00082  * \param[in] init_mem_ the initial memory allocation for each block. */
00083 container_periodic_poly::container_periodic_poly(double bx_,double bxy_,double by_,double bxz_,double byz_,double bz_,
00084         int nx_,int ny_,int nz_,int init_mem_)
00085         : container_periodic_base(bx_,bxy_,by_,bxz_,byz_,bz_,nx_,ny_,nz_,init_mem_,4),
00086         max_radius(0), vc(*this,2*nx_+1,2*ey+1,2*ez+1) {}
00087 
00088 /** Put a particle into the correct region of the container.
00089  * \param[in] n the numerical ID of the inserted particle.
00090  * \param[in] (x,y,z) the position vector of the inserted particle. */
00091 void container_periodic::put(int n,double x,double y,double z) {
00092         int ijk;
00093         put_locate_block(ijk,x,y,z);
00094         id[ijk][co[ijk]]=n;
00095         double *pp=p[ijk]+3*co[ijk]++;
00096         *(pp++)=x;*(pp++)=y;*pp=z;
00097 }
00098 
00099 /** Put a particle into the correct region of the container.
00100  * \param[in] n the numerical ID of the inserted particle.
00101  * \param[in] (x,y,z) the position vector of the inserted particle.
00102  * \param[in] r the radius of the particle. */
00103 void container_periodic_poly::put(int n,double x,double y,double z,double r) {
00104         int ijk;
00105         put_locate_block(ijk,x,y,z);
00106         id[ijk][co[ijk]]=n;
00107         double *pp=p[ijk]+4*co[ijk]++;
00108         *(pp++)=x;*(pp++)=y;*(pp++)=z;*pp=r;
00109         if(max_radius<r) max_radius=r;
00110 }
00111 
00112 /** Put a particle into the correct region of the container.
00113  * \param[in] n the numerical ID of the inserted particle.
00114  * \param[in] (x,y,z) the position vector of the inserted particle.
00115  * \param[out] (ai,aj,ak) the periodic image displacement that the particle is
00116  *                        in, with (0,0,0) corresponding to the primary domain.
00117  */
00118 void container_periodic::put(int n,double x,double y,double z,int &ai,int &aj,int &ak) {
00119         int ijk;
00120         put_locate_block(ijk,x,y,z,ai,aj,ak);
00121         id[ijk][co[ijk]]=n;
00122         double *pp=p[ijk]+3*co[ijk]++;
00123         *(pp++)=x;*(pp++)=y;*pp=z;
00124 }
00125 
00126 /** Put a particle into the correct region of the container.
00127  * \param[in] n the numerical ID of the inserted particle.
00128  * \param[in] (x,y,z) the position vector of the inserted particle.
00129  * \param[in] r the radius of the particle.
00130  * \param[out] (ai,aj,ak) the periodic image displacement that the particle is
00131  *                        in, with (0,0,0) corresponding to the primary domain.
00132  */
00133 void container_periodic_poly::put(int n,double x,double y,double z,double r,int &ai,int &aj,int &ak) {
00134         int ijk;
00135         put_locate_block(ijk,x,y,z,ai,aj,ak);
00136         id[ijk][co[ijk]]=n;
00137         double *pp=p[ijk]+4*co[ijk]++;
00138         *(pp++)=x;*(pp++)=y;*(pp++)=z;*pp=r;
00139         if(max_radius<r) max_radius=r;
00140 }
00141 
00142 /** Put a particle into the correct region of the container, also recording
00143  * into which region it was stored.
00144  * \param[in] vo the ordering class in which to record the region.
00145  * \param[in] n the numerical ID of the inserted particle.
00146  * \param[in] (x,y,z) the position vector of the inserted particle. */
00147 void container_periodic::put(particle_order &vo,int n,double x,double y,double z) {
00148         int ijk;
00149         put_locate_block(ijk,x,y,z);
00150         id[ijk][co[ijk]]=n;
00151         vo.add(ijk,co[ijk]);
00152         double *pp=p[ijk]+3*co[ijk]++;
00153         *(pp++)=x;*(pp++)=y;*pp=z;
00154 }
00155 
00156 /** Put a particle into the correct region of the container, also recording
00157  * into which region it was stored.
00158  * \param[in] vo the ordering class in which to record the region.
00159  * \param[in] n the numerical ID of the inserted particle.
00160  * \param[in] (x,y,z) the position vector of the inserted particle.
00161  * \param[in] r the radius of the particle. */
00162 void container_periodic_poly::put(particle_order &vo,int n,double x,double y,double z,double r) {
00163         int ijk;
00164         put_locate_block(ijk,x,y,z);
00165         id[ijk][co[ijk]]=n;
00166         vo.add(ijk,co[ijk]);
00167         double *pp=p[ijk]+4*co[ijk]++;
00168         *(pp++)=x;*(pp++)=y;*(pp++)=z;*pp=r;
00169         if(max_radius<r) max_radius=r;
00170 }
00171 
00172 /** Takes a particle position vector and computes the region index into which
00173  * it should be stored. If the container is periodic, then the routine also
00174  * maps the particle position to ensure it is in the primary domain. If the
00175  * container is not periodic, the routine bails out.
00176  * \param[out] ijk the region index.
00177  * \param[in,out] (x,y,z) the particle position, remapped into the primary
00178  *                        domain if necessary.
00179  * \return True if the particle can be successfully placed into the container,
00180  * false otherwise. */
00181 void container_periodic_base::put_locate_block(int &ijk,double &x,double &y,double &z) {
00182 
00183         // Remap particle in the z direction if necessary
00184         int k=step_int(z*zsp);
00185         if(k<0||k>=nz) {
00186                 int ak=step_div(k,nz);
00187                 z-=ak*bz;y-=ak*byz;x-=ak*bxz;k-=ak*nz;
00188         }
00189 
00190         // Remap particle in the y direction if necessary
00191         int j=step_int(y*ysp);
00192         if(j<0||j>=ny) {
00193                 int aj=step_div(j,ny);
00194                 y-=aj*by;x-=aj*bxy;j-=aj*ny;
00195         }
00196 
00197         // Remap particle in the x direction if necessary
00198         ijk=step_int(x*xsp);
00199         if(ijk<0||ijk>=nx) {
00200                 int ai=step_div(ijk,nx);
00201                 x-=ai*bx;ijk-=ai*nx;
00202         }
00203 
00204         // Compute the block index and check memory allocation
00205         j+=ey;k+=ez;
00206         ijk+=nx*(j+oy*k);
00207         if(co[ijk]==mem[ijk]) add_particle_memory(ijk);
00208 }
00209 
00210 /** Takes a particle position vector and computes the region index into which
00211  * it should be stored. If the container is periodic, then the routine also
00212  * maps the particle position to ensure it is in the primary domain. If the
00213  * container is not periodic, the routine bails out.
00214  * \param[out] ijk the region index.
00215  * \param[in,out] (x,y,z) the particle position, remapped into the primary
00216  *                        domain if necessary.
00217  * \param[out] (ai,aj,ak) the periodic image displacement that the particle is
00218  *                        in, with (0,0,0) corresponding to the primary domain.
00219  * \return True if the particle can be successfully placed into the container,
00220  * false otherwise. */
00221 void container_periodic_base::put_locate_block(int &ijk,double &x,double &y,double &z,int &ai,int &aj,int &ak) {
00222 
00223         // Remap particle in the z direction if necessary
00224         int k=step_int(z*zsp);
00225         if(k<0||k>=nz) {
00226                 ak=step_div(k,nz);
00227                 z-=ak*bz;y-=ak*byz;x-=ak*bxz;k-=ak*nz;
00228         } else ak=0;
00229 
00230         // Remap particle in the y direction if necessary
00231         int j=step_int(y*ysp);
00232         if(j<0||j>=ny) {
00233                 aj=step_div(j,ny);
00234                 y-=aj*by;x-=aj*bxy;j-=aj*ny;
00235         } else aj=0;
00236 
00237         // Remap particle in the x direction if necessary
00238         ijk=step_int(x*xsp);
00239         if(ijk<0||ijk>=nx) {
00240                 ai=step_div(ijk,nx);
00241                 x-=ai*bx;ijk-=ai*nx;
00242         } else ai=0;
00243 
00244         // Compute the block index and check memory allocation
00245         j+=ey;k+=ez;
00246         ijk+=nx*(j+oy*k);
00247         if(co[ijk]==mem[ijk]) add_particle_memory(ijk);
00248 }
00249 
00250 /** Takes a position vector and remaps it into the primary domain.
00251  * \param[out] (ai,aj,ak) the periodic image displacement that the vector is in,
00252  *                        with (0,0,0) corresponding to the primary domain.
00253  * \param[out] (ci,cj,ck) the index of the block that the position vector is
00254  *                        within, once it has been remapped.
00255  * \param[in,out] (x,y,z) the position vector to consider, which is remapped
00256  *                        into the primary domain during the routine.
00257  * \param[out] ijk the block index that the vector is within. */
00258 inline void container_periodic_base::remap(int &ai,int &aj,int &ak,int &ci,int &cj,int &ck,double &x,double &y,double &z,int &ijk) {
00259 
00260         // Remap particle in the z direction if necessary
00261         ck=step_int(z*zsp);
00262         if(ck<0||ck>=nz) {
00263                 ak=step_div(ck,nz);
00264                 z-=ak*bz;y-=ak*byz;x-=ak*bxz;ck-=ak*nz;
00265         } else ak=0;
00266 
00267         // Remap particle in the y direction if necessary
00268         cj=step_int(y*ysp);
00269         if(cj<0||cj>=ny) {
00270                 aj=step_div(cj,ny);
00271                 y-=aj*by;x-=aj*bxy;cj-=aj*ny;
00272         } else aj=0;
00273 
00274         // Remap particle in the x direction if necessary
00275         ci=step_int(x*xsp);
00276         if(ci<0||ci>=nx) {
00277                 ai=step_div(ci,nx);
00278                 x-=ai*bx;ci-=ai*nx;
00279         } else ai=0;
00280 
00281         cj+=ey;ck+=ez;
00282         ijk=ci+nx*(cj+oy*ck);
00283 }
00284 
00285 /** Takes a vector and finds the particle whose Voronoi cell contains that
00286  * vector. This is equivalent to finding the particle which is nearest to the
00287  * vector.
00288  * \param[in] (x,y,z) the vector to test.
00289  * \param[out] (rx,ry,rz) the position of the particle whose Voronoi cell
00290  *                        contains the vector. This may point to a particle in
00291  *                        a periodic image of the primary domain.
00292  * \param[out] pid the ID of the particle.
00293  * \return True if a particle was found. If the container has no particles,
00294  * then the search will not find a Voronoi cell and false is returned. */
00295 bool container_periodic::find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid) {
00296         int ai,aj,ak,ci,cj,ck,ijk;
00297         particle_record w;
00298         double mrs;
00299 
00300         // Remap the vector into the primary domain and then search for the
00301         // Voronoi cell that it is within
00302         remap(ai,aj,ak,ci,cj,ck,x,y,z,ijk);
00303         vc.find_voronoi_cell(x,y,z,ci,cj,ck,ijk,w,mrs);
00304 
00305         if(w.ijk!=-1) {
00306 
00307                 // Assemble the position vector of the particle to be returned,
00308                 // applying a periodic remapping if necessary
00309                 ci+=w.di;if(ci<0||ci>=nx) ai+=step_div(ci,nx);
00310                 rx=p[w.ijk][3*w.l]+ak*bxz+aj*bxy+ai*bx;
00311                 ry=p[w.ijk][3*w.l+1]+ak*byz+aj*by;
00312                 rz=p[w.ijk][3*w.l+2]+ak*bz;
00313                 pid=id[w.ijk][w.l];
00314                 return true;
00315         }
00316         return false;
00317 }
00318 
00319 /** Takes a vector and finds the particle whose Voronoi cell contains that
00320  * vector. Additional wall classes are not considered by this routine.
00321  * \param[in] (x,y,z) the vector to test.
00322  * \param[out] (rx,ry,rz) the position of the particle whose Voronoi cell
00323  *                        contains the vector. If the container is periodic,
00324  *                        this may point to a particle in a periodic image of
00325  *                        the primary domain.
00326  * \param[out] pid the ID of the particle.
00327  * \return True if a particle was found. If the container has no particles,
00328  * then the search will not find a Voronoi cell and false is returned. */
00329 bool container_periodic_poly::find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid) {
00330         int ai,aj,ak,ci,cj,ck,ijk;
00331         particle_record w;
00332         double mrs;
00333 
00334         // Remap the vector into the primary domain and then search for the
00335         // Voronoi cell that it is within
00336         remap(ai,aj,ak,ci,cj,ck,x,y,z,ijk);
00337         vc.find_voronoi_cell(x,y,z,ci,cj,ck,ijk,w,mrs);
00338 
00339         if(w.ijk!=-1) {
00340 
00341                 // Assemble the position vector of the particle to be returned,
00342                 // applying a periodic remapping if necessary
00343                 ci+=w.di;if(ci<0||ci>=nx) ai+=step_div(ci,nx);
00344                 rx=p[w.ijk][4*w.l]+ak*bxz+aj*bxy+ai*bx;
00345                 ry=p[w.ijk][4*w.l+1]+ak*byz+aj*by;
00346                 rz=p[w.ijk][4*w.l+2]+ak*bz;
00347                 pid=id[w.ijk][w.l];
00348                 return true;
00349         }
00350         return false;
00351 }
00352 
00353 /** Increase memory for a particular region.
00354  * \param[in] i the index of the region to reallocate. */
00355 void container_periodic_base::add_particle_memory(int i) {
00356 
00357         // Handle the case when no memory has been allocated for this block
00358         if(mem[i]==0) {
00359                 mem[i]=init_mem;
00360                 id[i]=new int[init_mem];
00361                 p[i]=new double[ps*init_mem];
00362                 return;
00363         }
00364 
00365         // Otherwise, double the memory allocation for this block. Carry out a
00366         // check on the memory allocation size, and print a status message if
00367         // requested.
00368         int l,nmem(mem[i]<<1);
00369         if(nmem>max_particle_memory)
00370                 voro_fatal_error("Absolute maximum memory allocation exceeded",VOROPP_MEMORY_ERROR);
00371 #if VOROPP_VERBOSE >=3
00372         fprintf(stderr,"Particle memory in region %d scaled up to %d\n",i,nmem);
00373 #endif
00374 
00375         // Allocate new memory and copy in the contents of the old arrays
00376         int *idp=new int[nmem];
00377         for(l=0;l<co[i];l++) idp[l]=id[i][l];
00378         double *pp=new double[ps*nmem];
00379         for(l=0;l<ps*co[i];l++) pp[l]=p[i][l];
00380 
00381         // Update pointers and delete old arrays
00382         mem[i]=nmem;
00383         delete [] id[i];id[i]=idp;
00384         delete [] p[i];p[i]=pp;
00385 }
00386 
00387 /** Import a list of particles from an open file stream into the container.
00388  * Entries of four numbers (Particle ID, x position, y position, z position)
00389  * are searched for. If the file cannot be successfully read, then the routine
00390  * causes a fatal error.
00391  * \param[in] fp the file handle to read from. */
00392 void container_periodic::import(FILE *fp) {
00393         int i,j;
00394         double x,y,z;
00395         while((j=fscanf(fp,"%d %lg %lg %lg",&i,&x,&y,&z))==4) put(i,x,y,z);
00396         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00397 }
00398 
00399 /** Import a list of particles from an open file stream, also storing the order
00400  * of that the particles are read. Entries of four numbers (Particle ID, x
00401  * position, y position, z position) are searched for. If the file cannot be
00402  * successfully read, then the routine causes a fatal error.
00403  * \param[in,out] vo a reference to an ordering class to use.
00404  * \param[in] fp the file handle to read from. */
00405 void container_periodic::import(particle_order &vo,FILE *fp) {
00406         int i,j;
00407         double x,y,z;
00408         while((j=fscanf(fp,"%d %lg %lg %lg",&i,&x,&y,&z))==4) put(vo,i,x,y,z);
00409         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00410 }
00411 
00412 /** Import a list of particles from an open file stream into the container.
00413  * Entries of five numbers (Particle ID, x position, y position, z position,
00414  * radius) are searched for. If the file cannot be successfully read, then the
00415  * routine causes a fatal error.
00416  * \param[in] fp the file handle to read from. */
00417 void container_periodic_poly::import(FILE *fp) {
00418         int i,j;
00419         double x,y,z,r;
00420         while((j=fscanf(fp,"%d %lg %lg %lg %lg",&i,&x,&y,&z,&r))==5) put(i,x,y,z,r);
00421         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00422 }
00423 
00424 /** Import a list of particles from an open file stream, also storing the order
00425  * of that the particles are read. Entries of four numbers (Particle ID, x
00426  * position, y position, z position, radius) are searched for. If the file
00427  * cannot be successfully read, then the routine causes a fatal error.
00428  * \param[in,out] vo a reference to an ordering class to use.
00429  * \param[in] fp the file handle to read from. */
00430 void container_periodic_poly::import(particle_order &vo,FILE *fp) {
00431         int i,j;
00432         double x,y,z,r;
00433         while((j=fscanf(fp,"%d %lg %lg %lg %lg",&i,&x,&y,&z,&r))==5) put(vo,i,x,y,z,r);
00434         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00435 }
00436 
00437 /** Outputs the a list of all the container regions along with the number of
00438  * particles stored within each. */
00439 void container_periodic_base::region_count() {
00440         int i,j,k,*cop=co;
00441         for(k=0;k<nz;k++) for(j=0;j<ny;j++) for(i=0;i<nx;i++)
00442                 printf("Region (%d,%d,%d): %d particles\n",i,j,k,*(cop++));
00443 }
00444 
00445 /** Clears a container of particles. */
00446 void container_periodic::clear() {
00447         for(int *cop=co;cop<co+nxyz;cop++) *cop=0;
00448 }
00449 
00450 /** Clears a container of particles, also clearing resetting the maximum radius
00451  * to zero. */
00452 void container_periodic_poly::clear() {
00453         for(int *cop=co;cop<co+nxyz;cop++) *cop=0;
00454         max_radius=0;
00455 }
00456 
00457 /** Computes all the Voronoi cells and saves customized information about them.
00458  * \param[in] format the custom output string to use.
00459  * \param[in] fp a file handle to write to. */
00460 void container_periodic::print_custom(const char *format,FILE *fp) {
00461         c_loop_all vl(*this);
00462         print_custom(vl,format,fp);
00463 }
00464 
00465 /** Computes all the Voronoi cells and saves customized
00466  * information about them.
00467  * \param[in] format the custom output string to use.
00468  * \param[in] fp a file handle to write to. */
00469 void container_periodic_poly::print_custom(const char *format,FILE *fp) {
00470         c_loop_all vl(*this);
00471         print_custom(vl,format,fp);
00472 }
00473 
00474 /** Computes all the Voronoi cells and saves customized information about them.
00475  * \param[in] format the custom output string to use.
00476  * \param[in] filename the name of the file to write to. */
00477 void container_periodic::print_custom(const char *format,const char *filename) {
00478         FILE *fp=safe_fopen(filename,"w");
00479         print_custom(format,fp);
00480         fclose(fp);
00481 }
00482 
00483 /** Computes all the Voronoi cells and saves customized
00484  * information about them
00485  * \param[in] format the custom output string to use.
00486  * \param[in] filename the name of the file to write to. */
00487 void container_periodic_poly::print_custom(const char *format,const char *filename) {
00488         FILE *fp=safe_fopen(filename,"w");
00489         print_custom(format,fp);
00490         fclose(fp);
00491 }
00492 
00493 /** Computes all of the Voronoi cells in the container, but does nothing
00494  * with the output. It is useful for measuring the pure computation time
00495  * of the Voronoi algorithm, without any additional calculations such as
00496  * volume evaluation or cell output. */
00497 void container_periodic::compute_all_cells() {
00498         voronoicell c;
00499         c_loop_all vl(*this);
00500         if(vl.start()) do compute_cell(c,vl);
00501         while(vl.inc());
00502 }
00503 
00504 /** Computes all of the Voronoi cells in the container, but does nothing
00505  * with the output. It is useful for measuring the pure computation time
00506  * of the Voronoi algorithm, without any additional calculations such as
00507  * volume evaluation or cell output. */
00508 void container_periodic_poly::compute_all_cells() {
00509         voronoicell c;
00510         c_loop_all vl(*this);
00511         if(vl.start()) do compute_cell(c,vl);while(vl.inc());
00512 }
00513 
00514 /** Calculates all of the Voronoi cells and sums their volumes. In most cases
00515  * without walls, the sum of the Voronoi cell volumes should equal the volume
00516  * of the container to numerical precision.
00517  * \return The sum of all of the computed Voronoi volumes. */
00518 double container_periodic::sum_cell_volumes() {
00519         voronoicell c;
00520         double vol=0;
00521         c_loop_all vl(*this);
00522         if(vl.start()) do if(compute_cell(c,vl)) vol+=c.volume();while(vl.inc());
00523         return vol;
00524 }
00525 
00526 /** Calculates all of the Voronoi cells and sums their volumes. In most cases
00527  * without walls, the sum of the Voronoi cell volumes should equal the volume
00528  * of the container to numerical precision.
00529  * \return The sum of all of the computed Voronoi volumes. */
00530 double container_periodic_poly::sum_cell_volumes() {
00531         voronoicell c;
00532         double vol=0;
00533         c_loop_all vl(*this);
00534         if(vl.start()) do if(compute_cell(c,vl)) vol+=c.volume();while(vl.inc());
00535         return vol;
00536 }
00537 
00538 /** This routine creates all periodic images of the particles. It is meant for
00539  * diagnostic purposes only, since usually periodic images are dynamically
00540  * created in when they are referenced. */
00541 void container_periodic_base::create_all_images() {
00542         int i,j,k;
00543         for(k=0;k<oz;k++) for(j=0;j<oy;j++) for(i=0;i<nx;i++) create_periodic_image(i,j,k);
00544 }
00545 
00546 /** Checks that the particles within each block lie within that block's bounds.
00547  * This is useful for diagnosing problems with periodic image computation. */
00548 void container_periodic_base::check_compartmentalized() {
00549         int c,l,i,j,k;
00550         double mix,miy,miz,max,may,maz,*pp;
00551         for(k=l=0;k<oz;k++) for(j=0;j<oy;j++) for(i=0;i<nx;i++,l++) if(mem[l]>0) {
00552 
00553                 // Compute the block's bounds, adding in a small tolerance
00554                 mix=i*boxx-tolerance;max=mix+boxx+tolerance;
00555                 miy=(j-ey)*boxy-tolerance;may=miy+boxy+tolerance;
00556                 miz=(k-ez)*boxz-tolerance;maz=miz+boxz+tolerance;
00557 
00558                 // Print entries for any particles that lie outside the block's
00559                 // bounds
00560                 for(pp=p[l],c=0;c<co[l];c++,pp+=ps) if(*pp<mix||*pp>max||pp[1]<miy||pp[1]>may||pp[2]<miz||pp[2]>maz)
00561                         printf("%d %d %d %d %f %f %f %f %f %f %f %f %f\n",
00562                                id[l][c],i,j,k,*pp,pp[1],pp[2],mix,max,miy,may,miz,maz);
00563         }
00564 }
00565 
00566 /** Creates particles within an image block that is aligned with the primary
00567  * domain in the z axis. In this case, the image block may be comprised of
00568  * particles from two primary blocks. The routine considers these two primary
00569  * blocks, and adds the needed particles to the image. The remaining particles
00570  * from the primary blocks are also filled into the neighboring images.
00571  * \param[in] (di,dj,dk) the index of the block to consider. The z index must
00572  *                       satisfy ez<=dk<wz. */
00573 void container_periodic_base::create_side_image(int di,int dj,int dk) {
00574         int l,dijk=di+nx*(dj+oy*dk),odijk,ima=step_div(dj-ey,ny);
00575         int qua=di+step_int(-ima*bxy*xsp),quadiv=step_div(qua,nx);
00576         int fi=qua-quadiv*nx,fijk=fi+nx*(dj-ima*ny+oy*dk);
00577         double dis=ima*bxy+quadiv*bx,switchx=di*boxx-ima*bxy-quadiv*bx,adis;
00578 
00579         // Left image computation
00580         if((img[dijk]&1)==0) {
00581                 if(di>0) {
00582                         odijk=dijk-1;adis=dis;
00583                 } else {
00584                         odijk=dijk+nx-1;adis=dis+bx;
00585                 }
00586                 img[odijk]|=2;
00587                 for(l=0;l<co[fijk];l++) {
00588                         if(p[fijk][ps*l]>switchx) put_image(dijk,fijk,l,dis,by*ima,0);
00589                         else put_image(odijk,fijk,l,adis,by*ima,0);
00590                 }
00591         }
00592 
00593         // Right image computation
00594         if((img[dijk]&2)==0) {
00595                 if(fi==nx-1) {
00596                         fijk+=1-nx;switchx+=(1-nx)*boxx;dis+=bx;
00597                 } else {
00598                         fijk++;switchx+=boxx;
00599                 }
00600                 if(di==nx-1) {
00601                         odijk=dijk-nx+1;adis=dis-bx;
00602                 } else {
00603                         odijk=dijk+1;adis=dis;
00604                 }
00605                 img[odijk]|=1;
00606                 for(l=0;l<co[fijk];l++) {
00607                         if(p[fijk][ps*l]<switchx) put_image(dijk,fijk,l,dis,by*ima,0);
00608                         else put_image(odijk,fijk,l,adis,by*ima,0);
00609                 }
00610         }
00611 
00612         // All contributions to the block now added, so set both two bits of
00613         // the image information
00614         img[dijk]=3;
00615 }
00616 
00617 /** Creates particles within an image block that is not aligned with the
00618  * primary domain in the z axis. In this case, the image block may be comprised
00619  * of particles from four primary blocks. The routine considers these four
00620  * primary blocks, and adds the needed particles to the image. The remaining
00621  * particles from the primary blocks are also filled into the neighboring
00622  * images.
00623  * \param[in] (di,dj,dk) the index of the block to consider. The z index must
00624  *                       satisfy dk<ez or dk>=wz. */
00625 void container_periodic_base::create_vertical_image(int di,int dj,int dk) {
00626         int l,dijk=di+nx*(dj+oy*dk),dijkl,dijkr,ima=step_div(dk-ez,nz);
00627         int qj=dj+step_int(-ima*byz*ysp),qjdiv=step_div(qj-ey,ny);
00628         int qi=di+step_int((-ima*bxz-qjdiv*bxy)*xsp),qidiv=step_div(qi,nx);
00629         int fi=qi-qidiv*nx,fj=qj-qjdiv*ny,fijk=fi+nx*(fj+oy*(dk-ima*nz)),fijk2;
00630         double disy=ima*byz+qjdiv*by,switchy=(dj-ey)*boxy-ima*byz-qjdiv*by;
00631         double disx=ima*bxz+qjdiv*bxy+qidiv*bx,switchx=di*boxx-ima*bxz-qjdiv*bxy-qidiv*bx;
00632         double switchx2,disxl,disxr,disx2,disxr2;
00633 
00634         if(di==0) {dijkl=dijk+nx-1;disxl=disx+bx;}
00635         else {dijkl=dijk-1;disxl=disx;}
00636 
00637         if(di==nx-1) {dijkr=dijk-nx+1;disxr=disx-bx;}
00638         else {dijkr=dijk+1;disxr=disx;}
00639 
00640         // Down-left image computation
00641         bool y_exist=dj!=0;
00642         if((img[dijk]&1)==0) {
00643                 img[dijkl]|=2;
00644                 if(y_exist) {
00645                         img[dijkl-nx]|=8;
00646                         img[dijk-nx]|=4;
00647                 }
00648                 for(l=0;l<co[fijk];l++) {
00649                         if(p[fijk][ps*l+1]>switchy) {
00650                                 if(p[fijk][ps*l]>switchx) put_image(dijk,fijk,l,disx,disy,bz*ima);
00651                                 else put_image(dijkl,fijk,l,disxl,disy,bz*ima);
00652                         } else {
00653                                 if(!y_exist) continue;
00654                                 if(p[fijk][ps*l]>switchx) put_image(dijk-nx,fijk,l,disx,disy,bz*ima);
00655                                 else put_image(dijkl-nx,fijk,l,disxl,disy,bz*ima);
00656                         }
00657                 }
00658         }
00659 
00660         // Down-right image computation
00661         if((img[dijk]&2)==0) {
00662                 if(fi==nx-1) {
00663                         fijk2=fijk+1-nx;switchx2=switchx+(1-nx)*boxx;disx2=disx+bx;disxr2=disxr+bx;
00664                 } else {
00665                         fijk2=fijk+1;switchx2=switchx+boxx;disx2=disx;disxr2=disxr;
00666                 }
00667                 img[dijkr]|=1;
00668                 if(y_exist) {
00669                         img[dijkr-nx]|=4;
00670                         img[dijk-nx]|=8;
00671                 }
00672                 for(l=0;l<co[fijk2];l++) {
00673                         if(p[fijk2][ps*l+1]>switchy) {
00674                                 if(p[fijk2][ps*l]>switchx2) put_image(dijkr,fijk2,l,disxr2,disy,bz*ima);
00675                                 else put_image(dijk,fijk2,l,disx2,disy,bz*ima);
00676                         } else {
00677                                 if(!y_exist) continue;
00678                                 if(p[fijk2][ps*l]>switchx2) put_image(dijkr-nx,fijk2,l,disxr2,disy,bz*ima);
00679                                 else put_image(dijk-nx,fijk2,l,disx2,disy,bz*ima);
00680                         }
00681                 }
00682         }
00683 
00684         // Recomputation of some intermediate quantities for boundary cases
00685         if(fj==wy-1) {
00686                 fijk+=nx*(1-ny)-fi;
00687                 switchy+=(1-ny)*boxy;
00688                 disy+=by;
00689                 qi=di+step_int(-(ima*bxz+(qjdiv+1)*bxy)*xsp);
00690                 int dqidiv=step_div(qi,nx)-qidiv;qidiv+=dqidiv;
00691                 fi=qi-qidiv*nx;
00692                 fijk+=fi;
00693                 disx+=bxy+bx*dqidiv;
00694                 disxl+=bxy+bx*dqidiv;
00695                 disxr+=bxy+bx*dqidiv;
00696                 switchx-=bxy+bx*dqidiv;
00697         } else {
00698                 fijk+=nx;switchy+=boxy;
00699         }
00700 
00701         // Up-left image computation
00702         y_exist=dj!=oy-1;
00703         if((img[dijk]&4)==0) {
00704                 img[dijkl]|=8;
00705                 if(y_exist) {
00706                         img[dijkl+nx]|=2;
00707                         img[dijk+nx]|=1;
00708                 }
00709                 for(l=0;l<co[fijk];l++) {
00710                         if(p[fijk][ps*l+1]>switchy) {
00711                                 if(!y_exist) continue;
00712                                 if(p[fijk][ps*l]>switchx) put_image(dijk+nx,fijk,l,disx,disy,bz*ima);
00713                                 else put_image(dijkl+nx,fijk,l,disxl,disy,bz*ima);
00714                         } else {
00715                                 if(p[fijk][ps*l]>switchx) put_image(dijk,fijk,l,disx,disy,bz*ima);
00716                                 else put_image(dijkl,fijk,l,disxl,disy,bz*ima);
00717                         }
00718                 }
00719         }
00720 
00721         // Up-right image computation
00722         if((img[dijk]&8)==0) {
00723                 if(fi==nx-1) {
00724                         fijk2=fijk+1-nx;switchx2=switchx+(1-nx)*boxx;disx2=disx+bx;disxr2=disxr+bx;
00725                 } else {
00726                         fijk2=fijk+1;switchx2=switchx+boxx;disx2=disx;disxr2=disxr;
00727                 }
00728                 img[dijkr]|=4;
00729                 if(y_exist) {
00730                         img[dijkr+nx]|=1;
00731                         img[dijk+nx]|=2;
00732                 }
00733                 for(l=0;l<co[fijk2];l++) {
00734                         if(p[fijk2][ps*l+1]>switchy) {
00735                                 if(!y_exist) continue;
00736                                 if(p[fijk2][ps*l]>switchx2) put_image(dijkr+nx,fijk2,l,disxr2,disy,bz*ima);
00737                                 else put_image(dijk+nx,fijk2,l,disx2,disy,bz*ima);
00738                         } else {
00739                                 if(p[fijk2][ps*l]>switchx2) put_image(dijkr,fijk2,l,disxr2,disy,bz*ima);
00740                                 else put_image(dijk,fijk2,l,disx2,disy,bz*ima);
00741                         }
00742                 }
00743         }
00744 
00745         // All contributions to the block now added, so set all four bits of
00746         // the image information
00747         img[dijk]=15;
00748 }
00749 
00750 /** Copies a particle position from the primary domain into an image block.
00751  * \param[in] reg the block index within the primary domain that the particle
00752  *                is within.
00753  * \param[in] fijk the index of the image block.
00754  * \param[in] l the index of the particle entry within the primary block.
00755  * \param[in] (dx,dy,dz) the displacement vector to add to the particle. */
00756 void container_periodic_base::put_image(int reg,int fijk,int l,double dx,double dy,double dz) {
00757         if(co[reg]==mem[reg]) add_particle_memory(reg);
00758         double *p1=p[reg]+ps*co[reg],*p2=p[fijk]+ps*l;
00759         *(p1++)=*(p2++)+dx;
00760         *(p1++)=*(p2++)+dy;
00761         *p1=*p2+dz;
00762         if(ps==4) *(++p1)=*(++p2);
00763         id[reg][co[reg]++]=id[fijk][l];
00764 }
00765 
00766 }