Class SteadystateProblem

Class Documentation

class amici::SteadystateProblem

The SteadystateProblem class solves a steady-state problem using Newton’s method and falls back to integration on failure.

Public Functions

explicit SteadystateProblem(const Solver &solver, const Model &model)

constructor

Parameters
void workSteadyStateProblem(Solver *solver, Model *model, int it)

Handles steady state computation in the forward case: tries to determine the steady state of the ODE system and computes steady state sensitivities if requested.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

  • it – integer with the index of the current time step

void workSteadyStateBackwardProblem(Solver *solver, Model *model, const BackwardProblem *bwd)

Integrates over the adjoint state backward in time by solving a linear system of equations, which gives the analytical solution. Computes the gradient via adjoint steady state sensitivities

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

  • bwd – backward problem

void findSteadyState(Solver *solver, NewtonSolver *newtonSolver, Model *model, int it)

Handles the computation of the steady state, throws an AmiException, if no steady state was found.

Parameters
  • solver – pointer to the solver object

  • newtonSolver – pointer to the newtonSolver solver object

  • model – pointer to the model object

  • it – integer with the index of the current time step

void findSteadyStateByNewtonsMethod(NewtonSolver *newtonSolver, Model *model, bool newton_retry)

Tries to determine the steady state by using Newton’s method.

Parameters
  • newtonSolver – pointer to the newtonSolver solver object

  • model – pointer to the model object

  • newton_retry – bool flag indicating whether being relaunched

void findSteadyStateBySimulation(const Solver *solver, Model *model, int it)

Tries to determine the steady state by using forward simulation.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

  • it – integer with the index of the current time step

void computeSteadyStateQuadrature(NewtonSolver *newtonSolver, const Solver *solver, Model *model)

Handles the computation of quadratures in adjoint mode.

Parameters
  • newtonSolver – pointer to the newtonSolver solver object

  • solver – pointer to the solver object

  • model – pointer to the model object

void getQuadratureByLinSolve(NewtonSolver *newtonSolver, Model *model)

Computes the quadrature in steady state backward mode by solving the linear system defined by the backward Jacobian.

Parameters
  • newtonSolver – pointer to the newtonSolver solver object

  • model – pointer to the model object

void getQuadratureBySimulation(const Solver *solver, Model *model)

Computes the quadrature in steady state backward mode by numerical integration of xB forward in time.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

void handleSteadyStateFailure(const Solver *solver, Model *model)

Stores state and throws an exception if equilibration failed.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

void writeErrorString(std::string *errorString, SteadyStateStatus status) const

Assembles the error message to be thrown.

Parameters
  • errorString – const pointer to string with error message

  • status – Entry of steady_state_status to be processed

bool getSensitivityFlag(const Model *model, const Solver *solver, int it, SteadyStateContext context)

Checks depending on the status of the Newton solver, solver settings, and the model, whether state sensitivities still need to be computed via a linear system solve or stored.

Parameters
  • model – pointer to the model object

  • solver – pointer to the solver object

  • it – integer with the index of the current time step

  • context – SteadyStateContext giving the situation for the flag

Returns

flag telling how to process state sensitivities

realtype getWrmsNorm(AmiVector const &x, AmiVector const &xdot, realtype atol, realtype rtol, AmiVector &ewt) const

Computes the weighted root mean square of xdot the weights are computed according to x: w_i = 1 / ( rtol * x_i + atol )

Parameters
  • x – current state (sx[ip] for sensitivities)

  • xdot – current rhs (sxdot[ip] for sensitivities)

  • atol – absolute tolerance

  • rtol – relative tolerance

  • ewt – error weight vector

Returns

root-mean-square norm

bool checkConvergence(const Solver *solver, Model *model, SensitivityMethod checkSensitivities)

Checks convergence for state and respective sensitivities.

Parameters
  • solverSolver instance

  • model – instance

  • checkSensitivities – flag whether sensitivities should be checked

Returns

boolean indicating convergence

void applyNewtonsMethod(Model *model, NewtonSolver *newtonSolver, bool newton_retry)

Runs the Newton solver iterations and checks for convergence to steady state.

Parameters
  • model – pointer to the model object

  • newtonSolver – pointer to the NewtonSolver object

  • newton_retry – flag indicating if Newton solver is rerun

void runSteadystateSimulation(const Solver *solver, Model *model, bool backward)

Simulation is launched, if Newton solver or linear system solve fails.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

  • backward – flag indicating adjoint mode (including quadrature)

std::unique_ptr<Solver> createSteadystateSimSolver(const Solver *solver, Model *model, bool forwardSensis, bool backward) const

Initialize CVodeSolver instance for preequilibration simulation.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

  • forwardSensis – flag switching on integration with FSA

  • backward – flag switching on quadratures computation

Returns

solver instance

bool initializeBackwardProblem(Solver *solver, Model *model, const BackwardProblem *bwd)

Initialize backward computation by setting state, time, adjoint state and checking for preequilibration mode.

Parameters
  • solver – pointer to the solver object

  • model – pointer to the model object

  • bwd – pointer to backward problem

Returns

flag indicating whether backward computation to be carried out

void computeQBfromQ(Model *model, const AmiVector &yQ, AmiVector &yQB) const

Compute the backward quadratures, which contribute to the gradient (xQB) from the quadrature over the backward state itself (xQ)

Parameters
  • model – pointer to the model object

  • yQ – vector to be multiplied with dxdotdp

  • yQB – resulting vector after multiplication

void storeSimulationState(Model *model, bool storesensi)

Store carbon copy of current simulation state variables as SimulationState.

Parameters
  • model – model carrying the ModelState to be used

  • storesensi – flag to enable storage of sensitivities

inline const SimulationState &getFinalSimulationState() const

Returns the stored SimulationState.

Returns

stored SimulationState

inline const AmiVector &getEquilibrationQuadratures() const

Returns the quadratures from pre- or postequilibration.

Returns

xQB Vector with quadratures

inline const AmiVector &getState() const

Returns state at steadystate.

Returns

x

inline const AmiVectorArray &getStateSensitivity() const

Returns state sensitivity at steadystate.

Returns

sx

inline std::vector<realtype> const &getDJydx() const

Accessor for dJydx.

Returns

dJydx

inline double getCPUTime() const

Accessor for run_time of the forward problem.

Returns

run_time

inline double getCPUTimeB() const

Accessor for run_time of the backward problem.

Returns

run_time

inline std::vector<SteadyStateStatus> const &getSteadyStateStatus() const

Accessor for steady_state_status.

Returns

steady_state_status

inline realtype getSteadyStateTime() const

Accessor for t.

Returns

t

inline realtype getResidualNorm() const

Accessor for wrms.

Returns

wrms

inline const std::vector<int> &getNumSteps() const

Accessor for numsteps.

Returns

numsteps

inline int getNumStepsB() const

Accessor for numstepsB.

Returns

numstepsB

inline const std::vector<int> &getNumLinSteps() const

Accessor for numlinsteps.

Returns

numlinsteps

void getAdjointUpdates(Model &model, const ExpData &edata)

computes adjoint updates dJydx according to provided model and expdata

Parameters
  • modelModel instance

  • edata – experimental data

inline AmiVector const &getAdjointState() const

Return the adjoint state.

Returns

xB adjoint state

inline AmiVector const &getAdjointQuadrature() const

Accessor for xQB.

Returns

xQB

inline bool hasQuadrature() const

Accessor for hasQuadrature_.

Returns

hasQuadrature_

bool checkSteadyStateSuccess() const

computes adjoint updates dJydx according to provided model and expdata

Returns

convergence of steady state solver