Scine::Kinetx
3.0.0
Kinetic models for reaction networks.
|
#include <Kinetx/Integrator/CashKarp5.h>
#include <Kinetx/Integrator/Cvode.h>
#include <Kinetx/Integrator/ExplicitEuler.h>
#include <Kinetx/Integrator/ImplicitEuler.h>
#include <Kinetx/Integrator/Integrator.h>
#include <pybind11/eigen.h>
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
#include <iomanip>
#include <memory>
Enumerations | |
enum | IntegratorName { CashKarp5, ExplicitEuler, ImplicitEuler, CvodeBdf } |
Enum class for the different numerical integration schemes. | |
Functions | |
std::tuple< Eigen::MatrixXd, Eigen::MatrixXd, Eigen::MatrixXd, Eigen::MatrixXd > | integrateDifferentialEquations (Network &network, std::vector< double > yStart, const double tStart, const double dt, const IntegratorName integratorName, const unsigned int batchInterval=1000, const unsigned int nBatches=100000, const double convergenceConcentrationChange=1e-10, bool integrateByTime=false, double maxTime=1e+5) |
A helper function to integrate the rate equations of a given reaction network in a batch-wise fashion. This function keeps track of the maximum concentration for each species in between integration batches. More... | |
IntegratorName | resolveIntegratorName (std::string integrator) |
Resolve a string to its enum representation for the numerical integration schemes. More... | |
void | init_numerical_integration (pybind11::module &m) |
Python bindings for the numerical integration. More... | |
void init_numerical_integration | ( | pybind11::module & | m | ) |
Python bindings for the numerical integration.
m |
std::tuple<Eigen::MatrixXd, Eigen::MatrixXd, Eigen::MatrixXd, Eigen::MatrixXd> integrateDifferentialEquations | ( | Network & | network, |
std::vector< double > | yStart, | ||
const double | tStart, | ||
const double | dt, | ||
const IntegratorName | integratorName, | ||
const unsigned int | batchInterval = 1000 , |
||
const unsigned int | nBatches = 100000 , |
||
const double | convergenceConcentrationChange = 1e-10 , |
||
bool | integrateByTime = false , |
||
double | maxTime = 1e+5 |
||
) |
A helper function to integrate the rate equations of a given reaction network in a batch-wise fashion. This function keeps track of the maximum concentration for each species in between integration batches.
network | The reaction network. |
yStart | The starting concentrations. |
tStart | The starting time. |
dt | The time step (note that this may be ignored by some numerical integration schemes) |
integratorName | The enum for the integration scheme. |
batchInterval | The number of propagation steps per batch. |
nBatches | The maximum number of batches. |
convergenceConcentrationChange | If the maximum change is concentration is lower than the given threshold. The integration is considered to be complete. |
integrateByTime | If true, the numerical integration is done up to the specified maxTime. |
maxTime | The maximum time to integrate. This is only used if integrateByTime is true. |
IntegratorName resolveIntegratorName | ( | std::string | integrator | ) |
Resolve a string to its enum representation for the numerical integration schemes.
integrator | The string. |