Class ReturnData

Inheritance Relationships

Base Type

Class Documentation

class ReturnData : public amici::ModelDimensions

Stores all data to be returned by amici::runAmiciSimulation.

NOTE: multi-dimensional arrays are stored in row-major order (C-style)

Public Functions

ReturnData() = default

Default constructor.

ReturnData(std::vector<realtype> ts_, ModelDimensions const &model_dimensions_, int nmaxevent_, int newton_maxsteps_, std::vector<int> plist_, std::vector<ParameterScaling> pscale_, SecondOrderMode o2mode_, SensitivityOrder sensi_, SensitivityMethod sensi_meth_, RDataReporting rdrm_, bool quadratic_llh_, bool sigma_res_, realtype sigma_offset_)

Constructor.

Parameters:
  • ts_ – see amici::SimulationParameters::ts

  • model_dimensions_Model dimensions

  • nmaxevent_ – see amici::ModelDimensions::nmaxevent

  • newton_maxsteps_ – see amici::Solver::newton_maxsteps

  • plist_ – see amici::Model::getParameterList

  • pscale_ – see amici::SimulationParameters::pscale

  • o2mode_ – see amici::SimulationParameters::o2mode

  • sensi_ – see amici::Solver::sensi

  • sensi_meth_ – see amici::Solver::sensi_meth

  • rdrm_ – see amici::Solver::rdata_reporting

  • quadratic_llh_ – whether model defines a quadratic nllh and computing res, sres and FIM makes sense

  • sigma_res_ – indicates whether additional residuals are to be added for each sigma

  • sigma_offset_ – offset to ensure real-valuedness of sigma residuals

ReturnData(Solver const &solver, Model const &model)

constructor that uses information from model and solver to appropriately initialize fields

Parameters:
  • solver – solver instance

  • model – model instance

~ReturnData() = default
void process_simulation_objects(ForwardProblem const *fwd, BackwardProblem const *bwd, Model &model, Solver const &solver, ExpData const *edata)

constructor that uses information from model and solver to appropriately initialize fields

Parameters:
  • fwd – simulated forward problem, pass nullptr to ignore

  • bwd – simulated backward problem, pass nullptr to ignore

  • model – matching model instance

  • solver – matching solver instance

  • edata – matching experimental data

Public Members

std::string id

Arbitrary (not necessarily unique) identifier.

std::vector<realtype> ts

Output or measurement timepoints (shape nt)

std::vector<realtype> xdot

time derivative (shape nx_solver) evaluated at t_last.

std::vector<realtype> J

Jacobian of differential equation right hand side.

The Jacobian of differential equation right hand side (shape nx_solver x nx_solver, row-major) evaluated at t_last.

The corresponding state variable IDs can be obtained from Model::getStateIdsSolver().

std::vector<realtype> w

Model expression values.

Values of model expressions (recurring terms in xdot, for imported SBML models from Python, this contains, e.g., the flux vector) at timepoints ReturnData::ts (shape nt x nw, row major).

The corresponding expression IDs can be obtained from Model::getExpressionIds().

std::vector<realtype> z

Event output (shape nmaxevent x nz, row-major)

std::vector<realtype> sigmaz

Event output sigma standard deviation (shape nmaxevent x nz, row-major)

std::vector<realtype> sz

Parameter derivative of event output (shape nmaxevent x nplist x nz, row-major)

std::vector<realtype> ssigmaz

Parameter derivative of event output standard deviation (shape nmaxevent x nplist x nz, row-major)

std::vector<realtype> rz

Event trigger output (shape nmaxevent x nz, row-major)

std::vector<realtype> srz

Parameter derivative of event trigger output (shape nmaxevent x nplist x nz, row-major)

std::vector<realtype> s2rz

Second-order parameter derivative of event trigger output (shape nmaxevent x nztrue x nplist x nplist, row-major)

std::vector<realtype> x

Model state.

The model state at timepoints ReturnData::ts (shape nt x nx_rdata, row-major).

The corresponding state variable IDs can be obtained from Model::getStateIds().

std::vector<realtype> sx

State sensitivities.

The derivative of the model state with respect to the chosen parameters (see Model::getParameterList() or ExpData::plist) at timepoints ReturnData::ts (shape nt x nplist x nx_rdata, row-major).

The corresponding state variable IDs can be obtained from Model::getStateIds().

std::vector<realtype> y

Observables.

The values of the observables at timepoints ReturnData::ts (shape nt x ny, row-major).

The corresponding observable IDs can be obtained from Model::getObservableIds().

