Class CVodeSolver

Inheritance Relationships

Base Type

Class Documentation

class amici::CVodeSolver : public amici::Solver

The CVodeSolver class is a wrapper around the SUNDIALS CVODES solver.

Public Functions

~CVodeSolver() override = default
Solver *clone() const override

Clone this instance.

Return

The clone

void reInit(realtype t0, const AmiVector &yy0, const AmiVector &yp0) const override

Reinitializes the states in the solver after an event occurrence.

Parameters
  • t0: reinitialization timepoint

  • yy0: initial state variables

  • yp0: initial derivative state variables (DAE only)

void sensReInit(const AmiVectorArray &yyS0, const AmiVectorArray &ypS0) const override

Reinitializes the state sensitivities in the solver after an event occurrence.

Parameters
  • yyS0: new state sensitivity

  • ypS0: new derivative state sensitivities (DAE only)

void sensToggleOff() const override

Switches off computation of state sensitivities without deallocating the memory for sensitivities.

void reInitB(int which, realtype tB0, const AmiVector &yyB0, const AmiVector &ypB0) const override

Reinitializes the adjoint states after an event occurrence.

Parameters
  • which: identifier of the backwards problem

  • tB0: reinitialization timepoint

  • yyB0: new adjoint state

  • ypB0: new adjoint derivative state

void quadReInitB(int which, const AmiVector &yQB0) const override

Reinitialize the adjoint states after an event occurrence.

Parameters
  • which: identifier of the backwards problem

  • yQB0: new adjoint quadrature state

int solve(realtype tout, int itask) const override

Solves the forward problem until a predefined timepoint.

Return

status flag indicating success of execution

Parameters
  • tout: timepoint until which simulation should be performed

  • itask: task identifier, can be CV_NORMAL or CV_ONE_STEP

int solveF(realtype tout, int itask, int *ncheckPtr) const override

Solves the forward problem until a predefined timepoint (adjoint only)

Return

status flag indicating success of execution

Parameters
  • tout: timepoint until which simulation should be performed

  • itask: task identifier, can be CV_NORMAL or CV_ONE_STEP

  • ncheckPtr: pointer to a number that counts the internal checkpoints

void solveB(realtype tBout, int itaskB) const override

Solves the backward problem until a predefined timepoint (adjoint only)

Parameters
  • tBout: timepoint until which simulation should be performed

  • itaskB: task identifier, can be CV_NORMAL or CV_ONE_STEP

void getDky(realtype t, int k) const override

interpolates the (derivative of the) solution at the requested timepoint

Parameters
  • t: timepoint

  • k: derivative order

void getSensDky(realtype t, int k) const override

interpolates the (derivative of the) solution at the requested timepoint

Parameters
  • t: timepoint

  • k: derivative order

void getQuadDkyB(realtype t, int k, int which) const override

interpolates the (derivative of the) solution at the requested timepoint

Parameters
  • t: timepoint

  • k: derivative order

  • which: index of backward problem

void getDkyB(realtype t, int k, int which) const override

interpolates the (derivative of the) solution at the requested timepoint

Parameters
  • t: timepoint

  • k: derivative order

  • which: index of backward problem

void getRootInfo(int *rootsfound) const override

getRootInfo extracts information which event occurred

Parameters
  • rootsfound: array with flags indicating whether the respective event occurred

void setStopTime(realtype tstop) const override

Sets a timepoint at which the simulation will be stopped.

Parameters
  • tstop: timepoint until which simulation should be performed

void turnOffRootFinding() const override

Disable rootfinding.

const Model *getModel() const override

Accessor function to the model stored in the user data

Return

user data model

void setLinearSolver() const override

Sets the linear solver for the forward problem.

void setLinearSolverB(int which) const override

Sets the linear solver for the backward problem.

Parameters
  • which: index of the backward problem

void setNonLinearSolver() const override

Set the non-linear solver for the forward problem.

void setNonLinearSolverSens() const override

Set the non-linear solver for sensitivities.

void setNonLinearSolverB(int which) const override

Set the non-linear solver for the backward problem.

Parameters
  • which: index of the backward problem

Protected Functions

void calcIC(realtype tout1) const override

Calculates consistent initial conditions, assumes initial states to be correct (DAE only)

Parameters
  • tout1: next timepoint to be computed (sets timescale)

void calcICB(int which, realtype tout1) const override

Calculates consistent initial conditions for the backwards problem, assumes initial states to be correct (DAE only)

Parameters
  • which: identifier of the backwards problem

  • tout1: next timepoint to be computed (sets timescale)

void getB(int which) const override

extracts the adjoint state at the current timepoint from solver memory and writes it to the xB member variable

Parameters
  • which: index of the backwards problem

void getSens() const override

extracts the state sensitivity at the current timepoint from solver memory and writes it to the sx member variable

