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