std::vector<realtype> sigmay

Observable standard deviation (shape nt x ny, row-major)

std::vector<realtype> sy

Observable sensitivities.

The derivative of the observables with respect to the chosen parameters (see Model::getParameterList() or ExpData::plist) at timepoints ReturnData::ts (shape nt x nplist x ny, row-major).

The corresponding observable IDs can be obtained from Model::getObservableIds().

std::vector<realtype> ssigmay

Parameter derivative of observable standard deviation (shape nt x nplist x ny, row-major)

std::vector<realtype> res

Residuals (shape nt*ny, row-major)

std::vector<realtype> sres

Parameter derivative of residual (shape nt*ny x nplist, row-major)

std::vector<realtype> FIM

Fisher information matrix (shape nplist x nplist, row-major)

std::vector<int> numsteps

Number of solver steps for the forward problem.

Cumulative number of integration steps for the forward problem for each output timepoint in ReturnData::ts (shape nt).

std::vector<int> numsteps_b

Number of solver steps for the backward problem.

Cumulative number of integration steps for the backward problem for each output timepoint in ReturnData::ts (shape nt).

std::vector<int> num_rhs_evals

Number of right hand side evaluations forward problem (shape nt)

std::vector<int> num_rhs_evals_b

Number of right hand side evaluations backward problem (shape nt)

std::vector<int> num_err_test_fails

Number of error test failures forward problem (shape nt)

std::vector<int> num_err_test_fails_b

Number of error test failures backward problem (shape nt)

std::vector<int> num_non_lin_solv_conv_fails

Number of linear solver convergence failures forward problem (shape nt)

std::vector<int> num_non_lin_solv_conv_fails_b

Number of linear solver convergence failures backward problem (shape nt)

std::vector<int> order

Employed order forward problem (shape nt)

double cpu_time = 0.0

Computation time of forward solve [ms].

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

double cpu_time_b = 0.0

Computation time of backward solve [ms].

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

double cpu_time_total = 0.0

Total CPU time from entering runAmiciSimulation until exiting [ms].

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

std::vector<SteadyStateStatus> preeq_status

Flags indicating success of steady-state solver (preequilibration)

double preeq_cpu_time = 0.0

Computation time of the steady state solver [ms] (pre-equilibration)

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

double preeq_cpu_time_b = 0.0

Computation time of the steady state solver of the backward problem [ms] (pre-equilibration)

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

std::vector<SteadyStateStatus> posteq_status

Flags indicating success of steady-state solver (postequilibration)

double posteq_cpu_time = 0.0

Computation time of the steady-state solver [ms] (postequilibration)

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

double posteq_cpu_time_b = 0.0

Computation time of the steady-state solver of the backward problem [ms] (postequilibration)

.. warning:: If AMICI was built without boost, this tracks the CPU-time of the current process. Therefore, in a multi-threaded context, this value may be incorrect.

std::vector<int> preeq_numsteps

Number of Newton steps for pre-equilibration.

[newton, simulation, newton] (length = 3)

int preeq_numsteps_b = 0

Number of simulation steps for adjoint pre-equilibration problem [== 0 if analytical solution worked, > 0 otherwise]

std::vector<int> posteq_numsteps

Number of Newton steps for post-equilibration [newton, simulation, newton].

int posteq_numsteps_b = 0

Number of simulation steps for the post-equilibration backward simulation [== 0 if analytical solution worked, > 0 otherwise]

realtype preeq_t = NAN

Model time at which the pre-equilibration steady state was reached via simulation.

realtype preeq_wrms = NAN

Weighted root-mean-square of the rhs at pre-equilibration steady state.

realtype posteq_t = NAN

Model time at which the post-equilibration steady state was reached via simulation.

realtype posteq_wrms = NAN

Weighted root-mean-square of the rhs at post-equilibration steady state.

std::vector<realtype> x0

Initial state of the main simulation (shape nx_rdata).

The corresponding state variable IDs can be obtained from Model::getStateIds().

std::vector<realtype> x_ss

Pre-equilibration steady state.

The values of the state variables at the pre-equilibration steady state (shape nx_rdata). The corresponding state variable IDs can be obtained from Model::getStateIds().

std::vector<realtype> sx0

Initial state sensitivities for the main simulation.

(shape nplist x nx_rdata, row-major).

std::vector<realtype> sx_ss

Pre-equilibration steady state sensitivities.

Sensitivities of the pre-equilibration steady state with respect to the selected parameters. (shape nplist x nx_rdata, row-major)

realtype llh = 0.0

Log-likelihood value

