Class Model_DAE
Defined in File model_dae.h
Inheritance Relationships
Base Type
public amici::Model
(Class Model)
Class Documentation
-
class 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
-
inline Model_DAE(ModelDimensions const &model_dimensions, SimulationParameters simulation_parameters, SecondOrderMode const o2mode, std::vector<realtype> const &idlist, std::vector<int> const &z2event, std::map<realtype, std::vector<int>> state_independent_events = {})
Constructor with model dimensions.
- Parameters:
model_dimensions – Model dimensions
simulation_parameters – Simulation parameters
o2mode – second order sensitivity mode
idlist – indexes indicating algebraic components (DAE only)
z2event – mapping of event outputs to events
state_independent_events – Map of events with state-independent triggers functions, mapping trigger timepoints to event indices.
-
virtual void fJ(realtype t, realtype cj, AmiVector const &x, AmiVector const &dx, AmiVector const &xdot, SUNMatrix J) override
Dense Jacobian function.
- Parameters:
t – time
cj – scaling factor (inverse of timestep, DAE only)
x – state
dx – 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 – timepoint
cj – scaling factor, inverse of the step size
x – Vector with the states
dx – Vector with the derivative states
xdot – Vector with the right hand side
J – Matrix to which the Jacobian will be written
-
virtual void fJB(realtype const t, realtype cj, AmiVector const &x, AmiVector const &dx, AmiVector const &xB, AmiVector const &dxB, AmiVector const &xBdot, SUNMatrix JB) override
Dense Jacobian function.
- Parameters:
t – time
cj – scaling factor (inverse of timestep, DAE only)
x – state
dx – time derivative of state (DAE only)
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
xBdot – 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 – timepoint
cj – scaling factor, inverse of the step size
x – Vector with the states
dx – Vector with the derivative states
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
JB – Matrix to which the Jacobian will be written
-
virtual void fJSparse(realtype t, realtype cj, AmiVector const &x, AmiVector const &dx, AmiVector const &xdot, SUNMatrix J) override
Sparse Jacobian function.
- Parameters:
t – time
cj – scaling factor (inverse of timestep, DAE only)
x – state
dx – 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 – timepoint
cj – scalar in Jacobian (inverse stepsize)
x – Vector with the states
dx – Vector with the derivative states
J – Matrix to which the Jacobian will be written
-
virtual void fJSparseB(realtype const t, realtype cj, AmiVector const &x, AmiVector const &dx, AmiVector const &xB, AmiVector const &dxB, AmiVector const &xBdot, SUNMatrix JB) override
Sparse Jacobian function.
- Parameters:
t – time
cj – scaling factor (inverse of timestep, DAE only)
x – state
dx – time derivative of state (DAE only)
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
xBdot – 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 – timepoint
cj – scalar in Jacobian
x – Vector with the states
dx – Vector with the derivative states
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
JB – Matrix to which the Jacobian will be written
-
virtual void fJDiag(realtype t, AmiVector &JDiag, realtype cj, AmiVector const &x, AmiVector const &dx) override
Diagonal of the Jacobian (for preconditioning)
- Parameters:
t – timepoint
JDiag – Vector to which the Jacobian diagonal will be written
cj – scaling factor, inverse of the step size
x – Vector with the states
dx – Vector with the derivative states
-
virtual void fJv(realtype t, AmiVector const &x, AmiVector const &dx, AmiVector const &xdot, AmiVector const &v, AmiVector &nJv, realtype cj) override
Jacobian multiply function.
- Parameters:
t – time
x – state
dx – 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 written
cj – 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 – timepoint
cj – scaling factor, inverse of the step size
x – Vector with the states
dx – Vector with the derivative states
v – Vector with which the Jacobian is multiplied
Jv – 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 – timepoint
x – Vector with the states
dx – Vector with the derivative states
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
vB – Vector with which the Jacobian is multiplied
JvB – Vector to which the Jacobian vector product will be written
cj – scalar in Jacobian (inverse stepsize)
-
virtual void froot(realtype t, AmiVector const &x, AmiVector const &dx, gsl::span<realtype> root) override
Root function.
- Parameters:
t – time
x – state
dx – 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 – timepoint
x – Vector with the states
dx – Vector with the derivative states
root – array with root function values
-
virtual void fxdot(realtype t, AmiVector const &x, AmiVector const &dx, AmiVector &xdot) override
Residual function.
- Parameters:
t – time
x – state
dx – 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 – timepoint
x – Vector with the states
dx – Vector with the derivative states
xdot – 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 – timepoint
x – Vector with the states
dx – Vector with the derivative states
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
xBdot – 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 – timepoint
x – Vector with the states
dx – Vector with the derivative states
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
qBdot – Vector with the adjoint quadrature right hand side
-
virtual void fxBdot_ss(realtype const t, AmiVector const &xB, AmiVector const &dxB, AmiVector &xBdot) override
Residual function backward when running in steady state mode.
- Parameters:
t – time
xB – adjoint state
dxB – 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 – timepoint
xB – Vector with the adjoint state
dxB – Vector with the adjoint derivative states
xBdot – 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 – timepoint
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
qBdot – Vector with the adjoint quadrature right hand side
-
virtual 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
-
virtual void writeSteadystateJB(realtype const t, realtype cj, AmiVector const &x, AmiVector const &dx, AmiVector const &xB, AmiVector const &dxB, AmiVector const &xBdot) override
Computes the sparse backward Jacobian for steadystate integration and writes it to the model member.
- Parameters:
t – timepoint
cj – scalar in Jacobian
x – Vector with the states
dx – Vector with the derivative states
xB – Vector with the adjoint states
dxB – Vector with the adjoint derivative states
xBdot – Vector with the adjoint state right hand side
-
void fdxdotdp(realtype t, const_N_Vector const x, const_N_Vector const dx)
Sensitivity of dx/dt wrt model parameters p.
- Parameters:
t – timepoint
x – Vector with the states
dx – Vector with the derivative states
-
inline virtual void fdxdotdp(realtype const t, AmiVector const &x, AmiVector const &dx) override
Model-specific sparse implementation of explicit parameter derivative of right hand side.
- Parameters:
t – time
x – state
dx – time derivative of state (DAE only)
-
virtual void fsxdot(realtype t, AmiVector const &x, AmiVector const &dx, int ip, AmiVector const &sx, AmiVector const &sdx, AmiVector &sxdot) override
Sensitivity Residual function.
- Parameters:
t – time
x – state
dx – time derivative of state (DAE only)
ip – parameter index
sx – sensitivity state
sdx – 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 – timepoint
x – Vector with the states
dx – Vector with the derivative states
ip – parameter index
sx – Vector with the state sensitivities
sdx – Vector with the derivative state sensitivities
sxdot – Vector with the sensitivity right hand side
-
void fM(realtype t, const_N_Vector x)
Mass matrix for DAE systems.
- Parameters:
t – timepoint
x – Vector with the states
Protected Functions
-
virtual void fJSparse(SUNMatrixContent_Sparse JSparse, realtype t, realtype const *x, double const *p, double const *k, realtype const *h, realtype cj, realtype const *dx, realtype const *w, realtype const *dwdx)
Model specific implementation for fJSparse.
- Parameters:
JSparse – Matrix to which the Jacobian will be written
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
cj – scaling factor, inverse of the step size
dx – Vector with the derivative states
w – vector with helper variables
dwdx – derivative of w wrt x
-
virtual void froot(realtype *root, realtype t, realtype const *x, double const *p, double const *k, realtype const *h, realtype const *dx)
Model specific implementation for froot.
- Parameters:
root – values of the trigger function
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
dx – Vector with the derivative states
-
virtual void fxdot(realtype *xdot, realtype t, realtype const *x, double const *p, double const *k, realtype const *h, realtype const *dx, realtype const *w) = 0
Model specific implementation for fxdot.
- Parameters:
xdot – residual function
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
w – vector with helper variables
dx – Vector with the derivative states
-
virtual void fdxdotdp(realtype *dxdotdp, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, int ip, realtype const *dx, realtype const *w, realtype const *dwdp)
Model specific implementation of fdxdotdp.
- Parameters:
dxdotdp – partial derivative xdot wrt p
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
ip – parameter index
dx – Vector with the derivative states
w – vector with helper variables
dwdp – derivative of w wrt p
-
virtual void fdxdotdp_explicit(realtype *dxdotdp_explicit, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *dx, realtype const *w)
Model specific implementation of fdxdotdp_explicit, no w chainrule (Py)
- Parameters:
dxdotdp_explicit – partial derivative xdot wrt p
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
dx – Vector with the derivative states
w – vector with helper variables
-
virtual void fdxdotdp_explicit_colptrs(SUNMatrixWrapper &dxdotdp)
Model specific implementation of fdxdotdp_explicit, colptrs part.
- Parameters:
dxdotdp – sparse matrix to which colptrs will be written
-
virtual void fdxdotdp_explicit_rowvals(SUNMatrixWrapper &dxdotdp)
Model specific implementation of fdxdotdp_explicit, rowvals part.
- Parameters:
dxdotdp – sparse matrix to which rowvals will be written
-
virtual void fdxdotdx_explicit(realtype *dxdotdx_explicit, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *dx, realtype const *w)
Model specific implementation of fdxdotdx_explicit, no w chainrule (Py)
- Parameters:
dxdotdx_explicit – partial derivative xdot wrt x
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – heavyside vector
dx – Vector with the derivative states
w – vector with helper variables
-
virtual void fdxdotdx_explicit_colptrs(SUNMatrixWrapper &dxdotdx)
Model specific implementation of fdxdotdx_explicit, colptrs part.
- Parameters:
dxdotdx – sparse matrix to which colptrs will be written
-
virtual void fdxdotdx_explicit_rowvals(SUNMatrixWrapper &dxdotdx)
Model specific implementation of fdxdotdx_explicit, rowvals part.
- Parameters:
dxdotdx – sparse matrix to which rowvals will be written
-
virtual void fdxdotdw(realtype *dxdotdw, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *dx, realtype const *w)
Model specific implementation of fdxdotdw, data part.
- Parameters:
dxdotdw – partial derivative xdot wrt w
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
dx – Vector with the derivative states
w – vector with helper variables
-
virtual void fdxdotdw_colptrs(SUNMatrixWrapper &dxdotdw)
Model specific implementation of fdxdotdw, colptrs part.
- Parameters:
dxdotdw – sparse matrix to which colptrs will be written
-
virtual void fdxdotdw_rowvals(SUNMatrixWrapper &dxdotdw)
Model specific implementation of fdxdotdw, rowvals part.
- Parameters:
dxdotdw – sparse matrix to which rowvals will be written
-
void fdxdotdw(realtype t, const_N_Vector x, const_N_Vector dx)
Sensitivity of dx/dt wrt model parameters w.
- Parameters:
t – timepoint
x – Vector with the states
dx – Vector with the derivative states
-
Model_DAE() = default