Class Model_DAE¶
Defined in File model_dae.h
Inheritance Relationships¶
Base Type¶
public amici::Model(Class Model)
Class Documentation¶
-
class
amici::Model_DAE: public amici::Model¶ The Model class represents an AMICI DAE model.
The model does not contain any data, but represents the state of the model at a specific time t. The states must not always be in sync, but may be updated asynchronously.
Public Functions
-
Model_DAE() = default¶ default constructor
-
Model_DAE(const ModelDimensions &model_dimensions, SimulationParameters simulation_parameters, const SecondOrderMode o2mode, std::vector<realtype> const &idlist, std::vector<int> const &z2event, const bool pythonGenerated = false, const int ndxdotdp_explicit = 0)¶ Constructor with model dimensions.
- Parameters
model_dimensions: Model dimensionssimulation_parameters: Simulation parameterso2mode: second order sensitivity modeidlist: indexes indicating algebraic components (DAE only)z2event: mapping of event outputs to eventspythonGenerated: flag indicating matlab or python wrappingndxdotdp_explicit: number of nonzero elements dxdotdp_explicit
-
void
fJ(realtype t, realtype cj, const AmiVector &x, const AmiVector &dx, const AmiVector &xdot, SUNMatrix J) override¶ Dense Jacobian function.
- Parameters
t: timecj: scaling factor (inverse of timestep, DAE only)x: statedx: time derivative of state (DAE only)xdot: values of residual function (unused)J: dense matrix to which values of the jacobian will be written
-
void
fJ(realtype t, realtype cj, const_N_Vector x, const_N_Vector dx, const_N_Vector xdot, SUNMatrix J)¶ Jacobian of xdot with respect to states x.
- Parameters
t: timepointcj: scaling factor, inverse of the step sizex: Vector with the statesdx: Vector with the derivative statesxdot: Vector with the right hand sideJ: Matrix to which the Jacobian will be written
-
void
fJB(const realtype t, realtype cj, const AmiVector &x, const AmiVector &dx, const AmiVector &xB, const AmiVector &dxB, const AmiVector &xBdot, SUNMatrix JB) override¶ Dense Jacobian function.
- Parameters
t: timecj: scaling factor (inverse of timestep, DAE only)x: statedx: time derivative of state (DAE only)xB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesxBdot: Vector with the adjoint right hand side (unused)JB: dense matrix to which values of the jacobian will be written
-
void
fJB(realtype t, realtype cj, const_N_Vector x, const_N_Vector dx, const_N_Vector xB, const_N_Vector dxB, SUNMatrix JB)¶ Jacobian of xBdot with respect to adjoint state xB.
- Parameters
t: timepointcj: scaling factor, inverse of the step sizex: Vector with the statesdx: Vector with the derivative statesxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesJB: Matrix to which the Jacobian will be written
-
void
fJSparse(realtype t, realtype cj, const AmiVector &x, const AmiVector &dx, const AmiVector &xdot, SUNMatrix J) override¶ Sparse Jacobian function.
- Parameters
t: timecj: scaling factor (inverse of timestep, DAE only)x: statedx: time derivative of state (DAE only)xdot: values of residual function (unused)J: sparse matrix to which values of the Jacobian will be written
-
void
fJSparse(realtype t, realtype cj, const_N_Vector x, const_N_Vector dx, SUNMatrix J)¶ J in sparse form (for sparse solvers from the SuiteSparse Package)
- Parameters
t: timepointcj: scalar in Jacobian (inverse stepsize)x: Vector with the statesdx: Vector with the derivative statesJ: Matrix to which the Jacobian will be written
-
void
fJSparseB(const realtype t, realtype cj, const AmiVector &x, const AmiVector &dx, const AmiVector &xB, const AmiVector &dxB, const AmiVector &xBdot, SUNMatrix JB) override¶ Sparse Jacobian function.
- Parameters
t: timecj: scaling factor (inverse of timestep, DAE only)x: statedx: time derivative of state (DAE only)xB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesxBdot: Vector with the adjoint right hand side (unused)JB: dense matrix to which values of the jacobian will be written
-
void
fJSparseB(realtype t, realtype cj, const_N_Vector x, const_N_Vector dx, const_N_Vector xB, const_N_Vector dxB, SUNMatrix JB)¶ JB in sparse form (for sparse solvers from the SuiteSparse Package)
- Parameters
t: timepointcj: scalar in Jacobianx: Vector with the statesdx: Vector with the derivative statesxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesJB: Matrix to which the Jacobian will be written
-
void
fJDiag(realtype t, AmiVector &JDiag, realtype cj, const AmiVector &x, const AmiVector &dx) override¶ Diagonal of the Jacobian (for preconditioning)
- Parameters
t: timepointJDiag: Vector to which the Jacobian diagonal will be writtencj: scaling factor, inverse of the step sizex: Vector with the statesdx: Vector with the derivative states
-
void
fJv(realtype t, const AmiVector &x, const AmiVector &dx, const AmiVector &xdot, const AmiVector &v, AmiVector &nJv, realtype cj) override¶ Jacobian multiply function.
- Parameters
t: timex: statedx: time derivative of state (DAE only)xdot: values of residual function (unused)v: multiplication vector (unused)nJv: array to which result of multiplication will be writtencj: scaling factor (inverse of timestep, DAE only)
-
void
fJv(realtype t, const_N_Vector x, const_N_Vector dx, const_N_Vector v, N_Vector Jv, realtype cj)¶ Matrix vector product of J with a vector v (for iterative solvers)
- Parameters
t: timepointcj: scaling factor, inverse of the step sizex: Vector with the statesdx: Vector with the derivative statesv: Vector with which the Jacobian is multipliedJv: Vector to which the Jacobian vector product will be written
-
void
fJvB(realtype t, const_N_Vector x, const_N_Vector dx, const_N_Vector xB, const_N_Vector dxB, const_N_Vector vB, N_Vector JvB, realtype cj)¶ Matrix vector product of JB with a vector v (for iterative solvers)
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative statesxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesvB: Vector with which the Jacobian is multipliedJvB: Vector to which the Jacobian vector product will be writtencj: scalar in Jacobian (inverse stepsize)
-
void
froot(realtype t, const AmiVector &x, const AmiVector &dx, gsl::span<realtype> root) override¶ Root function.
- Parameters
t: timex: statedx: time derivative of state (DAE only)root: array to which values of the root function will be written
-
void
froot(realtype t, const_N_Vector x, const_N_Vector dx, gsl::span<realtype> root)¶ Event trigger function for events.
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative statesroot: array with root function values
-
void
fxdot(realtype t, const AmiVector &x, const AmiVector &dx, AmiVector &xdot) override¶ Residual function.
- Parameters
t: timex: statedx: time derivative of state (DAE only)xdot: array to which values of the residual function will be written
-
void
fxdot(realtype t, const_N_Vector x, const_N_Vector dx, N_Vector xdot)¶ Residual function of the DAE.
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative statesxdot: Vector with the right hand side
-
void
fxBdot(realtype t, const_N_Vector x, const_N_Vector dx, const_N_Vector xB, const_N_Vector dxB, N_Vector xBdot)¶ Right hand side of differential equation for adjoint state xB.
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative statesxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesxBdot: Vector with the adjoint right hand side
-
void
fqBdot(realtype t, const_N_Vector x, const_N_Vector dx, const_N_Vector xB, const_N_Vector dxB, N_Vector qBdot)¶ Right hand side of integral equation for quadrature states qB.
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative statesxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesqBdot: Vector with the adjoint quadrature right hand side
-
void
fxBdot_ss(const realtype t, const AmiVector &xB, const AmiVector &dxB, AmiVector &xBdot) override¶ Residual function backward when running in steady state mode.
- Parameters
t: timexB: adjoint statedxB: time derivative of state (DAE only)xBdot: array to which values of the residual function will be written
-
void
fxBdot_ss(realtype t, const_N_Vector xB, const_N_Vector dxB, N_Vector xBdot) const¶ Implementation of fxBdot for steady state case at the N_Vector level.
- Parameters
t: timepointxB: Vector with the adjoint statedxB: Vector with the adjoint derivative statesxBdot: Vector with the adjoint right hand side
-
void
fqBdot_ss(realtype t, const_N_Vector xB, const_N_Vector dxB, N_Vector qBdot) const¶ Implementation of fqBdot for steady state at the N_Vector level.
- Parameters
t: timepointxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesqBdot: Vector with the adjoint quadrature right hand side
-
void
fJSparseB_ss(SUNMatrix JB) override¶ Sparse Jacobian function backward, steady state case.
- Parameters
JB: sparse matrix to which values of the Jacobian will be written
-
void
writeSteadystateJB(const realtype t, realtype cj, const AmiVector &x, const AmiVector &dx, const AmiVector &xB, const AmiVector &dxB, const AmiVector &xBdot) override¶ Computes the sparse backward Jacobian for steadystate integration and writes it to the model member.
- Parameters
t: timepointcj: scalar in Jacobianx: Vector with the statesdx: Vector with the derivative statesxB: Vector with the adjoint statesdxB: Vector with the adjoint derivative statesxBdot: Vector with the adjoint state right hand side
-
void
fdxdotdp(realtype t, const const_N_Vector x, const const_N_Vector dx)¶ Sensitivity of dx/dt wrt model parameters p.
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative states
-
void
fdxdotdp(const realtype t, const AmiVector &x, const AmiVector &dx) override¶ Model specific sparse implementation of explicit parameter derivative of right hand side.
- Parameters
t: timex: statedx: time derivative of state (DAE only)
-
void
fsxdot(realtype t, const AmiVector &x, const AmiVector &dx, int ip, const AmiVector &sx, const AmiVector &sdx, AmiVector &sxdot) override¶ Sensitivity Residual function.
- Parameters
t: timex: statedx: time derivative of state (DAE only)ip: parameter indexsx: sensitivity statesdx: time derivative of sensitivity state (DAE only)sxdot: array to which values of the sensitivity residual function will be written
-
void
fsxdot(realtype t, const_N_Vector x, const_N_Vector dx, int ip, const_N_Vector sx, const_N_Vector sdx, N_Vector sxdot)¶ Right hand side of differential equation for state sensitivities sx.
- Parameters
t: timepointx: Vector with the statesdx: Vector with the derivative statesip: parameter indexsx: Vector with the state sensitivitiessdx: Vector with the derivative state sensitivitiessxdot: Vector with the sensitivity right hand side
-
void
fM(realtype t, const_N_Vector x)¶ Mass matrix for DAE systems.
- Parameters
t: timepointx: Vector with the states
Protected Functions
-
void
fJSparse(SUNMatrixContent_Sparse JSparse, realtype t, const realtype *x, const double *p, const double *k, const realtype *h, realtype cj, const realtype *dx, const realtype *w, const realtype *dwdx) = 0¶ Model specific implementation for fJSparse.
- Parameters
JSparse: Matrix to which the Jacobian will be writtent: timepointx: Vector with the statesp: parameter vectork: constants vectorh: Heaviside vectorcj: scaling factor, inverse of the step sizedx: Vector with the derivative statesw: vector with helper variablesdwdx: derivative of w wrt x
-
void
froot(realtype *root, realtype t, const realtype *x, const double *p, const double *k, const realtype *h, const realtype *dx)¶ Model specific implementation for froot.
- Parameters
root: values of the trigger functiont: timepointx: Vector with the statesp: parameter vectork: constants vectorh: Heaviside vectordx: Vector with the derivative states
-
void
fxdot(realtype *xdot, realtype t, const realtype *x, const double *p, const double *k, const realtype *h, const realtype *dx, const realtype *w) = 0¶ Model specific implementation for fxdot.
- Parameters
xdot: residual functiont: timepointx: Vector with the statesp: parameter vectork: constants vectorh: Heaviside vectorw: vector with helper variablesdx: Vector with the derivative states
-
void
fdxdotdp(realtype *dxdotdp, realtype t, const realtype *x, const realtype *p, const realtype *k, const realtype *h, int ip, const realtype *dx, const realtype *w, const realtype *dwdp)¶ Model specific implementation of fdxdotdp.
- Parameters
dxdotdp: partial derivative xdot wrt pt: timepointx: Vector with the statesp: parameter vectork: constants vectorh: Heaviside vectorip: parameter indexdx: Vector with the derivative statesw: vector with helper variablesdwdp: derivative of w wrt p
-