realtype chi2 = 0.0

\(\chi^2\) value

std::vector<realtype> sllh

Parameter derivative of log-likelihood (shape nplist)

std::vector<realtype> s2llh

Second-order parameter derivative of log-likelihood (shape nJ-1 x nplist, row-major)

int status = AMICI_NOT_RUN

Simulation status code.

One of:

  • AMICI_SUCCESS, indicating successful simulation

  • AMICI_MAX_TIME_EXCEEDED, indicating that the simulation did not finish within the allowed time (see Solver.{set,get}MaxTime)

  • AMICI_ERROR, indicating that some error occurred during simulation (a more detailed error message will have been printed).

  • AMICI_NOT_RUN, if no simulation was started

int nplist = {0}

Number of parameters w.r.t. which sensitivities were requested

int nmaxevent = {0}

Maximal number of occurring events (for every event type)

int nt = {0}

Number of output timepoints (length of ReturnData::ts).

int newton_maxsteps = {0}

maximal number of newton iterations for steady state calculation

std::vector<ParameterScaling> pscale

Scaling of model parameters.

SecondOrderMode o2mode = {SecondOrderMode::none}

Flag indicating whether second-order sensitivities were requested.

SensitivityOrder sensi = {SensitivityOrder::none}

Sensitivity order.

SensitivityMethod sensi_meth = {SensitivityMethod::none}

Sensitivity method.

RDataReporting rdata_reporting = {RDataReporting::full}

Reporting mode.

bool sigma_res = {false}

Boolean indicating whether residuals for standard deviations have been added.

std::vector<LogItem> messages

Log messages.

realtype t_last = {std::numeric_limits<realtype>::quiet_NaN()}

The final internal time of the solver.

std::vector<int> plist

Indices of the parameters w.r.t. which sensitivities were computed.

The indices refer to parameter IDs in Model::getParameterIds().

Protected Functions

void initialize_likelihood_reporting(bool quadratic_llh)

initializes storage for likelihood reporting mode

Parameters:

quadratic_llh – whether model defines a quadratic nllh and computing res, sres and FIM makes sense.

void initializeObservablesLikelihoodReporting(bool quadratic_llh)

initializes storage for observables + likelihood reporting mode

Parameters:

quadratic_llh – whether model defines a quadratic nllh and computing res, sres and FIM makes sense.

void initialize_residual_reporting(bool enable_res)

initializes storage for residual reporting mode

Parameters:

enable_res – whether residuals are to be computed

void initialize_full_reporting(bool enable_fim)

initializes storage for full reporting mode

Parameters:

enable_fim – whether FIM Hessian approximation is to be computed

void initialize_objective_function(bool enable_chi2)

initialize values for chi2 and llh and derivatives

Parameters:

enable_chi2 – whether chi2 values are to be computed

void process_pre_equilibration(SteadyStateProblem const &preeq, SteadyStateBackwardProblem const *preeq_bwd, Model &model)

extracts data from a pre-equilibration SteadystateProblem

Parameters:
  • preeq – SteadystateProblem for pre-equilibration

  • preeq_bwdSteadyStateBackwardProblem from pre-equilibration

  • modelModel instance to compute return values

void process_post_equilibration(SteadyStateProblem const &posteq, SteadyStateBackwardProblem const *posteq_bwd, Model &model, ExpData const *edata)

extracts data from a post-equilibration SteadystateProblem

Parameters:
  • posteq – SteadystateProblem for post-equilibration

  • posteq_bwdSteadyStateBackwardProblem from post-equilibration

  • modelModel instance to compute return values

  • edataExpData instance containing observable data

void process_forward_problem(ForwardProblem const &fwd, Model &model, ExpData const *edata)

extracts results from forward problem

Parameters:
  • fwd – forward problem

  • model – model that was used for forward simulation

  • edataExpData instance containing observable data

void process_backward_problem(ForwardProblem const &fwd, BackwardProblem const &bwd, SteadyStateProblem const *preeq, SteadyStateBackwardProblem const *preeq_bwd, Model &model)

extracts results from backward problem

Parameters:
  • fwd – forward problem

  • bwd – backward problem

  • preeqSteadyStateProblem for preequilibration

  • preeq_bwdSteadyStateBackwardProblem for preequilibration

  • model – model that was used for forward/backward simulation

void process_solver(Solver const &solver)

extracts results from solver

Parameters:

solver – solver that was used for forward/backward simulation

template<class T>
inline void store_jacobian_and_derivative(T const &problem, Model &model)

