Scine::Kinetx  3.0.0
Kinetic models for reaction networks.
 All Classes Files Functions Variables Enumerations Pages
IntegratorsPython.cpp File Reference
#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>
Include dependency graph for IntegratorsPython.cpp:

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...
 

Detailed Description

Function Documentation

void init_numerical_integration ( pybind11::module &  m)

Python bindings for the numerical integration.

Parameters
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.

Parameters
networkThe reaction network.
yStartThe starting concentrations.
tStartThe starting time.
dtThe time step (note that this may be ignored by some numerical integration schemes)
integratorNameThe enum for the integration scheme.
batchIntervalThe number of propagation steps per batch.
nBatchesThe maximum number of batches.
convergenceConcentrationChangeIf the maximum change is concentration is lower than the given threshold. The integration is considered to be complete.
integrateByTimeIf true, the numerical integration is done up to the specified maxTime.
maxTimeThe maximum time to integrate. This is only used if integrateByTime is true.
Returns
A matrix of which the values of the first column are the final concentrations, the values of the second column are the maximum concentrations reached by the species, and the third column are the concentration flux through the species during the integration.
IntegratorName resolveIntegratorName ( std::string  integrator)

Resolve a string to its enum representation for the numerical integration schemes.

Parameters
integratorThe string.
Returns
The enum.