Class Model
Defined in File model.h
Inheritance Relationships
Base Types
public amici::AbstractModel(Class AbstractModel)public amici::ModelDimensions(Struct ModelDimensions)
Derived Types
public amici::Model_DAE(Class Model_DAE)public amici::Model_ODE(Class Model_ODE)
Class Documentation
-
class Model : public amici::AbstractModel, public amici::ModelDimensions
The Model class represents an AMICI ODE/DAE model.
The model can compute various model related quantities based on symbolically generated code.
Subclassed by amici::Model_DAE, amici::Model_ODE
Public Functions
-
Model() = default
Default constructor
-
Model(ModelDimensions const &model_dimensions, SimulationParameters simulation_parameters, SecondOrderMode o2mode, std::vector<realtype> idlist, std::vector<int> z2event, std::vector<Event> 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
events – Vector of events
-
~Model() override = default
Destructor.
-
Model &operator=(Model const &other) = delete
Copy assignment is disabled until const members are removed.
- Parameters:
other – Object to copy from
- Returns:
-
void initialize(realtype t, AmiVector &x, AmiVector &dx, AmiVectorArray &sx, AmiVectorArray &sdx, bool computeSensitivities, std::vector<int> &roots_found)
Initialize model properties.
- Parameters:
t – Timepoint
x – Reference to state variables
dx – Reference to time derivative of states (DAE only)
sx – Reference to state variable sensitivities
sdx – Reference to time derivative of state sensitivities (DAE only)
computeSensitivities – Flag indicating whether sensitivities are to be computed
roots_found – boolean indicators indicating whether roots were found at t0 by this fun
-
void reinitialize(realtype t, AmiVector &x, AmiVectorArray &sx, bool computeSensitivities)
Re-initialize model properties after changing simulation context.
- Parameters:
t – Timepoint
x – Reference to state variables
sx – Reference to state variable sensitivities
computeSensitivities – Flag indicating whether sensitivities are to be computed
-
void initialize_b(AmiVector &xB, AmiVector &dxB, AmiVector &xQB, bool posteq) const
Initialize model properties.
- Parameters:
xB – Adjoint state variables
dxB – Time derivative of adjoint states (DAE only)
xQB – Adjoint quadratures
posteq – Flag indicating whether postequilibration was performed
-
void initialize_state(realtype t, AmiVector &x)
Initialize model state.
- Parameters:
t – Initial timepoint
x – State vector to be initialized (size: nx_solver).
-
void initialize_state_sensitivities(realtype t, AmiVectorArray &sx, AmiVector const &x)
Initialize initial state sensitivities.
- Parameters:
t – Initial timepoint.
sx – Reference to state variable sensitivities
x – Reference to state variables
-
void initialize_splines()
Initialization of spline functions.
-
void initialize_spline_sensitivities()
Initialization of spline sensitivity functions.
-
void initialize_events(realtype t, AmiVector const &x, AmiVector const &dx, std::vector<int> &roots_found)
Initialize the Heaviside variables
hat the initial timet0.Heaviside variables activate/deactivate on event occurrences.
- Parameters:
t – Timepoint
x – Reference to state variables
dx – Reference to time derivative of states (DAE only)
roots_found – boolean indicators indicating whether roots were found at t0 by this fun
-
void reinit_events(realtype t, AmiVector const &x, AmiVector const &dx, std::vector<realtype> const &h_old, std::vector<int> &roots_found)
Re-initialize the Heaviside variables
hat the beginning of the second or subsequent simulation.- Parameters:
t – Timepoint
x – Reference to state variables
dx – Reference to time derivative of states (DAE only)
h_old – Previous values of Heaviside variables
roots_found – boolean indicators indicating whether roots were found at t0 by this fun
-
void reinit_explicit_roots()
Re-compute the explicit roots.
Re-compute the explicit roots based on the current model parameters.
-
int nplist() const
Get number of parameters wrt to which sensitivities are computed.
- Returns:
Length of sensitivity index vector
-
int np() const
Get total number of model parameters.
- Returns:
Length of parameter vector
-
int nk() const
Get number of constants.
- Returns:
Length of constant vector
-
int ncl() const
Get number of conservation laws.
- Returns:
Number of conservation laws (i.e., difference between
nx_rdataandnx_solver).
-
int nx_reinit() const
Get number of solver states subject to reinitialization.
- Returns:
Model member
nx_solver_reinit
-
double const *k() const
Get fixed parameters.
- Returns:
Pointer to constants array
-
int n_max_event() const
Get maximum number of events that may occur for each type.
- Returns:
Maximum number of events that may occur for each type
-
void set_n_max_event(int nmaxevent)
Set maximum number of events that may occur for each type.
- Parameters:
nmaxevent – Maximum number of events that may occur for each type
-
int nt() const
Get number of timepoints.
- Returns:
Number of timepoints
-
std::vector<ParameterScaling> const &get_parameter_scale() const
Get parameter scale for each parameter.
- Returns:
Vector of parameter scales
-
void set_parameter_scale(ParameterScaling pscale)
Set parameter scale for each parameter.
NOTE: Resets initial state sensitivities.
- Parameters:
pscale – Scalar parameter scale to be set for all parameters
-
void set_parameter_scale(std::vector<ParameterScaling> const &pscaleVec)
Set parameter scale for each parameter.
NOTE: Resets initial state sensitivities.
- Parameters:
pscaleVec – Vector of parameter scales
-
std::vector<realtype> const &get_unscaled_parameters() const
Get parameters with transformation according to parameter scale applied.
- Returns:
Unscaled parameters
-
std::vector<realtype> const &get_parameters() const
Get parameter vector.
- Returns:
The user-set parameters (see also
Model::getUnscaledParameters)
-
realtype get_parameter_by_id(std::string const &par_id) const
Get value of first model parameter with the specified ID.
- Parameters:
par_id – Parameter ID
- Returns:
Parameter value
-
realtype get_parameter_by_name(std::string const &par_name) const
Get value of first model parameter with the specified name.
- Parameters:
par_name – Parameter name
- Returns:
Parameter value
-
void set_parameters(std::vector<realtype> const &p)
Set the parameter vector.
- Parameters:
p – Vector of parameters
-
void set_parameter_by_id(std::map<std::string, realtype> const &p, bool ignoreErrors = false)
Set model parameters according to the parameter IDs and mapped values.
- Parameters:
p – Map of parameters IDs and values
ignoreErrors – Ignore errors such as parameter IDs in p which are not model parameters
-
void set_parameter_by_id(std::string const &par_id, realtype value)
Set value of first model parameter with the specified ID.
- Parameters:
par_id – Parameter ID
value – Parameter value
-
int set_parameters_by_id_regex(std::string const &par_id_regex, realtype value)
Set all values of model parameters with IDs matching the specified regular expression.
- Parameters:
par_id_regex – Parameter ID regex
value – Parameter value
- Returns:
Number of parameter IDs that matched the regex
-
void set_parameter_by_name(std::string const &par_name, realtype value)
Set value of first model parameter with the specified name.
- Parameters:
par_name – Parameter name
value – Parameter value
-
void set_parameter_by_name(std::map<std::string, realtype> const &p, bool ignoreErrors = false)
Set model parameters according to the parameter name and mapped values.
- Parameters:
p – Map of parameters names and values
ignoreErrors – Ignore errors such as parameter names in p which are not model parameters
-
int set_parameters_by_name_regex(std::string const &par_name_regex, realtype value)
Set all values of all model parameters with names matching the specified regex.
- Parameters:
par_name_regex – Parameter name regex
value – Parameter value
- Returns:
Number of fixed parameter names that matched the regex
-
std::vector<realtype> const &get_fixed_parameters() const
Get values of fixed parameters.
- Returns:
Vector of fixed parameters with same ordering as in Model::getFixedParameterIds
-
realtype get_fixed_parameter_by_id(std::string const &par_id) const
Get value of fixed parameter with the specified ID.
- Parameters:
par_id – Parameter ID
- Returns:
Parameter value
-
realtype get_fixed_parameter_by_name(std::string const &par_name) const
Get value of fixed parameter with the specified name.
If multiple parameters have the same name, the first parameter with matching name is returned.
- Parameters:
par_name – Parameter name
- Returns:
Parameter value
-
void set_fixed_parameters(std::vector<realtype> const &k)
Set values for constants.
- Parameters:
k – Vector of fixed parameters
-
void set_fixed_parameter_by_id(std::string const &par_id, realtype value)
Set value of first fixed parameter with the specified ID.
- Parameters:
par_id – Fixed parameter id
value – Fixed parameter value
-
int set_fixed_parameters_by_id_regex(std::string const &par_id_regex, realtype value)
Set values of all fixed parameters with the ID matching the specified regex.
- Parameters:
par_id_regex – Fixed parameter name regex
value – Fixed parameter value
- Returns:
Number of fixed parameter IDs that matched the regex
-
void set_fixed_parameter_by_name(std::string const &par_name, realtype value)
Set value of first fixed parameter with the specified name.
- Parameters:
par_name – Fixed parameter ID
value – Fixed parameter value
-
int set_fixed_parameters_by_name_regex(std::string const &par_name_regex, realtype value)
Set value of all fixed parameters with name matching the specified regex.
- Parameters:
par_name_regex – Fixed parameter name regex
value – Fixed parameter value
- Returns:
Number of fixed parameter names that matched the regex
-
virtual bool has_parameter_names() const
Report whether the model has parameter names set.
- Returns:
Boolean indicating whether parameter names were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_parameter_names() const
Get names of the model parameters.
- Returns:
The parameter names
-
virtual bool has_state_names() const
Report whether the model has state names set.
- Returns:
Boolean indicating whether state names were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_state_names() const
Get names of the model states.
- Returns:
State names
-
virtual std::vector<std::string> get_state_names_solver() const
Get names of the solver states.
- Returns:
State names
-
virtual bool has_fixed_parameter_names() const
Report whether the model has fixed parameter names set.
- Returns:
Boolean indicating whether fixed parameter names were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_fixed_parameter_names() const
Get names of the fixed model parameters.
- Returns:
Fixed parameter names
-
virtual bool has_observable_names() const
Report whether the model has observable names set.
- Returns:
Boolean indicating whether observable names were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_observable_names() const
Get names of the observables.
- Returns:
Observable names
-
virtual bool has_expression_names() const
Report whether the model has expression names set.
- Returns:
Boolean indicating whether expression names were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_expression_names() const
Get names of the expressions.
- Returns:
Expression names
-
virtual bool has_parameter_ids() const
Report whether the model has parameter IDs set.
- Returns:
Boolean indicating whether parameter IDs were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_parameter_ids() const
Get IDs of the model parameters.
- Returns:
Parameter IDs
-
virtual bool has_state_ids() const
Report whether the model has state IDs set.
- Returns:
Boolean indicating whether state IDs were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_state_ids() const
Get IDs of the model states.
- Returns:
State IDs
-
virtual std::vector<std::string> get_state_ids_solver() const
Get IDs of the solver states.
- Returns:
State IDs
-
virtual bool has_fixed_parameter_ids() const
Report whether the model has fixed parameter IDs set.
- Returns:
Boolean indicating whether fixed parameter IDs were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_fixed_parameter_ids() const
Get IDs of the fixed model parameters.
- Returns:
Fixed parameter IDs
-
virtual bool has_observable_ids() const
Report whether the model has observable IDs set.
- Returns:
Boolean indicating whether observable ids were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_observable_ids() const
Get IDs of the observables.
- Returns:
Observable IDs
-
virtual bool has_expression_ids() const
Report whether the model has expression IDs set.
- Returns:
Boolean indicating whether expression ids were set. Also returns
trueif the number of corresponding variables is just zero.
-
virtual std::vector<std::string> get_expression_ids() const
Get IDs of the expression.
- Returns:
Expression IDs
-
virtual bool has_quadratic_llh() const
Checks whether the defined noise model is Gaussian, i.e., the nllh is quadratic.
- Returns:
boolean flag
-
std::vector<realtype> const &get_timepoints() const
Get the timepoint vector.
- Returns:
Timepoint vector
-
realtype get_timepoint(int it) const
Get simulation timepoint for time index
it.- Parameters:
it – Time index
- Returns:
Timepoint
-
void set_timepoints(std::vector<realtype> const &ts)
Set the timepoint vector.
- Parameters:
ts – New timepoint vector
-
double t0() const
Get simulation start time.
- Returns:
Simulation start time
-
void set_t0(double t0)
Set simulation start time.
Output timepoints are absolute timepoints, independent of \( t_{0} \). For output timepoints \( t < t_{0} \), the initial state will be returned.
- Parameters:
t0 – Simulation start time
-
double t0_preeq() const
Get the initial time to use for pre-equilibration.
- Returns:
Initial time, or NAN to use the model’s t0.
-
void set_t0_preeq(double t0_preeq)
Set the initial time to use for pre-equilibration.
- Parameters:
t0_preeq – The initial time for pre-equilibration or NAN to use the model’s t0.
-
std::vector<bool> const &get_state_is_non_negative() const
Get flags indicating whether states should be treated as non-negative.
- Returns:
Vector of flags
-
void set_state_is_non_negative(std::vector<bool> const &stateIsNonNegative)
Set flags indicating whether states should be treated as non-negative.
- Parameters:
stateIsNonNegative – Vector of flags
-
void set_all_states_non_negative()
Set flags indicating that all states should be treated as non-negative.
-
inline ModelState const &get_model_state() const
Get the current model state.
- Returns:
Current model state
-
inline void set_model_state(ModelState const &state)
Set the current model state.
- Parameters:
state – Model state
-
inline void set_minimum_sigma_residuals(double const min_sigma)
Sets the estimated lower boundary for sigma_y. When :meth:
setAddSigmaResidualsis activated, this lower boundary must ensure that log(sigma) + min_sigma > 0.- Parameters:
min_sigma – lower boundary
-
inline realtype get_minimum_sigma_residuals() const
Gets the specified estimated lower boundary for sigma_y.
- Returns:
lower boundary
-
inline void set_add_sigma_residuals(bool const sigma_res)
Specifies whether residuals should be added to account for parameter dependent sigma.
If set to true, additional residuals of the form \( \sqrt{\log(\sigma) + C} \) will be added. This enables least-squares optimization for variables with Gaussian noise assumption and parameter dependent standard deviation sigma. The constant \( C \) can be set via :meth:
setMinimumSigmaResiduals.- Parameters:
sigma_res – if true, additional residuals are added
-
inline bool get_add_sigma_residuals() const
Checks whether residuals should be added to account for parameter dependent sigma.
- Returns:
sigma_res
-
std::vector<int> const &get_parameter_list() const
Get the list of parameters for which sensitivities are computed.
- Returns:
List of parameter indices
-
int plist(int pos) const
Get entry in parameter list by index.
- Parameters:
pos – Index in sensitivity parameter list
- Returns:
Index in parameter list
-
void set_parameter_list(std::vector<int> const &plist)
Set the list of parameters for which sensitivities are to be computed.
NOTE: Resets initial state sensitivities.
- Parameters:
plist – List of parameter indices
-
std::vector<realtype> get_initial_state(realtype t0)
Get the initial state.
- Parameters:
t0 – Custom t0 for which to get initial states.
- Returns:
Initial state vector, before any events are executed.
-
inline std::vector<realtype> get_initial_state()
Get the initial state for Model::t0()`.
- Returns:
Initial state vector, before any events are executed.
-
void set_initial_state(std::vector<realtype> const &x0)
Set the pre-event initial state.
- Parameters:
x0 – Initial state vector
-
bool has_custom_initial_state() const
Return whether custom initial state have been set.
- Returns:
trueif has custom initial state, otherwisefalse
-
inline std::vector<realtype> get_initial_state_sensitivities()
Get the initial state sensitivities.
- Returns:
vector of initial state sensitivities
-
std::vector<realtype> get_initial_state_sensitivities(realtype t0)
Get the initial states sensitivities.
- Parameters:
t0 – Custom t0 for which to get initial states.
- Returns:
vector of initial state sensitivities
-
void set_initial_state_sensitivities(std::vector<realtype> const &sx0)
Set the initial state sensitivities.
- Parameters:
sx0 – vector of initial state sensitivities with chain rule applied. This could be a slice of ReturnData::sx or ReturnData::sx0
-
bool has_custom_initial_state_sensitivities() const
Return whether custom initial state sensitivities have been set.
- Returns:
trueif has custom initial state sensitivities, otherwisefalse.
-
void set_unscaled_initial_state_sensitivities(std::vector<realtype> const &sx0)
Set the initial state sensitivities.
- Parameters:
sx0 – Vector of initial state sensitivities without chain rule applied. This could be the read-in from a
model.sx0datasaved to HDF5.
-
void set_steady_state_computation_mode(SteadyStateComputationMode mode)
Set the mode how steady state is computed in the steady state simulation.
- Parameters:
mode – Steady state computation mode
-
SteadyStateComputationMode get_steady_state_computation_mode() const
Gets the mode how steady state is computed in the steadystate simulation.
- Returns:
Mode
-
void set_steady_state_sensitivity_mode(SteadyStateSensitivityMode mode)
Set the mode how sensitivities are computed in the steady-state simulation.
- Parameters:
mode – Steady-state sensitivity mode
-
SteadyStateSensitivityMode get_steady_state_sensitivity_mode() const
Gets the mode how sensitivities are computed in the steady-state simulation.
- Returns:
Mode
-
void set_reinitialize_fixed_parameter_initial_states(bool flag)
Set whether initial states depending on fixed parameters are to be reinitialized after preequilibration and presimulation.
- Parameters:
flag – Fixed parameters reinitialized?
-
bool get_reinitialize_fixed_parameter_initial_states() const
Get whether initial states depending on fixedParameters are to be reinitialized after preequilibration and presimulation.
- Returns:
flag
true/false
-
void require_sensitivities_for_all_parameters()
Require computation of sensitivities for all parameters p [0..np[ in natural order.
NOTE: Resets initial state sensitivities.
-
void get_expression(gsl::span<realtype> w, realtype t, AmiVector const &x)
Get time-resolved
w.- Parameters:
w – Buffer (shape
nw)t – Current timepoint
x – Current state
-
void get_observable(gsl::span<realtype> y, realtype t, AmiVector const &x)
Get time-resolved observables.
- Parameters:
y – Buffer (shape
ny)t – Current timepoint
x – Current state
-
virtual ObservableScaling get_observable_scaling(int iy) const
Get scaling type for observable.
- Parameters:
iy – observable index
- Returns:
scaling type
-
void get_observable_sensitivity(gsl::span<realtype> sy, realtype t, AmiVector const &x, AmiVectorArray const &sx)
Get sensitivity of time-resolved observables.
Total derivative \( sy = dydx * sx + dydp\) (only for forward sensitivities).
- Parameters:
sy – buffer (shape
nyxnplist, row-major)t – Timpoint
x – State variables
sx – State sensitivities
-
void get_observable_sigma(gsl::span<realtype> sigmay, int it, ExpData const *edata)
Get time-resolved observable standard deviations.
- Parameters:
sigmay – Buffer (shape
ny)it – Timepoint index
edata – Pointer to experimental data instance (optional, pass
nullptrto ignore)
-
void get_observable_sigma_sensitivity(gsl::span<realtype> ssigmay, gsl::span<realtype const> sy, int it, ExpData const *edata)
Sensitivity of time-resolved observable standard deviation.
Total derivative (can be used with both adjoint and forward sensitivity).
- Parameters:
ssigmay – Buffer (shape
nyxnplist, row-major)sy – Sensitivity of time-resolved observables for current timepoint
it – Timepoint index
edata – Pointer to experimental data instance (optional, pass
nullptrto ignore)
-
void add_observable_objective(realtype &Jy, int it, AmiVector const &x, ExpData const &edata)
Add time-resolved measurement negative log-likelihood \( Jy \).
- Parameters:
Jy – Buffer (shape 1)
it – Timepoint index
x – State variables
edata – Experimental data
-
void add_observable_objective_sensitivity(std::vector<realtype> &sllh, std::vector<realtype> &s2llh, int it, AmiVector const &x, AmiVectorArray const &sx, ExpData const &edata)
Add sensitivity of time-resolved measurement negative log-likelihood \( Jy \).
- Parameters:
sllh – First-order buffer (shape
nplist)s2llh – Second-order buffer (shape
nJ - 1xnplist, row-major)it – Timepoint index
x – State variables
sx – State sensitivities
edata – Experimental data
-
void add_partial_observable_objective_sensitivity(std::vector<realtype> &sllh, std::vector<realtype> &s2llh, int it, AmiVector const &x, ExpData const &edata)
Add sensitivity of time-resolved measurement negative log-likelihood \( Jy \).
Partial derivative (to be used with adjoint sensitivities).
- Parameters:
sllh – First order output buffer (shape
nplist)s2llh – Second order output buffer (shape
nJ - 1xnplist, row-major)it – Timepoint index
x – State variables
edata – Experimental data
-
void get_adjoint_state_observable_update(gsl::span<realtype> dJydx, int it, AmiVector const &x, ExpData const &edata)
Get state sensitivity of the negative log-likelihood \( Jy \), partial derivative (to be used with adjoint sensitivities).
- Parameters:
dJydx – Output buffer (shape
nJxnx_solver, row-major)it – Timepoint index
x – State variables
edata – Experimental data instance
-
void get_event(gsl::span<realtype> z, int ie, realtype t, AmiVector const &x)
Get event-resolved observables.
- Parameters:
z – Output buffer (shape
nz)ie – Event index
t – Timepoint
x – State variables
-
void get_event_sensitivity(gsl::span<realtype> sz, int ie, realtype t, AmiVector const &x, AmiVectorArray const &sx)
Get sensitivities of event-resolved observables.
Total derivative (only forward sensitivities).
- Parameters:
sz – Output buffer (shape
nz x nplist, row-major)ie – Event index
t – Timepoint
x – State variables
sx – State sensitivities
-
void get_unobserved_event_sensitivity(gsl::span<realtype> sz, int ie)
Get sensitivity of
zat final timepoint.Ignores sensitivity of timepoint. Total derivative.
- Parameters:
sz – Output buffer (shape
nz x nplist, row-major)ie – Event index
-
void get_event_regularization(gsl::span<realtype> rz, int ie, realtype t, AmiVector const &x)
Get regularization for event-resolved observables.
- Parameters:
rz – Output buffer (shape
nz)ie – Event index
t – Timepoint
x – State variables
-
void get_event_regularization_sensitivity(gsl::span<realtype> srz, int ie, realtype t, AmiVector const &x, AmiVectorArray const &sx)
Get sensitivities of regularization for event-resolved observables.
Total derivative. Only forward sensitivities.
- Parameters:
srz – Output buffer (shape
nz x nplist, row-major)ie – Event index
t – Timepoint
x – State variables
sx – State sensitivities
-
void get_event_sigma(gsl::span<realtype> sigmaz, int ie, int nroots, realtype t, ExpData const *edata)
Get event-resolved observable standard deviations.
-
void get_event_sigma_sensitivity(gsl::span<realtype> ssigmaz, int ie, int nroots, realtype t, ExpData const *edata)
Get sensitivities of event-resolved observable standard deviations.
Total derivative (only forward sensitivities).
-
void add_event_objective(realtype &Jz, int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Add event-resolved observable negative log-likelihood.
-
void add_event_objective_regularization(realtype &Jrz, int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Add event-resolved observable negative log-likelihood.
-
void add_event_objective_sensitivity(std::vector<realtype> &sllh, std::vector<realtype> &s2llh, int ie, int nroots, realtype t, AmiVector const &x, AmiVectorArray const &sx, ExpData const &edata)
Add sensitivity of time-resolved measurement negative log-likelihood \( Jy \).
Total derivative (to be used with forward sensitivities).
-
void add_partial_event_objective_sensitivity(std::vector<realtype> &sllh, std::vector<realtype> &s2llh, int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Add sensitivity of time-resolved measurement negative log-likelihood \( Jy \).
Partial derivative (to be used with adjoint sensitivities).
-
void get_adjoint_state_event_update(gsl::span<realtype> dJzdx, int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
State sensitivity of the negative log-likelihood \( Jz \).
Partial derivative (to be used with adjoint sensitivities).
-
void get_event_time_sensitivity(std::vector<realtype> &stau, realtype t, int ie, AmiVector const &x, AmiVectorArray const &sx, AmiVector const &dx)
Sensitivity of event timepoint, total derivative.
Only forward sensitivities.
- Parameters:
stau – Timepoint sensitivity (shape
nplist)t – Timepoint
ie – Event index
x – State variables
sx – State sensitivities
dx – Current derivative of state (DAE only)
-
void add_state_event_update(AmiVector &x, int ie, realtype t, AmiVector const &xdot, AmiVector const &xdot_old, AmiVector const &x_old, ModelState const &state)
Update state variables after event.
- Parameters:
x – Current state (will be overwritten)
ie – Event index
t – Current timepoint
xdot – Current residual function values
xdot_old – Value of residual function before event
x_old – Current or old state from which to compute the state update
state – The model state based on which to compute the update.
-
void add_state_sensitivity_event_update(AmiVectorArray &sx, int ie, realtype t, AmiVector const &x, AmiVector const &x_old, AmiVector const &xdot, AmiVector const &xdot_old, AmiVectorArray const &sx_old, std::vector<realtype> const &stau)
Update state sensitivity after event.
- Parameters:
sx – Current state sensitivity (will be overwritten)
ie – Event index
t – Current timepoint
x – Current state
x_old – Pre-event state
xdot – Current residual function values
xdot_old – Value of residual function before event
sx_old – Pre-event state sensitivity
stau – Timepoint sensitivity, to be computed with
Model::getEventTimeSensitivity
-
void add_adjoint_state_event_update(AmiVector &xB, int ie, realtype t, AmiVector const &x, AmiVector const &xdot, AmiVector const &xdot_old, AmiVector const &x_old, AmiVector const &dx)
Update adjoint state after event.
- Parameters:
xB – Current adjoint state (will be overwritten)
ie – Event index
t – Current timepoint
x – Current state
xdot – Current residual function values
xdot_old – Value of residual function before event
x_old – Pre-event state
dx – Current derivative of state (DAE only)
-
void add_adjoint_quadrature_event_update(AmiVector &xQB, int ie, realtype t, AmiVector const &x, AmiVector const &xB, AmiVector const &xdot, AmiVector const &xdot_old, AmiVector const &x_old, AmiVector const &dx)
Update adjoint quadratures after event.
- Parameters:
xQB – Current quadrature state (will be overwritten)
ie – Event index
t – Current timepoint
x – Current state
xB – Current adjoint state
xdot – Current residual function values
xdot_old – Value of residual function before event
x_old – Pre-event state
dx – Current derivative of state (DAE only)
-
void update_heaviside(std::vector<int> const &rootsfound)
Update the Heaviside variables
hon event occurrences.- Parameters:
rootsfound – Provides the direction of the zero-crossing, so adding it will give the right update to the Heaviside variables (zero if no root was found)
-
int check_finite(gsl::span<realtype const> array, ModelQuantity model_quantity, realtype t) const
Check if the given array has only finite elements.
For (1D) spans.
- Parameters:
array –
model_quantity – The model quantity
arraycorresponds tot – Current timepoint
- Returns:
-
int check_finite(gsl::span<realtype const> array, ModelQuantity model_quantity, size_t num_cols, realtype t) const
Check if the given array has only finite elements.
For flattened 2D arrays.
- Parameters:
array – Flattened matrix
model_quantity – The model quantity
arraycorresponds tonum_cols – Number of columns of the non-flattened matrix
t – Current timepoint
- Returns:
-
int check_finite(SUNMatrix m, ModelQuantity model_quantity, realtype t) const
Check if the given array has only finite elements.
For SUNMatrix.
- Parameters:
m – Matrix to check
model_quantity – The model quantity
mcorresponds tot – current timepoint
- Returns:
-
void set_always_check_finite(bool alwaysCheck)
Set whether the result of every call to
Model::f*should be checked for finiteness.- Parameters:
alwaysCheck –
-
bool get_always_check_finite() const
Get setting of whether the result of every call to
Model::f*should be checked for finiteness.- Returns:
that
-
void fx0(realtype t, AmiVector &x)
Compute/get pre-event initial state.
- Parameters:
t – Timepoint.
x – Output buffer (size: nx_solver).
-
void fx0_fixedParameters(realtype t, AmiVector &x)
Set only those initial states that are specified via fixed parameters.
- Parameters:
t – Timepoint.
x – Output buffer.
-
void fsx0(realtype t, AmiVectorArray &sx, AmiVector const &x)
Compute/get initial value for initial state sensitivities.
- Parameters:
t – Timepoint.
sx – Output buffer for state sensitivities
x – State variables
-
void fsx0_fixedParameters(realtype t, AmiVectorArray &sx, AmiVector const &x)
Get only those initial states sensitivities that are affected from
amici::Model::fx0_fixedParameters.- Parameters:
t – Timepoint
sx – Output buffer for state sensitivities
x – State variables
-
virtual void fsdx0()
Compute sensitivity of derivative initial states sensitivities
sdx0.Only necessary for DAEs.
-
void fx_rdata(gsl::span<realtype> x_rdata, AmiVector const &x_solver)
Expand conservation law for states.
- Parameters:
x_rdata – Output buffer for state variables with conservation laws expanded (stored in
amici::ReturnData).x_solver – State variables with conservation laws applied (solver returns this)
-
void fsx_rdata(gsl::span<realtype> sx_rdata, AmiVectorArray const &sx_solver, AmiVector const &x_solver)
Expand conservation law for state sensitivities.
- Parameters:
sx_rdata – Output buffer for state variables sensitivities with conservation laws expanded (stored in
amici::ReturnDatashapenplistxnx, row-major).sx_solver – State variables sensitivities with conservation laws applied (solver returns this)
x_solver – State variables with conservation laws applied (solver returns this)
-
void set_reinitialization_state_idxs(std::vector<int> const &idxs)
Set indices of states to be reinitialized based on provided constants / fixed parameters.
- Parameters:
idxs – Array of state indices
-
std::vector<int> const &get_reinitialization_state_idxs() const
Return indices of states to be reinitialized based on provided constants / fixed parameters.
- Returns:
Those indices.
-
SUNMatrixWrapper const &get_dxdotdp_full() const
getter for dxdotdp
- Returns:
dxdotdp
-
virtual std::vector<double> get_trigger_timepoints() const
Get trigger times for events that don’t require root-finding.
To be called only after Model::initialize.
- Returns:
List of unique trigger points for events that don’t require root-finding (i.e. that trigger at predetermined timepoints), in ascending order.
-
inline std::vector<realtype> get_steadystate_mask() const
Get steady-state mask as std::vector.
See
set_steadystate_maskfor details.- Returns:
Steady-state mask
-
void set_steadystate_mask(std::vector<realtype> const &mask)
Set steady-state mask.
The mask is used to exclude certain state variables from the steady-state convergence check. Positive values indicate that the corresponding state variable should be included in the convergence check, while non-positive values indicate that the corresponding state variable should be excluded. An empty mask is interpreted as including all state variables.
- Parameters:
mask – Mask of length
nx_solver.
-
inline Event const &get_event(int const ie) const
Get event object for event index.
- Parameters:
ie – event index
- Returns:
The corresponding Event object.
-
inline bool get_any_state_nonnegative() const
Whether there is at least one state variable for which non-negativity is to be enforced.
- Returns:
Vector of all events.
-
inline virtual std::vector<std::vector<realtype>> fexplicit_roots(realtype const *p, realtype const *k) override
Compute explicit roots of the model.
- Parameters:
p – parameter vector
k – constant vector
- Returns:
A vector of length ne_solver, each containing a vector of explicit roots for the corresponding event.
-
inline std::vector<realtype> const &get_id_list() const
Get flag array for DAE equations.
- Returns:
Flag array for DAE equations
-
inline SecondOrderMode get_second_order_mode() const
Get second-order sensitivity mode.
- Returns:
Flag indicating whether for
SensitivityOrder::seconddirectional or full second order derivative will be computed
-
inline void set_logger(Logger *logger)
Set the logger.
- Parameters:
logger – Logger pointer or nullptr.
-
inline std::map<realtype, std::vector<int>> const &get_explicit_roots() const
Get map of trigger timepoints to event indices for events that don’t require root-finding.
- Returns:
Map of trigger timepoints to event indices
-
virtual void fdeltaqB(realtype *deltaqB, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *dx, int ip, int ie, realtype const *xdot, realtype const *xdot_old, realtype const *x_old, realtype const *xB)
Model-specific implementation of fdeltaqB.
- Parameters:
deltaqB – sensitivity update
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
dx – time derivative of state (DAE only)
ip – sensitivity index
ie – event index
xdot – new model right hand side
xdot_old – previous model right hand side
x_old – pre-event state
xB – adjoint state
-
virtual void fdeltasx(realtype *deltasx, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *w, int ip, int ie, realtype const *xdot, realtype const *xdot_old, realtype const *sx, realtype const *stau, realtype const *tcl, realtype const *x_old)
Model-specific implementation of fdeltasx.
- Parameters:
deltasx – sensitivity update
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
w – repeating elements vector
ip – sensitivity index
ie – event index
xdot – new model right hand side
xdot_old – previous model right hand side
sx – state sensitivity
stau – event-time sensitivity
tcl – total abundances for conservation laws
x_old – pre-event state
-
virtual void fdeltax(realtype *deltax, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, int ie, realtype const *xdot, realtype const *xdot_old, realtype const *x_old)
Model-specific implementation of fdeltax.
- Parameters:
deltax – state update
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
ie – event index
xdot – new model right hand side
xdot_old – previous model right hand side
x_old – pre-event state
-
virtual void fdeltaxB(realtype *deltaxB, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *dx, int ie, realtype const *xdot, realtype const *xdot_old, realtype const *x_old, realtype const *xB, realtype const *tcl)
Model-specific implementation of fdeltaxB.
- Parameters:
deltaxB – adjoint state update
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
dx – time derivative of state (DAE only)
ie – event index
xdot – new model right hand side
xdot_old – previous model right hand side
x_old – pre-event state
xB – current adjoint state
tcl – total abundances for conservation laws
-
virtual void fdJrzdsigma(realtype *dJrzdsigma, int iz, realtype const *p, realtype const *k, realtype const *rz, realtype const *sigmaz)
Model-specific implementation of fdJrzdsigma.
- Parameters:
dJrzdsigma – Sensitivity of event penalization Jrz w.r.t. standard deviation sigmaz
iz – event output index
p – parameter vector
k – constant vector
rz – model root output at timepoint
sigmaz – event measurement standard deviation at timepoint
-
virtual void fdJrzdz(realtype *dJrzdz, int iz, realtype const *p, realtype const *k, realtype const *rz, realtype const *sigmaz)
Model-specific implementation of fdJrzdz.
- Parameters:
dJrzdz – partial derivative of event penalization Jrz
iz – event output index
p – parameter vector
k – constant vector
rz – model root output at timepoint
sigmaz – event measurement standard deviation at timepoint
-
virtual void fdJydsigma(realtype *dJydsigma, int iy, realtype const *p, realtype const *k, realtype const *y, realtype const *sigmay, realtype const *my)
Model-specific implementation of fdJydsigma.
- Parameters:
dJydsigma – Sensitivity of time-resolved measurement negative log-likelihood Jy w.r.t. standard deviation sigmay
iy – output index
p – parameter vector
k – constant vector
y – model output at timepoint
sigmay – measurement standard deviation at timepoint
my – measurement at timepoint
-
virtual void fdJydy(realtype *dJydy, int iy, realtype const *p, realtype const *k, realtype const *y, realtype const *sigmay, realtype const *my)
Model-specific implementation of fdJydy.
- Parameters:
dJydy – partial derivative of time-resolved measurement negative log-likelihood Jy
iy – output index
p – parameter vector
k – constant vector
y – model output at timepoint
sigmay – measurement standard deviation at timepoint
my – measurement at timepoint
-
virtual void fdJydy_colptrs(SUNMatrixWrapper &dJydy, int index)
Model-specific implementation of fdJydy colptrs.
- Parameters:
dJydy – sparse matrix to which colptrs will be written
index – ytrue index
-
virtual void fdJydy_rowvals(SUNMatrixWrapper &dJydy, int index)
Model-specific implementation of fdJydy rowvals.
- Parameters:
dJydy – sparse matrix to which rowvals will be written
index –
ytrueindex
-
virtual void fdJzdsigma(realtype *dJzdsigma, int iz, realtype const *p, realtype const *k, realtype const *z, realtype const *sigmaz, realtype const *mz)
Model-specific implementation of fdJzdsigma.
- Parameters:
dJzdsigma – Sensitivity of event measurement negative log-likelihood Jz w.r.t. standard deviation sigmaz
iz – event output index
p – parameter vector
k – constant vector
z – model event output at timepoint
sigmaz – event measurement standard deviation at timepoint
mz – event measurement at timepoint
-
virtual void fdJzdz(realtype *dJzdz, int iz, realtype const *p, realtype const *k, realtype const *z, realtype const *sigmaz, realtype const *mz)
Model-specific implementation of fdJzdz.
- Parameters:
dJzdz – partial derivative of event measurement negative log-likelihood Jz
iz – event output index
p – parameter vector
k – constant vector
z – model event output at timepoint
sigmaz – event measurement standard deviation at timepoint
mz – event measurement at timepoint
-
virtual void fdrzdp(realtype *drzdp, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, int ip)
Model-specific implementation of fdrzdp.
- Parameters:
drzdp – partial derivative of root output rz w.r.t. model parameters p
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
ip – parameter index w.r.t. which the derivative is requested
-
virtual void fdrzdx(realtype *drzdx, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h)
Model-specific implementation of fdrzdx.
- Parameters:
drzdx – partial derivative of root output rz w.r.t. model states x
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
-
virtual void fdsigmaydp(realtype *dsigmaydp, realtype t, realtype const *p, realtype const *k, realtype const *y, int ip)
Model-specific implementation of fdsigmaydp.
- Parameters:
dsigmaydp – partial derivative of standard deviation of measurements
t – current time
p – parameter vector
k – constant vector
y – model output at timepoint t
ip – sensitivity index
-
virtual void fdsigmaydy(realtype *dsigmaydy, realtype t, realtype const *p, realtype const *k, realtype const *y)
Model-specific implementation of fsigmay.
- Parameters:
dsigmaydy – partial derivative of standard deviation of measurements w.r.t. model outputs
t – current time
p – parameter vector
k – constant vector
y – model output at timepoint t
-
virtual void fdsigmazdp(realtype *dsigmazdp, realtype t, realtype const *p, realtype const *k, int ip)
Model-specific implementation of fsigmaz.
- Parameters:
dsigmazdp – partial derivative of standard deviation of event measurements
t – current time
p – parameter vector
k – constant vector
ip – sensitivity index
-
virtual void fdtotal_cldp(realtype *dtotal_cldp, realtype const *x_rdata, realtype const *p, realtype const *k, int ip)
Compute dtotal_cl / dp.
- Parameters:
dtotal_cldp – dtotal_cl / dp
x_rdata – State variables with conservation laws applied
p – parameter vector
k – constant vector
ip – Sensitivity index
-
virtual void fdtotal_cldx_rdata(realtype *dtotal_cldx_rdata, realtype const *x_rdata, realtype const *p, realtype const *k, realtype const *tcl)
Compute dtotal_cl / dx_rdata.
- Parameters:
dtotal_cldx_rdata – dtotal_cl / dx_rdata
x_rdata – State variables with conservation laws applied
p – parameter vector
k – constant vector
tcl – Total abundances for conservation laws
-
virtual void fdtotal_cldx_rdata_colptrs(SUNMatrixWrapper &dtotal_cldx_rdata)
Model-specific implementation of fdtotal_cldx_rdata, colptrs part.
- Parameters:
dtotal_cldx_rdata – sparse matrix to which colptrs will be written
-
virtual void fdtotal_cldx_rdata_rowvals(SUNMatrixWrapper &dtotal_cldx_rdata)
Model-specific implementation of fdtotal_cldx_rdata, rowvals part.
- Parameters:
dtotal_cldx_rdata – sparse matrix to which rowvals will be written
-
virtual void fdwdp(realtype *dwdp, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *w, realtype const *tcl, realtype const *stcl, realtype const *spl, realtype const *sspl, bool include_static = true)
Model-specific sparse implementation of dwdp.
- Parameters:
dwdp – Recurring terms in xdot, parameter derivative
t – timepoint
x – vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
w – vector with helper variables
tcl – total abundances for conservation laws
stcl – sensitivities of total abundances for conservation laws
spl – spline value vector
sspl – sensitivities of spline values vector w.r.t. parameters \( p \)
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
virtual void fdwdp_colptrs(SUNMatrixWrapper &dwdp)
Model-specific implementation for dwdp, column pointers.
- Parameters:
dwdp – sparse matrix to which colptrs will be written
-
virtual void fdwdp_rowvals(SUNMatrixWrapper &dwdp)
Model-specific implementation for dwdp, row values.
- Parameters:
dwdp – sparse matrix to which rowvals will be written
-
virtual void fdwdw(realtype *dwdw, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *w, realtype const *tcl, bool include_static = true)
Model-specific implementation of fdwdw, no w chainrule.
- Parameters:
dwdw – partial derivative w wrt w
t – timepoint
x – Vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
w – vector with helper variables
tcl – Total abundances for conservation laws
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
virtual void fdwdw_colptrs(SUNMatrixWrapper &dwdw)
Model-specific implementation of fdwdw, colptrs part.
- Parameters:
dwdw – sparse matrix to which colptrs will be written
-
virtual void fdwdw_rowvals(SUNMatrixWrapper &dwdw)
Model-specific implementation of fdwdw, rowvals part.
- Parameters:
dwdw – sparse matrix to which rowvals will be written
-
virtual void fdwdx(realtype *dwdx, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *w, realtype const *tcl, realtype const *spl, bool include_static = true)
Model-specific implementation of dwdx, data part.
- Parameters:
dwdx – Recurring terms in xdot, state derivative
t – timepoint
x – vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
w – vector with helper variables
tcl – total abundances for conservation laws
spl – spline value vector
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
virtual void fdwdx_colptrs(SUNMatrixWrapper &dwdx)
Model-specific implementation for dwdx, column pointers.
- Parameters:
dwdx – sparse matrix to which colptrs will be written
-
virtual void fdwdx_rowvals(SUNMatrixWrapper &dwdx)
Model-specific implementation for dwdx, row values.
- Parameters:
dwdx – sparse matrix to which rowvals will be written
-
virtual void fdx_rdatadp(realtype *dx_rdatadp, realtype const *x, realtype const *tcl, realtype const *p, realtype const *k, int ip)
Compute dx_rdata / dp.
- Parameters:
dx_rdatadp – dx_rdata / dp
p – parameter vector
k – constant vector
x – State variables with conservation laws applied
tcl – Total abundances for conservation laws
ip – Sensitivity index
-
virtual void fdx_rdatadtcl(realtype *dx_rdatadtcl, realtype const *x, realtype const *tcl, realtype const *p, realtype const *k)
Compute dx_rdata / dtcl.
- Parameters:
dx_rdatadtcl – dx_rdata / dtcl
p – parameter vector
k – constant vector
x – State variables with conservation laws applied
tcl – Total abundances for conservation laws
-
virtual void fdx_rdatadtcl_colptrs(SUNMatrixWrapper &dx_rdatadtcl)
Model-specific implementation of fdx_rdatadtcl, colptrs part.
- Parameters:
dx_rdatadtcl – sparse matrix to which colptrs will be written
-
virtual void fdx_rdatadtcl_rowvals(SUNMatrixWrapper &dx_rdatadtcl)
Model-specific implementation of fdx_rdatadtcl, rowvals part.
- Parameters:
dx_rdatadtcl – sparse matrix to which rowvals will be written
-
virtual void fdx_rdatadx_solver(realtype *dx_rdatadx_solver, realtype const *x, realtype const *tcl, realtype const *p, realtype const *k)
Compute dx_rdata / dx_solver.
- Parameters:
dx_rdatadx_solver – dx_rdata / dx_solver
p – parameter vector
k – constant vector
x – State variables with conservation laws applied
tcl – Total abundances for conservation laws
-
virtual void fdx_rdatadx_solver_colptrs(SUNMatrixWrapper &dxrdatadxsolver)
Model-specific implementation of fdx_rdatadx_solver, colptrs part.
- Parameters:
dxrdatadxsolver – sparse matrix to which colptrs will be written
-
virtual void fdx_rdatadx_solver_rowvals(SUNMatrixWrapper &dxrdatadxsolver)
Model-specific implementation of fdx_rdatadx_solver, rowvals part.
- Parameters:
dxrdatadxsolver – sparse matrix to which rowvals will be written
-
virtual void fdydp(realtype *dydp, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, int ip, realtype const *w, realtype const *tcl, realtype const *dtcldp, realtype const *spl, realtype const *sspl)
Model-specific implementation of fdydp.
- Parameters:
dydp – partial derivative of observables y w.r.t. model parameters p
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
ip – parameter index w.r.t. which the derivative is requested
w – repeating elements vector
tcl – total abundances for conservation laws
dtcldp – Sensitivities of total abundances for conservation laws
spl – spline value vector
sspl – sensitivities of spline values vector w.r.t. parameters \( p \)
-
virtual void fdydx(realtype *dydx, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *w)
Model-specific implementation of fdydx.
- Parameters:
dydx – partial derivative of observables y w.r.t. model states x
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
w – repeating elements vector
-
virtual void fdzdp(realtype *dzdp, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, int ip)
Model-specific implementation of fdzdp.
- Parameters:
dzdp – partial derivative of event-resolved output z w.r.t. model parameters p
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
ip – parameter index w.r.t. which the derivative is requested
-
virtual void fdzdx(realtype *dzdx, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h)
Model-specific implementation of fdzdx.
- Parameters:
dzdx – partial derivative of event-resolved output z w.r.t. model states x
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
-
virtual void fJrz(realtype *nllh, int iz, realtype const *p, realtype const *k, realtype const *z, realtype const *sigmaz)
Model-specific implementation of fJrz.
- Parameters:
nllh – regularization for event measurements z
iz – event output index
p – parameter vector
k – constant vector
z – model event output at timepoint
sigmaz – event measurement standard deviation at timepoint
-
virtual void fJy(realtype *nllh, int iy, realtype const *p, realtype const *k, realtype const *y, realtype const *sigmay, realtype const *my)
Model-specific implementation of fJy.
- Parameters:
nllh – negative log-likelihood for measurements y
iy – output index
p – parameter vector
k – constant vector
y – model output at timepoint
sigmay – measurement standard deviation at timepoint
my – measurements at timepoint
-
virtual void fJz(realtype *nllh, int iz, realtype const *p, realtype const *k, realtype const *z, realtype const *sigmaz, realtype const *mz)
Model-specific implementation of fJz.
- Parameters:
nllh – negative log-likelihood for event measurements z
iz – event output index
p – parameter vector
k – constant vector
z – model event output at timepoint
sigmaz – event measurement standard deviation at timepoint
mz – event measurements at timepoint
-
virtual void frz(realtype *rz, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h)
Model-specific implementation of frz.
- Parameters:
rz – value of root function at current timepoint (non-output events not included)
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
-
virtual void fsigmay(realtype *sigmay, realtype t, realtype const *p, realtype const *k, realtype const *y)
Model-specific implementation of fsigmay.
- Parameters:
sigmay – standard deviation of measurements
t – current time
p – parameter vector
k – constant vector
y – model output at timepoint t
-
virtual void fsigmaz(realtype *sigmaz, realtype t, realtype const *p, realtype const *k)
Model-specific implementation of fsigmaz.
- Parameters:
sigmaz – standard deviation of event measurements
t – current time
p – parameter vector
k – constant vector
-
virtual void fsrz(realtype *srz, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *sx, int ip)
Model-specific implementation of fsrz.
- Parameters:
srz – Sensitivity of rz, total derivative
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
sx – current state sensitivity
h – Heaviside vector
ip – sensitivity index
-
virtual void fstau(realtype *stau, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *dx, realtype const *tcl, realtype const *sx, int ip, int ie)
Model-specific implementation of fstau.
- Parameters:
stau – total derivative of event timepoint
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
dx – time derivative of state (DAE only)
tcl – total abundances for conservation laws
sx – current state sensitivity
ip – sensitivity index
ie – event index
-
virtual void fsx0(realtype *sx0, realtype t, realtype const *x0, realtype const *p, realtype const *k, int ip)
Model-specific implementation of fsx0.
- Parameters:
sx0 – initial state sensitivities
t – initial time
x0 – initial state
p – parameter vector
k – constant vector
ip – sensitivity index
-
virtual void fsx0_fixedParameters(realtype *sx0, realtype t, realtype const *x0, realtype const *p, realtype const *k, int ip, gsl::span<int const> reinitialization_state_idxs)
Model-specific implementation of fsx0_fixedParameters.
- Parameters:
sx0 – initial state sensitivities
t – initial time
x0 – initial state
p – parameter vector
k – constant vector
ip – sensitivity index
reinitialization_state_idxs – Indices of states to be reinitialized based on provided constants / fixed parameters.
-
virtual void fsz(realtype *sz, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *sx, int ip)
Model-specific implementation of fsz.
- Parameters:
sz – Sensitivity of rz, total derivative
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
sx – current state sensitivity
ip – sensitivity index
-
virtual void fw(realtype *w, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *tcl, realtype const *spl, bool include_static = true)
Model-specific implementation of fw.
- Parameters:
w – Recurring terms in xdot
t – timepoint
x – vector with the states
p – parameter vector
k – constants vector
h – Heaviside vector
tcl – total abundances for conservation laws
spl – spline value vector
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
virtual void fx0(realtype *x0, realtype t, realtype const *p, realtype const *k)
Model-specific implementation of fx0.
- Parameters:
x0 – initial state
t – initial time
p – parameter vector
k – constant vector
-
virtual void fx0_fixedParameters(realtype *x0, realtype t, realtype const *p, realtype const *k, gsl::span<int const> reinitialization_state_idxs)
Model-specific implementation of fx0_fixedParameters.
- Parameters:
x0 – initial state
t – initial time
p – parameter vector
k – constant vector
reinitialization_state_idxs – Indices of states to be reinitialized based on provided constants / fixed parameters.
-
virtual void fy(realtype *y, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h, realtype const *w)
Model-specific implementation of fy.
- Parameters:
y – model output at current timepoint
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
w – repeating elements vector
-
virtual void fz(realtype *z, int ie, realtype t, realtype const *x, realtype const *p, realtype const *k, realtype const *h)
Model-specific implementation of fz.
- Parameters:
z – value of event output
ie – event index
t – current time
x – current state
p – parameter vector
k – constant vector
h – Heaviside vector
Protected Functions
-
void write_slice_event(gsl::span<realtype const> slice, gsl::span<realtype> buffer, int ie)
Write part of a slice to a buffer according to indices specified in z2event.
- Parameters:
slice – Input data slice
buffer – Output data slice
ie – Event index
-
void write_sensitivity_slice_event(gsl::span<realtype const> slice, gsl::span<realtype> buffer, int ie)
Write part of a sensitivity slice to a buffer according to indices specified in z2event.
- Parameters:
slice – source data slice
buffer – output data slice
ie – event index
-
void write_llh_sensitivity_slice(std::vector<realtype> const &dLLhdp, std::vector<realtype> &sllh, std::vector<realtype> &s2llh)
Separate first and second order objective sensitivity information and write them into the respective buffers.
- Parameters:
dLLhdp – Data with mangled first- and second-order information
sllh – First order buffer
s2llh – Second order buffer
-
void check_llh_buffer_size(std::vector<realtype> const &sllh, std::vector<realtype> const &s2llh) const
Verify that the provided buffers have the expected size.
- Parameters:
sllh – first order buffer
s2llh – second order buffer
-
void initialize_vectors()
Set the nplist-dependent vectors to their proper sizes.
-
void fy(realtype t, AmiVector const &x)
Compute observables / measurements.
- Parameters:
t – Current timepoint
x – Current state
-
void fdydp(realtype t, AmiVector const &x)
Compute partial derivative of observables \( y \) w.r.t. model parameters
p.- Parameters:
t – Current timepoint
x – Current state
-
void fdydx(realtype t, AmiVector const &x)
Compute partial derivative of observables \( y \) w.r.t. state variables
x.- Parameters:
t – Current timepoint
x – Current state
-
void fsigmay(int it, ExpData const *edata)
Compute standard deviation of measurements.
- Parameters:
it – Timepoint index
edata – Experimental data
-
void fdsigmaydp(int it, ExpData const *edata)
Compute partial derivative of standard deviation of measurements w.r.t. model parameters.
- Parameters:
it – Timepoint index
edata – pointer to
amici::ExpDatadata instance holding sigma values
-
void fdsigmaydy(int it, ExpData const *edata)
Compute partial derivative of standard deviation of measurements w.r.t. model outputs.
- Parameters:
it – Timepoint index
edata – pointer to
amici::ExpDatadata instance holding sigma values
-
void fdJydy(int it, AmiVector const &x, ExpData const &edata)
Compute partial derivative of time-resolved measurement negative log-likelihood \( Jy \).
- Parameters:
it – timepoint index
x – state variables
edata – Pointer to experimental data
-
void fdJydsigma(int it, AmiVector const &x, ExpData const &edata)
Sensitivity of time-resolved measurement negative log-likelihood Jy w.r.t. standard deviation sigma.
- Parameters:
it – timepoint index
x – state variables
edata – pointer to experimental data instance
-
void fdJydp(int it, AmiVector const &x, ExpData const &edata)
Compute sensitivity of time-resolved measurement negative log-likelihood \( Jy \) w.r.t. parameters for the given timepoint.
- Parameters:
it – timepoint index
x – state variables
edata – pointer to experimental data instance
-
void fdJydx(int it, AmiVector const &x, ExpData const &edata)
Sensitivity of time-resolved measurement negative log-likelihood \( Jy \) w.r.t. state variables.
- Parameters:
it – Timepoint index
x – State variables
edata – Pointer to experimental data instance
-
void fz(int ie, realtype t, AmiVector const &x)
Compute event-resolved output.
- Parameters:
ie – Event index
t – Current timepoint
x – Current state
-
void fdzdp(int ie, realtype t, AmiVector const &x)
Compute partial derivative of event-resolved output
zw.r.t. model parametersp- Parameters:
ie – event index
t – current timepoint
x – current state
-
void fdzdx(int ie, realtype t, AmiVector const &x)
Compute partial derivative of event-resolved output
zw.r.t. model statesx.- Parameters:
ie – Event index
t – Current timepoint
x – Current state
-
void frz(int ie, realtype t, AmiVector const &x)
Compute event root function of events.
Equal to
Model::frootbut does not include non-output events.- Parameters:
ie – Event index
t – Current timepoint
x – Current state
-
void fdrzdp(int ie, realtype t, AmiVector const &x)
Compute sensitivity of event-resolved root output w.r.t. model parameters
p.- Parameters:
ie – Event index
t – Current timepoint
x – Current state
-
void fdrzdx(int ie, realtype t, AmiVector const &x)
Compute sensitivity of event-resolved measurements \( rz \) w.r.t. model states
x.- Parameters:
ie – Event index
t – Current timepoint
x – Current state
-
void fsigmaz(int ie, int nroots, realtype t, ExpData const *edata)
Compute standard deviation of events.
-
void fdsigmazdp(int ie, int nroots, realtype t, ExpData const *edata)
Compute sensitivity of standard deviation of events measurements w.r.t. model parameters
p.
-
void fdJzdz(int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Compute partial derivative of event measurement negative log-likelihood \( Jz \).
-
void fdJzdsigma(int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Compute sensitivity of event measurement negative log-likelihood \( Jz \) w.r.t. standard deviation sigmaz.
-
void fdJzdp(int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Compute sensitivity of event-resolved measurement negative log-likelihood Jz w.r.t. parameters.
-
void fdJzdx(int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Compute sensitivity of event-resolved measurement negative log-likelihood Jz w.r.t. state variables.
-
void fdJrzdz(int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Compute partial derivative of event measurement negative log-likelihood J.
-
void fdJrzdsigma(int ie, int nroots, realtype t, AmiVector const &x, ExpData const &edata)
Compute sensitivity of event measurement negative log-likelihood Jz w.r.t. standard deviation sigmaz.
- Parameters:
ie – event index
nroots – event index
t – current timepoint
x – state variables
edata – pointer to experimental data instance
-
void fw(realtype t, realtype const *x, bool include_static = true)
Compute recurring terms in xdot.
- Parameters:
t – Timepoint
x – Array with the states
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
void fdwdp(realtype t, realtype const *x, bool include_static = true)
Compute parameter derivative for recurring terms in xdot.
- Parameters:
t – Timepoint
x – Array with the states
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
void fdwdx(realtype t, realtype const *x, bool include_static = true)
Compute state derivative for recurring terms in xdot.
- Parameters:
t – Timepoint
x – Array with the states
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
void fdwdw(realtype t, realtype const *x, bool include_static = true)
Compute self derivative for recurring terms in xdot.
- Parameters:
t – Timepoint
x – Array with the states
include_static – Whether to (re-)evaluate only dynamic expressions (false) or also static expressions (true). Dynamic expressions are those that depend directly or indirectly on time, static expressions are those that don’t.
-
virtual void fx_rdata(realtype *x_rdata, realtype const *x_solver, realtype const *tcl, realtype const *p, realtype const *k)
Compute fx_rdata.
To be implemented by derived class if applicable.
- Parameters:
x_rdata – State variables with conservation laws expanded
x_solver – State variables with conservation laws applied
tcl – Total abundances for conservation laws
p – parameter vector
k – constant vector
-
virtual void fsx_rdata(realtype *sx_rdata, realtype const *sx_solver, realtype const *stcl, realtype const *p, realtype const *k, realtype const *x_solver, realtype const *tcl, int ip)
Compute fsx_solver.
To be implemented by derived class if applicable.
- Parameters:
sx_rdata – State sensitivity variables with conservation laws expanded
sx_solver – State sensitivity variables with conservation laws applied
stcl – Sensitivities of total abundances for conservation laws
p – parameter vector
k – constant vector
x_solver – State variables with conservation laws applied
tcl – Total abundances for conservation laws
ip – Sensitivity index
-
virtual void fx_solver(realtype *x_solver, realtype const *x_rdata)
Compute fx_solver.
To be implemented by derived class if applicable.
- Parameters:
x_solver – State variables with conservation laws applied
x_rdata – State variables with conservation laws expanded
-
virtual void fsx_solver(realtype *sx_solver, realtype const *sx_rdata)
Compute fsx_solver.
To be implemented by derived class if applicable.
- Parameters:
sx_rdata – State sensitivity variables with conservation laws expanded
sx_solver – State sensitivity variables with conservation laws applied
-
virtual void ftotal_cl(realtype *total_cl, realtype const *x_rdata, realtype const *p, realtype const *k)
Compute ftotal_cl.
To be implemented by derived class if applicable.
- Parameters:
total_cl – Total abundances of conservation laws
x_rdata – State variables with conservation laws expanded
p – parameter vector
k – constant vector
-
virtual void fstotal_cl(realtype *stotal_cl, realtype const *sx_rdata, int ip, realtype const *x_rdata, realtype const *p, realtype const *k, realtype const *tcl)
Compute fstotal_cl.
To be implemented by derived class if applicable.
- Parameters:
stotal_cl – Sensitivities for the total abundances of conservation laws
sx_rdata – State sensitivity variables with conservation laws expanded
ip – Sensitivity index
x_rdata – State variables with conservation laws expanded
p – parameter vector
k – constant vector
tcl – Total abundances for conservation laws
-
const_N_Vector compute_x_pos(const_N_Vector x)
Compute non-negative state vector.
Compute non-negative state vector according to stateIsNonNegative. If anyStateNonNegative is set to
false, i.e., all entries in stateIsNonNegative arefalse, this function directly returnsx, otherwise all entries of x are copied in toamici::Model::x_pos_tmp_and negative values are replaced by0where applicable.- Parameters:
x – State vector possibly containing negative values
- Returns:
State vector with negative values replaced by
0according to stateIsNonNegative
-
realtype const *compute_x_pos(AmiVector const &x)
Compute non-negative state vector.
Compute non-negative state vector according to stateIsNonNegative. If anyStateNonNegative is set to
false, i.e., all entries in stateIsNonNegative arefalse, this function directly returnsx, otherwise all entries of x are copied in toamici::Model::x_pos_tmp_and negative values are replaced by0where applicable.- Parameters:
x – State vector possibly containing negative values
- Returns:
State vector with negative values replaced by
0according to stateIsNonNegative
Protected Attributes
-
ModelState state_
All variables necessary for function evaluation
-
ModelStateDerived derived_state_
Storage for model quantities beyond ModelState for the current timepoint
-
std::vector<HermiteSpline> splines_
Storage for splines of the model
-
std::vector<int> z2event_
index indicating to which event an event output belongs (size nz)
-
std::vector<bool> state_is_non_negative_
vector of bools indicating whether state variables are to be assumed to be positive
-
bool any_state_non_negative_ = {false}
boolean indicating whether any entry in stateIsNonNegative is
true
-
int nmaxevent_ = {10}
maximal number of events to track
-
SteadyStateComputationMode steadystate_computation_mode_{SteadyStateComputationMode::integrationOnly}
method for steady-state computation
-
SteadyStateSensitivityMode steadystate_sensitivity_mode_{SteadyStateSensitivityMode::integrationOnly}
method for steadystate sensitivities computation
-
bool always_check_finite_ = {true}
Indicates whether the result of every call to
Model::f*should be checked for finiteness
-
bool sigma_res_ = {false}
indicates whether sigma residuals are to be added for every datapoint
-
Model() = default