ppap4lmp  0.7.2
ProRadialDistributionFunction Class Reference

ProRadialDistributionFunction computes radial distribution function (RDF). More...

#include <pro_radial_distribution_function.h>

Inheritance diagram for ProRadialDistributionFunction:
Collaboration diagram for ProRadialDistributionFunction:

Public Member Functions

 ProRadialDistributionFunction (const ElPtr &targets, const ElPtr &box)
 Constructor of ProRadialDistributionFunction class for a snapshot of simulation. More...
 
 ProRadialDistributionFunction (const Vec< std::pair< ElPtr, ElPtr >> &pairs)
 Constructor of ProThicknessProfile class for multiple snapshots of simulation. More...
 
virtual void prepare () override
 Resize number_traj, volume_traj and counts_traj. More...
 
virtual void finish () override
 Compute rdf and rdf_traj from number_traj, volume_traj and counts_traj. More...
 
void set_bin (double bin_width_, int n_bins_)
 Set width and number of bins in the distance axis. Each bin is a spherical shell of which center is positioned at a reference particle. More...
 
void bin_from_r_to_r_plus_dr (bool bin_from_r_=true)
 By default, the bins are [0, 0.5*dr), [0.5*dr, 1.5*dr), ... . If this method is called without a parameter, the bins are set as [0, dr), [dr, 2*dr), ... . More...
 
void beyond_half_box_length (bool beyond_half_=true)
 Use a sample including particles of which distance from a reference particle are larger than half the simulation box length. More...
 
const ArrayXd get_r_axis ()
 Get a one-dimensional array of distance, [0, dr, 2*dr, ... ], where dr is width of a bin. More...
 
const ArrayXdget_rdf ()
 Get the radial distribution function as a one-dimensional array. More...
 
const Vec< ArrayXd > & get_rdf_traj ()
 Get a sequence of radial distribution functions as list of one-dimensional arrays. More...
 
- Public Member Functions inherited from Processor
 Processor ()=default
 Constructor of Processor class (default).
 
virtual bool run ()
 Run analysis using i th Generator object in generators, where i = i_generator. More...
 
void startup ()
 Startup this Processor object. More...
 

Protected Member Functions

virtual void run_impl (const int index) override
 Implementation of analysis using an element of generators. More...
 
- Protected Member Functions inherited from Processor
template<class GEN >
void register_generator (const ShPtr< GEN > &gen)
 
template<class GEN >
void register_generators (const Vec< ShPtr< GEN >> &gens)
 
virtual void use_generator_at (const int i)
 Call Generator::generate_data of i th Generator object in generators. More...
 
virtual void finish_using_generator_at (const int i)
 Call Generator::finish_using_generated_data of i* th Generator object in generators. More...
 

Protected Attributes

int n_bins = 1
 
double bin_width = 1.0
 
bool bin_from_r = false
 
bool beyond_half = false
 
ArrayXd rdf
 
Vec< ArrayXdrdf_traj
 
Vec< int > number_traj
 
Vec< double > volume_traj
 
Vec< ArrayXicounts_traj
 
- Protected Attributes inherited from Processor
int n_generators
 
Vec< ShPtr< Generator > > generators
 

Detailed Description

ProRadialDistributionFunction computes radial distribution function (RDF).

An object of this class makes a one-dimensional array of radial distribution function of particles (can be atoms, beads, molecules ...). Each element of the array corresponds to a ratio of local density at distance r to global density, where r is distance from a reference particle. The local density at r defined as a division of the number of particles in a spherical shell with radius of r by volume of the spherical shell.

List of radial distribution functions is also computed for a multiple snapshots of a simulation.

Note that a particle which is one of special bonds of a reference particle is excluded from a sample used for computing radial distribution function. About special bonds, please see AddSpecialBonds class.

About usage in Python, please see pybind::py_pro_radial_distribution_function.

Definition at line 37 of file pro_radial_distribution_function.h.

Constructor & Destructor Documentation

◆ ProRadialDistributionFunction() [1/2]

ProRadialDistributionFunction::ProRadialDistributionFunction ( const ElPtr targets,
const ElPtr box 
)

Constructor of ProRadialDistributionFunction class for a snapshot of simulation.

Parameters
targets

An Element object for particles (can be atoms, beads, molecules...).

Required property
  • id : integer
  • x : float
  • y : float
  • z : float
  • mass : float (for deformation)
  • I_xx : float (for deformation)
  • I_yy : float (for deformation)
  • I_zz : float (for deformation)
  • I_xy : float (for deformation)
  • I_xz : float (for deformation)
  • I_yz : float (for deformation)
  • special-bonds : array of integers (optional)

boxAn Element object for a simulation box.
Required property
  • lo_x : float
  • lo_y : float
  • lo_z : float
  • hi_x : float
  • hi_y : float
  • hi_z : float

A GenDict object is constructed taking the targets and box, and then put into generators by register_generator.

Definition at line 20 of file pro_radial_distribution_function.cpp.

Here is the call graph for this function:

◆ ProRadialDistributionFunction() [2/2]

ProRadialDistributionFunction::ProRadialDistributionFunction ( const Vec< std::pair< ElPtr, ElPtr >> &  pairs)

Constructor of ProThicknessProfile class for multiple snapshots of simulation.

Parameters
pairsList of pairs of two Element objects: the first object is for particles and the second object is for a simulation box. If each Element object has an array data, all the array should have the same length.
Required property (first)
  • id : integer
  • x : float
  • y : float
  • z : float
  • mass : float (for deformation)
  • I_xx : float (for deformation)
  • I_yy : float (for deformation)
  • I_zz : float (for deformation)
  • I_xy : float (for deformation)
  • I_xz : float (for deformation)
  • I_yz : float (for deformation)
  • special-bonds : array of integers (optional)
Required property (second)
  • lo_x : float
  • lo_y : float
  • lo_z : float
  • hi_x : float
  • hi_y : float
  • hi_z : float

GenDict objects are constructed taking each element of the pairs, and then put into generators by register_generators.

Definition at line 30 of file pro_radial_distribution_function.cpp.

Here is the call graph for this function:

Member Function Documentation

◆ beyond_half_box_length()

void ProRadialDistributionFunction::beyond_half_box_length ( bool  beyond_half_ = true)

Use a sample including particles of which distance from a reference particle are larger than half the simulation box length.

Parameters
beyond_half_A boolean, whether to use particles beyond half the simulation box from a reference particle (default is true). This parameter is assigned to beyond_half.
Returns
None.

Definition at line 246 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ bin_from_r_to_r_plus_dr()

void ProRadialDistributionFunction::bin_from_r_to_r_plus_dr ( bool  bin_from_r_ = true)

By default, the bins are [0, 0.5*dr), [0.5*dr, 1.5*dr), ... . If this method is called without a parameter, the bins are set as [0, dr), [dr, 2*dr), ... .

Parameters
bin_from_r_A boolean, whether to use bin of which inner radius is r and outer radius is r+dr (default is true). This parameter is assigned to bin_from_r.
Returns
None.

Definition at line 238 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ finish()

void ProRadialDistributionFunction::finish ( )
overridevirtual

Compute rdf and rdf_traj from number_traj, volume_traj and counts_traj.

Returns
None.

Note that rdf and an average over elements in rdf_traj are identical only when volume of the simulation box and the number of particles are unchanged during the simulation.

Reimplemented from Processor.

Reimplemented in ProRadialDistributionFunctionWithDeformation.

Definition at line 172 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ get_r_axis()

const ArrayXd ProRadialDistributionFunction::get_r_axis ( )

Get a one-dimensional array of distance, [0, dr, 2*dr, ... ], where dr is width of a bin.

Returns
ArrayXd (Numpy-Array in Python).

Definition at line 254 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ get_rdf()

const ArrayXd & ProRadialDistributionFunction::get_rdf ( )

Get the radial distribution function as a one-dimensional array.

Returns
rdf.

Definition at line 268 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ get_rdf_traj()

const Vec< ArrayXd > & ProRadialDistributionFunction::get_rdf_traj ( )

Get a sequence of radial distribution functions as list of one-dimensional arrays.

Returns
rdf_traj.

Definition at line 275 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ prepare()

void ProRadialDistributionFunction::prepare ( )
overridevirtual

Resize number_traj, volume_traj and counts_traj.

Returns
None.

Reimplemented from Processor.

Reimplemented in ProRadialDistributionFunctionWithDeformation.

Definition at line 163 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

◆ run_impl()

void ProRadialDistributionFunction::run_impl ( const int  index)
overrideprotectedvirtual

Implementation of analysis using an element of generators.

Parameters
indexAn index in generators for a Generator object to be analyzed.
Returns
None.

To be compatible with multithreading computation, class members should not be modified in this method. If modification of some members is necessary, please consider defining them as Vec, and modifying only their i th element, where i = index.

I am sorry to say that the best documentation for this method is its source code...

Implements Processor.

Reimplemented in ProRadialDistributionFunctionWithDeformation.

Definition at line 46 of file pro_radial_distribution_function.cpp.

Here is the call graph for this function:

◆ set_bin()

void ProRadialDistributionFunction::set_bin ( double  bin_width_,
int  n_bins_ 
)

Set width and number of bins in the distance axis. Each bin is a spherical shell of which center is positioned at a reference particle.

Parameters
bin_width_Width of a bin. This parameter is assigned to bin_width.
n_bins_Number of bins. This parameter is assigned to n_bins.
Returns
None.

Definition at line 228 of file pro_radial_distribution_function.cpp.

Here is the caller graph for this function:

Member Data Documentation

◆ beyond_half

bool ProRadialDistributionFunction::beyond_half = false
protected

If this member is false (default), a particle of which distance from a reference particle are larger than half the simulation box length is excluded from a sample used for computing radial distribution function. Such particles can be included in the sample for a case of small simulation box.

Definition at line 68 of file pro_radial_distribution_function.h.

◆ bin_from_r

bool ProRadialDistributionFunction::bin_from_r = false
protected

If this member is true, the bins are [0, dr), [dr, 2*dr), ..., [(n_bins-1)*dr, n_bins*dr). If this member is false, the bins are [0, 0.5*dr), [0.5*dr, 1.5*dr), ..., [(n_bins-1.5)*dr, (n_bins-0.5)*dr). Default is false.

Definition at line 59 of file pro_radial_distribution_function.h.

◆ bin_width

double ProRadialDistributionFunction::bin_width = 1.0
protected

Width of a bin. The maximum distance in a domain, where the radial distribution function is computed, is given as a product of n_bins (if bin_from_r is true) and bin_width.

Definition at line 52 of file pro_radial_distribution_function.h.

◆ counts_traj

Vec<ArrayXi> ProRadialDistributionFunction::counts_traj
protected

List consisting of temporary numbers of particles in a spherical shell with radius of r, where r is distance from a reference particle. Indices in this list corresponds those in generators.

Definition at line 94 of file pro_radial_distribution_function.h.

◆ n_bins

int ProRadialDistributionFunction::n_bins = 1
protected

Number of bins in the distance axis. The maximum distance in a domain, where the radial distribution function is computed, is given as a product of n_bins (if bin_from_r is true) and bin_width.

Definition at line 45 of file pro_radial_distribution_function.h.

◆ number_traj

Vec<int> ProRadialDistributionFunction::number_traj
protected

List consisting of temporary numbers of particles. Indices in this list corresponds those in generators.

Definition at line 82 of file pro_radial_distribution_function.h.

◆ rdf

ArrayXd ProRadialDistributionFunction::rdf
protected

A one-dimensional array of the radial distribution function.

Definition at line 72 of file pro_radial_distribution_function.h.

◆ rdf_traj

Vec<ArrayXd> ProRadialDistributionFunction::rdf_traj
protected

List consisting of temporary radial distribution functions. Indices in this list corresponds those in generators.

Definition at line 77 of file pro_radial_distribution_function.h.

◆ volume_traj

Vec<double> ProRadialDistributionFunction::volume_traj
protected

List consisting of temporary volumes of the simulation box. Indices in this list corresponds those in generators.

Definition at line 87 of file pro_radial_distribution_function.h.


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