Voro++
|
Extension of the container_base class for computing radical Voronoi tessellations. More...
#include <container.hh>
Public Member Functions | |
container_poly (double ax_, double bx_, double ay_, double by_, double az_, double bz_, int nx_, int ny_, int nz_, bool xperiodic_, bool yperiodic_, bool zperiodic_, int init_mem) | |
void | clear () |
void | put (int n, double x, double y, double z, double r) |
void | put (particle_order &vo, int n, double x, double y, double z, double r) |
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) |
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, double r) |
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) |
![]() | |
container_base (double ax_, double bx_, double ay_, double by_, double az_, double bz_, int nx_, int ny_, int nz_, bool xperiodic_, bool yperiodic_, bool zperiodic_, int init_mem, int ps_) | |
~container_base () | |
bool | point_inside (double x, double y, double z) |
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 | draw_domain_gnuplot (FILE *fp=stdout) |
void | draw_domain_gnuplot (const char *filename) |
void | draw_domain_pov (FILE *fp=stdout) |
void | draw_domain_pov (const char *filename) |
int | total_particles () |
![]() | |
bool | contains_neighbor (const char *format) |
voro_base (int nx_, int ny_, int nz_, double boxx_, double boxy_, double boxz_) | |
![]() | |
wall_list () | |
~wall_list () | |
void | add_wall (wall *w) |
void | add_wall (wall &w) |
void | add_wall (wall_list &wl) |
bool | point_inside_walls (double x, double y, double z) |
template<class c_class > | |
bool | apply_walls (c_class &c, double x, double y, double z) |
void | deallocate () |
![]() | |
radius_poly () | |
Friends | |
class | voro_compute< container_poly > |
Additional Inherited Members | |
![]() | |
const double | ax |
const double | bx |
const double | ay |
const double | by |
const double | az |
const double | bz |
const bool | xperiodic |
const bool | yperiodic |
const bool | zperiodic |
int ** | id |
double ** | p |
int * | co |
int * | mem |
const int | ps |
![]() | |
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 |
![]() | |
wall ** | walls |
wall ** | wep |
![]() | |
double ** | ppr |
double | max_radius |
![]() | |
static const unsigned int | wl [wl_seq_length *wl_hgridcu] |
![]() | |
void | add_particle_memory (int i) |
bool | put_locate_block (int &ijk, double &x, double &y, double &z) |
bool | put_remap (int &ijk, double &x, double &y, double &z) |
bool | remap (int &ai, int &aj, int &ak, int &ci, int &cj, int &ck, double &x, double &y, double &z, int &ijk) |
![]() | |
int | step_int (double a) |
int | step_mod (int a, int b) |
int | step_div (int a, int b) |
![]() | |
void | increase_wall_memory () |
![]() | |
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) |
![]() | |
wall ** | wel |
int | current_wall_size |
This class is an extension of container_base class that has routines specifically for computing the radical Voronoi tessellation that depends on the particle radii.
Definition at line 507 of file container.hh.
voro::container_poly::container_poly | ( | double | ax_, |
double | bx_, | ||
double | ay_, | ||
double | by_, | ||
double | az_, | ||
double | bz_, | ||
int | nx_, | ||
int | ny_, | ||
int | nz_, | ||
bool | xperiodic_, | ||
bool | yperiodic_, | ||
bool | zperiodic_, | ||
int | init_mem | ||
) |
The class constructor sets up the geometry of container.
[in] | (ax_,bx_) | the minimum and maximum x coordinates. |
[in] | (ay_,by_) | the minimum and maximum y coordinates. |
[in] | (az_,bz_) | the minimum and maximum z coordinates. |
[in] | (nx_,ny_,nz_) | the number of grid blocks in each of the three coordinate directions. |
[in] | (xperiodic_,yperiodic_,zperiodic_) | flags setting whether the container is periodic in each coordinate direction. |
[in] | init_mem | the initial memory allocation for each block. |
Definition at line 79 of file container.cc.
void voro::container_poly::clear | ( | ) |
Clears a container of particles, also clearing resetting the maximum radius to zero.
Definition at line 390 of file container.cc.
void voro::container_poly::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 446 of file container.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 681 of file container.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 693 of file container.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. |
[in] | r | the radius of the ghost particle. |
Definition at line 707 of file container.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 600 of file container.hh.
|
inline |
Compute all Voronoi cells and saves the output in gnuplot format.
[in] | fp | a file handle to write to. |
Definition at line 610 of file container.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 617 of file container.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 627 of file container.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 638 of file container.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 645 of file container.hh.
|
inline |
Dumps particle IDs, positions and radii to a file.
[in] | vl | the loop class to use. |
[in] | fp | a file handle to write to. |
Definition at line 548 of file container.hh.
|
inline |
Dumps all of the particle IDs, positions and radii to a file.
[in] | fp | a file handle to write to. |
Definition at line 558 of file container.hh.
|
inline |
Dumps all of the particle IDs, positions and radii to a file.
[in] | filename | the name of the file to write to. |
Definition at line 565 of file container.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 574 of file container.hh.
|
inline |
Dumps all the particle positions in POV-Ray format.
[in] | fp | a file handle to write to. |
Definition at line 584 of file container.hh.
|
inline |
Dumps all the particle positions in POV-Ray format.
[in] | filename | the name of the file to write to. |
Definition at line 590 of file container.hh.
bool voro::container_poly::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. Additional wall classes are not considered by this routine.
[in] | (x,y,z) | the vector to test. |
[out] | (rx,ry,rz) | the position of the particle whose Voronoi cell contains the vector. If the container is periodic, this may point to a particle in a periodic image of the primary domain. |
[out] | pid | the ID of the particle. |
Definition at line 272 of file container.cc.
void voro::container_poly::import | ( | FILE * | fp = stdin | ) |
Import a list of particles from an open file stream into the container. Entries of five numbers (Particle ID, x position, y position, z position, radius) 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 355 of file container.cc.
void voro::container_poly::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, radius) 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 368 of file container.cc.
|
inline |
Imports a list of particles from an open file stream into the container_poly class. Entries of five numbers (Particle ID, x position, y position, z position, radius) 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 523 of file container.hh.
|
inline |
Imports a list of particles from an open file stream into the container_poly class. Entries of five numbers (Particle ID, x position, y position, z position, radius) 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 537 of file container.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 656 of file container.hh.
void voro::container_poly::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 407 of file container.cc.
void voro::container_poly::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 425 of file container.cc.
void voro::container_poly::put | ( | int | n, |
double | x, | ||
double | y, | ||
double | z, | ||
double | r | ||
) |
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. |
[in] | r | the radius of the particle. |
Definition at line 100 of file container.cc.
void voro::container_poly::put | ( | particle_order & | vo, |
int | n, | ||
double | x, | ||
double | y, | ||
double | z, | ||
double | r | ||
) |
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. |
[in] | r | the radius of the particle. |
Definition at line 131 of file container.cc.
double voro::container_poly::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 468 of file container.cc.