ppap4lmp  0.7.2
pro_radial_distribution_function.h
Go to the documentation of this file.
1 
9 #ifndef PRO_RADIAL_DISTRIBUTION_FUNCTION_H
10 #define PRO_RADIAL_DISTRIBUTION_FUNCTION_H
11 
12 #include <processors/processor.h>
13 
38  protected:
45  int n_bins = 1;
52  double bin_width = 1.0;
59  bool bin_from_r = false;
68  bool beyond_half = false;
101  virtual void run_impl(
102  const int index) override;
103  public:
153  const ElPtr &targets,
154  const ElPtr &box);
202  const Vec<std::pair<ElPtr,ElPtr>> &pairs);
203  virtual ~ProRadialDistributionFunction() = default;
208  virtual void prepare() override;
219  virtual void finish() override;
239  void set_bin(
240  double bin_width_,
241  int n_bins_);
258  bool bin_from_r_ = true);
274  bool beyond_half_ = true);
281  const ArrayXd get_r_axis();
288  const ArrayXd &get_rdf();
295  const Vec<ArrayXd> &get_rdf_traj();
296 };
297 
298 #endif
const Vec< ArrayXd > & get_rdf_traj()
Get a sequence of radial distribution functions as list of one-dimensional arrays.
ShPtr< Element > ElPtr
An alias for a shared pointer of Element class.
Definition: element.h:378
virtual void run_impl(const int index) override
Implementation of analysis using an element of generators.
ProRadialDistributionFunction(const ElPtr &targets, const ElPtr &box)
Constructor of ProRadialDistributionFunction class for a snapshot of simulation.
virtual void finish() override
Compute rdf and rdf_traj from number_traj, volume_traj and counts_traj.
const ArrayXd & get_rdf()
Get the radial distribution function as a one-dimensional array.
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 para...
std::vector< T > Vec
Vec is an alias for vector (same as list in Python).
Definition: std.h:27
virtual void prepare() override
Resize number_traj, volume_traj and counts_traj.
const ArrayXd get_r_axis()
Get a one-dimensional array of distance, [0, dr, 2*dr, ... ], where dr is width of a bin...
This file has a definition of Processor class, where an analysis process is programmed.
Processor analyzes data contained in one or more Generator objects.
Definition: processor.h:18
ProRadialDistributionFunction computes radial distribution function (RDF).
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...
Eigen::Array< double, Eigen::Dynamic, 1 > ArrayXd
ArrayXd is an alias for a column array of float numbers.
Definition: eigen.h:42
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 p...