Voro++
Public Member Functions | Friends
voro::container_poly Class Reference

Extension of the container_base class for computing radical Voronoi tessellations. More...

#include <container.hh>

+ Inheritance diagram for voro::container_poly:

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)
 
- Public Member Functions inherited from voro::container_base
 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 ()
 
- 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_)
 
- Public Member Functions inherited from voro::wall_list
 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 ()
 
- Public Member Functions inherited from voro::radius_poly
 radius_poly ()
 

Friends

class voro_compute< container_poly >
 

Additional Inherited Members

- Data Fields inherited from voro::container_base
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
 
- 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
 
- Data Fields inherited from voro::wall_list
wall ** walls
 
wall ** wep
 
- Data Fields inherited from voro::radius_poly
double ** ppr
 
double max_radius
 
- 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_base
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)
 
- 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::wall_list
void increase_wall_memory ()
 
- Protected Member Functions inherited from voro::radius_poly
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::wall_list
wall ** wel
 
int current_wall_size
 

Detailed Description

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.

Constructor & Destructor Documentation

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.

Parameters
[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_memthe initial memory allocation for each block.

Definition at line 79 of file container.cc.

Member Function Documentation

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.

template<class v_cell , class c_loop >
bool voro::container_poly::compute_cell ( v_cell &  c,
c_loop &  vl 
)
inline

Computes the Voronoi cell for a particle currently being referenced by a loop class.

Parameters
[out]ca Voronoi cell class in which to store the computed cell.
[in]vlthe loop class to use.
Returns
True if the cell was computed. If the cell cannot be computed, if it is removed entirely by a wall or boundary condition, then the routine returns false.

Definition at line 681 of file container.hh.

template<class v_cell >
bool voro::container_poly::compute_cell ( v_cell &  c,
int  ijk,
int  q 
)
inline

Computes the Voronoi cell for given particle.

Parameters
[out]ca Voronoi cell class in which to store the computed cell.
[in]ijkthe block that the particle is within.
[in]qthe index of the particle within the block.
Returns
True if the cell was computed. If the cell cannot be computed, if it is removed entirely by a wall or boundary condition, then the routine returns false.

Definition at line 693 of file container.hh.

template<class v_cell >
bool voro::container_poly::compute_ghost_cell ( v_cell &  c,
double  x,
double  y,
double  z,
double  r 
)
inline

Computes the Voronoi cell for a ghost particle at a given location.

Parameters
[out]ca Voronoi cell class in which to store the computed cell.
[in](x,y,z)the location of the ghost particle.
[in]rthe radius of the ghost particle.
Returns
True if the cell was computed. If the cell cannot be computed, if it is removed entirely by a wall or boundary condition, then the routine returns false.

Definition at line 707 of file container.hh.

template<class c_loop >
void voro::container_poly::draw_cells_gnuplot ( c_loop &  vl,
FILE *  fp 
)
inline

Computes Voronoi cells and saves the output in gnuplot format.

Parameters
[in]vlthe loop class to use.
[in]fpa file handle to write to.

Definition at line 600 of file container.hh.

void voro::container_poly::draw_cells_gnuplot ( FILE *  fp = stdout)
inline

Compute all Voronoi cells and saves the output in gnuplot format.

Parameters
[in]fpa file handle to write to.

Definition at line 610 of file container.hh.

void voro::container_poly::draw_cells_gnuplot ( const char *  filename)
inline

Compute all Voronoi cells and saves the output in gnuplot format.

Parameters
[in]filenamethe name of the file to write to.

Definition at line 617 of file container.hh.

template<class c_loop >
void voro::container_poly::draw_cells_pov ( c_loop &  vl,
FILE *  fp 
)
inline

Computes Voronoi cells and saves the output in POV-Ray format.

Parameters
[in]vlthe loop class to use.
[in]fpa file handle to write to.

Definition at line 627 of file container.hh.

void voro::container_poly::draw_cells_pov ( FILE *  fp = stdout)
inline

Computes all Voronoi cells and saves the output in POV-Ray format.

Parameters
[in]fpa file handle to write to.

Definition at line 638 of file container.hh.

void voro::container_poly::draw_cells_pov ( const char *  filename)
inline

Computes all Voronoi cells and saves the output in POV-Ray format.

Parameters
[in]filenamethe name of the file to write to.

Definition at line 645 of file container.hh.

template<class c_loop >
void voro::container_poly::draw_particles ( c_loop &  vl,
FILE *  fp 
)
inline

Dumps particle IDs, positions and radii to a file.

Parameters
[in]vlthe loop class to use.
[in]fpa file handle to write to.

Definition at line 548 of file container.hh.

void voro::container_poly::draw_particles ( FILE *  fp = stdout)
inline

Dumps all of the particle IDs, positions and radii to a file.

Parameters
[in]fpa file handle to write to.

Definition at line 558 of file container.hh.

void voro::container_poly::draw_particles ( const char *  filename)
inline

Dumps all of the particle IDs, positions and radii to a file.

Parameters
[in]filenamethe name of the file to write to.

Definition at line 565 of file container.hh.

template<class c_loop >
void voro::container_poly::draw_particles_pov ( c_loop &  vl,
FILE *  fp 
)
inline

Dumps particle positions in POV-Ray format.

Parameters
[in]vlthe loop class to use.
[in]fpa file handle to write to.

Definition at line 574 of file container.hh.

void voro::container_poly::draw_particles_pov ( FILE *  fp = stdout)
inline

Dumps all the particle positions in POV-Ray format.

Parameters
[in]fpa file handle to write to.

Definition at line 584 of file container.hh.

void voro::container_poly::draw_particles_pov ( const char *  filename)
inline

Dumps all the particle positions in POV-Ray format.

Parameters
[in]filenamethe 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.

Parameters
[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]pidthe ID of the particle.
Returns
True if a particle was found. If the container has no particles, then the search will not find a Voronoi cell and false is returned.

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.

Parameters
[in]fpthe 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.

Parameters
[in,out]voa reference to an ordering class to use.
[in]fpthe file handle to read from.

Definition at line 368 of file container.cc.

void voro::container_poly::import ( const char *  filename)
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.

Parameters
[in]filenamethe name of the file to open and read from.

Definition at line 523 of file container.hh.

void voro::container_poly::import ( particle_order vo,
const char *  filename 
)
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.

Parameters
[in,out]vothe ordering class to use.
[in]filenamethe name of the file to open and read from.

Definition at line 537 of file container.hh.

template<class c_loop >
void voro::container_poly::print_custom ( c_loop &  vl,
const char *  format,
FILE *  fp 
)
inline

Computes the Voronoi cells and saves customized information about them.

Parameters
[in]vlthe loop class to use.
[in]formatthe custom output string to use.
[in]fpa 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.

Parameters
[in]formatthe custom output string to use.
[in]fpa 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

Parameters
[in]formatthe custom output string to use.
[in]filenamethe 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.

Parameters
[in]nthe numerical ID of the inserted particle.
[in](x,y,z)the position vector of the inserted particle.
[in]rthe 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.

Parameters
[in]vothe ordering class in which to record the region.
[in]nthe numerical ID of the inserted particle.
[in](x,y,z)the position vector of the inserted particle.
[in]rthe 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.

Returns
The sum of all of the computed Voronoi volumes.

Definition at line 468 of file container.cc.


The documentation for this class was generated from the following files: