Strip extra spaces from code.
[voro++.git] / branches / 2d / src / ctr_boundary_2d.hh
blob2299f4383619d52a404b5d9b4b57412021c1782d
1 // Voro++, a cell-based Voronoi library
2 //
3 // Authors : Chris H. Rycroft (LBL / UC Berkeley)
4 // Cody Robert Dance (UC Berkeley)
5 // Email : chr@alum.mit.edu
6 // Date : August 30th 2011
8 /** \file ctr_boundary_2d.hh
9 * \brief Header file for the container_boundary_2d and related classes. */
11 #ifndef VOROPP_CTR_BOUNDARY_2D_HH
12 #define VOROPP_CTR_BOUNDARY_2D_HH
14 #include <cstdio>
15 #include <cstdlib>
16 #include <cmath>
17 #include <vector>
18 using namespace std;
20 #include "config.hh"
21 #include "common.hh"
22 #include "v_base_2d.hh"
23 #include "cell_2d.hh"
24 #include "c_loops_2d.hh"
25 #include "rad_option.hh"
26 #include "v_compute_2d.hh"
28 namespace voro {
30 /** \brief Class for representing a particle system in a three-dimensional
31 * rectangular box.
33 * This class represents a system of particles in a three-dimensional
34 * rectangular box. Any combination of non-periodic and periodic coordinates
35 * can be used in the three coordinate directions. The class is not intended
36 * for direct use, but instead forms the base of the container and
37 * container_poly classes that add specialized routines for computing the
38 * regular and radical Voronoi tessellations respectively. It contains routines
39 * that are commonly between these two classes, such as those for drawing the
40 * domain, and placing particles within the internal data structure.
42 * The class is derived from the wall_list class, which encapsulates routines
43 * for associating walls with the container, and the voro_base class, which
44 * encapsulates routines about the underlying computational grid. */
45 class container_boundary_2d : public voro_base_2d, public radius_mono {
46 public:
48 /** The minimum x coordinate of the container. */
49 const double ax;
50 /** The maximum x coordinate of the container. */
51 const double bx;
52 /** The minimum y coordinate of the container. */
53 const double ay;
54 /** The maximum y coordinate of the container. */
55 const double by;
56 /** A boolean value that determines if the x coordinate in
57 * periodic or not. */
58 const bool xperiodic;
59 /** A boolean value that determines if the y coordinate in
60 * periodic or not. */
61 const bool yperiodic;
62 /** This array holds the numerical IDs of each particle in each
63 * computational box. */
64 int **id;
65 /** A two dimensional array holding particle positions. For the
66 * derived container_poly class, this also holds particle
67 * radii. */
68 double **p;
69 /** This array holds the number of particles within each
70 * computational box of the container. */
71 int *co;
72 /** This array holds the maximum amount of particle memory for
73 * each computational box of the container. If the number of
74 * particles in a particular box ever approaches this limit,
75 * more is allocated using the add_particle_memory() function.
77 int *mem;
78 int **wid;
79 int **nlab;
80 int ***plab;
81 int **bndpts;
82 int boundary_track;
83 int edbc;
84 int edbm;
85 int *edb;
86 double *bnds;
88 /** The amount of memory in the array structure for each
89 * particle. This is set to 2 when the basic class is
90 * initialized, so that the array holds (x,y) positions. If the
91 * 2D container class is initialized as part of the derived class
92 * container_poly_2d, then this is set to 3, to also hold the
93 * particle radii. */
94 const int ps;
95 container_boundary_2d(double ax_,double bx_,double ay_,double by_,
96 int nx_,int ny_,bool xperiodic_,bool yperiodic_,int init_mem);
97 ~container_boundary_2d();
98 void region_count();
99 /** Initializes the Voronoi cell prior to a compute_cell
100 * operation for a specific particle being carried out by a
101 * voro_compute class. The cell is initialized to fill the
102 * entire container. For non-periodic coordinates, this is set
103 * by the position of the walls. For periodic coordinates, the
104 * space is equally divided in either direction from the
105 * particle's initial position. Plane cuts made by any walls
106 * that have been added are then applied to the cell.
107 * \param[in,out] c a reference to a voronoicell_nonconvex_2d object.
108 * \param[in] ij the block that the particle is within.
109 * \param[in] q the index of the particle within its block.
110 * \param[in] (ci,cj) the coordinates of the block in the
111 * container coordinate system.
112 * \param[out] (i,j) the coordinates of the test block relative
113 * to the voro_compute coordinate system.
114 * \param[out] (x,y) the position of the particle.
115 * \param[out] disp a block displacement used internally by the
116 * compute_cell routine.
117 * \return False if the plane cuts applied by walls completely
118 * removed the cell, true otherwise. */
119 template<class v_cell_2d>
120 inline bool initialize_voronoicell(v_cell_2d &c,int ij,int q,int ci,int cj,
121 int &i,int &j,double &x,double &y,int &disp) {
122 double x1,x2,y1,y2,*pp=p[ij]+ps*q;
123 x=*(pp++);y=*(pp++);
124 if(xperiodic) {x1=-(x2=0.5*(bx-ax));i=nx;} else {x1=ax-x;x2=bx-x;i=ci;}
125 if(yperiodic) {y1=-(y2=0.5*(by-ay));j=ny;} else {y1=ay-y;y2=by-y;j=cj;}
126 if(bndpts[ij][q]==-1) c.init(x1,x2,y1,y2);
127 else {
128 int &bid=bndpts[ij][q];
129 double cx=bnds[2*bid],cy=bnds[2*bid+1];
130 int nwid=edb[2*bid],lwid=edb[2*bid+1];
131 double lx=bnds[2*lwid],ly=bnds[2*lwid+1];
132 double nx=bnds[2*nwid],ny=bnds[2*nwid+1];
133 c.init_nonconvex(x1,x2,y1,y2,nx-cx,ny-cy,lx-cx,ly-cy);
135 disp=ij-i-nx*j;
136 return true;
138 bool point_inside(double x,double y);
139 template<class v_cell_2d>
140 bool boundary_cuts(v_cell_2d &c,int ij,double x,double y);
141 /** Initializes parameters for a find_voronoi_cell call within
142 * the voro_compute template.
143 * \param[in] (ci,cj) the coordinates of the test block in
144 * the container coordinate system.
145 * \param[in] ij the index of the test block
146 * \param[out] (i,j) the coordinates of the test block relative
147 * to the voro_compute coordinate system.
148 * \param[out] disp a block displacement used internally by the
149 * find_voronoi_cell routine. */
150 inline void initialize_search(int ci,int cj,int ij,int &i,int &j,int &disp) {
151 i=xperiodic?nx:ci;
152 j=yperiodic?ny:cj;
153 disp=ij-i-nx*j;
155 /** Returns the position of a particle currently being computed
156 * relative to the computational block that it is within. It is
157 * used to select the optimal worklist entry to use.
158 * \param[in] (x,y) the position of the particle.
159 * \param[in] (ci,cj) the block that the particle is within.
160 * \param[out] (fx,fy) the position relative to the block.
162 inline void frac_pos(double x,double y,double ci,double cj,
163 double &fx,double &fy) {
164 fx=x-ax-boxx*ci;
165 fy=y-ay-boxy*cj;
167 /** Calculates the index of block in the container structure
168 * corresponding to given coordinates.
169 * \param[in] (ci,cj) the coordinates of the original block in
170 * the current computation, relative to the
171 * container coordinate system.
172 * \param[in] (ei,ej) the displacement of the current block
173 * from the original block.
174 * \param[in,out] (qx,qy) the periodic displacement that must
175 * be added to the particles within the
176 * computed block.
177 * \param[in] disp a block displacement used internally by the
178 * find_voronoi_cell and compute_cell routines.
179 * \return The block index. */
180 inline int region_index(int ci,int cj,int ei,int ej,double &qx,double &qy,int &disp) {
181 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;}
182 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;}
183 return disp+ei+nx*ej;
185 void draw_domain_gnuplot(FILE *fp=stdout);
186 /** Draws an outline of the domain in Gnuplot format.
187 * \param[in] filename the filename to write to. */
188 inline void draw_domain_gnuplot(const char* filename) {
189 FILE *fp=safe_fopen(filename,"w");
190 draw_domain_gnuplot(fp);
191 fclose(fp);
193 void draw_domain_pov(FILE *fp=stdout);
194 /** Draws an outline of the domain in Gnuplot format.
195 * \param[in] filename the filename to write to. */
196 inline void draw_domain_pov(const char* filename) {
197 FILE *fp=safe_fopen(filename,"w");
198 draw_domain_pov(fp);
199 fclose(fp);
201 void draw_boundary_gnuplot(FILE *fp=stdout);
202 inline void draw_boundary_gnuplot(const char* filename) {
203 FILE *fp=safe_fopen(filename,"w");
204 draw_boundary_gnuplot(fp);
205 fclose(fp);
207 /** Sums up the total number of stored particles.
208 * \return The number of particles. */
209 inline int total_particles() {
210 int tp=*co;
211 for(int *cop=co+1;cop<co+nxy;cop++) tp+=*cop;
212 return tp;
214 inline void start_boundary() {boundary_track=edbc;}
215 void end_boundary();
216 void register_boundary(double x,double y);
217 void clear();
218 void put(int n,double x,double y);
219 void put(particle_order &vo,int n,double x,double y);
227 void import(FILE *fp=stdin);
228 /** Imports a list of particles from an open file stream into
229 * the container. Entries of three numbers (Particle ID, x
230 * position, y position) are searched for. If the file cannot
231 * be successfully read, then the routine causes a fatal error.
232 * \param[in] filename the name of the file to open and read
233 * from. */
234 inline void import(const char* filename) {
235 FILE *fp=safe_fopen(filename,"r");
236 import(fp);
237 fclose(fp);
239 void compute_all_cells();
240 double sum_cell_areas();
241 /** Dumps particle IDs and positions to a file.
242 * \param[in] vl the loop class to use.
243 * \param[in] fp a file handle to write to. */
244 template<class c_loop_2d>
245 void draw_particles(c_loop_2d &vl,FILE *fp) {
246 double *pp;
247 if(vl.start()) do {
248 pp=p[vl.ij]+2*vl.q;
249 fprintf(fp,"%d %g %g\n",id[vl.ij][vl.q],*pp,pp[1]);
250 } while(vl.inc());
252 /** Dumps all of the particle IDs and positions to a file.
253 * \param[in] fp a file handle to write to. */
254 inline void draw_particles(FILE *fp=stdout) {
255 c_loop_all_2d vl(*this);
256 draw_particles(vl,fp);
258 /** Dumps all of the particle IDs and positions to a file.
259 * \param[in] filename the name of the file to write to. */
260 inline void draw_particles(const char *filename) {
261 FILE *fp=safe_fopen(filename,"w");
262 draw_particles(fp);
263 fclose(fp);
265 /** Dumps particle positions in POV-Ray format.
266 * \param[in] vl the loop class to use.
267 * \param[in] fp a file handle to write to. */
268 template<class c_loop_2d>
269 void draw_particles_pov(c_loop_2d &vl,FILE *fp) {
270 double *pp;
271 if(vl.start()) do {
272 pp=p[vl.ij]+2*vl.q;
273 fprintf(fp,"// id %d\nsphere{<%g,%g,0>,s}\n",
274 id[vl.ij][vl.q],*pp,pp[1]);
275 } while(vl.inc());
277 /** Dumps all particle positions in POV-Ray format.
278 * \param[in] fp a file handle to write to. */
279 inline void draw_particles_pov(FILE *fp=stdout) {
280 c_loop_all_2d vl(*this);
281 draw_particles_pov(vl,fp);
283 /** Dumps all particle positions in POV-Ray format.
284 * \param[in] filename the name of the file to write to. */
285 inline void draw_particles_pov(const char *filename) {
286 FILE *fp=safe_fopen(filename,"w");
287 draw_particles_pov(fp);
288 fclose(fp);
290 /** Computes Voronoi cells and saves the output in gnuplot
291 * format.
292 * \param[in] vl the loop class to use.
293 * \param[in] fp a file handle to write to. */
294 template<class c_loop_2d>
295 void draw_cells_gnuplot(c_loop_2d &vl,FILE *fp) {
296 voronoicell_nonconvex_2d c;double *pp;
297 if(vl.start()) do if(compute_cell(c,vl)) {
298 pp=p[vl.ij]+ps*vl.q;
299 fprintf(fp,"# [%d]\n",id[vl.ij][vl.q]);
300 c.draw_gnuplot(*pp,pp[1],fp);
301 fputs("\n",fp);
302 } while(vl.inc());
304 /** Computes all Voronoi cells and saves the output in gnuplot
305 * format.
306 * \param[in] fp a file handle to write to. */
307 inline void draw_cells_gnuplot(FILE *fp=stdout) {
308 c_loop_all_2d vl(*this);
309 draw_cells_gnuplot(vl,fp);
311 /** Computes all Voronoi cells and saves the output in gnuplot
312 * format.
313 * \param[in] filename the name of the file to write to. */
314 inline void draw_cells_gnuplot(const char *filename) {
315 FILE *fp=safe_fopen(filename,"w");
316 draw_cells_gnuplot(fp);
317 fclose(fp);
319 /** Computes Voronoi cells and saves the output in POV-Ray
320 * format.
321 * \param[in] vl the loop class to use.
322 * \param[in] fp a file handle to write to. */
323 template<class c_loop_2d>
324 void draw_cells_pov(c_loop_2d &vl,FILE *fp) {
325 voronoicell_nonconvex_2d c;double *pp;
326 if(vl.start()) do if(compute_cell(c,vl)) {
327 fprintf(fp,"// cell %d\n",id[vl.ij][vl.q]);
328 pp=p[vl.ij]+ps*vl.q;
329 c.draw_pov(*pp,pp[1],fp);
330 } while(vl.inc());
332 /** Computes all Voronoi cells and saves the output in POV-Ray
333 * format.
334 * \param[in] fp a file handle to write to. */
335 inline void draw_cells_pov(FILE *fp=stdout) {
336 c_loop_all_2d vl(*this);
337 draw_cells_pov(vl,fp);
339 /** Computes all Voronoi cells and saves the output in POV-Ray
340 * format.
341 * \param[in] filename the name of the file to write to. */
342 inline void draw_cells_pov(const char *filename) {
343 FILE *fp=safe_fopen(filename,"w");
344 draw_cells_pov(fp);
345 fclose(fp);
347 /** Computes the Voronoi cells and saves customized information
348 * about them.
349 * \param[in] vl the loop class to use.
350 * \param[in] format the custom output string to use.
351 * \param[in] fp a file handle to write to. */
352 template<class c_loop_2d>
353 void print_custom(c_loop_2d &vl,const char *format,FILE *fp) {
354 int ij,q;double *pp;
355 // bool glob=false, loc=false;
356 // if(contains_neighbor_global(format)){
357 // init_globne();
358 // glob=true;
359 // }
360 if(contains_neighbor(format)){
361 // loc=true;
362 // }
363 // if(glob || loc) {
364 voronoicell_nonconvex_neighbor_2d c;
365 if(vl.start()) do if(compute_cell(c,vl)) {
367 ij=vl.ij;q=vl.q;pp=p[ij]+ps*q;
368 // if(glob) add_globne_info(id[ij][q], c.ne, c.p);
369 c.output_custom(format,id[ij][q],*pp,pp[1],default_radius_2d,fp);
370 } while(vl.inc());
371 // if(glob) print_globne(fp);
372 } else {
373 voronoicell_nonconvex_2d c;
374 if(vl.start()) do if(compute_cell(c,vl)) {
375 ij=vl.ij;q=vl.q;pp=p[ij]+ps*q;
376 c.output_custom(format,id[ij][q],*pp,pp[1],default_radius_2d,fp);
377 } while(vl.inc());
380 void print_custom(const char *format,FILE *fp=stdout);
381 void print_custom(const char *format,const char *filename);
382 //bool find_voronoi_cell(double x,double y,double &rx,double &ry,int &pid);
383 /** Computes the Voronoi cell for a particle currently being
384 * referenced by a loop class.
385 * \param[out] c a Voronoi cell class in which to store the
386 * computed cell.
387 * \param[in] vl the loop class to use.
388 * \return True if the cell was computed. If the cell cannot be
389 * computed, if it is removed entirely by a wall or boundary
390 * condition, then the routine returns false. */
391 template<class v_cell_2d,class c_loop_2d>
392 inline bool compute_cell(v_cell_2d &c,c_loop_2d &vl) {
393 return vc.compute_cell(c,vl.ij,vl.q,vl.i,vl.j);
395 /** Computes the Voronoi cell for given particle.
396 * \param[out] c a Voronoi cell class in which to store the
397 * computed cell.
398 * \param[in] ij the block that the particle is within.
399 * \param[in] q the index of the particle within the block.
400 * \return True if the cell was computed. If the cell cannot be
401 * computed, if it is removed entirely by a wall or boundary
402 * condition, then the routine returns false. */
403 template<class v_cell_2d>
404 inline bool compute_cell(v_cell_2d &c,int ij,int q) {
405 int j=ij/nx,i=ij-j*nx;
406 return vc.compute_cell(c,ij,q,i,j);
408 void setup();
409 bool skip(int ij,int l,double x,double y);
410 private:
411 inline void draw_block(int ij) {
412 int i=ij%nx,j=ij/nx;
413 double lx=ax+i*boxx,ly=ay+j*boxy,ux=lx+boxx,uy=ly+boxy;
414 printf("%g %g\n%g %g\n%g %g\n%g %g\n%g %g\n\n\n",lx,ly,ux,ly,ux,uy,lx,uy,lx,ly);
416 int *soi;
417 int *tmp;
418 int *tmpp;
419 int *tmpe;
420 void create_label_table();
421 void add_particle_memory(int i);
422 inline bool put_locate_block(int &ij,double &x,double &y);
423 inline bool put_remap(int &ij,double &x,double &y);
424 inline bool remap(int &ai,int &aj,int &ci,int &cj,double &x,double &y,int &ij);
425 void add_temporary_label_memory();
426 void add_boundary_memory();
427 void tag_line(int &ij,int ije,int wid);
428 inline void tag(int ij,int wid);
429 void tag_walls(double x1,double y1,double x2,double y2,int wid);
430 inline bool cross_product(double x1,double y1,double x2,double y2) {
431 return x1*y2>x2*y1;
433 void semi_circle_labeling(double x1,double y1,double x2,double y2,int bid);
434 voro_compute_2d<container_boundary_2d> vc;
435 friend class voro_compute_2d<container_boundary_2d>;
440 #endif