Langevin C++ source
Loading...
Searching...
No Matches
sim_dplangevin.hpp
Go to the documentation of this file.
1
5
6#ifndef SIMDP_HPP
7#define SIMDP_HPP
8
9#include "dplangevin.hpp"
10
18class SimDP
19{
20private:
31
55 bool did_integrate = false;
57 bool is_initialized = false;
59 bool do_snapshot_grid = true;
61 bool do_verbose = false;
62
64 int count_epochs() const;
66 bool choose_integrator();
68 bool integrate(const int n_next_epochs);
69
71 bool pyprep_t_epochs();
76
77public:
79 SimDP(
80 const double linear, const double quadratic,
81 const double diffusion, const double noise,
82 const double t_final,
83 const double dx, const double dt,
84 const int random_seed,
85 const GridDimension grid_dimension,
86 const int_vec_t grid_size,
87 const gt_vec_t grid_topologies,
88 const bc_vec_t boundary_conditions,
89 const dbl_vec_t bc_values,
90 const InitialCondition initial_condition,
91 const dbl_vec_t ic_values,
92 const IntegrationMethod integration_method,
93 const bool do_snapshot_grid,
94 const bool do_verbose
95 );
97 bool initialize(int n_decimals);
99 bool run(const int n_next_epochs);
101 bool postprocess();
102
103 // Utilities provided to Python via the wrapper
104
106 int get_n_epochs() const;
108 int get_i_current_epoch() const;
110 int get_i_next_epoch() const;
112 double get_t_current_epoch() const;
114 double get_t_next_epoch() const;
116 py_array_t get_t_epochs() const;
120 py_array_t get_density() const;
121};
122
123
124#endif
DPLangevin model application of BaseLangevin class integrator.
bool initialize(int n_decimals)
Initialize the model simulation.
int get_i_current_epoch() const
Fetch the index of the current epoch of the simulation.
int i_current_epoch
Index of current epoch aka time step.
bool run(const int n_next_epochs)
Execute the model simulation for n_next_epochs
bool postprocess()
Process the model results data if available.
bool do_verbose
Flag whether to report sim parameters etc.
int i_next_epoch
Index of next epoch aka time step.
py_array_t get_density() const
Fetch the current Langevin density field grid as a Python array.
bool do_snapshot_grid
Flag whether to take density grid snapshots.
Parameters p
Model simulation parameters.
py_array_t pyarray_t_epochs
Python-compatible array of epochs time-series.
dbl_vec_t mean_densities
Vector time-series of grid-averaged field density values.
SimDP(const double linear, const double quadratic, const double diffusion, const double noise, const double t_final, const double dx, const double dt, const int random_seed, const GridDimension grid_dimension, const int_vec_t grid_size, const gt_vec_t grid_topologies, const bc_vec_t boundary_conditions, const dbl_vec_t bc_values, const InitialCondition initial_condition, const dbl_vec_t ic_values, const IntegrationMethod integration_method, const bool do_snapshot_grid, const bool do_verbose)
Constructor.
int count_epochs() const
Count upcoming number of epochs by running a dummy time-stepping loop.
Coefficients coefficients
Langevin equation coefficients.
bool choose_integrator()
Chooses function implementing either Runge-Kutta or Euler integration methods.
double t_next_epoch
Time of next epoch.
py_array_t get_t_epochs() const
Fetch a times-series vector of the simulation epochs as a Python array.
bool is_initialized
Flag whether simulation has been initialized or not.
void(DPLangevin::* integrator)(rng_t &)
Integrator: either a Runge-Kutta or an Euler method.
double get_t_next_epoch() const
Fetch the next epoch (time) of the simulation.
bool pyprep_mean_densities()
Generate a Python-compatible version of the mean densities time-series vector.
int n_epochs
Total number of simulation time steps aka "epochs".
dbl_vec_t t_epochs
Vector time-series of epochs.
py_array_t pyarray_density
Python-compatible array of current density grid.
py_array_t pyarray_mean_densities
Python-compatible array of mean density time-series.
bool pyprep_t_epochs()
Generate a Python-compatible version of the epochs time-series vector.
int n_decimals
Truncation number of decimal places when summing Δt.
double get_t_current_epoch() const
Fetch the current epoch (time) of the simulation.
DPLangevin * dpLangevin
Instance of DP Langevin integrator class (pointer to instance)
py_array_t get_mean_densities() const
Fetch a times-series vector of the grid-averaged density field over time as a Python array.
bool did_integrate
Flag whether integration step was successful or not.
int get_n_epochs() const
Fetch the total number of simulation epochs.
int get_i_next_epoch() const
Fetch the index of the next epoch of the simulation.
bool pyprep_density_grid()
Generate a Python-compatible version of the current density grid.
rng_t * rng
Random number generation function (Mersenne prime) (pointer to RNG)
bool integrate(const int n_next_epochs)
Perform Dornic-type integration of the DP Langevin equation for n_next_epochs
double t_current_epoch
Time of current epoch.
DPLangevin model application of BaseLangevin class integrator.
GridDimension
Density field grid dimension: only 1D or 2D grids implemented so far.
IntegrationMethod
Deterministic integration method: default is 4th-order Runge-Kutta; can be explicit Euler.
InitialCondition
Grid density field initial condition: random uniform or Gaussian variates; constant everywhere; or ce...
py::array_t< double, py::array::c_style > py_array_t
Type for Python arrays of doubles.
std::vector< GridTopology > gt_vec_t
Type for specifying grid topology in each direction x, y, z...
std::vector< int > int_vec_t
Type for vectors of integers.
std::vector< double > dbl_vec_t
Type for vectors of doubles.
std::mt19937 rng_t
Use Mersenne Twister random number generator.
Container for nonlinear Langevin equation coefficients.
Container for BaseLangevin integrator parameters.