void getQuadB(int which) const override

extracts the adjoint quadrature state at the current timepoint from solver memory and writes it to the xQB member variable

Parameters
  • which: index of the backwards problem

void getQuad(realtype &t) const override

extracts the quadrature at the current timepoint from solver memory and writes it to the xQ member variable

Parameters
  • t: timepoint for quadrature extraction

void getQuadDky(realtype t, int k) const override

interpolates the (derivative of the) solution at the requested timepoint

Parameters
  • t: timepoint

  • k: derivative order

void reInitPostProcessF(realtype tnext) const override

reInitPostProcessF postprocessing of the solver memory after a discontinuity in the forward problem

Parameters
  • tnext: next timepoint (defines integration direction)

void reInitPostProcessB(realtype tnext) const override

reInitPostProcessB postprocessing of the solver memory after a discontinuity in the backward problem

Parameters
  • tnext: next timepoint (defines integration direction)

void reInitPostProcess(void *cv_mem, realtype *t, AmiVector *yout, realtype tout) const

Postprocessing of the solver memory after a discontinuity.

Parameters
  • cv_mem: pointer to CVODES solver memory object

  • t: pointer to integration time

  • yout: new state vector

  • tout: anticipated next integration timepoint.

void allocateSolver() const override

Create specifies solver method and initializes solver memory for the forward problem.

void setSStolerances(double rtol, double atol) const override

sets scalar relative and absolute tolerances for the forward problem

Parameters
  • rtol: relative tolerances

  • atol: absolute tolerances

void setSensSStolerances(double rtol, const double *atol) const override

activates sets scalar relative and absolute tolerances for the sensitivity variables

Parameters
  • rtol: relative tolerances

  • atol: array of absolute tolerances for every sensitivity variable

void setSensErrCon(bool error_corr) const override

SetSensErrCon specifies whether error control is also enforced for sensitivities for the forward problem

Parameters
  • error_corr: activation flag

void setQuadErrConB(int which, bool flag) const override

Specifies whether error control is also enforced for the backward quadrature problem.

Parameters
  • which: identifier of the backwards problem

  • flag: activation flag

void setQuadErrCon(bool flag) const override

Specifies whether error control is also enforced for the forward quadrature problem.

Parameters
  • flag: activation flag

void setErrHandlerFn() const override

Attaches the error handler function (errMsgIdAndTxt) to the solver.

void setUserData(Model *model) const override

Attaches the user data instance (here this is a Model) to the forward problem.

Parameters

void setUserDataB(int which, Model *model) const override

attaches the user data instance (here this is a Model) to the backward problem

Parameters
  • which: identifier of the backwards problem

  • model: Model instance

void setMaxNumSteps(long int mxsteps) const override

specifies the maximum number of steps for the forward problem

Note

in contrast to the SUNDIALS method, this sets the overall maximum, not the maximum between output times.

Parameters
  • mxsteps: number of steps

void setStabLimDet(int stldet) const override

activates stability limit detection for the forward problem

Parameters
  • stldet: flag for stability limit detection (TRUE or FALSE)

void setStabLimDetB(int which, int stldet) const override

activates stability limit detection for the backward problem

Parameters
  • which: identifier of the backwards problem

  • stldet: flag for stability limit detection (TRUE or FALSE)

void setId(const Model *model) const override

specify algebraic/differential components (DAE only)

Parameters
  • model: model specification

void setSuppressAlg(bool flag) const override

deactivates error control for algebraic components (DAE only)

Parameters
  • flag: deactivation flag

void resetState(void *cv_mem, const_N_Vector y0) const

resetState reset the CVODES solver to restart integration after a rhs discontinuity.

Parameters
  • cv_mem: pointer to CVODES solver memory object

  • y0: new state vector

void setSensParams(const realtype *p, const realtype *pbar, const int *plist) const override

specifies the scaling and indexes for sensitivity computation

Parameters
  • p: parameters

  • pbar: parameter scaling constants

  • plist: parameter index list

void adjInit() const override

initializes the adjoint problem

void quadInit(const AmiVector &xQ0) const override

initializes the quadratures

Parameters
  • xQ0: vector with initial values for xQ

void allocateSolverB(int *which) const override

Specifies solver method and initializes solver memory for the backward problem.

Parameters
  • which: identifier of the backwards problem

void setSStolerancesB(int which, realtype relTolB, realtype absTolB) const override

sets relative and absolute tolerances for the backward problem

Parameters
  • which: identifier of the backwards problem

  • relTolB: relative tolerances

  • absTolB: absolute tolerances

void quadSStolerancesB(int which, realtype reltolQB, realtype abstolQB) const override

sets relative and absolute tolerances for the quadrature backward problem

Parameters
  • which: identifier of the backwards problem

  • reltolQB: relative tolerances

  • abstolQB: absolute tolerances

