Class ReturnData
Defined in File rdata.h
Inheritance Relationships
Base Type
public amici::ModelDimensions
(Struct ModelDimensions)
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 nplist, int nmaxevent, int nt, int newton_maxsteps, 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
nplist – see amici::ModelDimensions::nplist
nmaxevent – see amici::ModelDimensions::nmaxevent
nt – see amici::ModelDimensions::nt
newton_maxsteps – see amici::Solver::newton_maxsteps
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 processSimulationObjects(SteadystateProblem const *preeq, ForwardProblem const *fwd, BackwardProblem const *bwd, SteadystateProblem const *posteq, Model &model, Solver const &solver, ExpData const *edata)
constructor that uses information from model and solver to appropriately initialize fields
- Parameters:
preeq – simulated preequilibration problem, pass
nullptr
to ignorefwd – simulated forward problem, pass
nullptr
to ignorebwd – simulated backward problem, pass
nullptr
to ignoreposteq – simulated postequilibration problem, pass
nullptr
to ignoremodel – matching model instance
solver – matching solver instance
edata – matching experimental data
Public Members
-
std::string id
Arbitrary (not necessarily unique) identifier.
-
std::vector<realtype> J
Jacobian of differential equation right hand side (shape
nx_solver
xnx_solver
, row-major) evaluated att_last
.
-
std::vector<realtype> w
w data from the model (recurring terms in xdot, for imported SBML models from python, this contains the flux vector) (shape
nt
xnw
, row major)
-
std::vector<realtype> sigmaz
event output sigma standard deviation (shape
nmaxevent
xnz
, row-major)
-
std::vector<realtype> sz
parameter derivative of event output (shape
nmaxevent
xnplist
xnz
, row-major)
-
std::vector<realtype> ssigmaz
parameter derivative of event output standard deviation (shape
nmaxevent
xnplist
xnz
, row-major)
-
std::vector<realtype> srz
parameter derivative of event trigger output (shape
nmaxevent
xnplist
xnz
, row-major)
-
std::vector<realtype> s2rz
second-order parameter derivative of event trigger output (shape
nmaxevent
xnztrue
xnplist
xnplist
, row-major)
-
std::vector<realtype> ssigmay
parameter derivative of observable standard deviation (shape
nt
xnplist
xny
, row-major)
-
std::vector<int> numsteps
number of integration steps forward problem (shape
nt
)
-
std::vector<int> numstepsB
number of integration steps backward problem (shape
nt
)
-
std::vector<int> numrhsevals
number of right hand side evaluations forward problem (shape
nt
)
-
std::vector<int> numrhsevalsB
number of right hand side evaluations backward problem (shape
nt
)
-
std::vector<int> numerrtestfails
number of error test failures forward problem (shape
nt
)
-
std::vector<int> numerrtestfailsB
number of error test failures backward problem (shape
nt
)
-
std::vector<int> numnonlinsolvconvfails
number of linear solver convergence failures forward problem (shape
nt
)
-
std::vector<int> numnonlinsolvconvfailsB
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_timeB = 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] (preequilibration)
.. 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_timeB = 0.0
computation time of the steady state solver of the backward problem [ms] (preequilibration)
.. 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_timeB = 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 steady state problem (preequilibration) [newton, simulation, newton] (length = 3)
-
int preeq_numstepsB = 0
number of simulation steps for adjoint steady state problem (preequilibration) [== 0 if analytical solution worked, > 0 otherwise]
-
std::vector<int> posteq_numsteps
number of Newton steps for steady state problem (preequilibration) [newton, simulation, newton] (shape
3
) (postequilibration)
-
int posteq_numstepsB = 0
number of simulation steps for adjoint steady state problem (postequilibration) [== 0 if analytical solution worked, > 0 otherwise]
-
realtype preeq_wrms = NAN
weighted root-mean-square of the rhs when steadystate was reached (preequilibration)
-
realtype posteq_wrms = NAN
weighted root-mean-square of the rhs when steadystate was reached (postequilibration)
-
std::vector<realtype> s2llh
second-order parameter derivative of log-likelihood (shape
nJ-1
xnplist
, 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 nx = {0}
number of states (alias
nx_rdata
, kept for backward compatibility)
-
int nxtrue = {0}
number of states in the unaugmented system (alias nxtrue_rdata, kept for backward compatibility)
-
int nplist = {0}
number of parameter for which sensitivities were requested
-
int nmaxevent = {0}
maximal number of occurring events (for every event type)
-
int nt = {0}
number of considered timepoints
-
int newton_maxsteps = {0}
maximal number of newton iterations for steady state calculation
-
std::vector<ParameterScaling> pscale
scaling of parameterization
-
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
Protected Functions
-
void initializeLikelihoodReporting(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 initializeResidualReporting(bool enable_res)
initializes storage for residual reporting mode
- Parameters:
enable_res – whether residuals are to be computed
-
void initializeFullReporting(bool enable_fim)
initializes storage for full reporting mode
- Parameters:
enable_fim – whether FIM Hessian approximation is to be computed
-
void initializeObjectiveFunction(bool enable_chi2)
initialize values for chi2 and llh and derivatives
- Parameters:
enable_chi2 – whether chi2 values are to be computed
-
void processPreEquilibration(SteadystateProblem const &preeq, Model &model)
extracts data from a preequilibration SteadystateProblem
- Parameters:
preeq – SteadystateProblem for preequilibration
model – Model instance to compute return values
-
void processPostEquilibration(SteadystateProblem const &posteq, Model &model, ExpData const *edata)
extracts data from a preequilibration SteadystateProblem
- Parameters:
posteq – SteadystateProblem for postequilibration
model – Model instance to compute return values
edata – ExpData instance containing observable data
-
void processForwardProblem(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
edata – ExpData instance containing observable data
-
void processBackwardProblem(ForwardProblem const &fwd, BackwardProblem const &bwd, SteadystateProblem const *preeq, Model &model)
extracts results from backward problem
- Parameters:
fwd – forward problem
bwd – backward problem
preeq – SteadystateProblem for preequilibration
model – model that was used for forward/backward simulation
-
void processSolver(Solver const &solver)
extracts results from solver
- Parameters:
solver – solver that was used for forward/backward simulation
-
template<class T>
inline void storeJacobianAndDerivativeInReturnData(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, SimulationState const &simulation_state, ExpData const &edata)
Residual function.
- Parameters:
it – time index
model – model that was used for forward/backward simulation
simulation_state – simulation state the timepoint
it
edata – ExpData instance containing observable data
-
void fchi2(int it, ExpData const &edata)
Chi-squared function.
- Parameters:
it – time index
edata – ExpData instance containing observable data
-
void fsres(int it, Model &model, SimulationState const &simulation_state, ExpData const &edata)
Residual sensitivity function.
- Parameters:
it – time index
model – model that was used for forward/backward simulation
simulation_state – simulation state the timepoint
it
edata – ExpData instance containing observable data
-
void fFIM(int it, Model &model, SimulationState const &simulation_state, ExpData const &edata)
Fisher information matrix function.
- Parameters:
it – time index
model – model that was used for forward/backward simulation
simulation_state – simulation state the timepoint
it
edata – ExpData 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 invalidateLLH()
Set likelihood and chi2 to NaN (typically after integration failure)
-
void invalidateSLLH()
Set likelihood sensitivities to NaN (typically after integration failure)
-
void applyChainRuleFactorToSimulationResults(Model const &model)
applies the chain rule to account for parameter transformation in the sensitivities of simulation results
- Parameters:
model – Model from which the ReturnData was obtained
-
inline bool computingFSA() const
Checks whether forward sensitivity analysis is performed.
- Returns:
boolean indicator
-
void getDataOutput(int it, Model &model, SimulationState const &simulation_state, 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
simulation_state – simulation state the timepoint
it
edata – ExpData instance carrying experimental data
-
void getDataSensisFSA(int it, Model &model, SimulationState const &simulation_state, 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
simulation_state – simulation state the timepoint
it
edata – ExpData instance carrying experimental data
-
void getEventOutput(realtype t, std::vector<int> const &rootidx, Model &model, SimulationState const &simulation_state, ExpData const *edata)
Extracts output information for events, expects that the model state was set appropriately.
- Parameters:
t – event timepoint
rootidx – information about which roots fired (1 indicating fired, 0/-1 for not)
model – model that was used in forward solve
simulation_state – simulation state the timepoint
it
edata – ExpData instance carrying experimental data
-
void getEventSensisFSA(int ie, realtype t, Model &model, SimulationState const &simulation_state, ExpData const *edata)
Extracts event information for forward sensitivity analysis, expects the model state was set appropriately.
- Parameters:
ie – index of event type
t – event timepoint
model – model that was used in forward solve
simulation_state – simulation state the timepoint
it
edata – ExpData instance carrying experimental data
-
void handleSx0Backward(Model const &model, SteadystateProblem const &preeq, std::vector<realtype> &llhS0, AmiVector &xQB) 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
preeq – SteadystateProblem for preequilibration
llhS0 – contribution to likelihood for initial state sensitivities of preequilibration
xQB – vector with quadratures from adjoint computation
-
void handleSx0Forward(Model const &model, SimulationState const &simulation_state, std::vector<realtype> &llhS0, AmiVector &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
simulation_state – simulation state the timepoint
it
llhS0 – contribution to likelihood for initial state sensitivities
xB – vector with final adjoint state (excluding conservation laws)
Protected Attributes
-
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
-
ReturnData() = default