DP Langevin
Loading...
Searching...
No Matches
sim_dplangevin_private.cpp
Go to the documentation of this file.
1
5
6#include "sim_dplangevin.hpp"
7
10{
11 int n_epochs;
12 double t;
13 for (n_epochs=0, t=0; t<=p.t_final+p.dt; t+=p.dt, n_epochs++) {}
14 return n_epochs;
15}
16
18{
19 switch (p.integration_method)
20 {
21 case (IntegrationMethod::RUNGE_KUTTA):
23 return true;
24 case (IntegrationMethod::EULER):
26 return true;
27 default:
28 return false;
29 }
30}
31
32bool SimDP::integrate(const int n_next_epochs)
33{
34 // Check a further n_next_epochs won't exceed total permitted steps
35 if (t_epochs.size() < i_next_epoch+n_next_epochs)
36 {
37 std::cout << "Too many epochs: "
38 << t_epochs.size()
39 << " < "
40 << i_next_epoch+n_next_epochs
41 << std::endl;
42 return false;
43 }
44
45 // Perform (possibly another another) n_next_epochs integration steps
46 int i;
47 double t;
48 // For the very first epoch, record mean density right now
49 if (i_next_epoch==1) {
50 dpLangevin->apply_boundary_conditions(p, 0);
51 mean_densities[0] = dpLangevin->get_mean_density();
54 }
55 // Loop over integration steps.
56 // Effectively increment epoch counter and add to Δt to time counter
57 // so that both point the state *after* each integration step is complete.
58 // In so doing, we will record t_epochs.size() + 1 total integration steps.
59 for (
61 i<i_next_epoch+n_next_epochs;
62 t+=p.dt, i++
63 )
64 {
65 // Reapply boundary conditions prior to integrating
66 dpLangevin->apply_boundary_conditions(p, i);
67 // Perform a single integration over Δt
69 // Record this epoch
70 t_epochs[i] = t;
71 mean_densities[i] = dpLangevin->get_mean_density();
74 };
75 // Set epoch and time counters to point to *after* the last integration step
76 i_next_epoch = i;
77 t_next_epoch = t;
78 return true;
79}
80
82{
83 py_array_t epochs_array(n_epochs);
84 auto epochs_proxy = epochs_array.mutable_unchecked();
85 for (auto i=0; i<n_epochs; i++)
86 {
87 epochs_proxy(i) = t_epochs[i];
88 };
89 pyarray_t_epochs = epochs_array;
90 return true;
91}
92
94{
95 py_array_t mean_densities_array(n_epochs);
96 auto mean_densities_proxy = mean_densities_array.mutable_unchecked();
97 for (auto i=0; i<n_epochs; i++)
98 {
99 mean_densities_proxy(i) = mean_densities[i];
100 };
101 pyarray_mean_densities = mean_densities_array;
102 return true;
103}
104
106{
107 if (not (p.n_cells == p.n_x * p.n_y * p.n_z))
108 {
109 std::cout << "prep_density: grid size problem" << std::endl;
110 return false;
111 }
112 // Assume we're working with a <=2d grid for now
113 py_array_t density_array({p.n_x, p.n_y});
114 auto density_proxy = density_array.mutable_unchecked();
115 for (auto i=0; i<p.n_cells; i++)
116 {
117 int x = i % p.n_x;
118 int y = i / p.n_x;
119 density_proxy(x, y) = dpLangevin->get_density_grid_value(i);
120 };
121 pyarray_density = density_array;
122 return true;
123}
void integrate_rungekutta(rng_t &rng)
Runge-Kutta + stochastic integration + grid update.
void integrate_euler(rng_t &rng)
Explicit Euler + stochastic integration + grid update.
int i_current_epoch
Index of current epoch aka time step.
int i_next_epoch
Index of next epoch aka time step.
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.
int count_epochs() const
Count upcoming number of epochs by running a dummy time-stepping loop.
bool choose_integrator()
Chooses function implementing either Runge-Kutta or Euler integration methods.
double t_next_epoch
Time of next epoch.
void(DPLangevin::* integrator)(rng_t &)
Integrator: either a Runge-Kutta or an Euler method.
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.
DPLangevin * dpLangevin
Instance of DP Langevin integrator class (pointer to instance)
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.
py::array_t< double, py::array::c_style > py_array_t
Type for Python arrays of doubles.
Class that manages simulation of DPLangevin equation.