Scine::Kinetx  3.0.0
Kinetic models for reaction networks.
 All Classes Files Functions Variables Enumerations Pages
Scine::Kinetx::IntegratorBase Class Referenceabstract
Inheritance diagram for Scine::Kinetx::IntegratorBase:
Inheritance graph
Collaboration diagram for Scine::Kinetx::IntegratorBase:
Collaboration graph

Public Member Functions

 IntegratorBase (Network &net)
 Constructor. More...
 
virtual void propagate (Eigen::VectorXd &concentrations, Eigen::VectorXd &yFlux, Eigen::VectorXd &rFlux, Eigen::VectorXd &rForwardFlux, Eigen::VectorXd &rBackwardFlux, double &t, double &dt) const =0
 Propagate the numerical integration by one time step. More...
 
virtual Eigen::MatrixXd runIntegration (Eigen::VectorXd y, double tStart, double dt, Eigen::VectorXd &rFlux, Eigen::VectorXd &rForwardFlux, Eigen::VectorXd &rBackwardFlux, const unsigned int batchInterval=1000, const unsigned int nBatches=100000, const double convergenceConcentrationChange=1e-10)
 Run the numerical integration until a maximum number of steps or convergence is reached. More...
 
virtual Eigen::MatrixXd runIntegrationByTime (Eigen::VectorXd y, double tStart, double dt, Eigen::VectorXd &rFlux, Eigen::VectorXd &rForwardFlux, Eigen::VectorXd &rBackwardFlux, const double tMax, const unsigned int batchInterval=1000, const double convergenceConcentrationChange=1e-10)
 Run the numerical integration until a maximum time or convergence is reached. More...
 

Protected Member Functions

void trackVertexAndEdgeFluxes (const Eigen::VectorXd &y, const Eigen::VectorXd &yInitial, Eigen::VectorXd &yFlux, Eigen::VectorXd &rFlux, Eigen::VectorXd &rForwardFlux, Eigen::VectorXd &rBackwardFlux, const double &dt) const
 Keep track on the vertex and edge fluxes. More...
 
Eigen::VectorXd g (const Eigen::VectorXd &concentrations) const
 Calculate the concentration gradient. More...
 
Eigen::VectorXd gFlux (const Eigen::VectorXd &concentrationsT0, const Eigen::VectorXd &concentrationsT1, Eigen::VectorXd &flux, Eigen::VectorXd &forwardFlux, Eigen::VectorXd &backwardFlux) const
 Calculate the concentration flux gradient. More...
 
Eigen::MatrixXd fFlux (const Eigen::VectorXd &concentrationsT0, const Eigen::VectorXd &concentrationsT1) const
 Calculate the total edge flux gradient. More...
 
std::tuple< Eigen::MatrixXd,
Eigen::MatrixXd,
Eigen::MatrixXd > 
fFluxDirected (const Eigen::VectorXd &concentrationsT0, const Eigen::VectorXd &concentrationsT1) const
 Calculate the edge flux while keeping track on the direction. More...
 
Eigen::MatrixXd f (const Eigen::VectorXd &concentrations) const
 
std::pair< Eigen::MatrixXd,
Eigen::MatrixXd > 
fDirected (const Eigen::VectorXd &concentrations) const
 The same as f(..) but keeping track on the directions. More...
 
Eigen::SparseMatrix< double > jacobi (const Eigen::VectorXd &concentrations) const
 
void printHeader (const Eigen::VectorXd &y, const double dt, const double t)
 
bool printTimeAndCheckConvergenceStep (const Eigen::VectorXd &y, const Eigen::VectorXd &yOld, Eigen::VectorXd &yMax, const unsigned int batchInterval, const unsigned int iBatch, const double dt, const double t, const double convergenceConcentrationChange)
 

Static Protected Member Functions

static double fast_pow (const double &x, const int &n)
 

Protected Attributes

Network_net
 

Constructor & Destructor Documentation

Scine::Kinetx::IntegratorBase::IntegratorBase ( Network net)

Constructor.

Parameters
netThe network of reactions.

Member Function Documentation

Eigen::MatrixXd Scine::Kinetx::IntegratorBase::f ( const Eigen::VectorXd &  concentrations) const
protected

Calculate the reaction rates f for each reaction using a mass-weight ansatz. This returns a matrix in case several elementary steps are used for at least one reaction. Otherwise, this matrix will have only one column. Hence, it will generally be densely populated.

Parameters
concentrationsThe current concentrations.
Returns
The reaction rates as a vector/matrix.
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > Scine::Kinetx::IntegratorBase::fDirected ( const Eigen::VectorXd &  concentrations) const
protected

The same as f(..) but keeping track on the directions.

Parameters
concentrationsThe concentration.
Returns
The forward and backward reaction rates.
Eigen::MatrixXd Scine::Kinetx::IntegratorBase::fFlux ( const Eigen::VectorXd &  concentrationsT0,
const Eigen::VectorXd &  concentrationsT1 
) const
protected

Calculate the total edge flux gradient.