void quadSStolerances(realtype reltolQ, realtype abstolQ) const override

sets relative and absolute tolerances for the quadrature problem

Parameters
  • reltolQB: relative tolerances

  • abstolQB: absolute tolerances

void setMaxNumStepsB(int which, long int mxstepsB) const override

specifies the maximum number of steps for the forward problem

Note

in contrast to the SUNDIALS method, this sets the overall maximum, not the maximum between output times.

Parameters
  • which: identifier of the backwards problem

  • mxstepsB: number of steps

void diag() const override

attaches a diagonal linear solver to the forward problem

void diagB(int which) const override

attaches a diagonal linear solver to the backward problem

Parameters
  • which: identifier of the backwards problem

void getNumSteps(const void *ami_mem, long int *numsteps) const override

reports the number of solver steps

Parameters
  • ami_mem: pointer to the solver memory instance (can be from forward or backward problem)

  • numsteps: output array

void getNumRhsEvals(const void *ami_mem, long int *numrhsevals) const override

reports the number of right hand evaluations

Parameters
  • ami_mem: pointer to the solver memory instance (can be from forward or backward problem)

  • numrhsevals: output array

void getNumErrTestFails(const void *ami_mem, long int *numerrtestfails) const override

reports the number of local error test failures

Parameters
  • ami_mem: pointer to the solver memory instance (can be from forward or backward problem)

  • numerrtestfails: output array

void getNumNonlinSolvConvFails(const void *ami_mem, long int *numnonlinsolvconvfails) const override

reports the number of nonlinear convergence failures

Parameters
  • ami_mem: pointer to the solver memory instance (can be from forward or backward problem)

  • numnonlinsolvconvfails: output array

void getLastOrder(const void *ami_ami_mem, int *order) const override

Reports the order of the integration method during the last internal step.

Parameters
  • ami_mem: pointer to the solver memory instance (can be from forward or backward problem)

  • order: output array

void *getAdjBmem(void *ami_mem, int which) const override

Retrieves the solver memory instance for the backward problem.

Return

A (void *) pointer to the CVODES memory allocated for the backward problem.

Parameters
  • which: identifier of the backwards problem

  • ami_mem: pointer to the forward solver memory instance

void init(realtype t0, const AmiVector &x0, const AmiVector &dx0) const override

Initializes the states at the specified initial timepoint.

Parameters
  • t0: initial timepoint

  • x0: initial states

  • dx0: initial derivative states

void initSteadystate(const realtype t0, const AmiVector &x0, const AmiVector &dx0) const override

Initializes the states at the specified initial timepoint.

Parameters
  • t0: initial timepoint

  • x0: initial states

  • dx0: initial derivative states

void sensInit1(const AmiVectorArray &sx0, const AmiVectorArray &sdx0) const override

Initializes the forward sensitivities.

Parameters
  • sx0: initial states sensitivities

  • sdx0: initial derivative states sensitivities

void binit(int which, realtype tf, const AmiVector &xB0, const AmiVector &dxB0) const override

Initialize the adjoint states at the specified final timepoint.

Parameters
  • which: identifier of the backwards problem

  • tf: final timepoint

  • xB0: initial adjoint state

  • dxB0: initial adjoint derivative state

void qbinit(int which, const AmiVector &xQB0) const override

Initialize the quadrature states at the specified final timepoint.

Parameters
  • which: identifier of the backwards problem

  • xQB0: initial adjoint quadrature state

void rootInit(int ne) const override

Initializes the rootfinding for events.

Parameters
  • ne: number of different events

void setDenseJacFn() const override

Set the dense Jacobian function.

void setSparseJacFn() const override

sets the sparse Jacobian function

void setBandJacFn() const override

sets the banded Jacobian function

void setJacTimesVecFn() const override

sets the Jacobian vector multiplication function

void setDenseJacFnB(int which) const override

sets the dense Jacobian function

Parameters
  • which: identifier of the backwards problem

void setSparseJacFnB(int which) const override

sets the sparse Jacobian function

Parameters
  • which: identifier of the backwards problem

void setBandJacFnB(int which) const override

sets the banded Jacobian function

Parameters
  • which: identifier of the backwards problem

void setJacTimesVecFnB(int which) const override

sets the Jacobian vector multiplication function

Parameters
  • which: identifier of the backwards problem

void setSparseJacFn_ss() const override

sets the sparse Jacobian function for backward steady state case

Friends

template<class Archive>
friend void serialize(Archive &ar, CVodeSolver &s, unsigned int)

Serialize amici::CVodeSolver to boost archive.

Parameters
  • ar: Archive

  • s: Solver instance to serialize

friend bool operator==(const CVodeSolver &a, const CVodeSolver &b)

Equality operator.

Return

Whether a and b are equal

Parameters
  • a:

  • b: