Voro++
container.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.cc
00008  * \brief Function implementations for the container and related classes. */
00009 
00010 #include "container.hh"
00011 
00012 namespace voro {
00013 
00014 /** The class constructor sets up the geometry of container, initializing the
00015  * minimum and maximum coordinates in each direction, and setting whether each
00016  * direction is periodic or not. It divides the container into a rectangular
00017  * grid of blocks, and allocates memory for each of these for storing particle
00018  * positions and IDs.
00019  * \param[in] (ax_,bx_) the minimum and maximum x coordinates.
00020  * \param[in] (ay_,by_) the minimum and maximum y coordinates.
00021  * \param[in] (az_,bz_) the minimum and maximum z coordinates.
00022  * \param[in] (nx_,ny_,nz_) the number of grid blocks in each of the three
00023  *                       coordinate directions.
00024  * \param[in] (xperiodic_,yperiodic_,zperiodic_ ) flags setting whether the
00025  *                                                container is periodic in each
00026  *                                                coordinate direction.
00027  * \param[in] init_mem the initial memory allocation for each block.
00028  * \param[in] ps_ the number of floating point entries to store for each
00029  *                particle. */
00030 container_base::container_base(double ax_,double bx_,double ay_,double by_,double az_,double bz_,
00031                 int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem,int ps_)
00032         : voro_base(nx_,ny_,nz_,(bx_-ax_)/nx_,(by_-ay_)/ny_,(bz_-az_)/nz_),
00033         ax(ax_), bx(bx_), ay(ay_), by(by_), az(az_), bz(bz_),
00034         xperiodic(xperiodic_), yperiodic(yperiodic_), zperiodic(zperiodic_),
00035         id(new int*[nxyz]), p(new double*[nxyz]), co(new int[nxyz]), mem(new int[nxyz]), ps(ps_) {
00036         int l;
00037         for(l=0;l<nxyz;l++) co[l]=0;
00038         for(l=0;l<nxyz;l++) mem[l]=init_mem;
00039         for(l=0;l<nxyz;l++) id[l]=new int[init_mem];
00040         for(l=0;l<nxyz;l++) p[l]=new double[ps*init_mem];
00041 }
00042 
00043 /** The container destructor frees the dynamically allocated memory. */
00044 container_base::~container_base() {
00045         int l;
00046         for(l=0;l<nxyz;l++) delete [] p[l];
00047         for(l=0;l<nxyz;l++) delete [] id[l];
00048         delete [] id;
00049         delete [] p;
00050         delete [] co;
00051         delete [] mem;
00052 }
00053 
00054 /** The class constructor sets up the geometry of container.
00055  * \param[in] (ax_,bx_) the minimum and maximum x coordinates.
00056  * \param[in] (ay_,by_) the minimum and maximum y coordinates.
00057  * \param[in] (az_,bz_) the minimum and maximum z coordinates.
00058  * \param[in] (nx_,ny_,nz_) the number of grid blocks in each of the three
00059  *                       coordinate directions.
00060  * \param[in] (xperiodic_,yperiodic_,zperiodic_ ) flags setting whether the
00061  *                                                container is periodic in each
00062  *                                                coordinate direction.
00063  * \param[in] init_mem the initial memory allocation for each block. */
00064 container::container(double ax_,double bx_,double ay_,double by_,double az_,double bz_,
00065         int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem)
00066         : container_base(ax_,bx_,ay_,by_,az_,bz_,nx_,ny_,nz_,xperiodic_,yperiodic_,zperiodic_,init_mem,3),
00067         vc(*this,xperiodic_?2*nx_+1:nx_,yperiodic_?2*ny_+1:ny_,zperiodic_?2*nz_+1:nz_) {}
00068 
00069 /** The class constructor sets up the geometry of container.
00070  * \param[in] (ax_,bx_) the minimum and maximum x coordinates.
00071  * \param[in] (ay_,by_) the minimum and maximum y coordinates.
00072  * \param[in] (az_,bz_) the minimum and maximum z coordinates.
00073  * \param[in] (nx_,ny_,nz_) the number of grid blocks in each of the three
00074  *                       coordinate directions.
00075  * \param[in] (xperiodic_,yperiodic_,zperiodic_ ) flags setting whether the
00076  *                                                container is periodic in each
00077  *                                                coordinate direction.
00078  * \param[in] init_mem the initial memory allocation for each block. */
00079 container_poly::container_poly(double ax_,double bx_,double ay_,double by_,double az_,double bz_,
00080         int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem)
00081         : container_base(ax_,bx_,ay_,by_,az_,bz_,nx_,ny_,nz_,xperiodic_,yperiodic_,zperiodic_,init_mem,4),
00082         max_radius(0), vc(*this,xperiodic_?2*nx_+1:nx_,yperiodic_?2*ny_+1:ny_,zperiodic_?2*nz_+1:nz_) {}
00083 
00084 /** Put a particle into the correct region of the container.
00085  * \param[in] n the numerical ID of the inserted particle.
00086  * \param[in] (x,y,z) the position vector of the inserted particle. */
00087 void container::put(int n,double x,double y,double z) {
00088         int ijk;
00089         if(put_locate_block(ijk,x,y,z)) {
00090                 id[ijk][co[ijk]]=n;
00091                 double *pp=p[ijk]+3*co[ijk]++;
00092                 *(pp++)=x;*(pp++)=y;*pp=z;
00093         }
00094 }
00095 
00096 /** Put a particle into the correct region of the container.
00097  * \param[in] n the numerical ID of the inserted particle.
00098  * \param[in] (x,y,z) the position vector of the inserted particle.
00099  * \param[in] r the radius of the particle. */
00100 void container_poly::put(int n,double x,double y,double z,double r) {
00101         int ijk;
00102         if(put_locate_block(ijk,x,y,z)) {
00103                 id[ijk][co[ijk]]=n;
00104                 double *pp=p[ijk]+4*co[ijk]++;
00105                 *(pp++)=x;*(pp++)=y;*(pp++)=z;*pp=r;
00106                 if(max_radius<r) max_radius=r;
00107         }
00108 }
00109 
00110 /** Put a particle into the correct region of the container, also recording
00111  * into which region it was stored.
00112  * \param[in] vo the ordering class in which to record the region.
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 void container::put(particle_order &vo,int n,double x,double y,double z) {
00116         int ijk;
00117         if(put_locate_block(ijk,x,y,z)) {
00118                 id[ijk][co[ijk]]=n;
00119                 vo.add(ijk,co[ijk]);
00120                 double *pp=p[ijk]+3*co[ijk]++;
00121                 *(pp++)=x;*(pp++)=y;*pp=z;
00122         }
00123 }
00124 
00125 /** Put a particle into the correct region of the container, also recording
00126  * into which region it was stored.
00127  * \param[in] vo the ordering class in which to record the region.
00128  * \param[in] n the numerical ID of the inserted particle.
00129  * \param[in] (x,y,z) the position vector of the inserted particle.
00130  * \param[in] r the radius of the particle. */
00131 void container_poly::put(particle_order &vo,int n,double x,double y,double z,double r) {
00132         int ijk;
00133         if(put_locate_block(ijk,x,y,z)) {
00134                 id[ijk][co[ijk]]=n;
00135                 vo.add(ijk,co[ijk]);
00136                 double *pp=p[ijk]+4*co[ijk]++;
00137                 *(pp++)=x;*(pp++)=y;*(pp++)=z;*pp=r;
00138                 if(max_radius<r) max_radius=r;
00139         }
00140 }
00141 
00142 /** This routine takes a particle position vector, tries to remap it into the
00143  * primary domain. If successful, it computes the region into which it can be
00144  * stored and checks that there is enough memory within this region to store
00145  * it.
00146  * \param[out] ijk the region index.
00147  * \param[in,out] (x,y,z) the particle position, remapped into the primary
00148  *                        domain if necessary.
00149  * \return True if the particle can be successfully placed into the container,
00150  * false otherwise. */
00151 inline bool container_base::put_locate_block(int &ijk,double &x,double &y,double &z) {
00152         if(put_remap(ijk,x,y,z)) {
00153                 if(co[ijk]==mem[ijk]) add_particle_memory(ijk);
00154                 return true;
00155         }
00156 #if VOROPP_REPORT_OUT_OF_BOUNDS ==1
00157         fprintf(stderr,"Out of bounds: (x,y,z)=(%g,%g,%g)\n",x,y,z);
00158 #endif
00159         return false;
00160 }
00161 
00162 /** Takes a particle position vector and computes the region index into which
00163  * it should be stored. If the container is periodic, then the routine also
00164  * maps the particle position to ensure it is in the primary domain. If the
00165  * container is not periodic, the routine bails out.
00166  * \param[out] ijk the region index.
00167  * \param[in,out] (x,y,z) the particle position, remapped into the primary
00168  *                        domain if necessary.
00169  * \return True if the particle can be successfully placed into the container,
00170  * false otherwise. */
00171 inline bool container_base::put_remap(int &ijk,double &x,double &y,double &z) {
00172         int l;
00173 
00174         ijk=step_int((x-ax)*xsp);
00175         if(xperiodic) {l=step_mod(ijk,nx);x+=boxx*(l-ijk);ijk=l;}
00176         else if(ijk<0||ijk>=nx) return false;
00177 
00178         int j=step_int((y-ay)*ysp);
00179         if(yperiodic) {l=step_mod(j,ny);y+=boxy*(l-j);j=l;}
00180         else if(j<0||j>=ny) return false;
00181 
00182         int k=step_int((z-az)*zsp);
00183         if(xperiodic) {l=step_mod(k,nz);z+=boxz*(l-k);k=l;}
00184         else if(k<0||k>=nz) return false;
00185 
00186         ijk+=nx*j+nxy*k;
00187         return true;
00188 }
00189 
00190 /** Takes a position vector and attempts to remap it into the primary domain.
00191  * \param[out] (ai,aj,ak) the periodic image displacement that the vector is in,
00192  *                       with (0,0,0) corresponding to the primary domain.
00193  * \param[out] (ci,cj,ck) the index of the block that the position vector is
00194  *                        within, once it has been remapped.
00195  * \param[in,out] (x,y,z) the position vector to consider, which is remapped
00196  *                        into the primary domain during the routine.
00197  * \param[out] ijk the block index that the vector is within.
00198  * \return True if the particle is within the container or can be remapped into
00199  * it, false if it lies outside of the container bounds. */
00200 inline bool container_base::remap(int &ai,int &aj,int &ak,int &ci,int &cj,int &ck,double &x,double &y,double &z,int &ijk) {
00201         ci=step_int((x-ax)*xsp);
00202         if(ci<0||ci>=nx) {
00203                 if(xperiodic) {ai=step_div(ci,nx);x-=ai*(bx-ax);ci-=ai*nx;}
00204                 else return false;
00205         } else ai=0;
00206 
00207         cj=step_int((y-ay)*ysp);
00208         if(cj<0||cj>=ny) {
00209                 if(yperiodic) {aj=step_div(cj,ny);y-=aj*(by-ay);cj-=aj*ny;}
00210                 else return false;
00211         } else aj=0;
00212 
00213         ck=step_int((z-az)*zsp);
00214         if(ck<0||ck>=nz) {
00215                 if(zperiodic) {ak=step_div(ck,nz);z-=ak*(bz-az);ck-=ak*nz;}
00216                 else return false;
00217         } else ak=0;
00218 
00219         ijk=ci+nx*cj+nxy*ck;
00220         return true;
00221 }
00222 
00223 /** Takes a vector and finds the particle whose Voronoi cell contains that
00224  * vector. This is equivalent to finding the particle which is nearest to the
00225  * vector. Additional wall classes are not considered by this routine.
00226  * \param[in] (x,y,z) the vector to test.
00227  * \param[out] (rx,ry,rz) the position of the particle whose Voronoi cell
00228  *                        contains the vector. If the container is periodic,
00229  *                        this may point to a particle in a periodic image of
00230  *                        the primary domain.
00231  * \param[out] pid the ID of the particle.
00232  * \return True if a particle was found. If the container has no particles,
00233  * then the search will not find a Voronoi cell and false is returned. */
00234 bool container::find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid) {
00235         int ai,aj,ak,ci,cj,ck,ijk;
00236         particle_record w;
00237         double mrs;
00238 
00239         // If the given vector lies outside the domain, but the container
00240         // is periodic, then remap it back into the domain
00241         if(!remap(ai,aj,ak,ci,cj,ck,x,y,z,ijk)) return false;
00242         vc.find_voronoi_cell(x,y,z,ci,cj,ck,ijk,w,mrs);
00243 
00244         if(w.ijk!=-1) {
00245 
00246                 // Assemble the position vector of the particle to be returned,
00247                 // applying a periodic remapping if necessary
00248                 if(xperiodic) {ci+=w.di;if(ci<0||ci>=nx) ai+=step_div(ci,nx);}
00249                 if(yperiodic) {cj+=w.dj;if(cj<0||cj>=ny) aj+=step_div(cj,ny);}
00250                 if(zperiodic) {ck+=w.dk;if(ck<0||ck>=nz) ak+=step_div(ck,nz);}
00251                 rx=p[w.ijk][3*w.l]+ai*(bx-ax);
00252                 ry=p[w.ijk][3*w.l+1]+aj*(by-ay);
00253                 rz=p[w.ijk][3*w.l+2]+ak*(bz-az);
00254                 pid=id[w.ijk][w.l];
00255                 return true;
00256         }
00257 
00258         // If no particle if found then just return false
00259         return false;
00260 }
00261 
00262 /** Takes a vector and finds the particle whose Voronoi cell contains that
00263  * vector. Additional wall classes are not considered by this routine.
00264  * \param[in] (x,y,z) the vector to test.
00265  * \param[out] (rx,ry,rz) the position of the particle whose Voronoi cell
00266  *                        contains the vector. If the container is periodic,
00267  *                        this may point to a particle in a periodic image of
00268  *                        the primary domain.
00269  * \param[out] pid the ID of the particle.
00270  * \return True if a particle was found. If the container has no particles,
00271  * then the search will not find a Voronoi cell and false is returned. */
00272 bool container_poly::find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid) {
00273         int ai,aj,ak,ci,cj,ck,ijk;
00274         particle_record w;
00275         double mrs;
00276 
00277         // If the given vector lies outside the domain, but the container
00278         // is periodic, then remap it back into the domain
00279         if(!remap(ai,aj,ak,ci,cj,ck,x,y,z,ijk)) return false;
00280         vc.find_voronoi_cell(x,y,z,ci,cj,ck,ijk,w,mrs);
00281 
00282         if(w.ijk!=-1) {
00283 
00284                 // Assemble the position vector of the particle to be returned,
00285                 // applying a periodic remapping if necessary
00286                 if(xperiodic) {ci+=w.di;if(ci<0||ci>=nx) ai+=step_div(ci,nx);}
00287                 if(yperiodic) {cj+=w.dj;if(cj<0||cj>=ny) aj+=step_div(cj,ny);}
00288                 if(zperiodic) {ck+=w.dk;if(ck<0||ck>=nz) ak+=step_div(ck,nz);}
00289                 rx=p[w.ijk][4*w.l]+ai*(bx-ax);
00290                 ry=p[w.ijk][4*w.l+1]+aj*(by-ay);
00291                 rz=p[w.ijk][4*w.l+2]+ak*(bz-az);
00292                 pid=id[w.ijk][w.l];
00293                 return true;
00294         }
00295 
00296         // If no particle if found then just return false
00297         return false;
00298 }
00299 
00300 /** Increase memory for a particular region.
00301  * \param[in] i the index of the region to reallocate. */
00302 void container_base::add_particle_memory(int i) {
00303         int l,nmem=mem[i]<<1;
00304 
00305         // Carry out a check on the memory allocation size, and
00306         // print a status message if requested
00307         if(nmem>max_particle_memory)
00308                 voro_fatal_error("Absolute maximum memory allocation exceeded",VOROPP_MEMORY_ERROR);
00309 #if VOROPP_VERBOSE >=3
00310         fprintf(stderr,"Particle memory in region %d scaled up to %d\n",i,nmem);
00311 #endif
00312 
00313         // Allocate new memory and copy in the contents of the old arrays
00314         int *idp=new int[nmem];
00315         for(l=0;l<co[i];l++) idp[l]=id[i][l];
00316         double *pp=new double[ps*nmem];
00317         for(l=0;l<ps*co[i];l++) pp[l]=p[i][l];
00318 
00319         // Update pointers and delete old arrays
00320         mem[i]=nmem;
00321         delete [] id[i];id[i]=idp;
00322         delete [] p[i];p[i]=pp;
00323 }
00324 
00325 /** Import a list of particles from an open file stream into the container.
00326  * Entries of four numbers (Particle ID, x position, y position, z position)
00327  * are searched for. If the file cannot be successfully read, then the routine
00328  * causes a fatal error.
00329  * \param[in] fp the file handle to read from. */
00330 void container::import(FILE *fp) {
00331         int i,j;
00332         double x,y,z;
00333         while((j=fscanf(fp,"%d %lg %lg %lg",&i,&x,&y,&z))==4) put(i,x,y,z);
00334         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00335 }
00336 
00337 /** Import a list of particles from an open file stream, also storing the order
00338  * of that the particles are read. Entries of four numbers (Particle ID, x
00339  * position, y position, z position) are searched for. If the file cannot be
00340  * successfully read, then the routine causes a fatal error.
00341  * \param[in,out] vo a reference to an ordering class to use.
00342  * \param[in] fp the file handle to read from. */
00343 void container::import(particle_order &vo,FILE *fp) {
00344         int i,j;
00345         double x,y,z;
00346         while((j=fscanf(fp,"%d %lg %lg %lg",&i,&x,&y,&z))==4) put(vo,i,x,y,z);
00347         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00348 }
00349 
00350 /** Import a list of particles from an open file stream into the container.
00351  * Entries of five numbers (Particle ID, x position, y position, z position,
00352  * radius) are searched for. If the file cannot be successfully read, then the
00353  * routine causes a fatal error.
00354  * \param[in] fp the file handle to read from. */
00355 void container_poly::import(FILE *fp) {
00356         int i,j;
00357         double x,y,z,r;
00358         while((j=fscanf(fp,"%d %lg %lg %lg %lg",&i,&x,&y,&z,&r))==5) put(i,x,y,z,r);
00359         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00360 }
00361 
00362 /** Import a list of particles from an open file stream, also storing the order
00363  * of that the particles are read. Entries of four numbers (Particle ID, x
00364  * position, y position, z position, radius) are searched for. If the file
00365  * cannot be successfully read, then the routine causes a fatal error.
00366  * \param[in,out] vo a reference to an ordering class to use.
00367  * \param[in] fp the file handle to read from. */
00368 void container_poly::import(particle_order &vo,FILE *fp) {
00369         int i,j;
00370         double x,y,z,r;
00371         while((j=fscanf(fp,"%d %lg %lg %lg %lg",&i,&x,&y,&z,&r))==5) put(vo,i,x,y,z,r);
00372         if(j!=EOF) voro_fatal_error("File import error",VOROPP_FILE_ERROR);
00373 }
00374 
00375 /** Outputs the a list of all the container regions along with the number of
00376  * particles stored within each. */
00377 void container_base::region_count() {
00378         int i,j,k,*cop=co;
00379         for(k=0;k<nz;k++) for(j=0;j<ny;j++) for(i=0;i<nx;i++)
00380                 printf("Region (%d,%d,%d): %d particles\n",i,j,k,*(cop++));
00381 }
00382 
00383 /** Clears a container of particles. */
00384 void container::clear() {
00385         for(int *cop=co;cop<co+nxyz;cop++) *cop=0;
00386 }
00387 
00388 /** Clears a container of particles, also clearing resetting the maximum radius
00389  * to zero. */
00390 void container_poly::clear() {
00391         for(int *cop=co;cop<co+nxyz;cop++) *cop=0;
00392         max_radius=0;
00393 }
00394 
00395 /** Computes all the Voronoi cells and saves customized information about them.
00396  * \param[in] format the custom output string to use.
00397  * \param[in] fp a file handle to write to. */
00398 void container::print_custom(const char *format,FILE *fp) {
00399         c_loop_all vl(*this);
00400         print_custom(vl,format,fp);
00401 }
00402 
00403 /** Computes all the Voronoi cells and saves customized
00404  * information about them.
00405  * \param[in] format the custom output string to use.
00406  * \param[in] fp a file handle to write to. */
00407 void container_poly::print_custom(const char *format,FILE *fp) {
00408         c_loop_all vl(*this);
00409         print_custom(vl,format,fp);
00410 }
00411 
00412 /** Computes all the Voronoi cells and saves customized information about them.
00413  * \param[in] format the custom output string to use.
00414  * \param[in] filename the name of the file to write to. */
00415 void container::print_custom(const char *format,const char *filename) {
00416         FILE *fp=safe_fopen(filename,"w");
00417         print_custom(format,fp);
00418         fclose(fp);
00419 }
00420 
00421 /** Computes all the Voronoi cells and saves customized
00422  * information about them
00423  * \param[in] format the custom output string to use.
00424  * \param[in] filename the name of the file to write to. */
00425 void container_poly::print_custom(const char *format,const char *filename) {
00426         FILE *fp=safe_fopen(filename,"w");
00427         print_custom(format,fp);
00428         fclose(fp);
00429 }
00430 
00431 /** Computes all of the Voronoi cells in the container, but does nothing
00432  * with the output. It is useful for measuring the pure computation time
00433  * of the Voronoi algorithm, without any additional calculations such as
00434  * volume evaluation or cell output. */
00435 void container::compute_all_cells() {
00436         voronoicell c;
00437         c_loop_all vl(*this);
00438         if(vl.start()) do compute_cell(c,vl);
00439         while(vl.inc());
00440 }
00441 
00442 /** Computes all of the Voronoi cells in the container, but does nothing
00443  * with the output. It is useful for measuring the pure computation time
00444  * of the Voronoi algorithm, without any additional calculations such as
00445  * volume evaluation or cell output. */
00446 void container_poly::compute_all_cells() {
00447         voronoicell c;
00448         c_loop_all vl(*this);
00449         if(vl.start()) do compute_cell(c,vl);while(vl.inc());
00450 }
00451 
00452 /** Calculates all of the Voronoi cells and sums their volumes. In most cases
00453  * without walls, the sum of the Voronoi cell volumes should equal the volume
00454  * of the container to numerical precision.
00455  * \return The sum of all of the computed Voronoi volumes. */
00456 double container::sum_cell_volumes() {
00457         voronoicell c;
00458         double vol=0;
00459         c_loop_all vl(*this);
00460         if(vl.start()) do if(compute_cell(c,vl)) vol+=c.volume();while(vl.inc());
00461         return vol;
00462 }
00463 
00464 /** Calculates all of the Voronoi cells and sums their volumes. In most cases
00465  * without walls, the sum of the Voronoi cell volumes should equal the volume
00466  * of the container to numerical precision.
00467  * \return The sum of all of the computed Voronoi volumes. */
00468 double container_poly::sum_cell_volumes() {
00469         voronoicell c;
00470         double vol=0;
00471         c_loop_all vl(*this);
00472         if(vl.start()) do if(compute_cell(c,vl)) vol+=c.volume();while(vl.inc());
00473         return vol;
00474 }
00475 
00476 /** This function tests to see if a given vector lies within the container
00477  * bounds and any walls.
00478  * \param[in] (x,y,z) the position vector to be tested.
00479  * \return True if the point is inside the container, false if the point is
00480  *         outside. */
00481 bool container_base::point_inside(double x,double y,double z) {
00482         if(x<ax||x>bx||y<ay||y>by||z<az||z>bz) return false;
00483         return point_inside_walls(x,y,z);
00484 }
00485 
00486 /** Draws an outline of the domain in gnuplot format.
00487  * \param[in] fp the file handle to write to. */
00488 void container_base::draw_domain_gnuplot(FILE *fp) {
00489         fprintf(fp,"%g %g %g\n%g %g %g\n%g %g %g\n%g %g %g\n",ax,ay,az,bx,ay,az,bx,by,az,ax,by,az);
00490         fprintf(fp,"%g %g %g\n%g %g %g\n%g %g %g\n%g %g %g\n",ax,by,bz,bx,by,bz,bx,ay,bz,ax,ay,bz);
00491         fprintf(fp,"%g %g %g\n\n%g %g %g\n%g %g %g\n\n",ax,by,bz,ax,ay,az,ax,ay,bz);
00492         fprintf(fp,"%g %g %g\n%g %g %g\n\n%g %g %g\n%g %g %g\n\n",bx,ay,az,bx,ay,bz,bx,by,az,bx,by,bz);
00493 }
00494 
00495 /** Draws an outline of the domain in POV-Ray format.
00496  * \param[in] fp the file handle to write to. */
00497 void container_base::draw_domain_pov(FILE *fp) {
00498         fprintf(fp,"cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n"
00499                    "cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n",ax,ay,az,bx,ay,az,ax,by,az,bx,by,az);
00500         fprintf(fp,"cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n"
00501                    "cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n",ax,by,bz,bx,by,bz,ax,ay,bz,bx,ay,bz);
00502         fprintf(fp,"cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n"
00503                    "cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n",ax,ay,az,ax,by,az,bx,ay,az,bx,by,az);
00504         fprintf(fp,"cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n"
00505                    "cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n",bx,ay,bz,bx,by,bz,ax,ay,bz,ax,by,bz);
00506         fprintf(fp,"cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n"
00507                    "cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n",ax,ay,az,ax,ay,bz,bx,ay,az,bx,ay,bz);
00508         fprintf(fp,"cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n"
00509                    "cylinder{<%g,%g,%g>,<%g,%g,%g>,rr}\n",bx,by,az,bx,by,bz,ax,by,az,ax,by,bz);
00510         fprintf(fp,"sphere{<%g,%g,%g>,rr}\nsphere{<%g,%g,%g>,rr}\n"
00511                    "sphere{<%g,%g,%g>,rr}\nsphere{<%g,%g,%g>,rr}\n",ax,ay,az,bx,ay,az,ax,by,az,bx,by,az);
00512         fprintf(fp,"sphere{<%g,%g,%g>,rr}\nsphere{<%g,%g,%g>,rr}\n"
00513                    "sphere{<%g,%g,%g>,rr}\nsphere{<%g,%g,%g>,rr}\n",ax,ay,bz,bx,ay,bz,ax,by,bz,bx,by,bz);
00514 }
00515 
00516 
00517 /** The wall_list constructor sets up an array of pointers to wall classes. */
00518 wall_list::wall_list() : walls(new wall*[init_wall_size]), wep(walls), wel(walls+init_wall_size),
00519         current_wall_size(init_wall_size) {}
00520 
00521 /** The wall_list destructor frees the array of pointers to the wall classes.
00522  */
00523 wall_list::~wall_list() {
00524         delete [] walls;
00525 }
00526 
00527 /** Adds all of the walls on another wall_list to this class.
00528  * \param[in] wl a reference to the wall class. */
00529 void wall_list::add_wall(wall_list &wl) {
00530         for(wall **wp=wl.walls;wp<wl.wep;wp++) add_wall(*wp);
00531 }
00532 
00533 /** Deallocates all of the wall classes pointed to by the wall_list. */
00534 void wall_list::deallocate() {
00535         for(wall **wp=walls;wp<wep;wp++) delete *wp;
00536 }
00537 
00538 /** Increases the memory allocation for the walls array. */
00539 void wall_list::increase_wall_memory() {
00540         current_wall_size<<=1;
00541         if(current_wall_size>max_wall_size)
00542                 voro_fatal_error("Wall memory allocation exceeded absolute maximum",VOROPP_MEMORY_ERROR);
00543         wall **nwalls=new wall*[current_wall_size],**nwp=nwalls,**wp=walls;
00544         while(wp<wep) *(nwp++)=*(wp++);
00545         delete [] walls;
00546         walls=nwalls;wel=walls+current_wall_size;wep=nwp;
00547 }
00548 
00549 }