Parameters
concentrationsT0The original concentration.
concentrationsT1The updated concentration.
Returns
The total edge flux.
std::tuple< Eigen::MatrixXd, Eigen::MatrixXd, Eigen::MatrixXd > Scine::Kinetx::IntegratorBase::fFluxDirected ( const Eigen::VectorXd &  concentrationsT0,
const Eigen::VectorXd &  concentrationsT1 
) const
protected

Calculate the edge flux while keeping track on the direction.

Parameters
concentrationsT0The original concentration.
concentrationsT1The updated concentration.
Returns
The total, forward, and backward flux (in this order) as a tuple.
Eigen::VectorXd Scine::Kinetx::IntegratorBase::g ( const Eigen::VectorXd &  concentrations) const
protected

Calculate the concentration gradient.

Parameters
concentrationsThe current concentration.
Returns
The concentration gradient.
Eigen::VectorXd Scine::Kinetx::IntegratorBase::gFlux ( const Eigen::VectorXd &  concentrationsT0,
const Eigen::VectorXd &  concentrationsT1,
Eigen::VectorXd &  flux,
Eigen::VectorXd &  forwardFlux,
Eigen::VectorXd &  backwardFlux 
) const
protected

Calculate the concentration flux gradient.

Parameters
concentrationsT0The original concentraion.
concentrationsT1The updated concentation.
fluxThe total reaction edge flux.
forwardFluxThe reaction edge forward flux.
backwardFluxThe reaction edge backward flux.
Returns
The vertex flux gradient.
virtual void Scine::Kinetx::IntegratorBase::propagate ( Eigen::VectorXd &  concentrations,
Eigen::VectorXd &  yFlux,
Eigen::VectorXd &  rFlux,
Eigen::VectorXd &  rForwardFlux,
Eigen::VectorXd &  rBackwardFlux,
double &  t,
double &  dt 
) const
pure virtual

Propagate the numerical integration by one time step.

Parameters
concentrationsThe current concentrations. Updated inplace.
yFluxThe current vertex (aggregate) flux. Updated inplace.
rFluxThe current total edge (reaction) flux. Updated inplace.
rForwardFluxThe current forward edge flux. Updated inplace.
rBackwardFluxThe current backward edge flux. Updated inplace.
tThe current time. Updated inplace.
dtThe time increment. This may be updated inplace depending on the integration algorithm.

Implements Scine::Kinetx::Integrator.

Implemented in Scine::Kinetx::Cvode::Impl, and Scine::Kinetx::RungeKutta.

Eigen::MatrixXd Scine::Kinetx::IntegratorBase::runIntegration ( Eigen::VectorXd  y,
double  tStart,
double  dt,
Eigen::VectorXd &  rFlux,
Eigen::VectorXd &  rForwardFlux,
Eigen::VectorXd &  rBackwardFlux,
const unsigned int  batchInterval = 1000,
const unsigned int  nBatches = 100000,
const double  convergenceConcentrationChange = 1e-10 
)
virtual

Run the numerical integration until a maximum number of steps or convergence is reached.

Parameters
yThe input concentration.
tStartThe start time.
dtThe time increment.
rFluxThe reaction edge flux (total).
rForwardFluxThe forward reaction edge flux.
rBackwardFluxThe backward reaction edge flux.
batchIntervalThe number of steps per batch.
nBatchesThe number of integration batches.
convergenceConcentrationChangeThe concentration convergence threshold.
Returns
The final concentrations, max. concentrations, and vertex fluxes as a matrix with columns in this order.

Implements Scine::Kinetx::Integrator.

Eigen::MatrixXd Scine::Kinetx::IntegratorBase::runIntegrationByTime ( Eigen::VectorXd  y,
double  tStart,
double  dt,
Eigen::VectorXd &  rFlux,
Eigen::VectorXd &  rForwardFlux,
Eigen::VectorXd &  rBackwardFlux,
const double  tMax,
const unsigned int  batchInterval = 1000,
const double  convergenceConcentrationChange = 1e-10 
)
virtual

Run the numerical integration until a maximum time or convergence is reached.

Parameters
yThe input concentration.
tStartThe start time.
dtThe time increment.
rFluxThe reaction edge flux (total).
rForwardFluxThe forward reaction edge flux.
rBackwardFluxThe backward reaction edge flux.
tMaxThe final time to integrate to.
batchIntervalThe number of steps per batch.
convergenceConcentrationChangeThe concentration convergence threshold.
Returns
The final concentrations, max. concentrations, and vertex fluxes as a matrix with columns in this order.

Implements Scine::Kinetx::Integrator.

void Scine::Kinetx::IntegratorBase::trackVertexAndEdgeFluxes ( const Eigen::VectorXd &  y,
const Eigen::VectorXd &  yInitial,
Eigen::VectorXd &  yFlux,
Eigen::VectorXd &  rFlux,
Eigen::VectorXd &  rForwardFlux,
Eigen::VectorXd &  rBackwardFlux,
const double &  dt 
) const
protected

Keep track on the vertex and edge fluxes.

Parameters
yThe concentration after concentration propagation.
yInitialThe original concentration.
yFluxThe vertex flux. Updated inplace.
rFluxThe total edge flux.
rForwardFluxThe forward edge flux.
rBackwardFluxThe backward edge flux.
dtThe time increment.

The documentation for this class was generated from the following files: