Voro++
|
Extension of the container_periodic_base class for computing regular Voronoi tessellations. More...
#include <container_prd.hh>
Public Member Functions | |
container_periodic (double bx_, double bxy_, double by_, double bxz_, double byz_, double bz_, int nx_, int ny_, int nz_, int init_mem_) | |
void | clear () |
void | put (int n, double x, double y, double z) |
void | put (int n, double x, double y, double z, int &ai, int &aj, int &ak) |
void | put (particle_order &vo, int n, double x, double y, double z) |
void | import (FILE *fp=stdin) |
void | import (particle_order &vo, FILE *fp=stdin) |
void | import (const char *filename) |
void | import (particle_order &vo, const char *filename) |
void | compute_all_cells () |
double | sum_cell_volumes () |
template<class c_loop > | |
void | draw_particles (c_loop &vl, FILE *fp) |
void | draw_particles (FILE *fp=stdout) |
void | draw_particles (const char *filename) |
template<class c_loop > | |
void | draw_particles_pov (c_loop &vl, FILE *fp) |
void | draw_particles_pov (FILE *fp=stdout) |
void | draw_particles_pov (const char *filename) |
template<class c_loop > | |
void | draw_cells_gnuplot (c_loop &vl, FILE *fp) |
void | draw_cells_gnuplot (FILE *fp=stdout) |
void | draw_cells_gnuplot (const char *filename) |
template<class c_loop > | |
void | draw_cells_pov (c_loop &vl, FILE *fp) |
void | draw_cells_pov (FILE *fp=stdout) |
void | draw_cells_pov (const char *filename) |
template<class c_loop > | |
void | print_custom (c_loop &vl, const char *format, FILE *fp) |
void | print_custom (const char *format, FILE *fp=stdout) |
void | print_custom (const char *format, const char *filename) |
bool | find_voronoi_cell (double x, double y, double z, double &rx, double &ry, double &rz, int &pid) |
template<class v_cell , class c_loop > | |
bool | compute_cell (v_cell &c, c_loop &vl) |
template<class v_cell > | |
bool | compute_cell (v_cell &c, int ijk, int q) |
template<class v_cell > | |
bool | compute_ghost_cell (v_cell &c, double x, double y, double z) |
Public Member Functions inherited from voro::container_periodic_base | |
container_periodic_base (double bx_, double bxy_, double by_, double bxz_, double byz_, double bz_, int nx_, int ny_, int nz_, int init_mem_, int ps) | |
~container_periodic_base () | |
void | print_all_particles () |
void | region_count () |
template<class v_cell > | |
bool | initialize_voronoicell (v_cell &c, int ijk, int q, int ci, int cj, int ck, int &i, int &j, int &k, double &x, double &y, double &z, int &disp) |
void | initialize_search (int ci, int cj, int ck, int ijk, int &i, int &j, int &k, int &disp) |
void | frac_pos (double x, double y, double z, double ci, double cj, double ck, double &fx, double &fy, double &fz) |
int | region_index (int ci, int cj, int ck, int ei, int ej, int ek, double &qx, double &qy, double &qz, int &disp) |
void | create_all_images () |
void | check_compartmentalized () |
Public Member Functions inherited from voro::unitcell | |
unitcell (double bx_, double bxy_, double by_, double bxz_, double byz_, double bz_) | |
void | draw_domain_gnuplot (const char *filename) |
void | draw_domain_gnuplot (FILE *fp=stdout) |
void | draw_domain_pov (const char *filename) |
void | draw_domain_pov (FILE *fp=stdout) |
bool | intersects_image (double dx, double dy, double dz, double &vol) |
void | images (std::vector< int > &vi, std::vector< double > &vd) |
Public Member Functions inherited from voro::voro_base | |
bool | contains_neighbor (const char *format) |
voro_base (int nx_, int ny_, int nz_, double boxx_, double boxy_, double boxz_) | |
Friends | |
class | voro_compute< container_periodic > |
Additional Inherited Members | |
Data Fields inherited from voro::container_periodic_base | |
int | ey |
int | ez |
int | wy |
int | wz |
int | oy |
int | oz |
int | oxyz |
int ** | id |
double ** | p |
int * | co |
int * | mem |
char * | img |
const int | init_mem |
const int | ps |
Data Fields inherited from voro::unitcell | |
const double | bx |
const double | bxy |
const double | by |
const double | bxz |
const double | byz |
const double | bz |
voronoicell | unit_voro |
Data Fields inherited from voro::voro_base | |
const int | nx |
const int | ny |
const int | nz |
const int | nxy |
const int | nxyz |
const double | boxx |
const double | boxy |
const double | boxz |
const double | xsp |
const double | ysp |
const double | zsp |
double * | mrad |
Static Public Attributes inherited from voro::voro_base | |
static const unsigned int | wl [wl_seq_length *wl_hgridcu] |
Protected Member Functions inherited from voro::container_periodic_base | |
void | add_particle_memory (int i) |
void | put_locate_block (int &ijk, double &x, double &y, double &z) |
void | put_locate_block (int &ijk, double &x, double &y, double &z, int &ai, int &aj, int &ak) |
void | create_periodic_image (int di, int dj, int dk) |
void | create_side_image (int di, int dj, int dk) |
void | create_vertical_image (int di, int dj, int dk) |
void | put_image (int reg, int fijk, int l, double dx, double dy, double dz) |
void | remap (int &ai, int &aj, int &ak, int &ci, int &cj, int &ck, double &x, double &y, double &z, int &ijk) |
Protected Member Functions inherited from voro::voro_base | |
int | step_int (double a) |
int | step_mod (int a, int b) |
int | step_div (int a, int b) |
Protected Member Functions inherited from voro::radius_mono | |
void | r_init (int ijk, int s) |
void | r_prime (double rv) |
bool | r_ctest (double crs, double mrs) |
double | r_cutoff (double lrs) |
double | r_max_add (double rs) |
double | r_current_sub (double rs, int ijk, int q) |
double | r_scale (double rs, int ijk, int q) |
bool | r_scale_check (double &rs, double mrs, int ijk, int q) |
Protected Attributes inherited from voro::unitcell | |
double | max_uv_y |
double | max_uv_z |
This class is an extension of the container_periodic_base that has routines specifically for computing the regular Voronoi tessellation with no dependence on particle radii.
Definition at line 212 of file container_prd.hh.
voro::container_periodic::container_periodic | ( | double | bx_, |
double | bxy_, | ||
double | by_, | ||
double | bxz_, | ||
double | byz_, | ||
double | bz_, | ||
int | nx_, | ||
int | ny_, | ||
int | nz_, | ||
int | init_mem_ | ||
) |
The class constructor sets up the geometry of container.
[in] | (bx_) | The x coordinate of the first unit vector. |
[in] | (bxy_,by_) | The x and y coordinates of the second unit vector. |
[in] | (bxz_,byz_,bz_) | The x, y, and z coordinates of the third unit vector. |
[in] | (nx_,ny_,nz_) | the number of grid blocks in each of the three coordinate directions. |
[in] | init_mem_ | the initial memory allocation for each block. |
Definition at line 72 of file container_prd.cc.
void voro::container_periodic::clear | ( | ) |
Clears a container of particles.
Definition at line 448 of file container_prd.cc.
void voro::container_periodic::compute_all_cells | ( | ) |
Computes all of the Voronoi cells in the container, but does nothing with the output. It is useful for measuring the pure computation time of the Voronoi algorithm, without any additional calculations such as volume evaluation or cell output.
Definition at line 499 of file container_prd.cc.
|
inline |
Computes the Voronoi cell for a particle currently being referenced by a loop class.
[out] | c | a Voronoi cell class in which to store the computed cell. |
[in] | vl | the loop class to use. |
Definition at line 388 of file container_prd.hh.
|
inline |
Computes the Voronoi cell for given particle.
[out] | c | a Voronoi cell class in which to store the computed cell. |
[in] | ijk | the block that the particle is within. |
[in] | q | the index of the particle within the block. |
Definition at line 400 of file container_prd.hh.
|
inline |
Computes the Voronoi cell for a ghost particle at a given location.
[out] | c | a Voronoi cell class in which to store the computed cell. |
[in] | (x,y,z) | the location of the ghost particle. |
Definition at line 413 of file container_prd.hh.
|
inline |
Computes Voronoi cells and saves the output in gnuplot format.
[in] | vl | the loop class to use. |
[in] | fp | a file handle to write to. |
Definition at line 304 of file container_prd.hh.
|
inline |
Computes all Voronoi cells and saves the output in gnuplot format.
[in] | fp | a file handle to write to. |
Definition at line 314 of file container_prd.hh.
|
inline |
Compute all Voronoi cells and saves the output in gnuplot format.
[in] | filename | the name of the file to write to. |
Definition at line 321 of file container_prd.hh.
|
inline |
Computes Voronoi cells and saves the output in POV-Ray format.
[in] | vl | the loop class to use. |
[in] | fp | a file handle to write to. |
Definition at line 331 of file container_prd.hh.
|
inline |
Computes all Voronoi cells and saves the output in POV-Ray format.
[in] | fp | a file handle to write to. |
Definition at line 342 of file container_prd.hh.
|
inline |
Computes all Voronoi cells and saves the output in POV-Ray format.
[in] | filename | the name of the file to write to. |
Definition at line 349 of file container_prd.hh.
|
inline |
Dumps particle IDs and positions to a file.
[in] | vl | the loop class to use. |
[in] | fp | a file handle to write to. |
Definition at line 254 of file container_prd.hh.
|
inline |
Dumps all of the particle IDs and positions to a file.
[in] | fp | a file handle to write to. |
Definition at line 263 of file container_prd.hh.
|
inline |
Dumps all of the particle IDs and positions to a file.
[in] | filename | the name of the file to write to. |
Definition at line 269 of file container_prd.hh.
|
inline |
Dumps particle positions in POV-Ray format.
[in] | vl | the loop class to use. |
[in] | fp | a file handle to write to. |
Definition at line 278 of file container_prd.hh.
|
inline |
Dumps all particle positions in POV-Ray format.
[in] | fp | a file handle to write to. |
Definition at line 288 of file container_prd.hh.
|
inline |
Dumps all particle positions in POV-Ray format.
[in] | filename | the name of the file to write to. |
Definition at line 294 of file container_prd.hh.
bool voro::container_periodic::find_voronoi_cell | ( | double | x, |
double | y, | ||
double | z, | ||
double & | rx, | ||
double & | ry, | ||
double & | rz, | ||
int & | pid | ||
) |
Takes a vector and finds the particle whose Voronoi cell contains that vector. This is equivalent to finding the particle which is nearest to the vector.
[in] | (x,y,z) | the vector to test. |
[out] | (rx,ry,rz) | the position of the particle whose Voronoi cell contains the vector. This may point to a particle in a periodic image of the primary domain. |
[out] | pid | the ID of the particle. |
Definition at line 297 of file container_prd.cc.
void voro::container_periodic::import | ( | FILE * | fp = stdin | ) |
Import a list of particles from an open file stream into the container. Entries of four numbers (Particle ID, x position, y position, z position) are searched for. If the file cannot be successfully read, then the routine causes a fatal error.
[in] | fp | the file handle to read from. |
Definition at line 394 of file container_prd.cc.
void voro::container_periodic::import | ( | particle_order & | vo, |
FILE * | fp = stdin |
||
) |
Import a list of particles from an open file stream, also storing the order of that the particles are read. Entries of four numbers (Particle ID, x position, y position, z position) are searched for. If the file cannot be successfully read, then the routine causes a fatal error.
[in,out] | vo | a reference to an ordering class to use. |
[in] | fp | the file handle to read from. |
Definition at line 407 of file container_prd.cc.
|
inline |
Imports a list of particles from an open file stream into the container. Entries of four numbers (Particle ID, x position, y position, z position) are searched for. If the file cannot be successfully read, then the routine causes a fatal error.
[in] | filename | the name of the file to open and read from. |
Definition at line 229 of file container_prd.hh.
|
inline |
Imports a list of particles from an open file stream into the container. Entries of four numbers (Particle ID, x position, y position, z position) are searched for. In addition, the order in which particles are read is saved into an ordering class. If the file cannot be successfully read, then the routine causes a fatal error.
[in,out] | vo | the ordering class to use. |
[in] | filename | the name of the file to open and read from. |
Definition at line 243 of file container_prd.hh.
|
inline |
Computes the Voronoi cells and saves customized information about them.
[in] | vl | the loop class to use. |
[in] | format | the custom output string to use. |
[in] | fp | a file handle to write to. |
Definition at line 360 of file container_prd.hh.
void voro::container_periodic::print_custom | ( | const char * | format, |
FILE * | fp = stdout |
||
) |
Computes all the Voronoi cells and saves customized information about them.
[in] | format | the custom output string to use. |
[in] | fp | a file handle to write to. |
Definition at line 462 of file container_prd.cc.
void voro::container_periodic::print_custom | ( | const char * | format, |
const char * | filename | ||
) |
Computes all the Voronoi cells and saves customized information about them.
[in] | format | the custom output string to use. |
[in] | filename | the name of the file to write to. |
Definition at line 479 of file container_prd.cc.
void voro::container_periodic::put | ( | int | n, |
double | x, | ||
double | y, | ||
double | z | ||
) |
Put a particle into the correct region of the container.
[in] | n | the numerical ID of the inserted particle. |
[in] | (x,y,z) | the position vector of the inserted particle. |
Definition at line 93 of file container_prd.cc.
void voro::container_periodic::put | ( | int | n, |
double | x, | ||
double | y, | ||
double | z, | ||
int & | ai, | ||
int & | aj, | ||
int & | ak | ||
) |
Put a particle into the correct region of the container.
[in] | n | the numerical ID of the inserted particle. |
[in] | (x,y,z) | the position vector of the inserted particle. |
[out] | (ai,aj,ak) | the periodic image displacement that the particle is in, with (0,0,0) corresponding to the primary domain. |
Definition at line 120 of file container_prd.cc.
void voro::container_periodic::put | ( | particle_order & | vo, |
int | n, | ||
double | x, | ||
double | y, | ||
double | z | ||
) |
Put a particle into the correct region of the container, also recording into which region it was stored.
[in] | vo | the ordering class in which to record the region. |
[in] | n | the numerical ID of the inserted particle. |
[in] | (x,y,z) | the position vector of the inserted particle. |
Definition at line 149 of file container_prd.cc.
double voro::container_periodic::sum_cell_volumes | ( | ) |
Calculates all of the Voronoi cells and sums their volumes. In most cases without walls, the sum of the Voronoi cell volumes should equal the volume of the container to numerical precision.
Definition at line 520 of file container_prd.cc.