Scine::Kinetx
3.0.0
Kinetic models for reaction networks.
|
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 |
Scine::Kinetx::IntegratorBase::IntegratorBase | ( | Network & | net | ) |
Constructor.
net | The network of reactions. |
|
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.
concentrations | The current concentrations. |
|
protected |
The same as f(..) but keeping track on the directions.
concentrations | The concentration. |
|
protected |
Calculate the total edge flux gradient.
concentrationsT0 | The original concentration. |
concentrationsT1 | The updated concentration. |
|
protected |
Calculate the edge flux while keeping track on the direction.
concentrationsT0 | The original concentration. |
concentrationsT1 | The updated concentration. |
|
protected |
Calculate the concentration gradient.
concentrations | The current concentration. |
|
protected |
Calculate the concentration flux gradient.
concentrationsT0 | The original concentraion. |
concentrationsT1 | The updated concentation. |
flux | The total reaction edge flux. |
forwardFlux | The reaction edge forward flux. |
backwardFlux | The reaction edge backward flux. |
|
pure virtual |
Propagate the numerical integration by one time step.
concentrations | The current concentrations. Updated inplace. |
yFlux | The current vertex (aggregate) flux. Updated inplace. |
rFlux | The current total edge (reaction) flux. Updated inplace. |
rForwardFlux | The current forward edge flux. Updated inplace. |
rBackwardFlux | The current backward edge flux. Updated inplace. |
t | The current time. Updated inplace. |
dt | The 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.
|
virtual |
Run the numerical integration until a maximum number of steps or convergence is reached.
y | The input concentration. |
tStart | The start time. |
dt | The time increment. |
rFlux | The reaction edge flux (total). |
rForwardFlux | The forward reaction edge flux. |
rBackwardFlux | The backward reaction edge flux. |
batchInterval | The number of steps per batch. |
nBatches | The number of integration batches. |
convergenceConcentrationChange | The concentration convergence threshold. |
Implements Scine::Kinetx::Integrator.
|
virtual |
Run the numerical integration until a maximum time or convergence is reached.
y | The input concentration. |
tStart | The start time. |
dt | The time increment. |
rFlux | The reaction edge flux (total). |
rForwardFlux | The forward reaction edge flux. |
rBackwardFlux | The backward reaction edge flux. |
tMax | The final time to integrate to. |
batchInterval | The number of steps per batch. |
convergenceConcentrationChange | The concentration convergence threshold. |
Implements Scine::Kinetx::Integrator.
|
protected |
Keep track on the vertex and edge fluxes.
y | The concentration after concentration propagation. |
yInitial | The original concentration. |
yFlux | The vertex flux. Updated inplace. |
rFlux | The total edge flux. |
rForwardFlux | The forward edge flux. |
rBackwardFlux | The backward edge flux. |
dt | The time increment. |