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