Voro++
container_prd.hh
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.hh
00008  * \brief Header file for the container_periodic_base and related classes. */
00009 
00010 #ifndef VOROPP_CONTAINER_PRD_HH
00011 #define VOROPP_CONTAINER_PRD_HH
00012 
00013 #include <cstdio>
00014 #include <cstdlib>
00015 #include <cmath>
00016 #include <vector>
00017 using namespace std;
00018 
00019 #include "config.hh"
00020 #include "common.hh"
00021 #include "v_base.hh"
00022 #include "cell.hh"
00023 #include "c_loops.hh"
00024 #include "v_compute.hh"
00025 #include "unitcell.hh"
00026 
00027 namespace voro {
00028 
00029 /** \brief Class for representing a particle system in a 3D periodic
00030  * non-orthogonal periodic domain.
00031  *
00032  * This class represents a particle system in a three-dimensional
00033  * non-orthogonal periodic domain. The domain is defined by three periodicity
00034  * vectors (bx,0,0), (bxy,by,0), and (bxz,byz,bz) that represent a
00035  * parallelepiped. Internally, the class stores particles in the box 0<x<bx,
00036  * 0<y<by, 0<z<bz, and constructs periodic images of particles that displaced
00037  * by the three periodicity vectors when they are necessary for the
00038  * computation. The internal memory structure for this class is significantly
00039  * different from the container_base class in order to handle the dynamic
00040  * construction of these periodic images.
00041  *
00042  * The class is derived from the unitcell class, which encapsulates information
00043  * about the domain geometry, and the voro_base class, which encapsulates
00044  * information about the underlying computational grid. */
00045 class container_periodic_base : public unitcell, public voro_base {
00046         public:
00047                 /** The lower y index (inclusive) of the primary domain within
00048                  * the block structure. */
00049                 int ey;
00050                 /** The lower z index (inclusive) of the primary domain within
00051                  * the block structure. */
00052                 int ez;
00053                 /** The upper y index (exclusive) of the primary domain within
00054                  * the block structure. */
00055                 int wy;
00056                 /** The upper z index (exclusive) of the primary domain within
00057                  * the block structure. */
00058                 int wz;
00059                 /** The total size of the block structure (including images) in
00060                  * the y direction. */
00061                 int oy;
00062                 /** The total size of the block structure (including images) in
00063                  * the z direction. */
00064                 int oz;
00065                 /** The total number of blocks. */
00066                 int oxyz;
00067                 /** This array holds the numerical IDs of each particle in each
00068                  * computational box. */
00069                 int **id;
00070                 /** A two dimensional array holding particle positions. For the
00071                  * derived container_poly class, this also holds particle
00072                  * radii. */
00073                 double **p;
00074                 /** This array holds the number of particles within each
00075                  * computational box of the container. */
00076                 int *co;
00077                 /** This array holds the maximum amount of particle memory for
00078                  * each computational box of the container. If the number of
00079                  * particles in a particular box ever approaches this limit,
00080                  * more is allocated using the add_particle_memory() function.
00081                  */
00082                 int *mem;
00083                 /** An array holding information about periodic image
00084                  * construction at a given location. */
00085                 char *img;
00086                 /** The initial amount of memory to allocate for particles
00087                  * for each block. */
00088                 const int init_mem;
00089                 /** The amount of memory in the array structure for each
00090                  * particle. This is set to 3 when the basic class is
00091                  * initialized, so that the array holds (x,y,z) positions. If
00092                  * the container class is initialized as part of the derived
00093                  * class container_poly, then this is set to 4, to also hold
00094                  * the particle radii. */
00095                 const int ps;
00096                 container_periodic_base(double bx_,double bxy_,double by_,double bxz_,double byz_,double bz_,
00097                                 int nx_,int ny_,int nz_,int init_mem_,int ps);
00098                 ~container_periodic_base();
00099                 /** Prints all particles in the container, including those that
00100                  * have been constructed in image blocks. */
00101                 inline void print_all_particles() {
00102                         int ijk,q;
00103                         for(ijk=0;ijk<oxyz;ijk++) for(q=0;q<co[ijk];q++)
00104                                 printf("%d %g %g %g\n",id[ijk][q],p[ijk][ps*q],p[ijk][ps*q+1],p[ijk][ps*q+2]);
00105                 }
00106                 void region_count();
00107                 /** Initializes the Voronoi cell prior to a compute_cell
00108                  * operation for a specific particle being carried out by a
00109                  * voro_compute class. The cell is initialized to be the
00110                  * pre-computed unit Voronoi cell based on planes formed by
00111                  * periodic images of the particle.
00112                  * \param[in,out] c a reference to a voronoicell object.
00113                  * \param[in] ijk the block that the particle is within.
00114                  * \param[in] q the index of the particle within its block.
00115                  * \param[in] (ci,cj,ck) the coordinates of the block in the
00116                  *                       container coordinate system.
00117                  * \param[out] (i,j,k) the coordinates of the test block
00118                  *                     relative to the voro_compute
00119                  *                     coordinate system.
00120                  * \param[out] (x,y,z) the position of the particle.
00121                  * \param[out] disp a block displacement used internally by the
00122                  *                  compute_cell routine.
00123                  * \return False if the plane cuts applied by walls completely
00124                  * removed the cell, true otherwise. */
00125                 template<class v_cell>
00126                 inline bool initialize_voronoicell(v_cell &c,int ijk,int q,int ci,int cj,int ck,int &i,int &j,int &k,double &x,double &y,double &z,int &disp) {
00127                         c=unit_voro;
00128                         double *pp(p[ijk]+ps*q);
00129                         x=*(pp++);y=*(pp++);z=*pp;
00130                         i=nx;j=ey;k=ez;
00131                         return true;
00132                 }
00133                 /** Initializes parameters for a find_voronoi_cell call within
00134                  * the voro_compute template.
00135                  * \param[in] (ci,cj,ck) the coordinates of the test block in
00136                  *                       the container coordinate system.
00137                  * \param[in] ijk the index of the test block
00138                  * \param[out] (i,j,k) the coordinates of the test block
00139                  *                     relative to the voro_compute
00140                  *                     coordinate system.
00141                  * \param[out] disp a block displacement used internally by the
00142                  *                  find_voronoi_cell routine (but not needed
00143                  *                  in this instance.) */
00144                 inline void initialize_search(int ci,int cj,int ck,int ijk,int &i,int &j,int &k,int &disp) {
00145                         i=nx;j=ey;k=ez;
00146                 }
00147                 /** Returns the position of a particle currently being computed
00148                  * relative to the computational block that it is within. It is
00149                  * used to select the optimal worklist entry to use.
00150                  * \param[in] (x,y,z) the position of the particle.
00151                  * \param[in] (ci,cj,ck) the block that the particle is within.
00152                  * \param[out] (fx,fy,fz) the position relative to the block.
00153                  */
00154                 inline void frac_pos(double x,double y,double z,double ci,double cj,double ck,double &fx,double &fy,double &fz) {
00155                         fx=x-boxx*ci;
00156                         fy=y-boxy*(cj-ey);
00157                         fz=z-boxz*(ck-ez);
00158                 }
00159                 /** Calculates the index of block in the container structure
00160                  * corresponding to given coordinates.
00161                  * \param[in] (ci,cj,ck) the coordinates of the original block
00162                  *                       in the current computation, relative
00163                  *                       to the container coordinate system.
00164                  * \param[in] (ei,ej,ek) the displacement of the current block
00165                  *                       from the original block.
00166                  * \param[in,out] (qx,qy,qz) the periodic displacement that
00167                  *                           must be added to the particles
00168                  *                           within the computed block.
00169                  * \param[in] disp a block displacement used internally by the
00170                  *                  find_voronoi_cell and compute_cell routines
00171                  *                  (but not needed in this instance.)
00172                  * \return The block index. */
00173                 inline int region_index(int ci,int cj,int ck,int ei,int ej,int ek,double &qx,double &qy,double &qz,int &disp) {
00174                         int qi=ci+(ei-nx),qj=cj+(ej-ey),qk=ck+(ek-ez);
00175                         int iv(step_div(qi,nx));if(iv!=0) {qx=iv*bx;qi-=nx*iv;} else qx=0;
00176                         create_periodic_image(qi,qj,qk);
00177                         return qi+nx*(qj+oy*qk);
00178                 }
00179                 void create_all_images();
00180                 void check_compartmentalized();
00181         protected:
00182                 void add_particle_memory(int i);
00183                 void put_locate_block(int &ijk,double &x,double &y,double &z);
00184                 void put_locate_block(int &ijk,double &x,double &y,double &z,int &ai,int &aj,int &ak);
00185                 /** Creates particles within an image block by copying them
00186                  * from the primary domain and shifting them. If the given
00187                  * block is aligned with the primary domain in the z-direction,
00188                  * the routine calls the simpler create_side_image routine
00189                  * where the image block may comprise of particles from up to
00190                  * two primary blocks. Otherwise is calls the more complex
00191                  * create_vertical_image where the image block may comprise of
00192                  * particles from up to four primary blocks.
00193                  * \param[in] (di,dj,dk) the coordinates of the image block to
00194                  *                       create. */
00195                 inline void create_periodic_image(int di,int dj,int dk) {
00196                         if(di<0||di>=nx||dj<0||dj>=oy||dk<0||dk>=oz)
00197                                 voro_fatal_error("Constructing periodic image for nonexistent point",VOROPP_INTERNAL_ERROR);
00198                         if(dk>=ez&&dk<wz) {
00199                                 if(dj<ey||dj>=wy) create_side_image(di,dj,dk);
00200                         } else create_vertical_image(di,dj,dk);
00201                 }
00202                 void create_side_image(int di,int dj,int dk);
00203                 void create_vertical_image(int di,int dj,int dk);
00204                 void put_image(int reg,int fijk,int l,double dx,double dy,double dz);
00205                 inline void remap(int &ai,int &aj,int &ak,int &ci,int &cj,int &ck,double &x,double &y,double &z,int &ijk);
00206 };
00207 
00208 /** \brief Extension of the container_periodic_base class for computing regular
00209  * Voronoi tessellations.
00210  *
00211  * This class is an extension of the container_periodic_base that has routines
00212  * specifically for computing the regular Voronoi tessellation with no
00213  * dependence on particle radii. */
00214 class container_periodic : public container_periodic_base {
00215         public:
00216                 container_periodic(double bx_,double bxy_,double by_,double bxz_,double byz_,double bz_,
00217                                 int nx_,int ny_,int nz_,int init_mem_);
00218                 void clear();
00219                 void put(int n,double x,double y,double z);
00220                 void put(int n,double x,double y,double z,int &ai,int &aj,int &ak);
00221                 void put(particle_order &vo,int n,double x,double y,double z);
00222                 void import(FILE *fp=stdin);
00223                 void import(particle_order &vo,FILE *fp=stdin);
00224                 /** Imports a list of particles from an open file stream into
00225                  * the container. Entries of four numbers (Particle ID, x
00226                  * position, y position, z position) are searched for. If the
00227                  * file cannot be successfully read, then the routine causes a
00228                  * fatal error.
00229                  * \param[in] filename the name of the file to open and read
00230                  *                     from. */
00231                 inline void import(const char* filename) {
00232                         FILE *fp(safe_fopen(filename,"r"));
00233                         import(fp);
00234                         fclose(fp);
00235                 }
00236                 /** Imports a list of particles from an open file stream into
00237                  * the container. Entries of four numbers (Particle ID, x
00238                  * position, y position, z position) are searched for. In
00239                  * addition, the order in which particles are read is saved
00240                  * into an ordering class. If the file cannot be successfully
00241                  * read, then the routine causes a fatal error.
00242                  * \param[in,out] vo the ordering class to use.
00243                  * \param[in] filename the name of the file to open and read
00244                  *                     from. */
00245                 inline void import(particle_order &vo,const char* filename) {
00246                         FILE *fp(safe_fopen(filename,"r"));
00247                         import(vo,fp);
00248                         fclose(fp);
00249                 }
00250                 void compute_all_cells();
00251                 double sum_cell_volumes();
00252                 /** Dumps particle IDs and positions to a file.
00253                  * \param[in] vl the loop class to use.
00254                  * \param[in] fp a file handle to write to. */
00255                 template<class c_loop>
00256                 void draw_particles(c_loop &vl,FILE *fp) {
00257                         double *pp;
00258                         if(vl.start()) do {
00259                                 pp=p[vl.ijk]+3*vl.q;
00260                                 fprintf(fp,"%d %g %g %g\n",id[vl.ijk][vl.q],*pp,pp[1],pp[2]);
00261                         } while(vl.inc());
00262                 }
00263                 /** Dumps all of the particle IDs and positions to a file.
00264                  * \param[in] fp a file handle to write to. */
00265                 inline void draw_particles(FILE *fp=stdout) {
00266                         c_loop_all_periodic vl(*this);
00267                         draw_particles(vl,fp);
00268                 }
00269                 /** Dumps all of the particle IDs and positions to a file.
00270                  * \param[in] filename the name of the file to write to. */
00271                 inline void draw_particles(const char *filename) {
00272                         FILE *fp(safe_fopen(filename,"w"));
00273                         draw_particles(fp);
00274                         fclose(fp);
00275                 }
00276                 /** Dumps particle positions in POV-Ray format.
00277                  * \param[in] vl the loop class to use.
00278                  * \param[in] fp a file handle to write to. */
00279                 template<class c_loop>
00280                 void draw_particles_pov(c_loop &vl,FILE *fp) {
00281                         double *pp;
00282                         if(vl.start()) do {
00283                                 pp=p[vl.ijk]+3*vl.q;
00284                                 fprintf(fp,"// id %d\nsphere{<%g,%g,%g>,s}\n",
00285                                                 id[vl.ijk][vl.q],*pp,pp[1],pp[2]);
00286                         } while(vl.inc());
00287                 }
00288                 /** Dumps all particle positions in POV-Ray format.
00289                  * \param[in] fp a file handle to write to. */
00290                 inline void draw_particles_pov(FILE *fp=stdout) {
00291                         c_loop_all_periodic vl(*this);
00292                         draw_particles_pov(vl,fp);
00293                 }
00294                 /** Dumps all particle positions in POV-Ray format.
00295                  * \param[in] filename the name of the file to write to. */
00296                 inline void draw_particles_pov(const char *filename) {
00297                         FILE *fp(safe_fopen(filename,"w"));
00298                         draw_particles_pov(fp);
00299                         fclose(fp);
00300                 }
00301                 /** Computes Voronoi cells and saves the output in gnuplot
00302                  * format.
00303                  * \param[in] vl the loop class to use.
00304                  * \param[in] fp a file handle to write to. */
00305                 template<class c_loop>
00306                 void draw_cells_gnuplot(c_loop &vl,FILE *fp) {
00307                         voronoicell c;double *pp;
00308                         if(vl.start()) do if(compute_cell(c,vl)) {
00309                                 pp=p[vl.ijk]+ps*vl.q;
00310                                 c.draw_gnuplot(*pp,pp[1],pp[2],fp);
00311                         } while(vl.inc());
00312                 }
00313                 /** Computes all Voronoi cells and saves the output in gnuplot
00314                  * format.
00315                  * \param[in] fp a file handle to write to. */
00316                 inline void draw_cells_gnuplot(FILE *fp=stdout) {
00317                         c_loop_all_periodic vl(*this);
00318                         draw_cells_gnuplot(vl,fp);
00319                 }
00320                 /** Compute all Voronoi cells and saves the output in gnuplot
00321                  * format.
00322                  * \param[in] filename the name of the file to write to. */
00323                 inline void draw_cells_gnuplot(const char *filename) {
00324                         FILE *fp(safe_fopen(filename,"w"));
00325                         draw_cells_gnuplot(fp);
00326                         fclose(fp);
00327                 }
00328                 /** Computes Voronoi cells and saves the output in POV-Ray
00329                  * format.
00330                  * \param[in] vl the loop class to use.
00331                  * \param[in] fp a file handle to write to. */
00332                 template<class c_loop>
00333                 void draw_cells_pov(c_loop &vl,FILE *fp) {
00334                         voronoicell c;double *pp;
00335                         if(vl.start()) do if(compute_cell(c,vl)) {
00336                                 fprintf(fp,"// cell %d\n",id[vl.ijk][vl.q]);
00337                                 pp=p[vl.ijk]+ps*vl.q;
00338                                 c.draw_pov(*pp,pp[1],pp[2],fp);
00339                         } while(vl.inc());
00340                 }
00341                 /** Computes all Voronoi cells and saves the output in POV-Ray
00342                  * format.
00343                  * \param[in] fp a file handle to write to. */
00344                 inline void draw_cells_pov(FILE *fp=stdout) {
00345                         c_loop_all_periodic vl(*this);
00346                         draw_cells_pov(vl,fp);
00347                 }
00348                 /** Computes all Voronoi cells and saves the output in POV-Ray
00349                  * format.
00350                  * \param[in] filename the name of the file to write to. */
00351                 inline void draw_cells_pov(const char *filename) {
00352                         FILE *fp(safe_fopen(filename,"w"));
00353                         draw_cells_pov(fp);
00354                         fclose(fp);
00355                 }
00356                 /** Computes the Voronoi cells and saves customized information
00357                  * about them.
00358                  * \param[in] vl the loop class to use.
00359                  * \param[in] format the custom output string to use.
00360                  * \param[in] fp a file handle to write to. */
00361                 template<class c_loop>
00362                 void print_custom(c_loop &vl,const char *format,FILE *fp) {
00363                         int ijk,q;double *pp;
00364                         if(contains_neighbor(format)) {
00365                                 voronoicell_neighbor c;
00366                                 if(vl.start()) do if(compute_cell(c,vl)) {
00367                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00368                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],default_radius,fp);
00369                                 } while(vl.inc());
00370                         } else {
00371                                 voronoicell c;
00372                                 if(vl.start()) do if(compute_cell(c,vl)) {
00373                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00374                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],default_radius,fp);
00375                                 } while(vl.inc());
00376                         }
00377                 }
00378                 void print_custom(const char *format,FILE *fp=stdout);
00379                 void print_custom(const char *format,const char *filename);
00380                 bool find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid);
00381                 /** Computes the Voronoi cell for a particle currently being
00382                  * referenced by a loop class.
00383                  * \param[out] c a Voronoi cell class in which to store the
00384                  *               computed cell.
00385                  * \param[in] vl the loop class to use.
00386                  * \return True if the cell was computed. If the cell cannot be
00387                  * computed because it was removed entirely for some reason,
00388                  * then the routine returns false. */
00389                 template<class v_cell,class c_loop>
00390                 inline bool compute_cell(v_cell &c,c_loop &vl) {
00391                         return vc.compute_cell(c,vl.ijk,vl.q,vl.i,vl.j,vl.k);
00392                 }
00393                 /** Computes the Voronoi cell for given particle.
00394                  * \param[out] c a Voronoi cell class in which to store the
00395                  *               computed cell.
00396                  * \param[in] ijk the block that the particle is within.
00397                  * \param[in] q the index of the particle within the block.
00398                  * \return True if the cell was computed. If the cell cannot be
00399                  * computed because it was removed entirely for some reason,
00400                  * then the routine returns false. */
00401                 template<class v_cell>
00402                 inline bool compute_cell(v_cell &c,int ijk,int q) {
00403                         int k(ijk/(nx*oy)),ijkt(ijk-(nx*oy)*k),j(ijkt/nx),i(ijkt-j*nx);
00404                         return vc.compute_cell(c,ijk,q,i,j,k);
00405                 }
00406         private:
00407                 voro_compute<container_periodic> vc;
00408                 inline void r_init(int ijk,int s) {};
00409                 inline double r_cutoff(double lrs) {return lrs;}
00410                 inline double r_max_add(double rs) {return rs;}
00411                 inline double r_current_sub(double rs,int ijk,int q) {return rs;}
00412                 inline double r_scale(double rs,int ijk,int q) {return rs;}
00413                 friend class voro_compute<container_periodic>;
00414 };
00415 
00416 /** \brief Extension of the container_periodic_base class for computing radical
00417  * Voronoi tessellations.
00418  *
00419  * This class is an extension of container_periodic_base that has routines
00420  * specifically for computing the radical Voronoi tessellation that depends
00421  * on the particle radii. */
00422 class container_periodic_poly : public container_periodic_base {
00423         public:
00424                 /** The current maximum radius of any particle, used to
00425                  * determine when to cut off the radical Voronoi computation.
00426                  * */
00427                 double max_radius;
00428                 container_periodic_poly(double bx_,double bxy_,double by_,double bxz_,double byz_,double bz_,
00429                                 int nx_,int ny_,int nz_,int init_mem_);
00430                 void clear();
00431                 void put(int n,double x,double y,double z,double r);
00432                 void put(int n,double x,double y,double z,double r,int &ai,int &aj,int &ak);
00433                 void put(particle_order &vo,int n,double x,double y,double z,double r);
00434                 void import(FILE *fp=stdin);
00435                 void import(particle_order &vo,FILE *fp=stdin);
00436                 /** Imports a list of particles from an open file stream into
00437                  * the container_poly class. Entries of five numbers (Particle
00438                  * ID, x position, y position, z position, radius) are searched
00439                  * for. If the file cannot be successfully read, then the
00440                  * routine causes a fatal error.
00441                  * \param[in] filename the name of the file to open and read
00442                  *                     from. */
00443                 inline void import(const char* filename) {
00444                         FILE *fp(safe_fopen(filename,"r"));
00445                         import(fp);
00446                         fclose(fp);
00447                 }
00448                 /** Imports a list of particles from an open file stream into
00449                  * the container_poly class. Entries of five numbers (Particle
00450                  * ID, x position, y position, z position, radius) are searched
00451                  * for. In addition, the order in which particles are read is
00452                  * saved into an ordering class. If the file cannot be
00453                  * successfully read, then the routine causes a fatal error.
00454                  * \param[in,out] vo the ordering class to use.
00455                  * \param[in] filename the name of the file to open and read
00456                  *                     from. */
00457                 inline void import(particle_order &vo,const char* filename) {
00458                         FILE *fp(safe_fopen(filename,"r"));
00459                         import(vo,fp);
00460                         fclose(fp);
00461                 }
00462                 void compute_all_cells();
00463                 double sum_cell_volumes();
00464                 /** Dumps particle IDs, positions and radii to a file.
00465                  * \param[in] vl the loop class to use.
00466                  * \param[in] fp a file handle to write to. */
00467                 template<class c_loop>
00468                 void draw_particles(c_loop &vl,FILE *fp) {
00469                         double *pp;
00470                         if(vl.start()) do {
00471                                 pp=p[vl.ijk]+4*vl.q;
00472                                 fprintf(fp,"%d %g %g %g %g\n",id[vl.ijk][vl.q],*pp,pp[1],pp[2],pp[3]);
00473                         } while(vl.inc());
00474                 }
00475                 /** Dumps all of the particle IDs, positions and radii to a
00476                  * file.
00477                  * \param[in] fp a file handle to write to. */
00478                 inline void draw_particles(FILE *fp=stdout) {
00479                         c_loop_all_periodic vl(*this);
00480                         draw_particles(vl,fp);
00481                 }
00482                 /** Dumps all of the particle IDs, positions and radii to a
00483                  * file.
00484                  * \param[in] filename the name of the file to write to. */
00485                 inline void draw_particles(const char *filename) {
00486                         FILE *fp(safe_fopen(filename,"w"));
00487                         draw_particles(fp);
00488                         fclose(fp);
00489                 }
00490                 /** Dumps particle positions in POV-Ray format.
00491                  * \param[in] vl the loop class to use.
00492                  * \param[in] fp a file handle to write to. */
00493                 template<class c_loop>
00494                 void draw_particles_pov(c_loop &vl,FILE *fp) {
00495                         double *pp;
00496                         if(vl.start()) do {
00497                                 pp=p[vl.ijk]+4*vl.q;
00498                                 fprintf(fp,"// id %d\nsphere{<%g,%g,%g>,%g}\n",
00499                                                 id[vl.ijk][vl.q],*pp,pp[1],pp[2],pp[3]);
00500                         } while(vl.inc());
00501                 }
00502                 /** Dumps all the particle positions in POV-Ray format.
00503                  * \param[in] fp a file handle to write to. */
00504                 inline void draw_particles_pov(FILE *fp=stdout) {
00505                         c_loop_all_periodic vl(*this);
00506                         draw_particles_pov(vl,fp);
00507                 }
00508                 /** Dumps all the particle positions in POV-Ray format.
00509                  * \param[in] filename the name of the file to write to. */
00510                 inline void draw_particles_pov(const char *filename) {
00511                         FILE *fp(safe_fopen(filename,"w"));
00512                         draw_particles_pov(fp);
00513                         fclose(fp);
00514                 }
00515                 /** Computes Voronoi cells and saves the output in gnuplot
00516                  * format.
00517                  * \param[in] vl the loop class to use.
00518                  * \param[in] fp a file handle to write to. */
00519                 template<class c_loop>
00520                 void draw_cells_gnuplot(c_loop &vl,FILE *fp) {
00521                         voronoicell c;double *pp;
00522                         if(vl.start()) do if(compute_cell(c,vl)) {
00523                                 pp=p[vl.ijk]+ps*vl.q;
00524                                 c.draw_gnuplot(*pp,pp[1],pp[2],fp);
00525                         } while(vl.inc());
00526                 }
00527                 /** Compute all Voronoi cells and saves the output in gnuplot
00528                  * format.
00529                  * \param[in] fp a file handle to write to. */
00530                 inline void draw_cells_gnuplot(FILE *fp=stdout) {
00531                         c_loop_all_periodic vl(*this);
00532                         draw_cells_gnuplot(vl,fp);
00533                 }
00534                 /** Compute all Voronoi cells and saves the output in gnuplot
00535                  * format.
00536                  * \param[in] filename the name of the file to write to. */
00537                 inline void draw_cells_gnuplot(const char *filename) {
00538                         FILE *fp(safe_fopen(filename,"w"));
00539                         draw_cells_gnuplot(fp);
00540                         fclose(fp);
00541                 }
00542                 /** Computes Voronoi cells and saves the output in POV-Ray
00543                  * format.
00544                  * \param[in] vl the loop class to use.
00545                  * \param[in] fp a file handle to write to. */
00546                 template<class c_loop>
00547                 void draw_cells_pov(c_loop &vl,FILE *fp) {
00548                         voronoicell c;double *pp;
00549                         if(vl.start()) do if(compute_cell(c,vl)) {
00550                                 fprintf(fp,"// cell %d\n",id[vl.ijk][vl.q]);
00551                                 pp=p[vl.ijk]+ps*vl.q;
00552                                 c.draw_pov(*pp,pp[1],pp[2],fp);
00553                         } while(vl.inc());
00554                 }
00555                 /** Computes all Voronoi cells and saves the output in POV-Ray
00556                  * format.
00557                  * \param[in] fp a file handle to write to. */
00558                 inline void draw_cells_pov(FILE *fp=stdout) {
00559                         c_loop_all_periodic vl(*this);
00560                         draw_cells_pov(vl,fp);
00561                 }
00562                 /** Computes all Voronoi cells and saves the output in POV-Ray
00563                  * format.
00564                  * \param[in] filename the name of the file to write to. */
00565                 inline void draw_cells_pov(const char *filename) {
00566                         FILE *fp(safe_fopen(filename,"w"));
00567                         draw_cells_pov(fp);
00568                         fclose(fp);
00569                 }
00570                 /** Computes the Voronoi cells and saves customized information
00571                  * about them.
00572                  * \param[in] vl the loop class to use.
00573                  * \param[in] format the custom output string to use.
00574                  * \param[in] fp a file handle to write to. */
00575                 template<class c_loop>
00576                 void print_custom(c_loop &vl,const char *format,FILE *fp) {
00577                         int ijk,q;double *pp;
00578                         if(contains_neighbor(format)) {
00579                                 voronoicell_neighbor c;
00580                                 if(vl.start()) do if(compute_cell(c,vl)) {
00581                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00582                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],pp[3],fp);
00583                                 } while(vl.inc());
00584                         } else {
00585                                 voronoicell c;
00586                                 if(vl.start()) do if(compute_cell(c,vl)) {
00587                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00588                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],pp[3],fp);
00589                                 } while(vl.inc());
00590                         }
00591                 }
00592                 /** Computes the Voronoi cell for a particle currently being
00593                  * referenced by a loop class.
00594                  * \param[out] c a Voronoi cell class in which to store the
00595                  *               computed cell.
00596                  * \param[in] vl the loop class to use.
00597                  * \return True if the cell was computed. If the cell cannot be
00598                  * computed because it was removed entirely for some reason,
00599                  * then the routine returns false. */
00600                 template<class v_cell,class c_loop>
00601                 inline bool compute_cell(v_cell &c,c_loop &vl) {
00602                         return vc.compute_cell(c,vl.ijk,vl.q,vl.i,vl.j,vl.k);
00603                 }
00604                 /** Computes the Voronoi cell for given particle.
00605                  * \param[out] c a Voronoi cell class in which to store the
00606                  *               computed cell.
00607                  * \param[in] ijk the block that the particle is within.
00608                  * \param[in] q the index of the particle within the block.
00609                  * \return True if the cell was computed. If the cell cannot be
00610                  * computed because it was removed entirely for some reason,
00611                  * then the routine returns false. */
00612                 template<class v_cell>
00613                 inline bool compute_cell(v_cell &c,int ijk,int q) {
00614                         int k(ijk/(nx*oy)),ijkt(ijk-(nx*oy)*k),j(ijkt/nx),i(ijkt-j*nx);
00615                         return vc.compute_cell(c,ijk,q,i,j,k);
00616                 }
00617                 void print_custom(const char *format,FILE *fp=stdout);
00618                 void print_custom(const char *format,const char *filename);
00619                 bool find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid);
00620         private:
00621                 voro_compute<container_periodic_poly> vc;
00622                 double r_rad,r_mul;
00623                 inline void r_init(int ijk,int s) {
00624                         r_rad=p[ijk][4*s+3];
00625                         r_mul=1+(r_rad*r_rad-max_radius*max_radius)/((max_radius+r_rad)*(max_radius+r_rad));
00626                         r_rad*=r_rad;
00627                 }
00628                 inline double r_cutoff(double lrs) {return r_mul*lrs;}
00629                 inline double r_max_add(double rs) {return rs+max_radius*max_radius;}
00630                 inline double r_current_sub(double rs,int ijk,int q) {
00631                         return rs-p[ijk][4*q+3]*p[ijk][4*q+3];
00632                 }
00633                 inline double r_scale(double rs,int ijk,int q) {
00634                         return rs+r_rad-p[ijk][4*q+3]*p[ijk][4*q+3];
00635                 }
00636                 friend class voro_compute<container_periodic_poly>;
00637 };
00638 
00639 }
00640 
00641 #endif