1 |
francois |
979 |
// Voro++, a 3D cell-based Voronoi library |
2 |
|
|
// |
3 |
|
|
// Author : Chris H. Rycroft (LBL / UC Berkeley) |
4 |
|
|
// Email : chr@alum.mit.edu |
5 |
|
|
// Date : August 30th 2011 |
6 |
|
|
|
7 |
|
|
/** \file container.hh |
8 |
|
|
* \brief Header file for the container_base and related classes. */ |
9 |
|
|
|
10 |
|
|
#ifndef VOROPP_CONTAINER_HH |
11 |
|
|
#define VOROPP_CONTAINER_HH |
12 |
|
|
|
13 |
|
|
#include <cstdio> |
14 |
|
|
#include <vector> |
15 |
|
|
|
16 |
|
|
#include "config.hh" |
17 |
|
|
#include "common.hh" |
18 |
|
|
#include "v_base.hh" |
19 |
|
|
#include "cell.hh" |
20 |
|
|
#include "c_loops.hh" |
21 |
|
|
#include "v_compute.hh" |
22 |
|
|
#include "rad_option.hh" |
23 |
|
|
|
24 |
|
|
namespace voro { |
25 |
|
|
|
26 |
|
|
/** \brief Pure virtual class from which wall objects are derived. |
27 |
|
|
* |
28 |
|
|
* This is a pure virtual class for a generic wall object. A wall object |
29 |
|
|
* can be specified by deriving a new class from this and specifying the |
30 |
|
|
* functions.*/ |
31 |
|
|
class wall { |
32 |
|
|
public: |
33 |
|
|
virtual ~wall() {} |
34 |
|
|
/** A pure virtual function for testing whether a point is |
35 |
|
|
* inside the wall object. */ |
36 |
|
|
virtual bool point_inside(double x,double y,double z) = 0; |
37 |
|
|
/** A pure virtual function for cutting a cell without |
38 |
|
|
* neighbor-tracking with a wall. */ |
39 |
|
|
virtual bool cut_cell(voronoicell &c,double x,double y,double z) = 0; |
40 |
|
|
/** A pure virtual function for cutting a cell with |
41 |
|
|
* neighbor-tracking enabled with a wall. */ |
42 |
|
|
virtual bool cut_cell(voronoicell_neighbor &c,double x,double y,double z) = 0; |
43 |
|
|
}; |
44 |
|
|
|
45 |
|
|
/** \brief A class for storing a list of pointers to walls. |
46 |
|
|
* |
47 |
|
|
* This class stores a list of pointers to wall classes. It contains several |
48 |
|
|
* simple routines that make use of the wall classes (such as telling whether a |
49 |
|
|
* given position is inside all of the walls or not). It can be used by itself, |
50 |
|
|
* but also forms part of container_base, for associating walls with this |
51 |
|
|
* class. */ |
52 |
|
|
class wall_list { |
53 |
|
|
public: |
54 |
|
|
/** An array holding pointers to wall objects. */ |
55 |
|
|
wall **walls; |
56 |
|
|
/** A pointer to the next free position to add a wall pointer. |
57 |
|
|
*/ |
58 |
|
|
wall **wep; |
59 |
|
|
wall_list(); |
60 |
|
|
~wall_list(); |
61 |
|
|
/** Adds a wall to the list. |
62 |
|
|
* \param[in] w the wall to add. */ |
63 |
|
|
inline void add_wall(wall *w) { |
64 |
|
|
if(wep==wel) increase_wall_memory(); |
65 |
|
|
*(wep++)=w; |
66 |
|
|
} |
67 |
|
|
/** Adds a wall to the list. |
68 |
|
|
* \param[in] w a reference to the wall to add. */ |
69 |
|
|
inline void add_wall(wall &w) {add_wall(&w);} |
70 |
|
|
void add_wall(wall_list &wl); |
71 |
|
|
/** Determines whether a given position is inside all of the |
72 |
|
|
* walls on the list. |
73 |
|
|
* \param[in] (x,y,z) the position to test. |
74 |
|
|
* \return True if it is inside, false if it is outside. */ |
75 |
|
|
inline bool point_inside_walls(double x,double y,double z) { |
76 |
|
|
for(wall **wp=walls;wp<wep;wp++) if(!((*wp)->point_inside(x,y,z))) return false; |
77 |
|
|
return true; |
78 |
|
|
} |
79 |
|
|
/** Cuts a Voronoi cell by all of the walls currently on |
80 |
|
|
* the list. |
81 |
|
|
* \param[in] c a reference to the Voronoi cell class. |
82 |
|
|
* \param[in] (x,y,z) the position of the cell. |
83 |
|
|
* \return True if the cell still exists, false if the cell is |
84 |
|
|
* deleted. */ |
85 |
|
|
template<class c_class> |
86 |
|
|
bool apply_walls(c_class &c,double x,double y,double z) { |
87 |
|
|
for(wall **wp=walls;wp<wep;wp++) if(!((*wp)->cut_cell(c,x,y,z))) return false; |
88 |
|
|
return true; |
89 |
|
|
} |
90 |
|
|
void deallocate(); |
91 |
|
|
protected: |
92 |
|
|
void increase_wall_memory(); |
93 |
|
|
/** A pointer to the limit of the walls array, used to |
94 |
|
|
* determine when array is full. */ |
95 |
|
|
wall **wel; |
96 |
|
|
/** The current amount of memory allocated for walls. */ |
97 |
|
|
int current_wall_size; |
98 |
|
|
}; |
99 |
|
|
|
100 |
|
|
/** \brief Class for representing a particle system in a three-dimensional |
101 |
|
|
* rectangular box. |
102 |
|
|
* |
103 |
|
|
* This class represents a system of particles in a three-dimensional |
104 |
|
|
* rectangular box. Any combination of non-periodic and periodic coordinates |
105 |
|
|
* can be used in the three coordinate directions. The class is not intended |
106 |
|
|
* for direct use, but instead forms the base of the container and |
107 |
|
|
* container_poly classes that add specialized routines for computing the |
108 |
|
|
* regular and radical Voronoi tessellations respectively. It contains routines |
109 |
|
|
* that are commonly between these two classes, such as those for drawing the |
110 |
|
|
* domain, and placing particles within the internal data structure. |
111 |
|
|
* |
112 |
|
|
* The class is derived from the wall_list class, which encapsulates routines |
113 |
|
|
* for associating walls with the container, and the voro_base class, which |
114 |
|
|
* encapsulates routines about the underlying computational grid. */ |
115 |
|
|
class container_base : public voro_base, public wall_list { |
116 |
|
|
public: |
117 |
|
|
/** The minimum x coordinate of the container. */ |
118 |
|
|
const double ax; |
119 |
|
|
/** The maximum x coordinate of the container. */ |
120 |
|
|
const double bx; |
121 |
|
|
/** The minimum y coordinate of the container. */ |
122 |
|
|
const double ay; |
123 |
|
|
/** The maximum y coordinate of the container. */ |
124 |
|
|
const double by; |
125 |
|
|
/** The minimum z coordinate of the container. */ |
126 |
|
|
const double az; |
127 |
|
|
/** The maximum z coordinate of the container. */ |
128 |
|
|
const double bz; |
129 |
|
|
/** A boolean value that determines if the x coordinate in |
130 |
|
|
* periodic or not. */ |
131 |
|
|
const bool xperiodic; |
132 |
|
|
/** A boolean value that determines if the y coordinate in |
133 |
|
|
* periodic or not. */ |
134 |
|
|
const bool yperiodic; |
135 |
|
|
/** A boolean value that determines if the z coordinate in |
136 |
|
|
* periodic or not. */ |
137 |
|
|
const bool zperiodic; |
138 |
|
|
/** This array holds the numerical IDs of each particle in each |
139 |
|
|
* computational box. */ |
140 |
|
|
int **id; |
141 |
|
|
/** A two dimensional array holding particle positions. For the |
142 |
|
|
* derived container_poly class, this also holds particle |
143 |
|
|
* radii. */ |
144 |
|
|
double **p; |
145 |
|
|
/** This array holds the number of particles within each |
146 |
|
|
* computational box of the container. */ |
147 |
|
|
int *co; |
148 |
|
|
/** This array holds the maximum amount of particle memory for |
149 |
|
|
* each computational box of the container. If the number of |
150 |
|
|
* particles in a particular box ever approaches this limit, |
151 |
|
|
* more is allocated using the add_particle_memory() function. |
152 |
|
|
*/ |
153 |
|
|
int *mem; |
154 |
|
|
/** The amount of memory in the array structure for each |
155 |
|
|
* particle. This is set to 3 when the basic class is |
156 |
|
|
* initialized, so that the array holds (x,y,z) positions. If |
157 |
|
|
* the container class is initialized as part of the derived |
158 |
|
|
* class container_poly, then this is set to 4, to also hold |
159 |
|
|
* the particle radii. */ |
160 |
|
|
const int ps; |
161 |
|
|
container_base(double ax_,double bx_,double ay_,double by_,double az_,double bz_, |
162 |
|
|
int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_, |
163 |
|
|
int init_mem,int ps_); |
164 |
|
|
~container_base(); |
165 |
|
|
bool point_inside(double x,double y,double z); |
166 |
|
|
void region_count(); |
167 |
|
|
/** Initializes the Voronoi cell prior to a compute_cell |
168 |
|
|
* operation for a specific particle being carried out by a |
169 |
|
|
* voro_compute class. The cell is initialized to fill the |
170 |
|
|
* entire container. For non-periodic coordinates, this is set |
171 |
|
|
* by the position of the walls. For periodic coordinates, the |
172 |
|
|
* space is equally divided in either direction from the |
173 |
|
|
* particle's initial position. Plane cuts made by any walls |
174 |
|
|
* that have been added are then applied to the cell. |
175 |
|
|
* \param[in,out] c a reference to a voronoicell object. |
176 |
|
|
* \param[in] ijk the block that the particle is within. |
177 |
|
|
* \param[in] q the index of the particle within its block. |
178 |
|
|
* \param[in] (ci,cj,ck) the coordinates of the block in the |
179 |
|
|
* container coordinate system. |
180 |
|
|
* \param[out] (i,j,k) the coordinates of the test block |
181 |
|
|
* relative to the voro_compute |
182 |
|
|
* coordinate system. |
183 |
|
|
* \param[out] (x,y,z) the position of the particle. |
184 |
|
|
* \param[out] disp a block displacement used internally by the |
185 |
|
|
* compute_cell routine. |
186 |
|
|
* \return False if the plane cuts applied by walls completely |
187 |
|
|
* removed the cell, true otherwise. */ |
188 |
|
|
template<class v_cell> |
189 |
|
|
inline bool initialize_voronoicell(v_cell &c,int ijk,int q,int ci,int cj,int ck, |
190 |
|
|
int &i,int &j,int &k,double &x,double &y,double &z,int &disp) { |
191 |
|
|
double x1,x2,y1,y2,z1,z2,*pp=p[ijk]+ps*q; |
192 |
|
|
x=*(pp++);y=*(pp++);z=*pp; |
193 |
|
|
if(xperiodic) {x1=-(x2=0.5*(bx-ax));i=nx;} else {x1=ax-x;x2=bx-x;i=ci;} |
194 |
|
|
if(yperiodic) {y1=-(y2=0.5*(by-ay));j=ny;} else {y1=ay-y;y2=by-y;j=cj;} |
195 |
|
|
if(zperiodic) {z1=-(z2=0.5*(bz-az));k=nz;} else {z1=az-z;z2=bz-z;k=ck;} |
196 |
|
|
c.init(x1,x2,y1,y2,z1,z2); |
197 |
|
|
if(!apply_walls(c,x,y,z)) return false; |
198 |
|
|
disp=ijk-i-nx*(j+ny*k); |
199 |
|
|
return true; |
200 |
|
|
} |
201 |
|
|
/** Initializes parameters for a find_voronoi_cell call within |
202 |
|
|
* the voro_compute template. |
203 |
|
|
* \param[in] (ci,cj,ck) the coordinates of the test block in |
204 |
|
|
* the container coordinate system. |
205 |
|
|
* \param[in] ijk the index of the test block |
206 |
|
|
* \param[out] (i,j,k) the coordinates of the test block |
207 |
|
|
* relative to the voro_compute |
208 |
|
|
* coordinate system. |
209 |
|
|
* \param[out] disp a block displacement used internally by the |
210 |
|
|
* find_voronoi_cell routine. */ |
211 |
|
|
inline void initialize_search(int ci,int cj,int ck,int ijk,int &i,int &j,int &k,int &disp) { |
212 |
|
|
i=xperiodic?nx:ci; |
213 |
|
|
j=yperiodic?ny:cj; |
214 |
|
|
k=zperiodic?nz:ck; |
215 |
|
|
disp=ijk-i-nx*(j+ny*k); |
216 |
|
|
} |
217 |
|
|
/** Returns the position of a particle currently being computed |
218 |
|
|
* relative to the computational block that it is within. It is |
219 |
|
|
* used to select the optimal worklist entry to use. |
220 |
|
|
* \param[in] (x,y,z) the position of the particle. |
221 |
|
|
* \param[in] (ci,cj,ck) the block that the particle is within. |
222 |
|
|
* \param[out] (fx,fy,fz) the position relative to the block. |
223 |
|
|
*/ |
224 |
|
|
inline void frac_pos(double x,double y,double z,double ci,double cj,double ck, |
225 |
|
|
double &fx,double &fy,double &fz) { |
226 |
|
|
fx=x-ax-boxx*ci; |
227 |
|
|
fy=y-ay-boxy*cj; |
228 |
|
|
fz=z-az-boxz*ck; |
229 |
|
|
} |
230 |
|
|
/** Calculates the index of block in the container structure |
231 |
|
|
* corresponding to given coordinates. |
232 |
|
|
* \param[in] (ci,cj,ck) the coordinates of the original block |
233 |
|
|
* in the current computation, relative |
234 |
|
|
* to the container coordinate system. |
235 |
|
|
* \param[in] (ei,ej,ek) the displacement of the current block |
236 |
|
|
* from the original block. |
237 |
|
|
* \param[in,out] (qx,qy,qz) the periodic displacement that |
238 |
|
|
* must be added to the particles |
239 |
|
|
* within the computed block. |
240 |
|
|
* \param[in] disp a block displacement used internally by the |
241 |
|
|
* find_voronoi_cell and compute_cell routines. |
242 |
|
|
* \return The block index. */ |
243 |
|
|
inline int region_index(int ci,int cj,int ck,int ei,int ej,int ek,double &qx,double &qy,double &qz,int &disp) { |
244 |
|
|
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;} |
245 |
|
|
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;} |
246 |
|
|
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;} |
247 |
|
|
return disp+ei+nx*(ej+ny*ek); |
248 |
|
|
} |
249 |
|
|
void draw_domain_gnuplot(FILE *fp=stdout); |
250 |
|
|
/** Draws an outline of the domain in Gnuplot format. |
251 |
|
|
* \param[in] filename the filename to write to. */ |
252 |
|
|
inline void draw_domain_gnuplot(const char* filename) { |
253 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
254 |
|
|
draw_domain_gnuplot(fp); |
255 |
|
|
fclose(fp); |
256 |
|
|
} |
257 |
|
|
void draw_domain_pov(FILE *fp=stdout); |
258 |
|
|
/** Draws an outline of the domain in Gnuplot format. |
259 |
|
|
* \param[in] filename the filename to write to. */ |
260 |
|
|
inline void draw_domain_pov(const char* filename) { |
261 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
262 |
|
|
draw_domain_pov(fp); |
263 |
|
|
fclose(fp); |
264 |
|
|
} |
265 |
|
|
/** Sums up the total number of stored particles. |
266 |
|
|
* \return The number of particles. */ |
267 |
|
|
inline int total_particles() { |
268 |
|
|
int tp=*co; |
269 |
|
|
for(int *cop=co+1;cop<co+nxyz;cop++) tp+=*cop; |
270 |
|
|
return tp; |
271 |
|
|
} |
272 |
|
|
protected: |
273 |
|
|
void add_particle_memory(int i); |
274 |
|
|
bool put_locate_block(int &ijk,double &x,double &y,double &z); |
275 |
|
|
inline bool put_remap(int &ijk,double &x,double &y,double &z); |
276 |
|
|
inline bool remap(int &ai,int &aj,int &ak,int &ci,int &cj,int &ck,double &x,double &y,double &z,int &ijk); |
277 |
|
|
}; |
278 |
|
|
|
279 |
|
|
/** \brief Extension of the container_base class for computing regular Voronoi |
280 |
|
|
* tessellations. |
281 |
|
|
* |
282 |
|
|
* This class is an extension of the container_base class that has routines |
283 |
|
|
* specifically for computing the regular Voronoi tessellation with no |
284 |
|
|
* dependence on particle radii. */ |
285 |
|
|
class container : public container_base, public radius_mono { |
286 |
|
|
public: |
287 |
|
|
container(double ax_,double bx_,double ay_,double by_,double az_,double bz_, |
288 |
|
|
int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem); |
289 |
|
|
void clear(); |
290 |
|
|
void put(int n,double x,double y,double z); |
291 |
|
|
void put(particle_order &vo,int n,double x,double y,double z); |
292 |
|
|
void import(FILE *fp=stdin); |
293 |
|
|
void import(particle_order &vo,FILE *fp=stdin); |
294 |
|
|
/** Imports a list of particles from an open file stream into |
295 |
|
|
* the container. Entries of four numbers (Particle ID, x |
296 |
|
|
* position, y position, z position) are searched for. If the |
297 |
|
|
* file cannot be successfully read, then the routine causes a |
298 |
|
|
* fatal error. |
299 |
|
|
* \param[in] filename the name of the file to open and read |
300 |
|
|
* from. */ |
301 |
|
|
inline void import(const char* filename) { |
302 |
|
|
FILE *fp=safe_fopen(filename,"r"); |
303 |
|
|
import(fp); |
304 |
|
|
fclose(fp); |
305 |
|
|
} |
306 |
|
|
/** Imports a list of particles from an open file stream into |
307 |
|
|
* the container. Entries of four numbers (Particle ID, x |
308 |
|
|
* position, y position, z position) are searched for. In |
309 |
|
|
* addition, the order in which particles are read is saved |
310 |
|
|
* into an ordering class. If the file cannot be successfully |
311 |
|
|
* read, then the routine causes a fatal error. |
312 |
|
|
* \param[in,out] vo the ordering class to use. |
313 |
|
|
* \param[in] filename the name of the file to open and read |
314 |
|
|
* from. */ |
315 |
|
|
inline void import(particle_order &vo,const char* filename) { |
316 |
|
|
FILE *fp=safe_fopen(filename,"r"); |
317 |
|
|
import(vo,fp); |
318 |
|
|
fclose(fp); |
319 |
|
|
} |
320 |
|
|
void compute_all_cells(); |
321 |
|
|
double sum_cell_volumes(); |
322 |
|
|
/** Dumps particle IDs and positions to a file. |
323 |
|
|
* \param[in] vl the loop class to use. |
324 |
|
|
* \param[in] fp a file handle to write to. */ |
325 |
|
|
template<class c_loop> |
326 |
|
|
void draw_particles(c_loop &vl,FILE *fp) { |
327 |
|
|
double *pp; |
328 |
|
|
if(vl.start()) do { |
329 |
|
|
pp=p[vl.ijk]+3*vl.q; |
330 |
|
|
fprintf(fp,"%d %g %g %g\n",id[vl.ijk][vl.q],*pp,pp[1],pp[2]); |
331 |
|
|
} while(vl.inc()); |
332 |
|
|
} |
333 |
|
|
/** Dumps all of the particle IDs and positions to a file. |
334 |
|
|
* \param[in] fp a file handle to write to. */ |
335 |
|
|
inline void draw_particles(FILE *fp=stdout) { |
336 |
|
|
c_loop_all vl(*this); |
337 |
|
|
draw_particles(vl,fp); |
338 |
|
|
} |
339 |
|
|
/** Dumps all of the particle IDs and positions to a file. |
340 |
|
|
* \param[in] filename the name of the file to write to. */ |
341 |
|
|
inline void draw_particles(const char *filename) { |
342 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
343 |
|
|
draw_particles(fp); |
344 |
|
|
fclose(fp); |
345 |
|
|
} |
346 |
|
|
/** Dumps particle positions in POV-Ray format. |
347 |
|
|
* \param[in] vl the loop class to use. |
348 |
|
|
* \param[in] fp a file handle to write to. */ |
349 |
|
|
template<class c_loop> |
350 |
|
|
void draw_particles_pov(c_loop &vl,FILE *fp) { |
351 |
|
|
double *pp; |
352 |
|
|
if(vl.start()) do { |
353 |
|
|
pp=p[vl.ijk]+3*vl.q; |
354 |
|
|
fprintf(fp,"// id %d\nsphere{<%g,%g,%g>,s}\n", |
355 |
|
|
id[vl.ijk][vl.q],*pp,pp[1],pp[2]); |
356 |
|
|
} while(vl.inc()); |
357 |
|
|
} |
358 |
|
|
/** Dumps all particle positions in POV-Ray format. |
359 |
|
|
* \param[in] fp a file handle to write to. */ |
360 |
|
|
inline void draw_particles_pov(FILE *fp=stdout) { |
361 |
|
|
c_loop_all vl(*this); |
362 |
|
|
draw_particles_pov(vl,fp); |
363 |
|
|
} |
364 |
|
|
/** Dumps all particle positions in POV-Ray format. |
365 |
|
|
* \param[in] filename the name of the file to write to. */ |
366 |
|
|
inline void draw_particles_pov(const char *filename) { |
367 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
368 |
|
|
draw_particles_pov(fp); |
369 |
|
|
fclose(fp); |
370 |
|
|
} |
371 |
|
|
/** Computes Voronoi cells and saves the output in gnuplot |
372 |
|
|
* format. |
373 |
|
|
* \param[in] vl the loop class to use. |
374 |
|
|
* \param[in] fp a file handle to write to. */ |
375 |
|
|
template<class c_loop> |
376 |
|
|
void draw_cells_gnuplot(c_loop &vl,FILE *fp) { |
377 |
|
|
voronoicell c;double *pp; |
378 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
379 |
|
|
pp=p[vl.ijk]+ps*vl.q; |
380 |
|
|
c.draw_gnuplot(*pp,pp[1],pp[2],fp); |
381 |
|
|
} while(vl.inc()); |
382 |
|
|
} |
383 |
|
|
/** Computes all Voronoi cells and saves the output in gnuplot |
384 |
|
|
* format. |
385 |
|
|
* \param[in] fp a file handle to write to. */ |
386 |
|
|
inline void draw_cells_gnuplot(FILE *fp=stdout) { |
387 |
|
|
c_loop_all vl(*this); |
388 |
|
|
draw_cells_gnuplot(vl,fp); |
389 |
|
|
} |
390 |
|
|
/** Computes all Voronoi cells and saves the output in gnuplot |
391 |
|
|
* format. |
392 |
|
|
* \param[in] filename the name of the file to write to. */ |
393 |
|
|
inline void draw_cells_gnuplot(const char *filename) { |
394 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
395 |
|
|
draw_cells_gnuplot(fp); |
396 |
|
|
fclose(fp); |
397 |
|
|
} |
398 |
|
|
/** Computes Voronoi cells and saves the output in POV-Ray |
399 |
|
|
* format. |
400 |
|
|
* \param[in] vl the loop class to use. |
401 |
|
|
* \param[in] fp a file handle to write to. */ |
402 |
|
|
template<class c_loop> |
403 |
|
|
void draw_cells_pov(c_loop &vl,FILE *fp) { |
404 |
|
|
voronoicell c;double *pp; |
405 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
406 |
|
|
fprintf(fp,"// cell %d\n",id[vl.ijk][vl.q]); |
407 |
|
|
pp=p[vl.ijk]+ps*vl.q; |
408 |
|
|
c.draw_pov(*pp,pp[1],pp[2],fp); |
409 |
|
|
} while(vl.inc()); |
410 |
|
|
} |
411 |
|
|
/** Computes all Voronoi cells and saves the output in POV-Ray |
412 |
|
|
* format. |
413 |
|
|
* \param[in] fp a file handle to write to. */ |
414 |
|
|
inline void draw_cells_pov(FILE *fp=stdout) { |
415 |
|
|
c_loop_all vl(*this); |
416 |
|
|
draw_cells_pov(vl,fp); |
417 |
|
|
} |
418 |
|
|
/** Computes all Voronoi cells and saves the output in POV-Ray |
419 |
|
|
* format. |
420 |
|
|
* \param[in] filename the name of the file to write to. */ |
421 |
|
|
inline void draw_cells_pov(const char *filename) { |
422 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
423 |
|
|
draw_cells_pov(fp); |
424 |
|
|
fclose(fp); |
425 |
|
|
} |
426 |
|
|
/** Computes the Voronoi cells and saves customized information |
427 |
|
|
* about them. |
428 |
|
|
* \param[in] vl the loop class to use. |
429 |
|
|
* \param[in] format the custom output string to use. |
430 |
|
|
* \param[in] fp a file handle to write to. */ |
431 |
|
|
template<class c_loop> |
432 |
|
|
void print_custom(c_loop &vl,const char *format,FILE *fp) { |
433 |
|
|
int ijk,q;double *pp; |
434 |
|
|
if(contains_neighbor(format)) { |
435 |
|
|
voronoicell_neighbor c; |
436 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
437 |
|
|
ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q; |
438 |
|
|
c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],default_radius,fp); |
439 |
|
|
} while(vl.inc()); |
440 |
|
|
} else { |
441 |
|
|
voronoicell c; |
442 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
443 |
|
|
ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q; |
444 |
|
|
c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],default_radius,fp); |
445 |
|
|
} while(vl.inc()); |
446 |
|
|
} |
447 |
|
|
} |
448 |
|
|
void print_custom(const char *format,FILE *fp=stdout); |
449 |
|
|
void print_custom(const char *format,const char *filename); |
450 |
|
|
bool find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid); |
451 |
|
|
/** Computes the Voronoi cell for a particle currently being |
452 |
|
|
* referenced by a loop class. |
453 |
|
|
* \param[out] c a Voronoi cell class in which to store the |
454 |
|
|
* computed cell. |
455 |
|
|
* \param[in] vl the loop class to use. |
456 |
|
|
* \return True if the cell was computed. If the cell cannot be |
457 |
|
|
* computed, if it is removed entirely by a wall or boundary |
458 |
|
|
* condition, then the routine returns false. */ |
459 |
|
|
template<class v_cell,class c_loop> |
460 |
|
|
inline bool compute_cell(v_cell &c,c_loop &vl) { |
461 |
|
|
return vc.compute_cell(c,vl.ijk,vl.q,vl.i,vl.j,vl.k); |
462 |
|
|
} |
463 |
|
|
/** Computes the Voronoi cell for given particle. |
464 |
|
|
* \param[out] c a Voronoi cell class in which to store the |
465 |
|
|
* computed cell. |
466 |
|
|
* \param[in] ijk the block that the particle is within. |
467 |
|
|
* \param[in] q the index of the particle within the block. |
468 |
|
|
* \return True if the cell was computed. If the cell cannot be |
469 |
|
|
* computed, if it is removed entirely by a wall or boundary |
470 |
|
|
* condition, then the routine returns false. */ |
471 |
|
|
template<class v_cell> |
472 |
|
|
inline bool compute_cell(v_cell &c,int ijk,int q) { |
473 |
|
|
int k=ijk/nxy,ijkt=ijk-nxy*k,j=ijkt/nx,i=ijkt-j*nx; |
474 |
|
|
return vc.compute_cell(c,ijk,q,i,j,k); |
475 |
|
|
} |
476 |
|
|
/** Computes the Voronoi cell for a ghost particle at a given |
477 |
|
|
* location. |
478 |
|
|
* \param[out] c a Voronoi cell class in which to store the |
479 |
|
|
* computed cell. |
480 |
|
|
* \param[in] (x,y,z) the location of the ghost particle. |
481 |
|
|
* \return True if the cell was computed. If the cell cannot be |
482 |
|
|
* computed, if it is removed entirely by a wall or boundary |
483 |
|
|
* condition, then the routine returns false. */ |
484 |
|
|
template<class v_cell> |
485 |
|
|
inline bool compute_ghost_cell(v_cell &c,double x,double y,double z) { |
486 |
|
|
int ijk; |
487 |
|
|
if(put_locate_block(ijk,x,y,z)) { |
488 |
|
|
double *pp=p[ijk]+3*co[ijk]++; |
489 |
|
|
*(pp++)=x;*(pp++)=y;*pp=z; |
490 |
|
|
bool q=compute_cell(c,ijk,co[ijk]-1); |
491 |
|
|
co[ijk]--; |
492 |
|
|
return q; |
493 |
|
|
} |
494 |
|
|
return false; |
495 |
|
|
} |
496 |
|
|
private: |
497 |
|
|
voro_compute<container> vc; |
498 |
|
|
friend class voro_compute<container>; |
499 |
|
|
}; |
500 |
|
|
|
501 |
|
|
/** \brief Extension of the container_base class for computing radical Voronoi |
502 |
|
|
* tessellations. |
503 |
|
|
* |
504 |
|
|
* This class is an extension of container_base class that has routines |
505 |
|
|
* specifically for computing the radical Voronoi tessellation that depends on |
506 |
|
|
* the particle radii. */ |
507 |
|
|
class container_poly : public container_base, public radius_poly { |
508 |
|
|
public: |
509 |
|
|
container_poly(double ax_,double bx_,double ay_,double by_,double az_,double bz_, |
510 |
|
|
int nx_,int ny_,int nz_,bool xperiodic_,bool yperiodic_,bool zperiodic_,int init_mem); |
511 |
|
|
void clear(); |
512 |
|
|
void put(int n,double x,double y,double z,double r); |
513 |
|
|
void put(particle_order &vo,int n,double x,double y,double z,double r); |
514 |
|
|
void import(FILE *fp=stdin); |
515 |
|
|
void import(particle_order &vo,FILE *fp=stdin); |
516 |
|
|
/** Imports a list of particles from an open file stream into |
517 |
|
|
* the container_poly class. Entries of five numbers (Particle |
518 |
|
|
* ID, x position, y position, z position, radius) are searched |
519 |
|
|
* for. If the file cannot be successfully read, then the |
520 |
|
|
* routine causes a fatal error. |
521 |
|
|
* \param[in] filename the name of the file to open and read |
522 |
|
|
* from. */ |
523 |
|
|
inline void import(const char* filename) { |
524 |
|
|
FILE *fp=safe_fopen(filename,"r"); |
525 |
|
|
import(fp); |
526 |
|
|
fclose(fp); |
527 |
|
|
} |
528 |
|
|
/** Imports a list of particles from an open file stream into |
529 |
|
|
* the container_poly class. Entries of five numbers (Particle |
530 |
|
|
* ID, x position, y position, z position, radius) are searched |
531 |
|
|
* for. In addition, the order in which particles are read is |
532 |
|
|
* saved into an ordering class. If the file cannot be |
533 |
|
|
* successfully read, then the routine causes a fatal error. |
534 |
|
|
* \param[in,out] vo the ordering class to use. |
535 |
|
|
* \param[in] filename the name of the file to open and read |
536 |
|
|
* from. */ |
537 |
|
|
inline void import(particle_order &vo,const char* filename) { |
538 |
|
|
FILE *fp=safe_fopen(filename,"r"); |
539 |
|
|
import(vo,fp); |
540 |
|
|
fclose(fp); |
541 |
|
|
} |
542 |
|
|
void compute_all_cells(); |
543 |
|
|
double sum_cell_volumes(); |
544 |
|
|
/** Dumps particle IDs, positions and radii to a file. |
545 |
|
|
* \param[in] vl the loop class to use. |
546 |
|
|
* \param[in] fp a file handle to write to. */ |
547 |
|
|
template<class c_loop> |
548 |
|
|
void draw_particles(c_loop &vl,FILE *fp) { |
549 |
|
|
double *pp; |
550 |
|
|
if(vl.start()) do { |
551 |
|
|
pp=p[vl.ijk]+4*vl.q; |
552 |
|
|
fprintf(fp,"%d %g %g %g %g\n",id[vl.ijk][vl.q],*pp,pp[1],pp[2],pp[3]); |
553 |
|
|
} while(vl.inc()); |
554 |
|
|
} |
555 |
|
|
/** Dumps all of the particle IDs, positions and radii to a |
556 |
|
|
* file. |
557 |
|
|
* \param[in] fp a file handle to write to. */ |
558 |
|
|
inline void draw_particles(FILE *fp=stdout) { |
559 |
|
|
c_loop_all vl(*this); |
560 |
|
|
draw_particles(vl,fp); |
561 |
|
|
} |
562 |
|
|
/** Dumps all of the particle IDs, positions and radii to a |
563 |
|
|
* file. |
564 |
|
|
* \param[in] filename the name of the file to write to. */ |
565 |
|
|
inline void draw_particles(const char *filename) { |
566 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
567 |
|
|
draw_particles(fp); |
568 |
|
|
fclose(fp); |
569 |
|
|
} |
570 |
|
|
/** Dumps particle positions in POV-Ray format. |
571 |
|
|
* \param[in] vl the loop class to use. |
572 |
|
|
* \param[in] fp a file handle to write to. */ |
573 |
|
|
template<class c_loop> |
574 |
|
|
void draw_particles_pov(c_loop &vl,FILE *fp) { |
575 |
|
|
double *pp; |
576 |
|
|
if(vl.start()) do { |
577 |
|
|
pp=p[vl.ijk]+4*vl.q; |
578 |
|
|
fprintf(fp,"// id %d\nsphere{<%g,%g,%g>,%g}\n", |
579 |
|
|
id[vl.ijk][vl.q],*pp,pp[1],pp[2],pp[3]); |
580 |
|
|
} while(vl.inc()); |
581 |
|
|
} |
582 |
|
|
/** Dumps all the particle positions in POV-Ray format. |
583 |
|
|
* \param[in] fp a file handle to write to. */ |
584 |
|
|
inline void draw_particles_pov(FILE *fp=stdout) { |
585 |
|
|
c_loop_all vl(*this); |
586 |
|
|
draw_particles_pov(vl,fp); |
587 |
|
|
} |
588 |
|
|
/** Dumps all the particle positions in POV-Ray format. |
589 |
|
|
* \param[in] filename the name of the file to write to. */ |
590 |
|
|
inline void draw_particles_pov(const char *filename) { |
591 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
592 |
|
|
draw_particles_pov(fp); |
593 |
|
|
fclose(fp); |
594 |
|
|
} |
595 |
|
|
/** Computes Voronoi cells and saves the output in gnuplot |
596 |
|
|
* format. |
597 |
|
|
* \param[in] vl the loop class to use. |
598 |
|
|
* \param[in] fp a file handle to write to. */ |
599 |
|
|
template<class c_loop> |
600 |
|
|
void draw_cells_gnuplot(c_loop &vl,FILE *fp) { |
601 |
|
|
voronoicell c;double *pp; |
602 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
603 |
|
|
pp=p[vl.ijk]+ps*vl.q; |
604 |
|
|
c.draw_gnuplot(*pp,pp[1],pp[2],fp); |
605 |
|
|
} while(vl.inc()); |
606 |
|
|
} |
607 |
|
|
/** Compute all Voronoi cells and saves the output in gnuplot |
608 |
|
|
* format. |
609 |
|
|
* \param[in] fp a file handle to write to. */ |
610 |
|
|
inline void draw_cells_gnuplot(FILE *fp=stdout) { |
611 |
|
|
c_loop_all vl(*this); |
612 |
|
|
draw_cells_gnuplot(vl,fp); |
613 |
|
|
} |
614 |
|
|
/** Compute all Voronoi cells and saves the output in gnuplot |
615 |
|
|
* format. |
616 |
|
|
* \param[in] filename the name of the file to write to. */ |
617 |
|
|
inline void draw_cells_gnuplot(const char *filename) { |
618 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
619 |
|
|
draw_cells_gnuplot(fp); |
620 |
|
|
fclose(fp); |
621 |
|
|
} |
622 |
|
|
/** Computes Voronoi cells and saves the output in POV-Ray |
623 |
|
|
* format. |
624 |
|
|
* \param[in] vl the loop class to use. |
625 |
|
|
* \param[in] fp a file handle to write to. */ |
626 |
|
|
template<class c_loop> |
627 |
|
|
void draw_cells_pov(c_loop &vl,FILE *fp) { |
628 |
|
|
voronoicell c;double *pp; |
629 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
630 |
|
|
fprintf(fp,"// cell %d\n",id[vl.ijk][vl.q]); |
631 |
|
|
pp=p[vl.ijk]+ps*vl.q; |
632 |
|
|
c.draw_pov(*pp,pp[1],pp[2],fp); |
633 |
|
|
} while(vl.inc()); |
634 |
|
|
} |
635 |
|
|
/** Computes all Voronoi cells and saves the output in POV-Ray |
636 |
|
|
* format. |
637 |
|
|
* \param[in] fp a file handle to write to. */ |
638 |
|
|
inline void draw_cells_pov(FILE *fp=stdout) { |
639 |
|
|
c_loop_all vl(*this); |
640 |
|
|
draw_cells_pov(vl,fp); |
641 |
|
|
} |
642 |
|
|
/** Computes all Voronoi cells and saves the output in POV-Ray |
643 |
|
|
* format. |
644 |
|
|
* \param[in] filename the name of the file to write to. */ |
645 |
|
|
inline void draw_cells_pov(const char *filename) { |
646 |
|
|
FILE *fp=safe_fopen(filename,"w"); |
647 |
|
|
draw_cells_pov(fp); |
648 |
|
|
fclose(fp); |
649 |
|
|
} |
650 |
|
|
/** Computes the Voronoi cells and saves customized information |
651 |
|
|
* about them. |
652 |
|
|
* \param[in] vl the loop class to use. |
653 |
|
|
* \param[in] format the custom output string to use. |
654 |
|
|
* \param[in] fp a file handle to write to. */ |
655 |
|
|
template<class c_loop> |
656 |
|
|
void print_custom(c_loop &vl,const char *format,FILE *fp) { |
657 |
|
|
int ijk,q;double *pp; |
658 |
|
|
if(contains_neighbor(format)) { |
659 |
|
|
voronoicell_neighbor c; |
660 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
661 |
|
|
ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q; |
662 |
|
|
c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],pp[3],fp); |
663 |
|
|
} while(vl.inc()); |
664 |
|
|
} else { |
665 |
|
|
voronoicell c; |
666 |
|
|
if(vl.start()) do if(compute_cell(c,vl)) { |
667 |
|
|
ijk=vl.ijk;q=vl.q;pp=p[ijk]+ps*q; |
668 |
|
|
c.output_custom(format,id[ijk][q],*pp,pp[1],pp[2],pp[3],fp); |
669 |
|
|
} while(vl.inc()); |
670 |
|
|
} |
671 |
|
|
} |
672 |
|
|
/** Computes the Voronoi cell for a particle currently being |
673 |
|
|
* referenced by a loop class. |
674 |
|
|
* \param[out] c a Voronoi cell class in which to store the |
675 |
|
|
* computed cell. |
676 |
|
|
* \param[in] vl the loop class to use. |
677 |
|
|
* \return True if the cell was computed. If the cell cannot be |
678 |
|
|
* computed, if it is removed entirely by a wall or boundary |
679 |
|
|
* condition, then the routine returns false. */ |
680 |
|
|
template<class v_cell,class c_loop> |
681 |
|
|
inline bool compute_cell(v_cell &c,c_loop &vl) { |
682 |
|
|
return vc.compute_cell(c,vl.ijk,vl.q,vl.i,vl.j,vl.k); |
683 |
|
|
} |
684 |
|
|
/** Computes the Voronoi cell for given particle. |
685 |
|
|
* \param[out] c a Voronoi cell class in which to store the |
686 |
|
|
* computed cell. |
687 |
|
|
* \param[in] ijk the block that the particle is within. |
688 |
|
|
* \param[in] q the index of the particle within the block. |
689 |
|
|
* \return True if the cell was computed. If the cell cannot be |
690 |
|
|
* computed, if it is removed entirely by a wall or boundary |
691 |
|
|
* condition, then the routine returns false. */ |
692 |
|
|
template<class v_cell> |
693 |
|
|
inline bool compute_cell(v_cell &c,int ijk,int q) { |
694 |
|
|
int k=ijk/nxy,ijkt=ijk-nxy*k,j=ijkt/nx,i=ijkt-j*nx; |
695 |
|
|
return vc.compute_cell(c,ijk,q,i,j,k); |
696 |
|
|
} |
697 |
|
|
/** Computes the Voronoi cell for a ghost particle at a given |
698 |
|
|
* location. |
699 |
|
|
* \param[out] c a Voronoi cell class in which to store the |
700 |
|
|
* computed cell. |
701 |
|
|
* \param[in] (x,y,z) the location of the ghost particle. |
702 |
|
|
* \param[in] r the radius of the ghost particle. |
703 |
|
|
* \return True if the cell was computed. If the cell cannot be |
704 |
|
|
* computed, if it is removed entirely by a wall or boundary |
705 |
|
|
* condition, then the routine returns false. */ |
706 |
|
|
template<class v_cell> |
707 |
|
|
inline bool compute_ghost_cell(v_cell &c,double x,double y,double z,double r) { |
708 |
|
|
int ijk; |
709 |
|
|
if(put_locate_block(ijk,x,y,z)) { |
710 |
|
|
double *pp=p[ijk]+4*co[ijk]++,tm=max_radius; |
711 |
|
|
*(pp++)=x;*(pp++)=y;*(pp++)=z;*pp=r; |
712 |
|
|
if(r>max_radius) max_radius=r; |
713 |
|
|
bool q=compute_cell(c,ijk,co[ijk]-1); |
714 |
|
|
co[ijk]--;max_radius=tm; |
715 |
|
|
return q; |
716 |
|
|
} |
717 |
|
|
return false; |
718 |
|
|
} |
719 |
|
|
void print_custom(const char *format,FILE *fp=stdout); |
720 |
|
|
void print_custom(const char *format,const char *filename); |
721 |
|
|
bool find_voronoi_cell(double x,double y,double z,double &rx,double &ry,double &rz,int &pid); |
722 |
|
|
private: |
723 |
|
|
voro_compute<container_poly> vc; |
724 |
|
|
friend class voro_compute<container_poly>; |
725 |
|
|
}; |
726 |
|
|
|
727 |
|
|
} |
728 |
|
|
|
729 |
|
|
#endif |