00001 // Voro++, a 3D cell-based Voronoi library 00002 // 00003 // Author : Chris H. Rycroft (LBL / UC Berkeley) 00004 // Email : [email protected] 00005 // Date : July 1st 2008 00006 00007 /** \file voro++.hh 00008 * \brief A file that loads all of the Voro++ header files. */ 00009 00010 /** \mainpage Voro++ class reference manual 00011 * \section intro Introduction 00012 * Voro++ is a software library for carrying out 3D cell-based calculations of 00013 * the Voronoi tessellation. It is primarily designed for applications in 00014 * physics and materials science, where the Voronoi tessellation can be a 00015 * useful tool in analyzing particle systems. 00016 * 00017 * Voro++ is comprised of several C++ classes, and is designed to be 00018 * incorporated into other programs. This manual provides a reference for every 00019 * function in the class structure. For a general overview of the program, see 00020 * the Voro++ website at http://math.lbl.gov/voro++/ and in particular the 00021 * example programs at http://math.lbl.gov/voro++/examples/ that demonstrate 00022 * many of the library's features. 00023 * 00024 * \section class C++ class structure 00025 * The code is structured around two main C++ classes. The voronoicell class 00026 * contains all of the routines for constructing a single Voronoi cell. It 00027 * represents the cell as a collection of vertices that are connected by edges, 00028 * and there are routines for initializing, making, and outputting the cell. 00029 * The container class represents a three-dimensional simulation region into 00030 * which particles can be added. The class can then carry out a variety of 00031 * Voronoi calculations by computing cells using the voronoicell class. It also 00032 * has a general mechanism using virtual functions to implement walls, 00033 * discussed in more detail below. To implement the radical Voronoi 00034 * tessellation and the neighbor calculations, two class variants called 00035 * voronoicell_neighbor and container_poly are provided by making use of 00036 * templates, which is discussed below. 00037 * 00038 * \section voronoicell The voronoicell class 00039 * The voronoicell class represents a single Voronoi cell as a convex 00040 * polyhedron, with a set of vertices that are connected by edges. The class 00041 * contains a variety of functions that can be used to compute and output the 00042 * Voronoi cell corresponding to a particular particle. The command init() 00043 * can be used to initialize a cell as a large rectangular box. The Voronoi cell 00044 * can then be computed by repeatedly cutting it with planes that correspond to 00045 * the perpendicular bisectors between that particle and its neighbors. 00046 * 00047 * This is achieved by using the plane() routine, which will recompute the 00048 * cell's vertices and edges after cutting it with a single plane. This is the 00049 * key routine in voronoicell class. It begins by exploiting the convexity 00050 * of the underlying cell, tracing between edges to work out if the cell 00051 * intersects the cutting plane. If it does not intersect, then the routine 00052 * immediately exits. Otherwise, it finds an edge or vertex that intersects 00053 * the plane, and from there, traces out a new face on the cell, recomputing 00054 * the edge and vertex structure accordingly. 00055 * 00056 * Once the cell is computed, it can be drawn using commands such as 00057 * draw_gnuplot() and draw_pov(), or its volume can be evaluated using the 00058 * volume() function. Many more routines are available, and are described in 00059 * the online reference manual. 00060 * 00061 * \subsection internal Internal data representation 00062 * The voronoicell class has a public member p representing the 00063 * number of vertices. The polyhedral structure of the cell is stored 00064 * in the following arrays: 00065 * 00066 * - pts[]: an array of floating point numbers, that represent the position 00067 * vectors x_0, x_1, ..., x_{p-1} of the polyhedron vertices. 00068 * - nu[]: the order of each vertex n_0, n_1, ..., n_{p-1}, corresponding to 00069 * the number of other vertices to which each is connected. 00070 * - ed[][]: a table of edges and relations. For the ith vertex, ed[i] has 00071 * 2n_i+1 elements. The first n_i elements are the edges e(j,i), where e(j,i) 00072 * is the jth neighbor of vertex i. The edges are ordered according to a 00073 * right-hand rule with respect to an outward-pointing normal. The next n_i 00074 * elements are the relations l(j,i) which satisfy the property 00075 * e(l(j,i),e(j,i)) = i. The final element of the ed[i] list is a back 00076 * pointer used in memory allocation. 00077 * 00078 * In a very large number of cases, the values of n_i will be 3. This is because 00079 * the only way that a higher-order vertex can be created in the plane() 00080 * routine is if the cutting plane perfectly intersects an existing vertex. For 00081 * random particle arrangements with position vectors specified to double 00082 * precision this should happen very rarely. A preliminary version of this code 00083 * was quite successful with only making use of vertices of order 3. However, 00084 * when calculating millions of cells, it was found that this approach is not 00085 * robust, since a single floating point error can invalidate the computation. 00086 * This can also be a problem for cases featuring crystalline arrangements of 00087 * particles where the corresponding Voronoi cells may have high-order vertices 00088 * by construction. 00089 * 00090 * Because of this, Voro++ takes the approach that it if an existing vertex is 00091 * within a small numerical tolerance of the cutting plane, it is treated as 00092 * being exactly on the plane, and the polyhedral topology is recomputed 00093 * accordingly. However, while this improves robustness, it also adds the 00094 * complexity that n_i may no longer always be 3. This causes memory management 00095 * to be significantly more complicated, as different vertices require a 00096 * different number of elements in the ed[][] array. To accommodate this, the 00097 * voronoicell class allocated edge memory in a different array called mep[][], 00098 * in such a way that all vertices of order k are held in mep[k]. If vertex 00099 * i has order k, then ed[i] points to memory within mep[k]. The array ed[][] 00100 * is never directly initialized as a two-dimensional array itself, but points 00101 * at allocations within mep[][]. To the user, it appears as though each row of 00102 * ed[][] has a different number of elements. When vertices are added or 00103 * deleted, care must be taken to reorder and reassign elements in these 00104 * arrays. 00105 * 00106 * During the plane() routine, the code traces around the vertices of the cell, 00107 * and adds new vertices along edges which intersect the cutting plane to 00108 * create a new face. The values of l(j,i) are used in this computation, as 00109 * when the code is traversing from one vertex on the cell to another, this 00110 * information allows the code to immediately work out which edge of a vertex 00111 * points back to the one it came from. As new vertices are created, the 00112 * l(j,i) are also updated to ensure consistency. To ensure robustness, the 00113 * plane cutting algorithm should work with any possible combination of 00114 * vertices which are inside, outside, or exactly on the cutting plane. 00115 * 00116 * Vertices exactly on the cutting plane create some additional computational 00117 * difficulties. If there are two marginal vertices connected by an existing 00118 * edge, then it would be possible for duplicate edges to be created between 00119 * those two vertices, if the plane routine traces along both sides of this 00120 * edge while constructing the new face. The code recognizes these cases and 00121 * prevents the double edge from being formed. Another possibility is the 00122 * formation of vertices of order two or one. At the end of the plane cutting 00123 * routine, the code checks to see if any of these are present, removing the 00124 * order one vertices by just deleting them, and removing the order two 00125 * vertices by connecting the two neighbors of each vertex together. It is 00126 * possible that the removal of a single low-order vertex could result in the 00127 * creation of additional low-order vertices, so the process is applied 00128 * recursively until no more are left. 00129 * 00130 * \section container The container class 00131 * The container class represents a three-dimensional rectangular box of 00132 * particles. The constructor for this class sets up the coordinate ranges, 00133 * sets whether each direction is periodic or not, and divides the box into a 00134 * rectangular subgrid of regions. Particles can be added to the container 00135 * using the put() command, that adds a particle's position and an integer 00136 * numerical ID label to the corresponding region. Alternatively, the command 00137 * import() can be used to read large numbers of particles from a text file. 00138 * 00139 * The key routine in this class is compute_cell(), which makes use of the 00140 * voronoicell class to construct a Voronoi cell for a specific particle in the 00141 * container. The basic approach that this function takes is to repeatedly cut 00142 * the Voronoi cell by planes corresponding neighboring particles, and stop 00143 * when it recognizes that all the remaining particles in the container are too 00144 * far away to possibly influence cell's shape. The code makes use of two 00145 * possible methods for working out when a cell computation is complete: 00146 * 00147 * - Radius test: if the maximum distance of a Voronoi cell 00148 * vertex from the cell center is R, then no particles more than a distance 00149 * 2R away can possibly influence the cell. This a very fast computation to 00150 * do, but it has no directionality: if the cell extends a long way in one 00151 * direction then particles a long distance in other directions will still 00152 * need to be tested. 00153 * - Region test: it is possible to test whether a specific region can 00154 * possibly influence the cell by applying a series of plane tests at the 00155 * point on the region which is closest to the Voronoi cell center. This is a 00156 * slower computation to do, but it has directionality. 00157 * 00158 * Another useful observation is that the regions that need to be tested 00159 * are simply connected, meaning that if a particular region does not need 00160 * to be tested, then neighboring regions which are further away do not 00161 * need to be tested. 00162 * 00163 * For maximum efficiency, it was found that a hybrid approach making use of both 00164 * of the above tests worked well in practice. Radius tests work well for the 00165 * first few blocks, but switching to region tests after then prevent the code 00166 * from becoming extremely slow, due to testing over very large spherical shells of 00167 * particles. The compute_cell() routine therefore takes the following 00168 * approach: 00169 * 00170 * - Initialize the voronoicell class to fill the entire computational domain. 00171 * - Cut the cell by any wall objects that have been added to the container. 00172 * - Apply plane cuts to the cell corresponding to the other particles which 00173 * are within the current particle's region. 00174 * - Test over a pre-computed worklist of neighboring regions, that have been 00175 * ordered according to the minimum distance away from the particle's 00176 * position. Apply radius tests after every few regions to see if the 00177 * calculation can terminate. 00178 * - If the code reaches the end of the worklist, add all the neighboring 00179 * regions to a new list. 00180 * - Carry out a region test on the first item of the list. If the region needs 00181 * to be tested, apply the plane() routine for all of its particles, and then 00182 * add any neighboring regions to the end of the list that need to be tested. 00183 * Continue until the list has no elements left. 00184 * 00185 * The compute_cell() routine forms the basis of many other routines, such as 00186 * store_cell_volumes() and draw_cells_gnuplot() that can be used to calculate 00187 * and draw the cells in the entire container or in a subdomain. 00188 * 00189 * \section walls Wall computation 00190 * Wall computations are handled by making use of a pure virtual wall class. 00191 * Specific wall types are derived from this class, and require the 00192 * specification of two routines: point_inside() that tests to see if a point 00193 * is inside a wall or not, and cut_cell() that cuts a cell according to the 00194 * wall's position. The walls can be added to the container using the 00195 * add_wall() command, and these are called each time a compute_cell() command 00196 * is carried out. At present, wall types for planes, spheres, cylinders, and 00197 * cones are provided, although custom walls can be added by creating new 00198 * classes derived from the pure virtual class. Currently all wall types 00199 * approximate the wall surface with a single plane, which produces some small 00200 * errors, but generally gives good results for dense particle packings in 00201 * direct contact with a wall surface. It would be possible to create more 00202 * accurate walls by making cut_cell() routines that approximate the curved 00203 * surface with multiple plane cuts. 00204 * 00205 * The wall objects can used for periodic calculations, although to obtain 00206 * valid results, the walls should also be periodic as well. For example, in a 00207 * domain that is periodic in the x direction, a cylinder aligned along the x 00208 * axis could be added. At present, the interior of all wall objects are convex 00209 * domains, and consequently any superposition of them will be a convex domain 00210 * also. Carrying out computations in non-convex domains poses some problems, 00211 * since this could theoretically lead to non-convex Voronoi cells, which the 00212 * internal data representation of the voronoicell class does not support. For 00213 * non-convex cases where the wall surfaces feature just a small amount of 00214 * negative curvature (eg. a torus) approximating the curved surface with a 00215 * single plane cut may give an acceptable level of accuracy. For non-convex 00216 * cases that feature internal angles, the best strategy may be to decompose 00217 * the domain into several convex subdomains, carry out a calculation in each, 00218 * and then add the results together. The voronoicell class cannot be easily 00219 * modified to handle non-convex cells as this would fundamentally alter the 00220 * algorithms that it uses, and cases could arise where a single plane cut 00221 * could create several new faces as opposed to just one. 00222 * 00223 * \section templates Extra functionality via the use of templates 00224 * C++ templates are often presented as a mechanism for allowing functions to 00225 * be coded to work with several different data types. However, they also 00226 * provide an extremely powerful mechanism for achieving static polymorphism, 00227 * allowing several variations of a program to be compiled from a single source 00228 * code. Voro++ makes use of templates in order to handle the radical Voronoi 00229 * tessellation and the neighbor calculations, both of which require only 00230 * relatively minimal alterations to the main body of code. 00231 * 00232 * The main body of the voronoicell class is written as template named 00233 * voronoicell_base. Two additional small classes are then written: 00234 * neighbor_track, which contains small, inlined functions that encapsulate all 00235 * of the neighbor calculations, and neighbor_none, which contains the same 00236 * function names left blank. By making use of the typedef command, two classes 00237 * are then created from the template: 00238 * 00239 * - voronoicell: an instance of voronoicell_base with the neighbor_none class. 00240 * - voronoicell_neighbor: an instance of voronoicell_base with the 00241 * neighbor_track class. 00242 * 00243 * The two classes will be the same, except that the second will get all of the 00244 * additional neighbor-tracking functionality compiled into it through the 00245 * neighbor_track class. Since the two instances of the template are created 00246 * during the compilation, and since all of the functions in neighbor_none and 00247 * neighbor_track are inlined, there should be no speed overhead with this 00248 * construction; it should have the same efficiency as writing two completely 00249 * separate classes. C++ has other methods for achieving similar results, such 00250 * as virtual functions and class inheritance, but these are more focused on 00251 * dynamic polymorphism, switching between functionality at run-time, resulting 00252 * in a drop in performance. This would be particularly apparent in this case, 00253 * as the neighbor computation code, while small, is heavily integrated into 00254 * the low-level details of the plane() routine, and a virtual function 00255 * approach would require a very large number of function address look-ups. 00256 * 00257 * In a similar manner, two small classes called radius_mono and radius_poly 00258 * are provided. The first contains all routines suitable for calculate the 00259 * standard Voronoi tessellation associated with a monodisperse particle 00260 * packing, while the second incorporates variations to carry out the radical 00261 * Voronoi tessellation associated with a polydisperse particle packing. Two 00262 * classes are then created via typedef commands: 00263 * 00264 * - container: an instance of container_base with the radius_mono class. 00265 * - container_poly: an instance of container_base with the radius_poly class. 00266 * 00267 * The container_poly class accepts an additional variable in the put() command 00268 * for the particle's radius. These radii are then used to weight the plane 00269 * positions in the compute_cell() routine. 00270 * 00271 * It should be noted that the underlying template structure is largely hidden 00272 * from a typical user accessing the library's functionality, and as 00273 * demonstrated in the examples, the classes listed above behave like regular 00274 * C++ classes, and can be used in all the same ways. However, the template 00275 * structure may provide an additional method of customizing the code; for 00276 * example, an additional radius class could be written to implement a Voronoi 00277 * tessellation variant. */ 00278 00279 #ifndef VOROPP_HH 00280 #define VOROPP_HH 00281 00282 #include "cell.hh" 00283 #include "container.hh" 00284 #include "wall.hh" 00285 00286 #endif