Evaluates and stores the Jacobian and right hand side at final timepoint.

Parameters:
  • problem – forward problem or steadystate problem

  • model – model that was used for forward/backward simulation

void fres(int it, Model &model, SolutionState const &sol, ExpData const &edata)

Residual function.

Parameters:
  • it – time index

  • model – model that was used for forward/backward simulation

  • sol – Solution state the timepoint it

  • edataExpData instance containing observable data

void fchi2(int it, ExpData const &edata)

Chi-squared function.

Parameters:
  • it – time index

  • edataExpData instance containing observable data

void fsres(int it, Model &model, SolutionState const &sol, ExpData const &edata)

Residual sensitivity function.

Parameters:
  • it – time index

  • model – model that was used for forward/backward simulation

  • sol – solution state the timepoint it

  • edataExpData instance containing observable data

void fFIM(int it, Model &model, SolutionState const &sol, ExpData const &edata)

Fisher information matrix function.

Parameters:
  • it – time index

  • model – model that was used for forward/backward simulation

  • sol – Solution state the timepoint it

  • edataExpData instance containing observable data

void invalidate(int it_start)

Set likelihood, state variables, outputs and respective sensitivities to NaN (typically after integration failure)

Parameters:

it_start – time index at which to start invalidating

void invalidate_llh()

Set likelihood and chi2 to NaN (typically after integration failure)

void invalidate_sllh()

Set likelihood sensitivities to NaN (typically after integration failure)

void apply_chain_rule_factor_to_simulation_results(Model const &model)

applies the chain rule to account for parameter transformation in the sensitivities of simulation results

Parameters:

modelModel from which the ReturnData was obtained

inline bool computing_fsa() const

Checks whether forward sensitivity analysis is performed.

Returns:

boolean indicator

void get_data_output(int it, Model &model, SolutionState const &sol, ExpData const *edata)

Extracts output information for data-points, expects that the model state was set appropriately.

Parameters:
  • it – timepoint index

  • model – model that was used in forward solve

  • sol – solution state the timepoint it

  • edataExpData instance carrying experimental data

void get_data_sensis_fsa(int it, Model &model, SolutionState const &sol, ExpData const *edata)

Extracts data information for forward sensitivity analysis, expects that the model state was set appropriately.

Parameters:
  • it – index of current timepoint

  • model – model that was used in forward solve

  • sol – Solution state the timepoint it

  • edataExpData instance carrying experimental data

void get_event_output(std::vector<int> const &rootidx, Model &model, SolutionState const &sol, ExpData const *edata)

Extracts output information for events, expects that the model state was set appropriately.

Parameters:
  • rootidx – information about which roots fired (1 indicating fired, 0/-1 for not)

  • model – model that was used in forward solve

  • sol – Solution state the timepoint it

  • edataExpData instance carrying experimental data

void get_event_sensis_fsa(int ie, Model &model, SolutionState const &sol, ExpData const *edata)

Extracts event information for forward sensitivity analysis, expects the model state was set appropriately.

Parameters:
  • ie – index of event type

  • model – model that was used in forward solve

  • sol – Solution state the timepoint it

  • edataExpData instance carrying experimental data

void handle_sx0_backward(Model const &model, AmiVectorArray const &sx0, AmiVector const &xB, std::vector<realtype> &llhS0) const

Updates contribution to likelihood from quadratures (xQB), if preequilibration was run in adjoint mode.

Parameters:
  • model – model that was used for forward/backward simulation

  • sx0 – State sensitivities at pre-equilibration steady state

  • xB – Adjoint state from pre-equilibration.

  • llhS0 – contribution to likelihood for initial state sensitivities of preequilibration

void handle_sx0_forward(Model const &model, SolutionState const &sol, std::vector<realtype> &llhS0, AmiVector const &xB) const

Updates contribution to likelihood for initial state sensitivities (llhS0), if no preequilibration was run or if forward sensitivities were used.

Parameters:
  • model – model that was used for forward/backward simulation

  • sol – Solution state the timepoint it

  • llhS0 – contribution to likelihood for initial state sensitivities

  • xB – vector with final adjoint state (excluding conservation laws)

Protected Attributes

realtype sigma_offset = {0.0}

offset for sigma_residuals

std::vector<int> nroots_

array of number of found roots for a certain event type (shape ne)

Friends

template<class Archive>
friend void serialize(Archive &ar, ReturnData &r, unsigned int version)

Serialize ReturnData (see boost::serialization::serialize)

Parameters:
  • ar – Archive to serialize to

  • r – Data to serialize

  • version – Version number