Voro++
container.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.hh
00008  * \brief Header file for the container_base and related classes. */
00009 
00010 #ifndef VOROPP_CONTAINER_HH
00011 #define VOROPP_CONTAINER_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 "rad_option.hh"
00026 
00027 namespace voro {
00028 
00029 /** \brief Pure virtual class from which wall objects are derived.
00030  *
00031  * This is a pure virtual class for a generic wall object. A wall object
00032  * can be specified by deriving a new class from this and specifying the
00033  * functions.*/
00034 class wall {
00035         public:
00036                 virtual ~wall() {}
00037                 /** A pure virtual function for testing whether a point is
00038                  * inside the wall object. */
00039                 virtual bool point_inside(double x,double y,double z) = 0;
00040                 /** A pure virtual function for cutting a cell without
00041                  * neighbor-tracking with a wall. */
00042                 virtual bool cut_cell(voronoicell &c,double x,double y,double z) = 0;
00043                 /** A pure virtual function for cutting a cell with
00044                  * neighbor-tracking enabled with a wall. */
00045                 virtual bool cut_cell(voronoicell_neighbor &c,double x,double y,double z) = 0;
00046 };
00047 
00048 /** \brief A class for storing a list of pointers to walls.
00049  *
00050  * This class stores a list of pointers to wall classes. It contains several
00051  * simple routines that make use of the wall classes (such as telling whether a
00052  * given position is inside all of the walls or not). It can be used by itself,
00053  * but also forms part of container_base, for associating walls with this
00054  * class. */
00055 class wall_list {
00056         public:
00057                 /** An array holding pointers to wall objects. */
00058                 wall **walls;
00059                 /** A pointer to the next free position to add a wall pointer.
00060                  */
00061                 wall **wep;
00062                 wall_list();
00063                 ~wall_list();
00064                 /** Adds a wall to the list.
00065                  * \param[in] w the wall to add. */
00066                 inline void add_wall(wall *w) {
00067                         if(wep==wel) increase_wall_memory();
00068                         *(wep++)=w;
00069                 }
00070                 /** Adds a wall to the list.
00071                  * \param[in] w a reference to the wall to add. */
00072                 inline void add_wall(wall &w) {add_wall(&w);}
00073                 void add_wall(wall_list &wl);
00074                 /** Determines whether a given position is inside all of the
00075                  * walls on the list.
00076                  * \param[in] (x,y,z) the position to test.
00077                  * \return True if it is inside, false if it is outside. */
00078                 inline bool point_inside_walls(double x,double y,double z) {
00079                         for(wall **wp=walls;wp<wep;wp++) if(!((*wp)->point_inside(x,y,z))) return false;
00080                         return true;
00081                 }
00082                 /** Cuts a Voronoi cell by all of the walls currently on
00083                  * the list.
00084                  * \param[in] c a reference to the Voronoi cell class.
00085                  * \param[in] (x,y,z) the position of the cell.
00086                  * \return True if the cell still exists, false if the cell is
00087                  * deleted. */
00088                 template<class c_class>
00089                 bool apply_walls(c_class &c,double x,double y,double z) {
00090                         for(wall **wp=walls;wp<wep;wp++) if(!((*wp)->cut_cell(c,x,y,z))) return false;
00091                         return true;
00092                 }
00093                 void deallocate();
00094         protected:
00095                 void increase_wall_memory();
00096                 /** A pointer to the limit of the walls array, used to
00097                  * determine when array is full. */
00098                 wall **wel;
00099                 /** The current amount of memory allocated for walls. */
00100                 int current_wall_size;
00101 };
00102 
00103 /** \brief Class for representing a particle system in a three-dimensional
00104  * rectangular box.
00105  *
00106  * This class represents a system of particles in a three-dimensional
00107  * rectangular box. Any combination of non-periodic and periodic coordinates
00108  * can be used in the three coordinate directions. The class is not intended
00109  * for direct use, but instead forms the base of the container and
00110  * container_poly classes that add specialized routines for computing the
00111  * regular and radical Voronoi tessellations respectively. It contains routines
00112  * that are commonly between these two classes, such as those for drawing the
00113  * domain, and placing particles within the internal data structure.
00114  *
00115  * The class is derived from the wall_list class, which encapsulates routines
00116  * for associating walls with the container, and the voro_base class, which
00117  * encapsulates routines about the underlying computational grid. */
00118 class container_base : public voro_base, public wall_list {
00119         public:
00120                 /** The minimum x coordinate of the container. */
00121                 const double ax;
00122                 /** The maximum x coordinate of the container. */
00123                 const double bx;
00124                 /** The minimum y coordinate of the container. */
00125                 const double ay;
00126                 /** The maximum y coordinate of the container. */
00127                 const double by;
00128                 /** The minimum z coordinate of the container. */
00129                 const double az;
00130                 /** The maximum z coordinate of the container. */
00131                 const double bz;
00132                 /** A boolean value that determines if the x coordinate in
00133                  * periodic or not. */
00134                 const bool xperiodic;
00135                 /** A boolean value that determines if the y coordinate in
00136                  * periodic or not. */
00137                 const bool yperiodic;
00138                 /** A boolean value that determines if the z coordinate in
00139                  * periodic or not. */
00140                 const bool zperiodic;
00141                 /** This array holds the numerical IDs of each particle in each
00142                  * computational box. */
00143                 int **id;
00144                 /** A two dimensional array holding particle positions. For the
00145                  * derived container_poly class, this also holds particle
00146                  * radii. */
00147                 double **p;
00148                 /** This array holds the number of particles within each
00149                  * computational box of the container. */
00150                 int *co;
00151                 /** This array holds the maximum amount of particle memory for
00152                  * each computational box of the container. If the number of
00153                  * particles in a particular box ever approaches this limit,
00154                  * more is allocated using the add_particle_memory() function.
00155                  */
00156                 int *mem;
00157                 /** The amount of memory in the array structure for each
00158                  * particle. This is set to 3 when the basic class is
00159                  * initialized, so that the array holds (x,y,z) positions. If
00160                  * the container class is initialized as part of the derived
00161                  * class container_poly, then this is set to 4, to also hold
00162                  * the particle radii. */
00163                 const int ps;
00164                 container_base(double ax_,double bx_,double ay_,double by_,double az_,double bz_,
00165                                 int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,
00166                                 int init_mem,int ps_);
00167                 ~container_base();
00168                 bool point_inside(double x,double y,double z);
00169                 void region_count();
00170                 /** Initializes the Voronoi cell prior to a compute_cell
00171                  * operation for a specific particle being carried out by a
00172                  * voro_compute class. The cell is initialized to fill the
00173                  * entire container. For non-periodic coordinates, this is set
00174                  * by the position of the walls. For periodic coordinates, the
00175                  * space is equally divided in either direction from the
00176                  * particle's initial position. Plane cuts made by any walls
00177                  * that have been added are then applied to the cell.
00178                  * \param[in,out] c a reference to a voronoicell object.
00179                  * \param[in] ijk the block that the particle is within.
00180                  * \param[in] q the index of the particle within its block.
00181                  * \param[in] (ci,cj,ck) the coordinates of the block in the
00182                  *                       container coordinate system.
00183                  * \param[out] (i,j,k) the coordinates of the test block
00184                  *                     relative to the voro_compute
00185                  *                     coordinate system.
00186                  * \param[out] (x,y,z) the position of the particle.
00187                  * \param[out] disp a block displacement used internally by the
00188                  *                  compute_cell routine.
00189                  * \return False if the plane cuts applied by walls completely
00190                  * removed the cell, true otherwise. */
00191                 template<class v_cell>
00192                 inline bool initialize_voronoicell(v_cell &c,int ijk,int q,int ci,int cj,int ck,
00193                                 int &i,int &j,int &k,double &x,double &y,double &z,int &disp) {
00194                         double x1,x2,y1,y2,z1,z2,*pp=p[ijk]+ps*q;
00195                         x=*(pp++);y=*(pp++);z=*pp;
00196                         if(xperiodic) {x1=-(x2=0.5*(bx-ax));i=nx;} else {x1=ax-x;x2=bx-x;i=ci;}
00197                         if(yperiodic) {y1=-(y2=0.5*(by-ay));j=ny;} else {y1=ay-y;y2=by-y;j=cj;}
00198                         if(zperiodic) {z1=-(z2=0.5*(bz-az));k=nz;} else {z1=az-z;z2=bz-z;k=ck;}
00199                         c.init(x1,x2,y1,y2,z1,z2);
00200                         if(!apply_walls(c,x,y,z)) return false;
00201                         disp=ijk-i-nx*(j+ny*k);
00202                         return true;
00203                 }
00204                 /** Initializes parameters for a find_voronoi_cell call within
00205                  * the voro_compute template.
00206                  * \param[in] (ci,cj,ck) the coordinates of the test block in
00207                  *                       the container coordinate system.
00208                  * \param[in] ijk the index of the test block
00209                  * \param[out] (i,j,k) the coordinates of the test block
00210                  *                     relative to the voro_compute
00211                  *                     coordinate system.
00212                  * \param[out] disp a block displacement used internally by the
00213                  *                  find_voronoi_cell routine. */
00214                 inline void initialize_search(int ci,int cj,int ck,int ijk,int &i,int &j,int &k,int &disp) {
00215                         i=xperiodic?nx:ci;
00216                         j=yperiodic?ny:cj;
00217                         k=zperiodic?nz:ck;
00218                         disp=ijk-i-nx*(j+ny*k);
00219                 }
00220                 /** Returns the position of a particle currently being computed
00221                  * relative to the computational block that it is within. It is
00222                  * used to select the optimal worklist entry to use.
00223                  * \param[in] (x,y,z) the position of the particle.
00224                  * \param[in] (ci,cj,ck) the block that the particle is within.
00225                  * \param[out] (fx,fy,fz) the position relative to the block.
00226                  */
00227                 inline void frac_pos(double x,double y,double z,double ci,double cj,double ck,
00228                                 double &fx,double &fy,double &fz) {
00229                         fx=x-ax-boxx*ci;
00230                         fy=y-ay-boxy*cj;
00231                         fz=z-az-boxz*ck;
00232                 }
00233                 /** Calculates the index of block in the container structure
00234                  * corresponding to given coordinates.
00235                  * \param[in] (ci,cj,ck) the coordinates of the original block
00236                  *                       in the current computation, relative
00237                  *                       to the container coordinate system.
00238                  * \param[in] (ei,ej,ek) the displacement of the current block
00239                  *                       from the original block.
00240                  * \param[in,out] (qx,qy,qz) the periodic displacement that
00241                  *                           must be added to the particles
00242                  *                           within the computed block.
00243                  * \param[in] disp a block displacement used internally by the
00244                  *                  find_voronoi_cell and compute_cell routines.
00245                  * \return The block index. */
00246                 inline int region_index(int ci,int cj,int ck,int ei,int ej,int ek,double &qx,double &qy,double &qz,int &disp) {
00247                         if(xperiodic) {if(ci+ei<nx) {ei+=nx;qx=-(bx-ax);} else if(ci+ei>=(nx<<1)) {ei-=nx;qx=bx-ax;} else qx=0;}
00248                         if(yperiodic) {if(cj+ej<ny) {ej+=ny;qy=-(by-ay);} else if(cj+ej>=(ny<<1)) {ej-=ny;qy=by-ay;} else qy=0;}
00249                         if(zperiodic) {if(ck+ek<nz) {ek+=nz;qz=-(bz-az);} else if(ck+ek>=(nz<<1)) {ek-=nz;qz=bz-az;} else qz=0;}
00250                         return disp+ei+nx*(ej+ny*ek);
00251                 }
00252                 void draw_domain_gnuplot(FILE *fp=stdout);
00253                 /** Draws an outline of the domain in Gnuplot format.
00254                  * \param[in] filename the filename to write to. */
00255                 inline void draw_domain_gnuplot(const char* filename) {
00256                         FILE *fp=safe_fopen(filename,"w");
00257                         draw_domain_gnuplot(fp);
00258                         fclose(fp);
00259                 }
00260                 void draw_domain_pov(FILE *fp=stdout);
00261                 /** Draws an outline of the domain in Gnuplot format.
00262                  * \param[in] filename the filename to write to. */
00263                 inline void draw_domain_pov(const char* filename) {
00264                         FILE *fp=safe_fopen(filename,"w");
00265                         draw_domain_pov(fp);
00266                         fclose(fp);
00267                 }
00268                 /** Sums up the total number of stored particles.
00269                  * \return The number of particles. */
00270                 inline int total_particles() {
00271                         int tp=*co;
00272                         for(int *cop=co+1;cop<co+nxyz;cop++) tp+=*cop;
00273                         return tp;
00274                 }
00275         protected:
00276                 void add_particle_memory(int i);
00277                 inline bool put_locate_block(int &ijk,double &x,double &y,double &z);
00278                 inline bool put_remap(int &ijk,double &x,double &y,double &z);
00279                 inline bool remap(int &ai,int &aj,int &ak,int &ci,int &cj,int &ck,double &x,double &y,double &z,int &ijk);
00280 };
00281 
00282 /** \brief Extension of the container_base class for computing regular Voronoi
00283  * tessellations.
00284  *
00285  * This class is an extension of the container_base class that has routines
00286  * specifically for computing the regular Voronoi tessellation with no
00287  * dependence on particle radii. */
00288 class container : public container_base, public radius_mono {
00289         public:
00290                 container(double ax_,double bx_,double ay_,double by_,double az_,double bz_,
00291                                 int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem);
00292                 void clear();
00293                 void put(int n,double x,double y,double z);
00294                 void put(particle_order &vo,int n,double x,double y,double z);
00295                 void import(FILE *fp=stdin);
00296                 void import(particle_order &vo,FILE *fp=stdin);
00297                 /** Imports a list of particles from an open file stream into
00298                  * the container. Entries of four numbers (Particle ID, x
00299                  * position, y position, z position) are searched for. If the
00300                  * file cannot be successfully read, then the routine causes a
00301                  * fatal error.
00302                  * \param[in] filename the name of the file to open and read
00303                  *                     from. */
00304                 inline void import(const char* filename) {
00305                         FILE *fp=safe_fopen(filename,"r");
00306                         import(fp);
00307                         fclose(fp);
00308                 }
00309                 /** Imports a list of particles from an open file stream into
00310                  * the container. Entries of four numbers (Particle ID, x
00311                  * position, y position, z position) are searched for. In
00312                  * addition, the order in which particles are read is saved
00313                  * into an ordering class. If the file cannot be successfully
00314                  * read, then the routine causes a fatal error.
00315                  * \param[in,out] vo the ordering class to use.
00316                  * \param[in] filename the name of the file to open and read
00317                  *                     from. */
00318                 inline void import(particle_order &vo,const char* filename) {
00319                         FILE *fp=safe_fopen(filename,"r");
00320                         import(vo,fp);
00321                         fclose(fp);
00322                 }
00323                 void compute_all_cells();
00324                 double sum_cell_volumes();
00325                 /** Dumps particle IDs and positions to a file.
00326                  * \param[in] vl the loop class to use.
00327                  * \param[in] fp a file handle to write to. */
00328                 template<class c_loop>
00329                 void draw_particles(c_loop &vl,FILE *fp) {
00330                         double *pp;
00331                         if(vl.start()) do {
00332                                 pp=p[vl.ijk]+3*vl.q;
00333                                 fprintf(fp,"%d %g %g %g\n",id[vl.ijk][vl.q],*pp,pp[1],pp[2]);
00334                         } while(vl.inc());
00335                 }
00336                 /** Dumps all of the particle IDs and positions to a file.
00337                  * \param[in] fp a file handle to write to. */
00338                 inline void draw_particles(FILE *fp=stdout) {
00339                         c_loop_all vl(*this);
00340                         draw_particles(vl,fp);
00341                 }
00342                 /** Dumps all of the particle IDs and positions to a file.
00343                  * \param[in] filename the name of the file to write to. */
00344                 inline void draw_particles(const char *filename) {
00345                         FILE *fp=safe_fopen(filename,"w");
00346                         draw_particles(fp);
00347                         fclose(fp);
00348                 }
00349                 /** Dumps particle positions in POV-Ray format.
00350                  * \param[in] vl the loop class to use.
00351                  * \param[in] fp a file handle to write to. */
00352                 template<class c_loop>
00353                 void draw_particles_pov(c_loop &vl,FILE *fp) {
00354                         double *pp;
00355                         if(vl.start()) do {
00356                                 pp=p[vl.ijk]+3*vl.q;
00357                                 fprintf(fp,"// id %d\nsphere{<%g,%g,%g>,s}\n",
00358                                                 id[vl.ijk][vl.q],*pp,pp[1],pp[2]);
00359                         } while(vl.inc());
00360                 }
00361                 /** Dumps all particle positions in POV-Ray format.
00362                  * \param[in] fp a file handle to write to. */
00363                 inline void draw_particles_pov(FILE *fp=stdout) {
00364                         c_loop_all vl(*this);
00365                         draw_particles_pov(vl,fp);
00366                 }
00367                 /** Dumps all particle positions in POV-Ray format.
00368                  * \param[in] filename the name of the file to write to. */
00369                 inline void draw_particles_pov(const char *filename) {
00370                         FILE *fp=safe_fopen(filename,"w");
00371                         draw_particles_pov(fp);
00372                         fclose(fp);
00373                 }
00374                 /** Computes Voronoi cells and saves the output in gnuplot
00375                  * format.
00376                  * \param[in] vl the loop class to use.
00377                  * \param[in] fp a file handle to write to. */
00378                 template<class c_loop>
00379                 void draw_cells_gnuplot(c_loop &vl,FILE *fp) {
00380                         voronoicell c;double *pp;
00381                         if(vl.start()) do if(compute_cell(c,vl)) {
00382                                 pp=p[vl.ijk]+ps*vl.q;
00383                                 c.draw_gnuplot(*pp,pp[1],pp[2],fp);
00384                         } while(vl.inc());
00385                 }
00386                 /** Computes all Voronoi cells and saves the output in gnuplot
00387                  * format.
00388                  * \param[in] fp a file handle to write to. */
00389                 inline void draw_cells_gnuplot(FILE *fp=stdout) {
00390                         c_loop_all vl(*this);
00391                         draw_cells_gnuplot(vl,fp);
00392                 }
00393                 /** Computes all Voronoi cells and saves the output in gnuplot
00394                  * format.
00395                  * \param[in] filename the name of the file to write to. */
00396                 inline void draw_cells_gnuplot(const char *filename) {
00397                         FILE *fp=safe_fopen(filename,"w");
00398                         draw_cells_gnuplot(fp);
00399                         fclose(fp);
00400                 }
00401                 /** Computes Voronoi cells and saves the output in POV-Ray
00402                  * format.
00403                  * \param[in] vl the loop class to use.
00404                  * \param[in] fp a file handle to write to. */
00405                 template<class c_loop>
00406                 void draw_cells_pov(c_loop &vl,FILE *fp) {
00407                         voronoicell c;double *pp;
00408                         if(vl.start()) do if(compute_cell(c,vl)) {
00409                                 fprintf(fp,"// cell %d\n",id[vl.ijk][vl.q]);
00410                                 pp=p[vl.ijk]+ps*vl.q;
00411                                 c.draw_pov(*pp,pp[1],pp[2],fp);
00412                         } while(vl.inc());
00413                 }
00414                 /** Computes all Voronoi cells and saves the output in POV-Ray
00415                  * format.
00416                  * \param[in] fp a file handle to write to. */
00417                 inline void draw_cells_pov(FILE *fp=stdout) {
00418                         c_loop_all vl(*this);
00419                         draw_cells_pov(vl,fp);
00420                 }
00421                 /** Computes all Voronoi cells and saves the output in POV-Ray
00422                  * format.
00423                  * \param[in] filename the name of the file to write to. */
00424                 inline void draw_cells_pov(const char *filename) {
00425                         FILE *fp=safe_fopen(filename,"w");
00426                         draw_cells_pov(fp);
00427                         fclose(fp);
00428                 }
00429                 /** Computes the Voronoi cells and saves customized information
00430                  * about them.
00431                  * \param[in] vl the loop class to use.
00432                  * \param[in] format the custom output string to use.
00433                  * \param[in] fp a file handle to write to. */
00434                 template<class c_loop>
00435                 void print_custom(c_loop &vl,const char *format,FILE *fp) {
00436                         int ijk,q;double *pp;
00437                         if(contains_neighbor(format)) {
00438                                 voronoicell_neighbor c;
00439                                 if(vl.start()) do if(compute_cell(c,vl)) {
00440                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00441                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],default_radius,fp);
00442                                 } while(vl.inc());
00443                         } else {
00444                                 voronoicell c;
00445                                 if(vl.start()) do if(compute_cell(c,vl)) {
00446                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00447                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],default_radius,fp);
00448                                 } while(vl.inc());
00449                         }
00450                 }
00451                 void print_custom(const char *format,FILE *fp=stdout);
00452                 void print_custom(const char *format,const char *filename);
00453                 bool find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid);
00454                 /** Computes the Voronoi cell for a particle currently being
00455                  * referenced by a loop class.
00456                  * \param[out] c a Voronoi cell class in which to store the
00457                  *               computed cell.
00458                  * \param[in] vl the loop class to use.
00459                  * \return True if the cell was computed. If the cell cannot be
00460                  * computed, if it is removed entirely by a wall or boundary
00461                  * condition, then the routine returns false. */
00462                 template<class v_cell,class c_loop>
00463                 inline bool compute_cell(v_cell &c,c_loop &vl) {
00464                         return vc.compute_cell(c,vl.ijk,vl.q,vl.i,vl.j,vl.k);
00465                 }
00466                 /** Computes the Voronoi cell for given particle.
00467                  * \param[out] c a Voronoi cell class in which to store the
00468                  *               computed cell.
00469                  * \param[in] ijk the block that the particle is within.
00470                  * \param[in] q the index of the particle within the block.
00471                  * \return True if the cell was computed. If the cell cannot be
00472                  * computed, if it is removed entirely by a wall or boundary
00473                  * condition, then the routine returns false. */
00474                 template<class v_cell>
00475                 inline bool compute_cell(v_cell &c,int ijk,int q) {
00476                         int k=ijk/nxy,ijkt=ijk-nxy*k,j=ijkt/nx,i=ijkt-j*nx;
00477                         return vc.compute_cell(c,ijk,q,i,j,k);
00478                 }
00479         private:
00480                 voro_compute<container> vc;
00481                 friend class voro_compute<container>;
00482 };
00483 
00484 /** \brief Extension of the container_base class for computing radical Voronoi
00485  * tessellations.
00486  *
00487  * This class is an extension of container_base class that has routines
00488  * specifically for computing the radical Voronoi tessellation that depends on
00489  * the particle radii. */
00490 class container_poly : public container_base, public radius_poly {
00491         public:
00492                 container_poly(double ax_,double bx_,double ay_,double by_,double az_,double bz_,
00493                                 int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem);
00494                 void clear();
00495                 void put(int n,double x,double y,double z,double r);
00496                 void put(particle_order &vo,int n,double x,double y,double z,double r);
00497                 void import(FILE *fp=stdin);
00498                 void import(particle_order &vo,FILE *fp=stdin);
00499                 /** Imports a list of particles from an open file stream into
00500                  * the container_poly class. Entries of five numbers (Particle
00501                  * ID, x position, y position, z position, radius) are searched
00502                  * for. If the file cannot be successfully read, then the
00503                  * routine causes a fatal error.
00504                  * \param[in] filename the name of the file to open and read
00505                  *                     from. */
00506                 inline void import(const char* filename) {
00507                         FILE *fp=safe_fopen(filename,"r");
00508                         import(fp);
00509                         fclose(fp);
00510                 }
00511                 /** Imports a list of particles from an open file stream into
00512                  * the container_poly class. Entries of five numbers (Particle
00513                  * ID, x position, y position, z position, radius) are searched
00514                  * for. In addition, the order in which particles are read is
00515                  * saved into an ordering class. If the file cannot be
00516                  * successfully read, then the routine causes a fatal error.
00517                  * \param[in,out] vo the ordering class to use.
00518                  * \param[in] filename the name of the file to open and read
00519                  *                     from. */
00520                 inline void import(particle_order &vo,const char* filename) {
00521                         FILE *fp=safe_fopen(filename,"r");
00522                         import(vo,fp);
00523                         fclose(fp);
00524                 }
00525                 void compute_all_cells();
00526                 double sum_cell_volumes();
00527                 /** Dumps particle IDs, positions and radii to a file.
00528                  * \param[in] vl the loop class to use.
00529                  * \param[in] fp a file handle to write to. */
00530                 template<class c_loop>
00531                 void draw_particles(c_loop &vl,FILE *fp) {
00532                         double *pp;
00533                         if(vl.start()) do {
00534                                 pp=p[vl.ijk]+4*vl.q;
00535                                 fprintf(fp,"%d %g %g %g %g\n",id[vl.ijk][vl.q],*pp,pp[1],pp[2],pp[3]);
00536                         } while(vl.inc());
00537                 }
00538                 /** Dumps all of the particle IDs, positions and radii to a
00539                  * file.
00540                  * \param[in] fp a file handle to write to. */
00541                 inline void draw_particles(FILE *fp=stdout) {
00542                         c_loop_all vl(*this);
00543                         draw_particles(vl,fp);
00544                 }
00545                 /** Dumps all of the particle IDs, positions and radii to a
00546                  * file.
00547                  * \param[in] filename the name of the file to write to. */
00548                 inline void draw_particles(const char *filename) {
00549                         FILE *fp=safe_fopen(filename,"w");
00550                         draw_particles(fp);
00551                         fclose(fp);
00552                 }
00553                 /** Dumps particle positions in POV-Ray format.
00554                  * \param[in] vl the loop class to use.
00555                  * \param[in] fp a file handle to write to. */
00556                 template<class c_loop>
00557                 void draw_particles_pov(c_loop &vl,FILE *fp) {
00558                         double *pp;
00559                         if(vl.start()) do {
00560                                 pp=p[vl.ijk]+4*vl.q;
00561                                 fprintf(fp,"// id %d\nsphere{<%g,%g,%g>,%g}\n",
00562                                                 id[vl.ijk][vl.q],*pp,pp[1],pp[2],pp[3]);
00563                         } while(vl.inc());
00564                 }
00565                 /** Dumps all the particle positions in POV-Ray format.
00566                  * \param[in] fp a file handle to write to. */
00567                 inline void draw_particles_pov(FILE *fp=stdout) {
00568                         c_loop_all vl(*this);
00569                         draw_particles_pov(vl,fp);
00570                 }
00571                 /** Dumps all the particle positions in POV-Ray format.
00572                  * \param[in] filename the name of the file to write to. */
00573                 inline void draw_particles_pov(const char *filename) {
00574                         FILE *fp=safe_fopen(filename,"w");
00575                         draw_particles_pov(fp);
00576                         fclose(fp);
00577                 }
00578                 /** Computes Voronoi cells and saves the output in gnuplot
00579                  * format.
00580                  * \param[in] vl the loop class to use.
00581                  * \param[in] fp a file handle to write to. */
00582                 template<class c_loop>
00583                 void draw_cells_gnuplot(c_loop &vl,FILE *fp) {
00584                         voronoicell c;double *pp;
00585                         if(vl.start()) do if(compute_cell(c,vl)) {
00586                                 pp=p[vl.ijk]+ps*vl.q;
00587                                 c.draw_gnuplot(*pp,pp[1],pp[2],fp);
00588                         } while(vl.inc());
00589                 }
00590                 /** Compute all Voronoi cells and saves the output in gnuplot
00591                  * format.
00592                  * \param[in] fp a file handle to write to. */
00593                 inline void draw_cells_gnuplot(FILE *fp=stdout) {
00594                         c_loop_all vl(*this);
00595                         draw_cells_gnuplot(vl,fp);
00596                 }
00597                 /** Compute all Voronoi cells and saves the output in gnuplot
00598                  * format.
00599                  * \param[in] filename the name of the file to write to. */
00600                 inline void draw_cells_gnuplot(const char *filename) {
00601                         FILE *fp=safe_fopen(filename,"w");
00602                         draw_cells_gnuplot(fp);
00603                         fclose(fp);
00604                 }
00605                 /** Computes Voronoi cells and saves the output in POV-Ray
00606                  * format.
00607                  * \param[in] vl the loop class to use.
00608                  * \param[in] fp a file handle to write to. */
00609                 template<class c_loop>
00610                 void draw_cells_pov(c_loop &vl,FILE *fp) {
00611                         voronoicell c;double *pp;
00612                         if(vl.start()) do if(compute_cell(c,vl)) {
00613                                 fprintf(fp,"// cell %d\n",id[vl.ijk][vl.q]);
00614                                 pp=p[vl.ijk]+ps*vl.q;
00615                                 c.draw_pov(*pp,pp[1],pp[2],fp);
00616                         } while(vl.inc());
00617                 }
00618                 /** Computes all Voronoi cells and saves the output in POV-Ray
00619                  * format.
00620                  * \param[in] fp a file handle to write to. */
00621                 inline void draw_cells_pov(FILE *fp=stdout) {
00622                         c_loop_all vl(*this);
00623                         draw_cells_pov(vl,fp);
00624                 }
00625                 /** Computes all Voronoi cells and saves the output in POV-Ray
00626                  * format.
00627                  * \param[in] filename the name of the file to write to. */
00628                 inline void draw_cells_pov(const char *filename) {
00629                         FILE *fp=safe_fopen(filename,"w");
00630                         draw_cells_pov(fp);
00631                         fclose(fp);
00632                 }
00633                 /** Computes the Voronoi cells and saves customized information
00634                  * about them.
00635                  * \param[in] vl the loop class to use.
00636                  * \param[in] format the custom output string to use.
00637                  * \param[in] fp a file handle to write to. */
00638                 template<class c_loop>
00639                 void print_custom(c_loop &vl,const char *format,FILE *fp) {
00640                         int ijk,q;double *pp;
00641                         if(contains_neighbor(format)) {
00642                                 voronoicell_neighbor c;
00643                                 if(vl.start()) do if(compute_cell(c,vl)) {
00644                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00645                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],pp[3],fp);
00646                                 } while(vl.inc());
00647                         } else {
00648                                 voronoicell c;
00649                                 if(vl.start()) do if(compute_cell(c,vl)) {
00650                                         ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q;
00651                                         c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],pp[3],fp);
00652                                 } while(vl.inc());
00653                         }
00654                 }
00655                 /** Computes the Voronoi cell for a particle currently being
00656                  * referenced by a loop class.
00657                  * \param[out] c a Voronoi cell class in which to store the
00658                  *               computed cell.
00659                  * \param[in] vl the loop class to use.
00660                  * \return True if the cell was computed. If the cell cannot be
00661                  * computed, if it is removed entirely by a wall or boundary
00662                  * condition, then the routine returns false. */
00663                 template<class v_cell,class c_loop>
00664                 inline bool compute_cell(v_cell &c,c_loop &vl) {
00665                         return vc.compute_cell(c,vl.ijk,vl.q,vl.i,vl.j,vl.k);
00666                 }
00667                 /** Computes the Voronoi cell for given particle.
00668                  * \param[out] c a Voronoi cell class in which to store the
00669                  *               computed cell.
00670                  * \param[in] ijk the block that the particle is within.
00671                  * \param[in] q the index of the particle within the block.
00672                  * \return True if the cell was computed. If the cell cannot be
00673                  * computed, if it is removed entirely by a wall or boundary
00674                  * condition, then the routine returns false. */
00675                 template<class v_cell>
00676                 inline bool compute_cell(v_cell &c,int ijk,int q) {
00677                         int k=ijk/nxy,ijkt=ijk-nxy*k,j=ijkt/nx,i=ijkt-j*nx;
00678                         return vc.compute_cell(c,ijk,q,i,j,k);
00679                 }
00680                 void print_custom(const char *format,FILE *fp=stdout);
00681                 void print_custom(const char *format,const char *filename);
00682                 bool find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid);
00683         private:
00684                 voro_compute<container_poly> vc;
00685                 friend class voro_compute<container_poly>;
00686 };
00687 
00688 }
00689 
00690 #endif