Class ReturnData¶
Defined in File rdata.h
Class Documentation¶
-
class
amici::ReturnData¶ Stores all data to be returned by amici::runAmiciSimulation.
NOTE: multidimensional arrays are stored in row-major order (FORTRAN-style)
Public Functions
-
ReturnData() = default¶ default constructor
-
ReturnData(std::vector<realtype> ts, int np, int nk, int nx, int nx_solver, int nxtrue, int nx_solver_reinit, int ny, int nytrue, int nz, int nztrue, int ne, int nJ, int nplist, int nmaxevent, int nt, int newton_maxsteps, int nw, std::vector<ParameterScaling> pscale, SecondOrderMode o2mode, SensitivityOrder sensi, SensitivityMethod sensi_meth, RDataReporting rdrm, bool quadratic_llh)¶ -
- Parameters
ts: see amici::Model::tsnp: see amici::Model::npnk: see amici::Model::nknx: see amici::Model::nx_rdatanx_solver: see amici::Model::nx_solvernxtrue: see amici::Model::nxtrue_rdatanx_solver_reinit: see amici::Model::nx_solver_reinitny: see amici::Model::nynytrue: see amici::Model::nytruenz: see amici::Model::nznztrue: see amici::Model::nztruene: see amici::Model::nenJ: see amici::Model::nJnplist: see amici::Model::nplistnmaxevent: see amici::Model::nmaxeventnt: see amici::Model::ntnewton_maxsteps: see amici::Solver::newton_maxstepsnw: see amici::Model::nwpscale: see amici::Model::pscaleo2mode: see amici::Model::o2modesensi: see amici::Solver::sensisensi_meth: see amici::Solver::sensi_methrdrm: see amici::Solver::rdata_reportingquadratic_llh: whether model defines a quadratic nllh and computing res, sres and FIM makes sense
-
ReturnData(Solver const &solver, const Model &model)¶ constructor that uses information from model and solver to appropriately initialize fields
- Parameters
solver: solver instancemodel: 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 instancesolver: matching solver instanceedata: matching experimental data
Public Members
-
std::vector<realtype>
J¶ Jacobian of differential equation right hand side (dimension: nx x nx, row-major)
-
std::vector<realtype>
w¶ w data from the model (recurring terms in xdot, for imported SBML models from python, this contains the flux vector) (dimensions: nt x nw, row major)
-
std::vector<realtype>
sigmaz¶ event output sigma standard deviation (dimension: nmaxevent x nz, row-major)
-
std::vector<realtype>
sz¶ parameter derivative of event output (dimension: nmaxevent x nz, row-major)
-
std::vector<realtype>
ssigmaz¶ parameter derivative of event output standard deviation (dimension: nmaxevent x nz, row-major)
-
std::vector<realtype>
srz¶ parameter derivative of event trigger output (dimension: nmaxevent x nz x nplist, row-major)
-
std::vector<realtype>
s2rz¶ second order parameter derivative of event trigger output (dimension: nmaxevent x nztrue x nplist x nplist, row-major)
-
std::vector<realtype>
sy¶ parameter derivative of observable (dimension: nt x nplist x ny, row-major)
-
std::vector<realtype>
ssigmay¶ parameter derivative of observable standard deviation (dimension: nt x nplist x ny, row-major)
-
std::vector<int>
numsteps¶ number of integration steps forward problem (dimension: nt)
-
std::vector<int>
numstepsB¶ number of integration steps backward problem (dimension: nt)
-
std::vector<int>
numrhsevals¶ number of right hand side evaluations forward problem (dimension: nt)
-
std::vector<int>
numrhsevalsB¶ number of right hand side evaluations backward problem (dimension: nt)
-
std::vector<int>
numerrtestfails¶ number of error test failures forward problem (dimension: nt)
-
std::vector<int>
numerrtestfailsB¶ number of error test failures backward problem (dimension: nt)
-
std::vector<int>
numnonlinsolvconvfails¶ number of linear solver convergence failures forward problem (dimension: nt)
-
std::vector<int>
numnonlinsolvconvfailsB¶ number of linear solver convergence failures backward problem (dimension: nt)
-
std::vector<int>
order¶ employed order forward problem (dimension: nt)
-
double
cpu_time= 0.0¶ computation time of forward solve [ms]
-
double
cpu_timeB= 0.0¶ computation time of backward solve [ms]
-
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)
-
double
preeq_cpu_timeB= 0.0¶ computation time of the steady state solver of the backward problem [ms] (preequilibration)
-
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)
-
double
posteq_cpu_timeB= 0.0¶ computation time of the steady state solver of the backward problem [ms] (postequilibration)
-
std::vector<int>
preeq_numsteps¶ number of Newton steps for steady state problem (preequilibration) [newton, simulation, newton] (length = 3)
-
std::vector<int>
preeq_numlinsteps¶ number of linear steps by Newton step for steady state problem. this will only be filled for iterative solvers (preequilibration) (length = newton_maxsteps * 2)
-
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] (length = 3) (postequilibration)
-
std::vector<int>
posteq_numlinsteps¶ number of linear steps by Newton step for steady state problem. this will only be filled for iterative solvers (postequilibration) (length = newton_maxsteps * 2)
-
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>
sx_ss¶ preequilibration sensitivities found by Newton solver (dimension: nplist x nx, row-major)
-
std::vector<realtype>
s2llh¶ second order parameter derivative of loglikelihood (dimension: (nJ-1) x nplist, row-major)
-
int
status= 0¶ status code
-
int
np= {0}¶ total number of model parameters
-
int
nk= {0}¶ number of fixed parameters
-
int
nx= {0}¶ number of states
-
int
nx_solver= {0}¶ number of states with conservation laws applied
-
int
nxtrue= {0}¶ number of states in the unaugmented system
-
int
nx_solver_reinit= {0}¶ number of solver states to be reinitialized after preequilibration
-
int
ny= {0}¶ number of observables
-
int
nytrue= {0}¶ number of observables in the unaugmented system
-
int
nz= {0}¶ number of event outputs
-
int
nztrue= {0}¶ number of event outputs in the unaugmented system
-
int
ne= {0}¶ number of events
-
int
nJ= {0}¶ dimension of the augmented objective function for 2nd order ASA
-
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
nw= {0}¶ number of columns in w
-
int
newton_maxsteps= {0}¶ maximal number of newton iterations for steady state calculation
-
std::vector<ParameterScaling>
pscale¶ scaling of parameterization (lin,log,log10)
-
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
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: wheter residuals are to be computed
-
void
initializeFullReporting(bool enable_fim)¶ initializes storage for full reporting mode
- Parameters
enable_fim: wheter FIM Hessian approximation is to be computed
-
void
initializeObjectiveFunction(bool enable_chi2)¶ initialize values for chi2 and llh and derivatives
- Parameters
enable_chi2: wheter chi2 values are to be computed
-
void
processPreEquilibration(SteadystateProblem const &preeq, Model &model)¶ extracts data from a preequilibration steadystateproblem
- Parameters
preeq: Steadystateproblem for preequilibrationmodel: Model instance to compute return values
-
void
processPostEquilibration(SteadystateProblem const &posteq, Model &model, ExpData const *edata)¶ extracts data from a preequilibration steadystateproblem
-
void
processForwardProblem(ForwardProblem const &fwd, Model &model, ExpData const *edata)¶ extracts results from forward problem
- Parameters
fwd: forward problemmodel: model that was used for forward simulationedata: 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 problembwd: backward problempreeq: Steadystateproblem for preequilibrationmodel: 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>
voidstoreJacobianAndDerivativeInReturnData(T const &problem, Model &model)¶ Evaluates and stores the Jacobian and right hand side at final timepoint.
- Parameters
problem: forward problem or steadystate problemmodel: model that was used for forward/backward simulation
-
void
readSimulationState(SimulationState const &state, Model &model)¶ sets member variables and model state according to provided simulation state
- Parameters
state: simulation state provided by Problemmodel: model that was used for forward/backward simulation
-
void
fres(int it, Model &model, const ExpData &edata)¶ Residual function.
- Parameters
it: time indexmodel: model that was used for forward/backward simulationedata: ExpData instance containing observable data
-
void
fchi2(int it)¶ Chi-squared function.
- Parameters
it: time index
-
void
fsres(int it, Model &model, const ExpData &edata)¶ Residual sensitivity function.
- Parameters
it: time indexmodel: model that was used for forward/backward simulationedata: ExpData instance containing observable data
-
void
fFIM(int it, Model &model, const ExpData &edata)¶ Fisher information matrix function.
- Parameters
it: time indexmodel: model that was used for forward/backward simulationedata: 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(const Model &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
-
bool
computingFSA() const¶ Checks whether forward sensitivity analysis is performed.
- Return
boolean indicator
-
void
getDataOutput(int it, Model &model, ExpData const *edata)¶ Extracts output information for data-points, expects that x_solver and sx_solver were were set appropriately.
- Parameters
it: timepoint indexmodel: model that was used in forward solveedata: ExpData instance carrying experimental data
-
void
getDataSensisFSA(int it, Model &model, ExpData const *edata)¶ Extracts data information for forward sensitivity analysis, expects that x_solver and sx_solver were were set appropriately.
- Parameters
it: index of current timepointmodel: model that was used in forward solveedata: ExpData instance carrying experimental data
-
void
getEventOutput(int iroot, realtype t, const std::vector<int> rootidx, Model &model, ExpData const *edata)¶ Extracts output information for events, expects that x_solver and sx_solver were were set appropriately.
- Parameters
iroot: event indext: event timepointrootidx: information about which roots fired (1 indicating fired, 0/-1 for not)model: model that was used in forward solveedata: ExpData instance carrying experimental data
-
void
getEventSensisFSA(int iroot, int ie, realtype t, Model &model, ExpData const *edata)¶ Extracts event information for forward sensitivity analysis, expects that x_solver and sx_solver were set appropriately.
- Parameters
iroot: event indexie: index of event typet: event timepointmodel: model that was used in forward solveedata: ExpData instance carrying experimental data
-
void
handleSx0Backward(const Model &model, SteadystateProblem const &preeq, 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 simulationpreeq: Steadystateproblem for preequilibrationxQB: vector with quadratures from adjoint computation
-
void
handleSx0Forward(const Model &model, 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 simulationllhS0: contribution to likelihood for initial state sensitivitiesxB: vector with final adjoint state (exluding conservation laws)
Protected Attributes
-
AmiVector
dx_solver_¶ partial time derivative of state vector, excluding states eliminated from conservation laws
-
AmiVectorArray
sx_solver_¶ partial sensitivity state vector array, excluding states eliminated from conservation laws
-
AmiVectorArray
sx_rdata_¶ full sensitivity state vector array, including states eliminated from conservation laws
-
std::vector<int>
nroots_¶ array of number of found roots for a certain event type (dimension: ne)
Friends
-
template<class
Archive>
friend voidserialize(Archive &ar, ReturnData &r, unsigned int version)¶ Serialize ReturnData (see boost::serialization::serialize)
- Parameters
ar: Archive to serialize tor: Data to serializeversion: Version number
-