amici.amici

Core C++ bindings

This module encompasses the complete public C++ API of AMICI, which was exposed via swig. All functions listed here are directly accessible in the main amici package, i.e., amici.amici.ExpData is available as amici.ExpData. Usage of functions and classes from the base amici package is generally recommended as they often include convenience wrappers that avoid common pitfalls when accessing C++ types from python and implement some nonstandard type conversions.

Module Attributes

SensitivityOrder_none

Don't compute sensitivities.

SensitivityOrder_first

First-order sensitivities.

SensitivityOrder_second

Second-order sensitivities.

SensitivityMethod_none

Don't compute sensitivities.

SensitivityMethod_forward

Forward sensitivity analysis.

SensitivityMethod_adjoint

Adjoint sensitivity analysis.

NonlinearSolverIteration_fixedpoint

deprecated

Functions

compiledWithOpenMP()

AMICI extension was compiled with OpenMP?

enum(prefix)

getScaledParameter(unscaledParameter,Β scaling)

Apply parameter scaling according to scaling

getUnscaledParameter(scaledParameter,Β scaling)

Remove parameter scaling according to scaling

parameterScalingFromIntVector(intVec)

Swig-Generated class, which, in contrast to other Vector classes, does not allow for simple interoperability with common Python types, but must be created using amici.amici.parameterScalingFromIntVector()

runAmiciSimulation(solver,Β edata,Β model[,Β ...])

Core integration routine.

runAmiciSimulations(solver,Β edatas,Β model,Β ...)

Same as runAmiciSimulation, but for multiple ExpData instances.

simulation_status_to_str(status)

Get the string representation of the given simulation status code (see ReturnData::status).

Classes

BoolVector(*args)

Swig-Generated class templating common python types including Iterable [bool] and numpy.array [bool] to facilitate interfacing with C++ bindings.

Constraint(value[,Β names,Β module,Β qualname,Β ...])

CpuTimer()

Tracks elapsed CPU time using std::clock.

DoubleVector(*args)

Swig-Generated class templating common python types including Iterable [float] and numpy.array [float] to facilitate interfacing with C++ bindings.

ExpData(*args)

ExpData carries all information about experimental or condition-specific data.

ExpDataPtr(*args)

Swig-Generated class that implements smart pointers to ExpData as objects.

ExpDataPtrVector(*args)

Swig-Generated class templating common python types including Iterable [amici.amici.ExpData] and numpy.array [amici.amici.ExpData] to facilitate interfacing with C++ bindings.

FixedParameterContext(value[,Β names,Β ...])

IntVector(*args)

Swig-Generated class templating common python types including Iterable [int] and numpy.array [int] to facilitate interfacing with C++ bindings.

InternalSensitivityMethod(value[,Β names,Β ...])

InterpolationType(value[,Β names,Β module,Β ...])

LinearMultistepMethod(value[,Β names,Β ...])

LinearSolver(value[,Β names,Β module,Β ...])

LogItem(*args)

A log item.

LogItemVector(*args)

Model(*args,Β **kwargs)

The Model class represents an AMICI ODE/DAE model.

ModelDimensions(*args)

Container for model dimensions.

ModelPtr(*args)

Swig-Generated class that implements smart pointers to Model as objects.

NewtonDampingFactorMode(value[,Β names,Β ...])

NonlinearSolverIteration(value[,Β names,Β ...])

ObservableScaling(value[,Β names,Β module,Β ...])

ParameterScaling(value[,Β names,Β module,Β ...])

ParameterScalingVector(*args)

RDataReporting(value[,Β names,Β module,Β ...])

ReturnData(*args)

Stores all data to be returned by amici.amici.runAmiciSimulation().

ReturnDataPtr(*args)

Swig-Generated class that implements smart pointers to ReturnData as objects.

SecondOrderMode(value[,Β names,Β module,Β ...])

SensitivityMethod(value[,Β names,Β module,Β ...])

SensitivityOrder(value[,Β names,Β module,Β ...])

SimulationParameters(*args)

Container for various simulation parameters.

Solver(*args,Β **kwargs)

The Solver class provides a generic interface to CVODES and IDAS solvers, individual realizations are realized in the CVodeSolver and the IDASolver class.

SolverPtr(*args)

Swig-Generated class that implements smart pointers to Solver as objects.

SteadyStateComputationMode(value[,Β names,Β ...])

SteadyStateSensitivityMode(value[,Β names,Β ...])

SteadyStateStatus(value[,Β names,Β module,Β ...])

SteadyStateStatusVector(*args)

StringDoubleMap(*args)

Swig-Generated class templating Dict [str, float] to facilitate interfacing with C++ bindings.

StringVector(*args)

Swig-Generated class templating common python types including Iterable [str] and numpy.array [str] to facilitate interfacing with C++ bindings.

class amici.amici.BoolVector(*args)[source]

Swig-Generated class templating common python types including Iterable [bool] and numpy.array [bool] to facilitate interfacing with C++ bindings.

class amici.amici.Constraint(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
negative = -2
non_negative = 1
non_positive = -1
none = 0
positive = 2
class amici.amici.CpuTimer[source]

Tracks elapsed CPU time using std::clock.

__init__()[source]

Constructor

elapsed_milliseconds() float[source]

Get elapsed CPU time in milliseconds since initialization or last reset

Return type:

float

Returns:

CPU time in milliseconds

elapsed_seconds() float[source]

Get elapsed CPU time in seconds since initialization or last reset

Return type:

float

Returns:

CPU time in seconds

reset()[source]

Reset the timer

uses_thread_clock = False

Whether the timer uses a thread clock (i.e. provides proper, thread-specific CPU time).

class amici.amici.DoubleVector(*args)[source]

Swig-Generated class templating common python types including Iterable [float] and numpy.array [float] to facilitate interfacing with C++ bindings.

class amici.amici.ExpData(*args)[source]

ExpData carries all information about experimental or condition-specific data.

__init__(*args)[source]

Overload 1:

Default constructor.


Overload 2:

Copy constructor.


Overload 3:

Constructor that only initializes dimensions.

Parameters:
  • nytrue (int) – Number of observables

  • nztrue (int) – Number of event outputs

  • nmaxevent (int) – Maximal number of events to track


Overload 4:

constructor that initializes timepoints from vectors

Parameters:
  • nytrue (int) – Number of observables

  • nztrue (int) – Number of event outputs

  • nmaxevent (int) – Maximal number of events to track

  • ts (DoubleVector) – Timepoints (dimension: nt)


Overload 5:

constructor that initializes timepoints and fixed parameters from vectors

Parameters:
  • nytrue (int) – Number of observables

  • nztrue (int) – Number of event outputs

  • nmaxevent (int) – Maximal number of events to track

  • ts (DoubleVector) – Timepoints (dimension: nt)

  • fixedParameters (DoubleVector) – Model constants (dimension: nk)


Overload 6:

constructor that initializes timepoints and data from vectors

Parameters:
  • nytrue (int) – Number of observables

  • nztrue (int) – Number of event outputs

  • nmaxevent (int) – Maximal number of events to track

  • ts (DoubleVector) – Timepoints (dimension: nt)

  • observedData (DoubleVector) – observed data (dimension: nt x nytrue, row-major)

  • observedDataStdDev (DoubleVector) – standard deviation of observed data (dimension: nt x nytrue, row-major)

  • observedEvents (DoubleVector) – observed events (dimension: nmaxevents x nztrue, row-major)

  • observedEventsStdDev (DoubleVector) – standard deviation of observed events/roots (dimension: nmaxevents x nztrue, row-major)


Overload 7:

constructor that initializes with Model

Parameters:

model (Model) – pointer to model specification object


Overload 8:

constructor that initializes with returnData, adds noise according to specified sigmas

Parameters:
  • rdata (ReturnData) – return data pointer with stored simulation results

  • sigma_y (float) – scalar standard deviations for all observables

  • sigma_z (float) – scalar standard deviations for all event observables


Overload 9:

constructor that initializes with returnData, adds noise according to specified sigmas

Parameters:
  • rdata (ReturnData) – return data pointer with stored simulation results

  • sigma_y (DoubleVector) – vector of standard deviations for observables (dimension: nytrue or nt x nytrue, row-major)

  • sigma_z (DoubleVector) – vector of standard deviations for event observables (dimension: nztrue or nmaxevent x nztrue, row-major)

clear_observations()[source]

Set all observations and their standard deviations to NaN.

Useful, e.g., after calling ExpData::setTimepoints.

property fixedParameters

Model constants

Vector of size Model::nk() or empty

property fixedParametersPreequilibration

Model constants for pre-equilibration

Vector of size Model::nk() or empty.

property fixedParametersPresimulation

Model constants for pre-simulation

Vector of size Model::nk() or empty.

getObservedData() Sequence[float][source]

Get all measurements.

Return type:

typing.Sequence[float]

Returns:

observed data (dimension: nt x nytrue, row-major)

getObservedDataPtr(it: int) float[source]

Get measurements for a given timepoint index.

Parameters:

it (int) – timepoint index

Return type:

float

Returns:

pointer to observed data at index (dimension: nytrue)

getObservedDataStdDev() Sequence[float][source]

Get measurement standard deviations.

Return type:

typing.Sequence[float]

Returns:

standard deviation of observed data

getObservedDataStdDevPtr(it: int) float[source]

Get pointer to measurement standard deviations.

Parameters:

it (int) – timepoint index

Return type:

float

Returns:

pointer to standard deviation of observed data at index

getObservedEvents() Sequence[float][source]

Get observed event data.

Return type:

typing.Sequence[float]

Returns:

observed event data

getObservedEventsPtr(ie: int) float[source]

Get pointer to observed data at ie-th occurrence.

Parameters:

ie (int) – event occurrence

Return type:

float

Returns:

pointer to observed event data at ie-th occurrence

getObservedEventsStdDev() Sequence[float][source]

Get standard deviation of observed event data.

Return type:

typing.Sequence[float]

Returns:

standard deviation of observed event data

getObservedEventsStdDevPtr(ie: int) float[source]

Get pointer to standard deviation of observed event data at ie-th occurrence.

Parameters:

ie (int) – event occurrence

Return type:

float

Returns:

pointer to standard deviation of observed event data at ie-th occurrence

getTimepoint(it: int) float[source]

Get timepoint for the given index

Parameters:

it (int) – timepoint index

Return type:

float

Returns:

timepoint timepoint at index

getTimepoints() Sequence[float][source]

Get output timepoints.

Return type:

typing.Sequence[float]

Returns:

ExpData::ts

property id

Arbitrary (not necessarily unique) identifier.

isSetObservedData(it: int, iy: int) bool[source]

Whether there is a measurement for the given time- and observable- index.

Parameters:
  • it (int) – time index

  • iy (int) – observable index

Return type:

bool

Returns:

boolean specifying if data was set

isSetObservedDataStdDev(it: int, iy: int) bool[source]

Whether standard deviation for a measurement at specified timepoint- and observable index has been set.

Parameters:
  • it (int) – time index

  • iy (int) – observable index

Return type:

bool

Returns:

boolean specifying if standard deviation of data was set

isSetObservedEvents(ie: int, iz: int) bool[source]

Check whether event data at specified indices has been set.

Parameters:
  • ie (int) – event index

  • iz (int) – event observable index

Return type:

bool

Returns:

boolean specifying if data was set

isSetObservedEventsStdDev(ie: int, iz: int) bool[source]

Check whether standard deviation of event data at specified indices has been set.

Parameters:
  • ie (int) – event index

  • iz (int) – event observable index

Return type:

bool

Returns:

boolean specifying if standard deviation of event data was set

nmaxevent() int[source]

maximal number of events to track

Return type:

int

Returns:

maximal number of events to track

nt() int[source]

number of timepoints

Return type:

int

Returns:

number of timepoints

nytrue() int[source]

number of observables of the non-augmented model

Return type:

int

Returns:

number of observables of the non-augmented model

nztrue() int[source]

number of event observables of the non-augmented model

Return type:

int

Returns:

number of event observables of the non-augmented model

property parameters

Model parameters

Vector of size Model::np() or empty with parameter scaled according to SimulationParameter::pscale.

property plist

Parameter indices w.r.t. which to compute sensitivities

property pscale

Parameter scales

Vector of parameter scale of size Model::np(), indicating how/if each parameter is to be scaled.

property reinitialization_state_idxs_presim

Indices of states to be reinitialized based on provided presimulation constants / fixed parameters.

property reinitialization_state_idxs_sim

Indices of states to be reinitialized based on provided constants / fixed parameters.

reinitializeAllFixedParameterDependentInitialStates(nx_rdata: int)

Set reinitialization of all states based on model constants for all simulation phases.

Convenience function to populate reinitialization_state_idxs_presim and reinitialization_state_idxs_sim

Parameters:

nx_rdata (int) – Number of states (Model::nx_rdata)

reinitializeAllFixedParameterDependentInitialStatesForPresimulation(nx_rdata: int)

Set reinitialization of all states based on model constants for presimulation (only meaningful if preequilibration is performed).

Convenience function to populate reinitialization_state_idxs_presim and reinitialization_state_idxs_sim

Parameters:

nx_rdata (int) – Number of states (Model::nx_rdata)

reinitializeAllFixedParameterDependentInitialStatesForSimulation(nx_rdata: int)

Set reinitialization of all states based on model constants for the β€˜main’ simulation (only meaningful if presimulation or preequilibration is performed).

Convenience function to populate reinitialization_state_idxs_presim and reinitialization_state_idxs_sim

Parameters:

nx_rdata (int) – Number of states (Model::nx_rdata)

property reinitializeFixedParameterInitialStates

Flag indicating whether reinitialization of states depending on fixed parameters is activated

setObservedData(*args)[source]

Overload 1:

Set all measurements.

Parameters:

observedData (DoubleVector) – observed data (dimension: nt x nytrue, row-major)


Overload 2:

Set measurements for a given observable index

Parameters:
  • observedData (DoubleVector) – observed data (dimension: nt)

  • iy (int) – observed data index

setObservedDataStdDev(*args)[source]

Overload 1:

Set standard deviations for measurements.

Parameters:

observedDataStdDev (DoubleVector) – standard deviation of observed data (dimension: nt x nytrue, row-major)


Overload 2:

Set identical standard deviation for all measurements.

Parameters:

stdDev (float) – standard deviation (dimension: scalar)


Overload 3:

Set standard deviations of observed data for a specific observable index.

Parameters:
  • observedDataStdDev (DoubleVector) – standard deviation of observed data (dimension: nt)

  • iy (int) – observed data index


Overload 4:

Set all standard deviation for a given observable index to the input value.

Parameters:
  • stdDev (float) – standard deviation (dimension: scalar)

  • iy (int) – observed data index

setObservedEvents(*args)[source]

Overload 1:

Set observed event data.

Parameters:

observedEvents (DoubleVector) – observed data (dimension: nmaxevent x nztrue, row-major)


Overload 2:

Set observed event data for specific event observable.

Parameters:
  • observedEvents (DoubleVector) – observed data (dimension: nmaxevent)

  • iz (int) – observed event data index

setObservedEventsStdDev(*args)[source]

Overload 1:

Set standard deviation of observed event data.

Parameters:

observedEventsStdDev (DoubleVector) – standard deviation of observed event data


Overload 2:

Set standard deviation of observed event data.

Parameters:

stdDev (float) – standard deviation (dimension: scalar)


Overload 3:

Set standard deviation of observed data for a specific observable.

Parameters:
  • observedEventsStdDev (DoubleVector) – standard deviation of observed data (dimension: nmaxevent)

  • iz (int) – observed data index


Overload 4:

Set all standard deviations of a specific event-observable.

Parameters:
  • stdDev (float) – standard deviation (dimension: scalar)

  • iz (int) – observed data index

setTimepoints(ts: Sequence[float])[source]

Set output timepoints.

If the number of timepoint increases, this will grow the observation/sigma matrices and fill new entries with NaN. If the number of timepoints decreases, this will shrink the observation/sigma matrices.

Note that the mapping from timepoints to measurements will not be preserved. E.g., say there are measurements at t = 2, and this function is called with [1, 2], then the old measurements will belong to t = 1.

Parameters:

ts (typing.Sequence[float]) – timepoints

property sx0

Initial state sensitivities

Dimensions: Model::nx() * Model::nplist(), Model::nx() * ExpData::plist.size(), if ExpData::plist is not empty, or empty

property t_presim

Duration of pre-simulation.

If this is > 0, presimulation will be performed from (model->t0 - t_presim) to model->t0 using the fixedParameters in fixedParametersPresimulation

property ts_

Timepoints for which model state/outputs/… are requested

Vector of timepoints.

property tstart_

Starting time of the simulation.

Output timepoints are absolute timepoints, independent of \(t_{start}\). For output timepoints \(t < t_{start}\), the initial state will be returned.

property x0

Initial state

Vector of size Model::nx() or empty

class amici.amici.ExpDataPtr(*args)[source]

Swig-Generated class that implements smart pointers to ExpData as objects.

property fixedParameters

Model constants

Vector of size Model::nk() or empty

property fixedParametersPreequilibration

Model constants for pre-equilibration

Vector of size Model::nk() or empty.

property fixedParametersPresimulation

Model constants for pre-simulation

Vector of size Model::nk() or empty.

property id

Arbitrary (not necessarily unique) identifier.

property parameters

Model parameters

Vector of size Model::np() or empty with parameter scaled according to SimulationParameter::pscale.

property plist

Parameter indices w.r.t. which to compute sensitivities

property pscale

Parameter scales

Vector of parameter scale of size Model::np(), indicating how/if each parameter is to be scaled.

property reinitialization_state_idxs_presim

Indices of states to be reinitialized based on provided presimulation constants / fixed parameters.

property reinitialization_state_idxs_sim

Indices of states to be reinitialized based on provided constants / fixed parameters.

property reinitializeFixedParameterInitialStates

Flag indicating whether reinitialization of states depending on fixed parameters is activated

property sx0

Initial state sensitivities

Dimensions: Model::nx() * Model::nplist(), Model::nx() * ExpData::plist.size(), if ExpData::plist is not empty, or empty

property t_presim

Duration of pre-simulation.

If this is > 0, presimulation will be performed from (model->t0 - t_presim) to model->t0 using the fixedParameters in fixedParametersPresimulation

property ts_

Timepoints for which model state/outputs/… are requested

Vector of timepoints.

property tstart_

Starting time of the simulation.

Output timepoints are absolute timepoints, independent of \(t_{start}\). For output timepoints \(t < t_{start}\), the initial state will be returned.

property x0

Initial state

Vector of size Model::nx() or empty

class amici.amici.ExpDataPtrVector(*args)[source]

Swig-Generated class templating common python types including Iterable [amici.amici.ExpData] and numpy.array [amici.amici.ExpData] to facilitate interfacing with C++ bindings.

class amici.amici.FixedParameterContext(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
preequilibration = 1
presimulation = 2
simulation = 0
class amici.amici.IntVector(*args)[source]

Swig-Generated class templating common python types including Iterable [int] and numpy.array [int] to facilitate interfacing with C++ bindings.

class amici.amici.InternalSensitivityMethod(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
simultaneous = 1
staggered = 2
staggered1 = 3
class amici.amici.InterpolationType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
hermite = 1
polynomial = 2
class amici.amici.LinearMultistepMethod(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
BDF = 2
__init__(*args, **kwds)
adams = 1
class amici.amici.LinearSolver(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
KLU = 9
LAPACKBand = 4
LAPACKDense = 3
SPBCG = 7
SPGMR = 6
SPTFQMR = 8
SuperLUMT = 10
__init__(*args, **kwds)
band = 2
dense = 1
diag = 5
class amici.amici.LogItem(*args)[source]

A log item.

__init__(*args)[source]

Overload 1:

Default ctor.


Overload 2:

Construct a LogItem

Parameters:
  • severity (int)

  • identifier (str)

  • message (str)

property identifier

Short identifier for the logged event

property message

A more detailed and readable message

property severity

Severity level

class amici.amici.LogItemVector(*args)[source]
class amici.amici.Model(*args, **kwargs)[source]

The Model class represents an AMICI ODE/DAE model.

The model can compute various model related quantities based on symbolically generated code.

__init__(*args, **kwargs)[source]
Overload 1:

Default ctor


Overload 2:

Constructor with model dimensions

Parameters:
  • nx_rdata (int) – Number of state variables

  • nxtrue_rdata (int) – Number of state variables of the non-augmented model

  • nx_solver (int) – Number of state variables with conservation laws applied

  • nxtrue_solver (int) – Number of state variables of the non-augmented model with conservation laws applied

  • nx_solver_reinit (int) – Number of state variables with conservation laws subject to reinitialization

  • np (int) – Number of parameters

  • nk (int) – Number of constants

  • ny (int) – Number of observables

  • nytrue (int) – Number of observables of the non-augmented model

  • nz (int) – Number of event observables

  • nztrue (int) – Number of event observables of the non-augmented model

  • ne (int) – Number of events

  • ne_solver (int) – Number of events that require root-finding

  • nspl (int) – Number of splines

  • nJ (int) – Number of objective functions

  • nw (int) – Number of repeating elements

  • ndwdx (int) – Number of nonzero elements in the x derivative of the repeating elements

  • ndwdp (int) – Number of nonzero elements in the p derivative of the repeating elements

  • ndwdw (int) – Number of nonzero elements in the w derivative of the repeating elements

  • ndxdotdw (int) – Number of nonzero elements in the \(w\) derivative of \(xdot\)

  • ndJydy (IntVector) – Number of nonzero elements in the \(y\) derivative of \(dJy\) (shape nytrue)

  • ndxrdatadxsolver (int) – Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

  • ndxrdatadtcl (int) – Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

  • ndtotal_cldx_rdata (int) – Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

  • nnz (int) – Number of nonzero elements in Jacobian

  • ubw (int) – Upper matrix bandwidth in the Jacobian

  • lbw (int) – Lower matrix bandwidth in the Jacobian

clone() Model[source]

Clone this instance.

Return type:

amici.amici.Model

Returns:

The clone

fdsigmaydy(dsigmaydy: float, t: float, p: float, k: float, y: float)

Model-specific implementation of fsigmay

Parameters:
  • dsigmaydy (float) – partial derivative of standard deviation of measurements w.r.t. model outputs

  • t (float) – current time

  • p (float) – parameter vector

  • k (float) – constant vector

  • y (float) – model output at timepoint t

fdspline_slopesdp(dspline_slopesdp: float, p: float, k: float, ip: int)

Model-specific implementation the parametric derivatives of slopevalues at spline nodes

Parameters:
  • dspline_slopesdp (float) – vector to which derivatives will be written

  • p (float) – parameter vector

  • k (float) – constants vector

  • ip (int) – Sensitivity index

fdspline_valuesdp(dspline_valuesdp: float, p: float, k: float, ip: int)

Model-specific implementation the parametric derivatives of spline node values

Parameters:
  • dspline_valuesdp (float) – vector to which derivatives will be written

  • p (float) – parameter vector

  • k (float) – constants vector

  • ip (int) – Sensitivity index

fdtotal_cldp(dtotal_cldp: float, x_rdata: float, p: float, k: float, ip: int)

Compute dtotal_cl / dp

Parameters:
  • dtotal_cldp (float) – dtotal_cl / dp

  • x_rdata (float) – State variables with conservation laws applied

  • p (float) – parameter vector

  • k (float) – constant vector

  • ip (int) – Sensitivity index

fdtotal_cldx_rdata(dtotal_cldx_rdata: float, x_rdata: float, p: float, k: float, tcl: float)

Compute dtotal_cl / dx_rdata

Parameters:
  • dtotal_cldx_rdata (float) – dtotal_cl / dx_rdata

  • x_rdata (float) – State variables with conservation laws applied

  • p (float) – parameter vector

  • k (float) – constant vector

  • tcl (float) – Total abundances for conservation laws

fdx_rdatadp(dx_rdatadp: float, x: float, tcl: float, p: float, k: float, ip: int)

Compute dx_rdata / dp

Parameters:
  • dx_rdatadp (float) – dx_rdata / dp

  • p (float) – parameter vector

  • k (float) – constant vector

  • x (float) – State variables with conservation laws applied

  • tcl (float) – Total abundances for conservation laws

  • ip (int) – Sensitivity index

fdx_rdatadtcl(dx_rdatadtcl: float, x: float, tcl: float, p: float, k: float)

Compute dx_rdata / dtcl

Parameters:
  • dx_rdatadtcl (float) – dx_rdata / dtcl

  • p (float) – parameter vector

  • k (float) – constant vector

  • x (float) – State variables with conservation laws applied

  • tcl (float) – Total abundances for conservation laws

fdx_rdatadx_solver(dx_rdatadx_solver: float, x: float, tcl: float, p: float, k: float)

Compute dx_rdata / dx_solver

Parameters:
  • dx_rdatadx_solver (float) – dx_rdata / dx_solver

  • p (float) – parameter vector

  • k (float) – constant vector

  • x (float) – State variables with conservation laws applied

  • tcl (float) – Total abundances for conservation laws

getAddSigmaResiduals() bool[source]

Checks whether residuals should be added to account for parameter dependent sigma.

Return type:

bool

Returns:

sigma_res

getAlwaysCheckFinite() bool[source]

Get setting of whether the result of every call to Model::f* should be checked for finiteness.

Return type:

bool

Returns:

that

getAmiciCommit() str

Returns the AMICI commit that was used to generate the model

Return type:

str

Returns:

AMICI commit string

getAmiciVersion() str

Returns the AMICI version that was used to generate the model

Return type:

str

Returns:

AMICI version string

getExpressionIds() Sequence[str][source]

Get IDs of the expression.

Return type:

typing.Sequence[str]

Returns:

Expression IDs

getExpressionNames() Sequence[str][source]

Get names of the expressions.

Return type:

typing.Sequence[str]

Returns:

Expression names

getFixedParameterById(par_id: str) float[source]

Get value of fixed parameter with the specified ID.

Parameters:

par_id (str) – Parameter ID

Return type:

float

Returns:

Parameter value

getFixedParameterByName(par_name: str) float[source]

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 (str) – Parameter name

Return type:

float

Returns:

Parameter value

getFixedParameterIds() Sequence[str][source]

Get IDs of the fixed model parameters.

Return type:

typing.Sequence[str]

Returns:

Fixed parameter IDs

getFixedParameterNames() Sequence[str][source]

Get names of the fixed model parameters.

Return type:

typing.Sequence[str]

Returns:

Fixed parameter names

getFixedParameters() Sequence[float][source]

Get values of fixed parameters.

Return type:

typing.Sequence[float]

Returns:

Vector of fixed parameters with same ordering as in Model::getFixedParameterIds

getInitialStateSensitivities() Sequence[float][source]

Get the initial states sensitivities.

Return type:

typing.Sequence[float]

Returns:

vector of initial state sensitivities

getInitialStates() Sequence[float][source]

Get the initial states.

Return type:

typing.Sequence[float]

Returns:

Initial state vector

getMinimumSigmaResiduals() float[source]

Gets the specified estimated lower boundary for sigma_y.

Return type:

float

Returns:

lower boundary

getName() str[source]

Get the model name.

Return type:

str

Returns:

Model name

getObservableIds() Sequence[str][source]

Get IDs of the observables.

Return type:

typing.Sequence[str]

Returns:

Observable IDs

getObservableNames() Sequence[str][source]

Get names of the observables.

Return type:

typing.Sequence[str]

Returns:

Observable names

getObservableScaling(iy: int) int[source]

Get scaling type for observable

Parameters:

iy (int) – observable index

Return type:

int

Returns:

scaling type

getParameterById(par_id: str) float[source]

Get value of first model parameter with the specified ID.

Parameters:

par_id (str) – Parameter ID

Return type:

float

Returns:

Parameter value

getParameterByName(par_name: str) float[source]

Get value of first model parameter with the specified name.

Parameters:

par_name (str) – Parameter name

Return type:

float

Returns:

Parameter value

getParameterIds() Sequence[str][source]

Get IDs of the model parameters.

Return type:

typing.Sequence[str]

Returns:

Parameter IDs

getParameterList() Sequence[int][source]

Get the list of parameters for which sensitivities are computed.

Return type:

typing.Sequence[int]

Returns:

List of parameter indices

getParameterNames() Sequence[str][source]

Get names of the model parameters.

Return type:

typing.Sequence[str]

Returns:

The parameter names

getParameterScale() ParameterScalingVector[source]

Get parameter scale for each parameter.

Return type:

amici.amici.ParameterScalingVector

Returns:

Vector of parameter scales

getParameters() Sequence[float][source]

Get parameter vector.

Return type:

typing.Sequence[float]

Returns:

The user-set parameters (see also Model::getUnscaledParameters)

getReinitializationStateIdxs() Sequence[int][source]

Return indices of states to be reinitialized based on provided constants / fixed parameters

Return type:

typing.Sequence[int]

Returns:

Those indices.

getReinitializeFixedParameterInitialStates() bool[source]

Get whether initial states depending on fixedParameters are to be reinitialized after preequilibration and presimulation.

Return type:

bool

Returns:

flag true / false

getSolver() Solver

Retrieves the solver object

Return type:

amici.amici.Solver

Returns:

The Solver instance

getStateIds() Sequence[str][source]

Get IDs of the model states.

Return type:

typing.Sequence[str]

Returns:

State IDs

getStateIdsSolver() Sequence[str][source]

Get IDs of the solver states.

Return type:

typing.Sequence[str]

Returns:

State IDs

getStateIsNonNegative() Sequence[bool][source]

Get flags indicating whether states should be treated as non-negative.

Return type:

typing.Sequence[bool]

Returns:

Vector of flags

getStateNames() Sequence[str][source]

Get names of the model states.

Return type:

typing.Sequence[str]

Returns:

State names

getStateNamesSolver() Sequence[str][source]

Get names of the solver states.

Return type:

typing.Sequence[str]

Returns:

State names

getSteadyStateComputationMode() int[source]

Gets the mode how steady state is computed in the steadystate simulation.

Return type:

int

Returns:

Mode

getSteadyStateSensitivityMode() SteadyStateSensitivityMode[source]

Gets the mode how sensitivities are computed in the steadystate simulation.

Return type:

amici.amici.SteadyStateSensitivityMode

Returns:

Mode

getTimepoint(it: int) float[source]

Get simulation timepoint for time index it.

Parameters:

it (int) – Time index

Return type:

float

Returns:

Timepoint

getTimepoints() Sequence[float][source]

Get the timepoint vector.

Return type:

typing.Sequence[float]

Returns:

Timepoint vector

getUnscaledParameters() Sequence[float][source]

Get parameters with transformation according to parameter scale applied.

Return type:

typing.Sequence[float]

Returns:

Unscaled parameters

get_trigger_timepoints() Sequence[float][source]

Get trigger times for events that don’t require root-finding.

Return type:

typing.Sequence[float]

Returns:

List of unique trigger points for events that don’t require root-finding (i.e. that trigger at predetermined timepoints), in ascending order.

hasCustomInitialStateSensitivities() bool[source]

Return whether custom initial state sensitivities have been set.

Return type:

bool

Returns:

true if has custom initial state sensitivities, otherwise false.

hasCustomInitialStates() bool[source]

Return whether custom initial states have been set.

Return type:

bool

Returns:

true if has custom initial states, otherwise false

hasExpressionIds() bool[source]

Report whether the model has expression IDs set.

Return type:

bool

Returns:

Boolean indicating whether expression ids were set. Also returns true if the number of corresponding variables is just zero.

hasExpressionNames() bool[source]

Report whether the model has expression names set.

Return type:

bool

Returns:

Boolean indicating whether expression names were set. Also returns true if the number of corresponding variables is just zero.

hasFixedParameterIds() bool[source]

Report whether the model has fixed parameter IDs set.

Return type:

bool

Returns:

Boolean indicating whether fixed parameter IDs were set. Also returns true if the number of corresponding variables is just zero.

hasFixedParameterNames() bool[source]

Report whether the model has fixed parameter names set.

Return type:

bool

Returns:

Boolean indicating whether fixed parameter names were set. Also returns true if the number of corresponding variables is just zero.

hasObservableIds() bool[source]

Report whether the model has observable IDs set.

Return type:

bool

Returns:

Boolean indicating whether observable ids were set. Also returns true if the number of corresponding variables is just zero.

hasObservableNames() bool[source]

Report whether the model has observable names set.

Return type:

bool

Returns:

Boolean indicating whether observable names were set. Also returns true if the number of corresponding variables is just zero.

hasParameterIds() bool[source]

Report whether the model has parameter IDs set.

Return type:

bool

Returns:

Boolean indicating whether parameter IDs were set. Also returns true if the number of corresponding variables is just zero.

hasParameterNames() bool[source]

Report whether the model has parameter names set.

Return type:

bool

Returns:

Boolean indicating whether parameter names were set. Also returns true if the number of corresponding variables is just zero.

hasQuadraticLLH() bool[source]

Checks whether the defined noise model is gaussian, i.e., the nllh is quadratic

Return type:

bool

Returns:

boolean flag

hasStateIds() bool[source]

Report whether the model has state IDs set.

Return type:

bool

Returns:

Boolean indicating whether state IDs were set. Also returns true if the number of corresponding variables is just zero.

hasStateNames() bool[source]

Report whether the model has state names set.

Return type:

bool

Returns:

Boolean indicating whether state names were set. Also returns true if the number of corresponding variables is just zero.

property idlist

Flag array for DAE equations

initializeSplineSensitivities()[source]

Initialization of spline sensitivity functions

initializeSplines()[source]

Initialization of spline functions

isFixedParameterStateReinitializationAllowed() bool

Function indicating whether reinitialization of states depending on fixed parameters is permissible

Return type:

bool

Returns:

flag indicating whether reinitialization of states depending on fixed parameters is permissible

k() float[source]

Get fixed parameters.

Return type:

float

Returns:

Pointer to constants array

property lbw

Lower bandwidth of the Jacobian

property logger

Logger

property nJ

Dimension of the augmented objective function for 2nd order ASA

nMaxEvent() int[source]

Get maximum number of events that may occur for each type.

Return type:

int

Returns:

Maximum number of events that may occur for each type

ncl() int[source]

Get number of conservation laws.

Return type:

int

Returns:

Number of conservation laws (i.e., difference between nx_rdata and nx_solver).

property ndJydy

Number of nonzero elements in the \(y\) derivative of \(dJy\) (dimension nytrue)

property ndtotal_cldx_rdata

Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

property ndwdp

Number of nonzero elements in the p derivative of the repeating elements

property ndwdw

Number of nonzero elements in the w derivative of the repeating elements

property ndwdx

Number of nonzero elements in the x derivative of the repeating elements

property ndxdotdw

Number of nonzero elements in the \(w\) derivative of \(xdot\)

property ndxrdatadtcl

Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

property ndxrdatadxsolver

Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

property ne

Number of events

property ne_solver

Number of events that require root-finding

nk() int[source]

Get number of constants

Return type:

int

Returns:

Length of constant vector

property nnz

Number of nonzero entries in Jacobian

np() int[source]

Get total number of model parameters.

Return type:

int

Returns:

Length of parameter vector

nplist() int[source]

Get number of parameters wrt to which sensitivities are computed.

Return type:

int

Returns:

Length of sensitivity index vector

property nspl

Number of spline functions in the model

nt() int[source]

Get number of timepoints.

Return type:

int

Returns:

Number of timepoints

property nw

Number of common expressions

property nx_rdata

Number of states

nx_reinit() int[source]

Get number of solver states subject to reinitialization.

Return type:

int

Returns:

Model member nx_solver_reinit

property nx_solver

Number of states with conservation laws applied

property nx_solver_reinit

Number of solver states subject to reinitialization

property nxtrue_rdata

Number of states in the unaugmented system

property nxtrue_solver

Number of states in the unaugmented system with conservation laws applied

property ny

Number of observables

property nytrue

Number of observables in the unaugmented system

property nz

Number of event outputs

property nztrue

Number of event outputs in the unaugmented system

property o2mode

Flag indicating whether for amici::Solver::sensi_ == amici::SensitivityOrder::second directional or full second order derivative will be computed

plist(pos: int) int[source]

Get entry in parameter list by index.

Parameters:

pos (int) – Index in sensitivity parameter list

Return type:

int

Returns:

Index in parameter list

property pythonGenerated

Flag indicating Matlab- or Python-based model generation

requireSensitivitiesForAllParameters()[source]

Require computation of sensitivities for all parameters p [0..np[ in natural order.

NOTE: Resets initial state sensitivities.

setAddSigmaResiduals(sigma_res: bool)[source]

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 setMinimumSigmaResiduals().

Parameters:

sigma_res (bool) – if true, additional residuals are added

setAllStatesNonNegative()[source]

Set flags indicating that all states should be treated as non-negative.

setAlwaysCheckFinite(alwaysCheck: bool)[source]

Set whether the result of every call to Model::f* should be checked for finiteness.

Parameters:

alwaysCheck (bool)

setFixedParameterById(par_id: str, value: float)[source]

Set value of first fixed parameter with the specified ID.

Parameters:
  • par_id (str) – Fixed parameter id

  • value (float) – Fixed parameter value

setFixedParameterByName(par_name: str, value: float)[source]

Set value of first fixed parameter with the specified name.

Parameters:
  • par_name (str) – Fixed parameter ID

  • value (float) – Fixed parameter value

setFixedParameters(k: Sequence[float])[source]

Set values for constants.

Parameters:

k (typing.Sequence[float]) – Vector of fixed parameters

setFixedParametersByIdRegex(par_id_regex: str, value: float) int[source]

Set values of all fixed parameters with the ID matching the specified regex.

Parameters:
  • par_id_regex (str) – Fixed parameter name regex

  • value (float) – Fixed parameter value

Return type:

int

Returns:

Number of fixed parameter IDs that matched the regex

setFixedParametersByNameRegex(par_name_regex: str, value: float) int[source]

Set value of all fixed parameters with name matching the specified regex.

Parameters:
  • par_name_regex (str) – Fixed parameter name regex

  • value (float) – Fixed parameter value

Return type:

int

Returns:

Number of fixed parameter names that matched the regex

setInitialStateSensitivities(sx0: Sequence[float])[source]

Set the initial state sensitivities.

Parameters:

sx0 (typing.Sequence[float]) – vector of initial state sensitivities with chainrule applied. This could be a slice of ReturnData::sx or ReturnData::sx0

setInitialStates(x0: Sequence[float])[source]

Set the initial states.

Parameters:

x0 (typing.Sequence[float]) – Initial state vector

setMinimumSigmaResiduals(min_sigma: float)[source]

Sets the estimated lower boundary for sigma_y. When setAddSigmaResiduals() is activated, this lower boundary must ensure that log(sigma) + min_sigma > 0.

Parameters:

min_sigma (float) – lower boundary

setNMaxEvent(nmaxevent: int)[source]

Set maximum number of events that may occur for each type.

Parameters:

nmaxevent (int) – Maximum number of events that may occur for each type

setParameterById(*args)[source]

Overload 1:

Set model parameters according to the parameter IDs and mapped values.

Parameters:
  • p (StringDoubleMap) – Map of parameters IDs and values

  • ignoreErrors (boolean, optional) – Ignore errors such as parameter IDs in p which are not model parameters


Overload 2:

Set value of first model parameter with the specified ID.

Parameters:
  • par_id (str) – Parameter ID

  • value (float) – Parameter value

setParameterByName(*args)[source]

Overload 1:

Set value of first model parameter with the specified name.

Parameters:
  • par_name (str) – Parameter name

  • value (float) – Parameter value


Overload 2:

Set model parameters according to the parameter name and mapped values.

Parameters:
  • p (StringDoubleMap) – Map of parameters names and values

  • ignoreErrors (boolean, optional) – Ignore errors such as parameter names in p which are not model parameters


Overload 3:

Set model parameters according to the parameter name and mapped values.

Parameters:
  • p (StringDoubleMap) – Map of parameters names and values

  • ignoreErrors – Ignore errors such as parameter names in p which are not model parameters

setParameterList(plist: Sequence[int])[source]

Set the list of parameters for which sensitivities are to be computed.

NOTE: Resets initial state sensitivities.

Parameters:

plist (typing.Sequence[int]) – List of parameter indices

setParameterScale(*args)[source]
setParameters(p: Sequence[float])[source]

Set the parameter vector.

Parameters:

p (typing.Sequence[float]) – Vector of parameters

setParametersByIdRegex(par_id_regex: str, value: float) int[source]

Set all values of model parameters with IDs matching the specified regular expression.

Parameters:
  • par_id_regex (str) – Parameter ID regex

  • value (float) – Parameter value

Return type:

int

Returns:

Number of parameter IDs that matched the regex

setParametersByNameRegex(par_name_regex: str, value: float) int[source]

Set all values of all model parameters with names matching the specified regex.

Parameters:
  • par_name_regex (str) – Parameter name regex

  • value (float) – Parameter value

Return type:

int

Returns:

Number of fixed parameter names that matched the regex

setReinitializationStateIdxs(idxs: Sequence[int])[source]

Set indices of states to be reinitialized based on provided constants / fixed parameters

Parameters:

idxs (typing.Sequence[int]) – Array of state indices

setReinitializeFixedParameterInitialStates(flag: bool)[source]

Set whether initial states depending on fixed parameters are to be reinitialized after preequilibration and presimulation.

Parameters:

flag (bool) – Fixed parameters reinitialized?

setStateIsNonNegative(stateIsNonNegative: Sequence[bool])[source]

Set flags indicating whether states should be treated as non-negative.

Parameters:

stateIsNonNegative (typing.Sequence[bool]) – Vector of flags

setSteadyStateComputationMode(mode: int)[source]

Set the mode how steady state is computed in the steadystate simulation.

Parameters:

mode (int) – Steadystate computation mode

setSteadyStateSensitivityMode(mode: SteadyStateSensitivityMode)[source]

Set the mode how sensitivities are computed in the steadystate simulation.

Parameters:

mode (amici.amici.SteadyStateSensitivityMode) – Steadystate sensitivity mode

setT0(t0: float)[source]

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 (float) – Simulation start time

setTimepoints(ts: Sequence[float])[source]

Set the timepoint vector.

Parameters:

ts (typing.Sequence[float]) – New timepoint vector

setUnscaledInitialStateSensitivities(sx0: Sequence[float])[source]

Set the initial state sensitivities.

Parameters:

sx0 (typing.Sequence[float]) – Vector of initial state sensitivities without chainrule applied. This could be the readin from a model.sx0data saved to HDF5.

property state_independent_events_

Map of trigger timepoints to event indices for events that don’t require root-finding.

t0() float[source]

Get simulation start time.

Return type:

float

Returns:

Simulation start time

property ubw

Upper bandwidth of the Jacobian

class amici.amici.ModelDimensions(*args)[source]

Container for model dimensions.

Holds number of states, observables, etc.

__init__(*args)[source]
Overload 1:

Default ctor


Overload 2:

Constructor with model dimensions

Parameters:
  • nx_rdata (int) – Number of state variables

  • nxtrue_rdata (int) – Number of state variables of the non-augmented model

  • nx_solver (int) – Number of state variables with conservation laws applied

  • nxtrue_solver (int) – Number of state variables of the non-augmented model with conservation laws applied

  • nx_solver_reinit (int) – Number of state variables with conservation laws subject to reinitialization

  • np (int) – Number of parameters

  • nk (int) – Number of constants

  • ny (int) – Number of observables

  • nytrue (int) – Number of observables of the non-augmented model

  • nz (int) – Number of event observables

  • nztrue (int) – Number of event observables of the non-augmented model

  • ne (int) – Number of events

  • ne_solver (int) – Number of events that require root-finding

  • nspl (int) – Number of splines

  • nJ (int) – Number of objective functions

  • nw (int) – Number of repeating elements

  • ndwdx (int) – Number of nonzero elements in the x derivative of the repeating elements

  • ndwdp (int) – Number of nonzero elements in the p derivative of the repeating elements

  • ndwdw (int) – Number of nonzero elements in the w derivative of the repeating elements

  • ndxdotdw (int) – Number of nonzero elements in the \(w\) derivative of \(xdot\)

  • ndJydy (IntVector) – Number of nonzero elements in the \(y\) derivative of \(dJy\) (shape nytrue)

  • ndxrdatadxsolver (int) – Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

  • ndxrdatadtcl (int) – Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

  • ndtotal_cldx_rdata (int) – Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

  • nnz (int) – Number of nonzero elements in Jacobian

  • ubw (int) – Upper matrix bandwidth in the Jacobian

  • lbw (int) – Lower matrix bandwidth in the Jacobian

property lbw

Lower bandwidth of the Jacobian

property nJ

Dimension of the augmented objective function for 2nd order ASA

property ndJydy

Number of nonzero elements in the \(y\) derivative of \(dJy\) (dimension nytrue)

property ndtotal_cldx_rdata

Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

property ndwdp

Number of nonzero elements in the p derivative of the repeating elements

property ndwdw

Number of nonzero elements in the w derivative of the repeating elements

property ndwdx

Number of nonzero elements in the x derivative of the repeating elements

property ndxdotdw

Number of nonzero elements in the \(w\) derivative of \(xdot\)

property ndxrdatadtcl

Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

property ndxrdatadxsolver

Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

property ne

Number of events

property ne_solver

Number of events that require root-finding

property nk

Number of constants

property nnz

Number of nonzero entries in Jacobian

property np

Number of parameters

property nspl

Number of spline functions in the model

property nw

Number of common expressions

property nx_rdata

Number of states

property nx_solver

Number of states with conservation laws applied

property nx_solver_reinit

Number of solver states subject to reinitialization

property nxtrue_rdata

Number of states in the unaugmented system

property nxtrue_solver

Number of states in the unaugmented system with conservation laws applied

property ny

Number of observables

property nytrue

Number of observables in the unaugmented system

property nz

Number of event outputs

property nztrue

Number of event outputs in the unaugmented system

property ubw

Upper bandwidth of the Jacobian

class amici.amici.ModelPtr(*args)[source]

Swig-Generated class that implements smart pointers to Model as objects.

property idlist

Flag array for DAE equations

property lbw

Lower bandwidth of the Jacobian

property logger

Logger

property nJ

Dimension of the augmented objective function for 2nd order ASA

property ndJydy

Number of nonzero elements in the \(y\) derivative of \(dJy\) (dimension nytrue)

property ndtotal_cldx_rdata

Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

property ndwdp

Number of nonzero elements in the p derivative of the repeating elements

property ndwdw

Number of nonzero elements in the w derivative of the repeating elements

property ndwdx

Number of nonzero elements in the x derivative of the repeating elements

property ndxdotdw

Number of nonzero elements in the \(w\) derivative of \(xdot\)

property ndxrdatadtcl

Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

property ndxrdatadxsolver

Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

property ne

Number of events

property ne_solver

Number of events that require root-finding

property nnz

Number of nonzero entries in Jacobian

property nspl

Number of spline functions in the model

property nw

Number of common expressions

property nx_rdata

Number of states

property nx_solver

Number of states with conservation laws applied

property nx_solver_reinit

Number of solver states subject to reinitialization

property nxtrue_rdata

Number of states in the unaugmented system

property nxtrue_solver

Number of states in the unaugmented system with conservation laws applied

property ny

Number of observables

property nytrue

Number of observables in the unaugmented system

property nz

Number of event outputs

property nztrue

Number of event outputs in the unaugmented system

property o2mode

Flag indicating whether for amici::Solver::sensi_ == amici::SensitivityOrder::second directional or full second order derivative will be computed

property pythonGenerated

Flag indicating Matlab- or Python-based model generation

property state_independent_events_

Map of trigger timepoints to event indices for events that don’t require root-finding.

property ubw

Upper bandwidth of the Jacobian

class amici.amici.NewtonDampingFactorMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
off = 0
on = 1
class amici.amici.NonlinearSolverIteration(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
fixedpoint = 1
functional = 1
newton = 2
amici.amici.NonlinearSolverIteration_fixedpoint = 1

deprecated

class amici.amici.ObservableScaling(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
lin = 0
log = 1
log10 = 2
class amici.amici.ParameterScaling(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
ln = 1
log10 = 2
none = 0
class amici.amici.ParameterScalingVector(*args)[source]
class amici.amici.RDataReporting(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
full = 0
likelihood = 2
residuals = 1
class amici.amici.ReturnData(*args)[source]

Stores all data to be returned by amici.amici.runAmiciSimulation().

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

property FIM

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

property J

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

__init__(*args)[source]

Overload 1:

Default constructor


Overload 2:

Constructor

Parameters:
  • ts (DoubleVector) – see amici::SimulationParameters::ts

  • model_dimensions (ModelDimensions) – Model dimensions

  • nplist (int) – see amici::ModelDimensions::nplist

  • nmaxevent (int) – see amici::ModelDimensions::nmaxevent

  • nt (int) – see amici::ModelDimensions::nt

  • newton_maxsteps (int) – see amici::Solver::newton_maxsteps

  • pscale (ParameterScalingVector) – see amici::SimulationParameters::pscale

  • o2mode (int) – see amici::SimulationParameters::o2mode

  • sensi (SensitivityOrder) – see amici::Solver::sensi

  • sensi_meth (SensitivityMethod) – see amici::Solver::sensi_meth

  • rdrm (RDataReporting) – see amici::Solver::rdata_reporting

  • quadratic_llh (boolean) – whether model defines a quadratic nllh and computing res, sres and FIM makes sense

  • sigma_res (boolean) – indicates whether additional residuals are to be added for each sigma

  • sigma_offset (float) – offset to ensure real-valuedness of sigma residuals


Overload 3:

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

Parameters:
  • solver (Solver) – solver instance

  • model (Model) – model instance

property chi2

\(\chi^2\) value

property cpu_time

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.

property cpu_timeB

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.

property cpu_time_total

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.

property id

Arbitrary (not necessarily unique) identifier.

property lbw

Lower bandwidth of the Jacobian

property llh

log-likelihood value

property messages

log messages

property nJ

Dimension of the augmented objective function for 2nd order ASA

property ndJydy

Number of nonzero elements in the \(y\) derivative of \(dJy\) (dimension nytrue)

property ndtotal_cldx_rdata

Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

property ndwdp

Number of nonzero elements in the p derivative of the repeating elements

property ndwdw

Number of nonzero elements in the w derivative of the repeating elements

property ndwdx

Number of nonzero elements in the x derivative of the repeating elements

property ndxdotdw

Number of nonzero elements in the \(w\) derivative of \(xdot\)

property ndxrdatadtcl

Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

property ndxrdatadxsolver

Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

property ne

Number of events

property ne_solver

Number of events that require root-finding

property newton_maxsteps

maximal number of newton iterations for steady state calculation

property nk

Number of constants

property nmaxevent

maximal number of occurring events (for every event type)

property nnz

Number of nonzero entries in Jacobian

property np

Number of parameters

property nplist

number of parameter for which sensitivities were requested

property nspl

Number of spline functions in the model

property nt

number of considered timepoints

property numerrtestfails

number of error test failures forward problem (shape nt)

property numerrtestfailsB

number of error test failures backward problem (shape nt)

property numnonlinsolvconvfails

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

property numnonlinsolvconvfailsB

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

property numrhsevals

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

property numrhsevalsB

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

property numsteps

number of integration steps forward problem (shape nt)

property numstepsB

number of integration steps backward problem (shape nt)

property nw

Number of common expressions

property nx

number of states (alias nx_rdata, kept for backward compatibility)

property nx_rdata

Number of states

property nx_solver

Number of states with conservation laws applied

property nx_solver_reinit

Number of solver states subject to reinitialization

property nxtrue

number of states in the unaugmented system (alias nxtrue_rdata, kept for backward compatibility)

property nxtrue_rdata

Number of states in the unaugmented system

property nxtrue_solver

Number of states in the unaugmented system with conservation laws applied

property ny

Number of observables

property nytrue

Number of observables in the unaugmented system

property nz

Number of event outputs

property nztrue

Number of event outputs in the unaugmented system

property o2mode

flag indicating whether second-order sensitivities were requested

property order

employed order forward problem (shape nt)

property posteq_cpu_time

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.

property posteq_cpu_timeB

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.

property posteq_numsteps

number of Newton steps for steady state problem (preequilibration) [newton, simulation, newton] (shape 3) (postequilibration)

property posteq_numstepsB

number of simulation steps for adjoint steady state problem (postequilibration) [== 0 if analytical solution worked, > 0 otherwise]

property posteq_status

flags indicating success of steady state solver (postequilibration)

property posteq_t

time when steadystate was reached via simulation (postequilibration)

property posteq_wrms

weighted root-mean-square of the rhs when steadystate was reached (postequilibration)

property preeq_cpu_time

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.

property preeq_cpu_timeB

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.

property preeq_numsteps

number of Newton steps for steady state problem (preequilibration) [newton, simulation, newton] (length = 3)

property preeq_numstepsB

number of simulation steps for adjoint steady state problem (preequilibration) [== 0 if analytical solution worked, > 0 otherwise]

property preeq_status

flags indicating success of steady state solver (preequilibration)

property preeq_t

time when steadystate was reached via simulation (preequilibration)

property preeq_wrms

weighted root-mean-square of the rhs when steadystate was reached (preequilibration)

property pscale

scaling of parameterization

property rdata_reporting

reporting mode

property res

observable (shape nt*ny, row-major)

property rz

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

property s2llh

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

property s2rz

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

property sensi

sensitivity order

property sensi_meth

sensitivity method

property sigma_res

boolean indicating whether residuals for standard deviations have been added

property sigmay

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

property sigmaz

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

property sllh

parameter derivative of log-likelihood (shape nplist)

property sres

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

property srz

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

property ssigmay

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

property ssigmaz

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

property status

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

property sx

parameter derivative of state (shape nt x nplist x nx, row-major)

property sx0

initial sensitivities (shape nplist x nx, row-major)

property sx_ss

preequilibration sensitivities (shape nplist x nx, row-major)

property sy

parameter derivative of observable (shape nt x nplist x ny, row-major)

property sz

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

property t_last

The final internal time of the solver.

property ts

timepoints (shape nt)

property ubw

Upper bandwidth of the Jacobian

property w

w data from the model (recurring terms in xdot, for imported SBML models from python, this contains the flux vector) (shape nt x nw, row major)

property x

state (shape nt x nx, row-major)

property x0

initial state (shape nx)

property x_ss

preequilibration steady state (shape nx)

property xdot

time derivative (shape nx) evaluated at t_last.

property y

observable (shape nt x ny, row-major)

property z

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

class amici.amici.ReturnDataPtr(*args)[source]

Swig-Generated class that implements smart pointers to ReturnData as objects.

property FIM

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

property J

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

property chi2

\(\chi^2\) value

property cpu_time

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.

property cpu_timeB

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.

property cpu_time_total

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.

property id

Arbitrary (not necessarily unique) identifier.

property lbw

Lower bandwidth of the Jacobian

property llh

log-likelihood value

property messages

log messages

property nJ

Dimension of the augmented objective function for 2nd order ASA

property ndJydy

Number of nonzero elements in the \(y\) derivative of \(dJy\) (dimension nytrue)

property ndtotal_cldx_rdata

Number of nonzero elements in the \(x_rdata\) derivative of \(total_cl\)

property ndwdp

Number of nonzero elements in the p derivative of the repeating elements

property ndwdw

Number of nonzero elements in the w derivative of the repeating elements

property ndwdx

Number of nonzero elements in the x derivative of the repeating elements

property ndxdotdw

Number of nonzero elements in the \(w\) derivative of \(xdot\)

property ndxrdatadtcl

Number of nonzero elements in the \(tcl\) derivative of \(x_rdata\)

property ndxrdatadxsolver

Number of nonzero elements in the \(x\) derivative of \(x_rdata\)

property ne

Number of events

property ne_solver

Number of events that require root-finding

property newton_maxsteps

maximal number of newton iterations for steady state calculation

property nk

Number of constants

property nmaxevent

maximal number of occurring events (for every event type)

property nnz

Number of nonzero entries in Jacobian

property np

Number of parameters

property nplist

number of parameter for which sensitivities were requested

property nspl

Number of spline functions in the model

property nt

number of considered timepoints

property numerrtestfails

number of error test failures forward problem (shape nt)

property numerrtestfailsB

number of error test failures backward problem (shape nt)

property numnonlinsolvconvfails

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

property numnonlinsolvconvfailsB

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

property numrhsevals

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

property numrhsevalsB

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

property numsteps

number of integration steps forward problem (shape nt)

property numstepsB

number of integration steps backward problem (shape nt)

property nw

Number of common expressions

property nx

number of states (alias nx_rdata, kept for backward compatibility)

property nx_rdata

Number of states

property nx_solver

Number of states with conservation laws applied

property nx_solver_reinit

Number of solver states subject to reinitialization

property nxtrue

number of states in the unaugmented system (alias nxtrue_rdata, kept for backward compatibility)

property nxtrue_rdata

Number of states in the unaugmented system

property nxtrue_solver

Number of states in the unaugmented system with conservation laws applied

property ny

Number of observables

property nytrue

Number of observables in the unaugmented system

property nz

Number of event outputs

property nztrue

Number of event outputs in the unaugmented system

property o2mode

flag indicating whether second-order sensitivities were requested

property order

employed order forward problem (shape nt)

property posteq_cpu_time

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.

property posteq_cpu_timeB

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.

property posteq_numsteps

number of Newton steps for steady state problem (preequilibration) [newton, simulation, newton] (shape 3) (postequilibration)

property posteq_numstepsB

number of simulation steps for adjoint steady state problem (postequilibration) [== 0 if analytical solution worked, > 0 otherwise]

property posteq_status

flags indicating success of steady state solver (postequilibration)

property posteq_t

time when steadystate was reached via simulation (postequilibration)

property posteq_wrms

weighted root-mean-square of the rhs when steadystate was reached (postequilibration)

property preeq_cpu_time

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.

property preeq_cpu_timeB

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.

property preeq_numsteps

number of Newton steps for steady state problem (preequilibration) [newton, simulation, newton] (length = 3)

property preeq_numstepsB

number of simulation steps for adjoint steady state problem (preequilibration) [== 0 if analytical solution worked, > 0 otherwise]

property preeq_status

flags indicating success of steady state solver (preequilibration)

property preeq_t

time when steadystate was reached via simulation (preequilibration)

property preeq_wrms

weighted root-mean-square of the rhs when steadystate was reached (preequilibration)

property pscale

scaling of parameterization

property rdata_reporting

reporting mode

property res

observable (shape nt*ny, row-major)

property rz

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

property s2llh

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

property s2rz

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

property sensi

sensitivity order

property sensi_meth

sensitivity method

property sigma_res

boolean indicating whether residuals for standard deviations have been added

property sigmay

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

property sigmaz

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

property sllh

parameter derivative of log-likelihood (shape nplist)

property sres

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

property srz

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

property ssigmay

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

property ssigmaz

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

property status

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

property sx

parameter derivative of state (shape nt x nplist x nx, row-major)

property sx0

initial sensitivities (shape nplist x nx, row-major)

property sx_ss

preequilibration sensitivities (shape nplist x nx, row-major)

property sy

parameter derivative of observable (shape nt x nplist x ny, row-major)

property sz

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

property t_last

The final internal time of the solver.

property ts

timepoints (shape nt)

property ubw

Upper bandwidth of the Jacobian

property w

w data from the model (recurring terms in xdot, for imported SBML models from python, this contains the flux vector) (shape nt x nw, row major)

property x

state (shape nt x nx, row-major)

property x0

initial state (shape nx)

property x_ss

preequilibration steady state (shape nx)

property xdot

time derivative (shape nx) evaluated at t_last.

property y

observable (shape nt x ny, row-major)

property z

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

class amici.amici.SecondOrderMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
directional = 2
full = 1
none = 0
class amici.amici.SensitivityMethod(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
adjoint = 2
forward = 1
none = 0
amici.amici.SensitivityMethod_adjoint = 2

Adjoint sensitivity analysis.

amici.amici.SensitivityMethod_forward = 1

Forward sensitivity analysis.

amici.amici.SensitivityMethod_none = 0

Don’t compute sensitivities.

class amici.amici.SensitivityOrder(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
first = 1
none = 0
second = 2
amici.amici.SensitivityOrder_first = 1

First-order sensitivities.

amici.amici.SensitivityOrder_none = 0

Don’t compute sensitivities.

amici.amici.SensitivityOrder_second = 2

Second-order sensitivities.

class amici.amici.SimulationParameters(*args)[source]

Container for various simulation parameters.

__init__(*args)[source]

Overload 1:

Constructor

Parameters:

timepoints (DoubleVector) – Timepoints for which simulation results are requested


Overload 2:

Constructor

Parameters:
property fixedParameters

Model constants

Vector of size Model::nk() or empty

property fixedParametersPreequilibration

Model constants for pre-equilibration

Vector of size Model::nk() or empty.

property fixedParametersPresimulation

Model constants for pre-simulation

Vector of size Model::nk() or empty.

property parameters

Model parameters

Vector of size Model::np() or empty with parameter scaled according to SimulationParameter::pscale.

property plist

Parameter indices w.r.t. which to compute sensitivities

property pscale

Parameter scales

Vector of parameter scale of size Model::np(), indicating how/if each parameter is to be scaled.

property reinitialization_state_idxs_presim

Indices of states to be reinitialized based on provided presimulation constants / fixed parameters.

property reinitialization_state_idxs_sim

Indices of states to be reinitialized based on provided constants / fixed parameters.

reinitializeAllFixedParameterDependentInitialStates(nx_rdata: int)[source]

Set reinitialization of all states based on model constants for all simulation phases.

Convenience function to populate reinitialization_state_idxs_presim and reinitialization_state_idxs_sim

Parameters:

nx_rdata (int) – Number of states (Model::nx_rdata)

reinitializeAllFixedParameterDependentInitialStatesForPresimulation(nx_rdata: int)[source]

Set reinitialization of all states based on model constants for presimulation (only meaningful if preequilibration is performed).

Convenience function to populate reinitialization_state_idxs_presim and reinitialization_state_idxs_sim

Parameters:

nx_rdata (int) – Number of states (Model::nx_rdata)

reinitializeAllFixedParameterDependentInitialStatesForSimulation(nx_rdata: int)[source]

Set reinitialization of all states based on model constants for the β€˜main’ simulation (only meaningful if presimulation or preequilibration is performed).

Convenience function to populate reinitialization_state_idxs_presim and reinitialization_state_idxs_sim

Parameters:

nx_rdata (int) – Number of states (Model::nx_rdata)

property reinitializeFixedParameterInitialStates

Flag indicating whether reinitialization of states depending on fixed parameters is activated

property sx0

Initial state sensitivities

Dimensions: Model::nx() * Model::nplist(), Model::nx() * ExpData::plist.size(), if ExpData::plist is not empty, or empty

property t_presim

Duration of pre-simulation.

If this is > 0, presimulation will be performed from (model->t0 - t_presim) to model->t0 using the fixedParameters in fixedParametersPresimulation

property ts_

Timepoints for which model state/outputs/… are requested

Vector of timepoints.

property tstart_

Starting time of the simulation.

Output timepoints are absolute timepoints, independent of \(t_{start}\). For output timepoints \(t < t_{start}\), the initial state will be returned.

property x0

Initial state

Vector of size Model::nx() or empty

class amici.amici.Solver(*args, **kwargs)[source]

The Solver class provides a generic interface to CVODES and IDAS solvers, individual realizations are realized in the CVodeSolver and the IDASolver class. All transient private/protected members (CVODES/IDAS memory, interface variables and status flags) are specified as mutable and not included in serialization or equality checks. No solver setting parameter should be marked mutable.

NOTE: Any changes in data members here must be propagated to copy ctor, equality operator, serialization functions in serialization.h, and amici::hdf5::(read/write)SolverSettings(From/To)HDF5 in hdf5.cpp.

__init__(*args, **kwargs)[source]
clone() Solver[source]

Clone this instance

Return type:

amici.amici.Solver

Returns:

The clone

computingASA() bool[source]

check if ASA is being computed

Return type:

bool

Returns:

flag

computingFSA() bool[source]

check if FSA is being computed

Return type:

bool

Returns:

flag

getAbsoluteTolerance() float[source]

Get the absolute tolerances for the forward problem

Same tolerance is used for the backward problem if not specified differently via setAbsoluteToleranceASA.

Return type:

float

Returns:

absolute tolerances

getAbsoluteToleranceB() float[source]

Returns the absolute tolerances for the backward problem for adjoint sensitivity analysis

Return type:

float

Returns:

absolute tolerances

getAbsoluteToleranceFSA() float[source]

Returns the absolute tolerances for the forward sensitivity problem

Return type:

float

Returns:

absolute tolerances

getAbsoluteToleranceQuadratures() float[source]

returns the absolute tolerance for the quadrature problem

Return type:

float

Returns:

absolute tolerance

getAbsoluteToleranceSteadyState() float[source]

returns the absolute tolerance for the steady state problem

Return type:

float

Returns:

absolute tolerance

getAbsoluteToleranceSteadyStateSensi() float[source]

returns the absolute tolerance for the sensitivities of the steady state problem

Return type:

float

Returns:

absolute tolerance

getConstraints() Sequence[float][source]

Get constraints on the model state.

Return type:

typing.Sequence[float]

Returns:

constraints

getInternalSensitivityMethod() InternalSensitivityMethod[source]

returns the internal sensitivity method

Return type:

amici.amici.InternalSensitivityMethod

Returns:

internal sensitivity method

getInterpolationType() InterpolationType[source]
Return type:

amici.amici.InterpolationType

Returns:

getLinearMultistepMethod() LinearMultistepMethod[source]

returns the linear system multistep method

Return type:

amici.amici.LinearMultistepMethod

Returns:

linear system multistep method

getLinearSolver() LinearSolver[source]
Return type:

amici.amici.LinearSolver

Returns:

getMaxConvFails() int[source]

Get the maximum number of nonlinear solver convergence failures permitted per step.

Return type:

int

Returns:

maximum number of nonlinear solver convergence

getMaxNonlinIters() int[source]

Get the maximum number of nonlinear solver iterations permitted per step.

Return type:

int

Returns:

maximum number of nonlinear solver iterations

getMaxStepSize() float[source]

Get the maximum step size

Return type:

float

Returns:

maximum step size

getMaxSteps() int[source]

returns the maximum number of solver steps for the forward problem

Return type:

int

Returns:

maximum number of solver steps

getMaxStepsBackwardProblem() int[source]

returns the maximum number of solver steps for the backward problem

Return type:

int

Returns:

maximum number of solver steps

getMaxTime() float[source]

Returns the maximum time allowed for integration

Return type:

float

Returns:

Time in seconds

getNewtonDampingFactorLowerBound() float[source]

Get a lower bound of the damping factor used in the Newton solver

Return type:

float

Returns:

getNewtonDampingFactorMode() NewtonDampingFactorMode[source]

Get a state of the damping factor used in the Newton solver

Return type:

amici.amici.NewtonDampingFactorMode

Returns:

getNewtonMaxSteps() int[source]

Get maximum number of allowed Newton steps for steady state computation

Return type:

int

Returns:

getNewtonStepSteadyStateCheck() bool[source]

Returns how convergence checks for steadystate computation are performed. If activated, convergence checks are limited to every 25 steps in the simulation solver to limit performance impact.

Return type:

bool

Returns:

boolean flag indicating newton step (true) or the right hand side (false)

getNonlinearSolverIteration() NonlinearSolverIteration[source]

returns the nonlinear system solution method

Return type:

amici.amici.NonlinearSolverIteration

Returns:

getRelativeTolerance() float[source]

Get the relative tolerances for the forward problem

Same tolerance is used for the backward problem if not specified differently via setRelativeToleranceASA.

Return type:

float

Returns:

relative tolerances

getRelativeToleranceB() float[source]

Returns the relative tolerances for the adjoint sensitivity problem

Return type:

float

Returns:

relative tolerances

getRelativeToleranceFSA() float[source]

Returns the relative tolerances for the forward sensitivity problem

Return type:

float

Returns:

relative tolerances

getRelativeToleranceQuadratures() float[source]

Returns the relative tolerance for the quadrature problem

Return type:

float

Returns:

relative tolerance

getRelativeToleranceSteadyState() float[source]

returns the relative tolerance for the steady state problem

Return type:

float

Returns:

relative tolerance

getRelativeToleranceSteadyStateSensi() float[source]

returns the relative tolerance for the sensitivities of the steady state problem

Return type:

float

Returns:

relative tolerance

getReturnDataReportingMode() RDataReporting[source]

returns the ReturnData reporting mode

Return type:

amici.amici.RDataReporting

Returns:

ReturnData reporting mode

getSensiSteadyStateCheck() bool[source]

Returns how convergence checks for steadystate computation are performed.

Return type:

bool

Returns:

boolean flag indicating state and sensitivity equations (true) or only state variables (false).

getSensitivityMethod() SensitivityMethod[source]

Return current sensitivity method

Return type:

amici.amici.SensitivityMethod

Returns:

method enum

getSensitivityMethodPreequilibration() SensitivityMethod[source]

Return current sensitivity method during preequilibration

Return type:

amici.amici.SensitivityMethod

Returns:

method enum

getSensitivityOrder() SensitivityOrder[source]

Get sensitivity order

Return type:

amici.amici.SensitivityOrder

Returns:

sensitivity order

getStabilityLimitFlag() bool[source]

returns stability limit detection mode

Return type:

bool

Returns:

stldet can be false (deactivated) or true (activated)

getStateOrdering() int[source]

Gets KLU / SuperLUMT state ordering mode

Return type:

int

Returns:

State-ordering as integer according to SUNLinSolKLU::StateOrdering or SUNLinSolSuperLUMT::StateOrdering (which differ).

getSteadyStateSensiToleranceFactor() float[source]

returns the steady state sensitivity simulation tolerance factor.

Steady state sensitivity simulation tolerances are the product of the sensitivity simulation tolerances and this factor, unless manually set with set(Absolute/Relative)ToleranceSteadyStateSensi().

Return type:

float

Returns:

steady state simulation tolerance factor

getSteadyStateToleranceFactor() float[source]

returns the steady state simulation tolerance factor.

Steady state simulation tolerances are the product of the simulation tolerances and this factor, unless manually set with set(Absolute/Relative)ToleranceSteadyState().

Return type:

float

Returns:

steady state simulation tolerance factor

property logger
nplist() int[source]

number of parameters with which the solver was initialized

Return type:

int

Returns:

sx.getLength()

nquad() int[source]

number of quadratures with which the solver was initialized

Return type:

int

Returns:

xQB.getLength()

nx() int[source]

number of states with which the solver was initialized

Return type:

int

Returns:

x.getLength()

setAbsoluteTolerance(atol: float)[source]

Sets the absolute tolerances for the forward problem

Same tolerance is used for the backward problem if not specified differently via setAbsoluteToleranceASA.

Parameters:

atol (float) – absolute tolerance (non-negative number)

setAbsoluteToleranceB(atol: float)[source]

Sets the absolute tolerances for the backward problem for adjoint sensitivity analysis

Parameters:

atol (float) – absolute tolerance (non-negative number)

setAbsoluteToleranceFSA(atol: float)[source]

Sets the absolute tolerances for the forward sensitivity problem

Parameters:

atol (float) – absolute tolerance (non-negative number)

setAbsoluteToleranceQuadratures(atol: float)[source]

sets the absolute tolerance for the quadrature problem

Parameters:

atol (float) – absolute tolerance (non-negative number)

setAbsoluteToleranceSteadyState(atol: float)[source]

sets the absolute tolerance for the steady state problem

Parameters:

atol (float) – absolute tolerance (non-negative number)

setAbsoluteToleranceSteadyStateSensi(atol: float)[source]

sets the absolute tolerance for the sensitivities of the steady state problem

Parameters:

atol (float) – absolute tolerance (non-negative number)

setConstraints(constraints: Sequence[float])[source]

Set constraints on the model state.

See https://sundials.readthedocs.io/en/latest/cvode/Usage/index.html#c.CVodeSetConstraints.

Parameters:

constraints (typing.Sequence[float])

setInternalSensitivityMethod(ism: InternalSensitivityMethod)[source]

sets the internal sensitivity method

Parameters:

ism (amici.amici.InternalSensitivityMethod) – internal sensitivity method

setInterpolationType(interpType: InterpolationType)[source]

sets the interpolation of the forward solution that is used for the backwards problem

Parameters:

interpType (amici.amici.InterpolationType) – interpolation type

setLinearMultistepMethod(lmm: LinearMultistepMethod)[source]

sets the linear system multistep method

Parameters:

lmm (amici.amici.LinearMultistepMethod) – linear system multistep method

setLinearSolver(linsol: LinearSolver)[source]
Parameters:

linsol (amici.amici.LinearSolver)

setMaxConvFails(max_conv_fails: int)[source]

Set the maximum number of nonlinear solver convergence failures permitted per step.

Parameters:

max_conv_fails (int) – maximum number of nonlinear solver convergence

setMaxNonlinIters(max_nonlin_iters: int)[source]

Set the maximum number of nonlinear solver iterations permitted per step.

Parameters:

max_nonlin_iters (int) – maximum number of nonlinear solver iterations

setMaxStepSize(max_step_size: float)[source]

Set the maximum step size

Parameters:

max_step_size (float) – maximum step size. 0.0 means no limit.

setMaxSteps(maxsteps: int)[source]

sets the maximum number of solver steps for the forward problem

Parameters:

maxsteps (int) – maximum number of solver steps (positive number)

setMaxStepsBackwardProblem(maxsteps: int)[source]

sets the maximum number of solver steps for the backward problem

Parameters:

maxsteps (int) – maximum number of solver steps (non-negative number)

Notes: default behaviour (100 times the value for the forward problem) can be restored by passing maxsteps=0

setMaxTime(maxtime: float)[source]

Set the maximum CPU time allowed for integration

Parameters:

maxtime (float) – Time in seconds. Zero means infinite time.

setNewtonDampingFactorLowerBound(dampingFactorLowerBound: float)[source]

Set a lower bound of the damping factor in the Newton solver

Parameters:

dampingFactorLowerBound (float)

setNewtonDampingFactorMode(dampingFactorMode: NewtonDampingFactorMode)[source]

Turn on/off a damping factor in the Newton method

Parameters:

dampingFactorMode (amici.amici.NewtonDampingFactorMode)

setNewtonMaxSteps(newton_maxsteps: int)[source]

Set maximum number of allowed Newton steps for steady state computation

Parameters:

newton_maxsteps (int)

setNewtonStepSteadyStateCheck(flag: bool)[source]

Sets how convergence checks for steadystate computation are performed.

Parameters:

flag (bool) – boolean flag to pick newton step (true) or the right hand side (false, default)

setNonlinearSolverIteration(iter: NonlinearSolverIteration)[source]

sets the nonlinear system solution method

Parameters:

iter (amici.amici.NonlinearSolverIteration) – nonlinear system solution method

setRelativeTolerance(rtol: float)[source]

Sets the relative tolerances for the forward problem

Same tolerance is used for the backward problem if not specified differently via setRelativeToleranceASA.

Parameters:

rtol (float) – relative tolerance (non-negative number)

setRelativeToleranceB(rtol: float)[source]

Sets the relative tolerances for the adjoint sensitivity problem

Parameters:

rtol (float) – relative tolerance (non-negative number)

setRelativeToleranceFSA(rtol: float)[source]

Sets the relative tolerances for the forward sensitivity problem

Parameters:

rtol (float) – relative tolerance (non-negative number)

setRelativeToleranceQuadratures(rtol: float)[source]

sets the relative tolerance for the quadrature problem

Parameters:

rtol (float) – relative tolerance (non-negative number)

setRelativeToleranceSteadyState(rtol: float)[source]

sets the relative tolerance for the steady state problem

Parameters:

rtol (float) – relative tolerance (non-negative number)

setRelativeToleranceSteadyStateSensi(rtol: float)[source]

sets the relative tolerance for the sensitivities of the steady state problem

Parameters:

rtol (float) – relative tolerance (non-negative number)

setReturnDataReportingMode(rdrm: RDataReporting)[source]

sets the ReturnData reporting mode

Parameters:

rdrm (amici.amici.RDataReporting) – ReturnData reporting mode

setSensiSteadyStateCheck(flag: bool)[source]

Sets for which variables convergence checks for steadystate computation are performed.

Parameters:

flag (bool) – boolean flag to pick state and sensitivity equations (true, default) or only state variables (false).

setSensitivityMethod(sensi_meth: SensitivityMethod)[source]

Set sensitivity method

Parameters:

sensi_meth (amici.amici.SensitivityMethod)

setSensitivityMethodPreequilibration(sensi_meth_preeq: SensitivityMethod)[source]

Set sensitivity method for preequilibration

Parameters:

sensi_meth_preeq (amici.amici.SensitivityMethod)

setSensitivityOrder(sensi: SensitivityOrder)[source]

Set the sensitivity order

Parameters:

sensi (amici.amici.SensitivityOrder) – sensitivity order

setStabilityLimitFlag(stldet: bool)[source]

set stability limit detection mode

Parameters:

stldet (bool) – can be false (deactivated) or true (activated)

setStateOrdering(ordering: int)[source]

Sets KLU / SuperLUMT state ordering mode

This only applies when linsol is set to LinearSolver::KLU or LinearSolver::SuperLUMT. Mind the difference between SUNLinSolKLU::StateOrdering and SUNLinSolSuperLUMT::StateOrdering.

Parameters:

ordering (int) – state ordering

setSteadyStateSensiToleranceFactor(factor: float)[source]

set the steady state sensitivity simulation tolerance factor.

Steady state sensitivity simulation tolerances are the product of the sensitivity simulation tolerances and this factor, unless manually set with set(Absolute/Relative)ToleranceSteadyStateSensi().

Parameters:

factor (float) – tolerance factor (non-negative number)

setSteadyStateToleranceFactor(factor: float)[source]

set the steady state simulation tolerance factor.

Steady state simulation tolerances are the product of the simulation tolerances and this factor, unless manually set with set(Absolute/Relative)ToleranceSteadyState().

Parameters:

factor (float) – tolerance factor (non-negative number)

class amici.amici.SolverPtr(*args)[source]

Swig-Generated class that implements smart pointers to Solver as objects.

property logger
class amici.amici.SteadyStateComputationMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
integrateIfNewtonFails = 2
integrationOnly = 1
newtonOnly = 0
class amici.amici.SteadyStateSensitivityMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
integrateIfNewtonFails = 2
integrationOnly = 1
newtonOnly = 0
class amici.amici.SteadyStateStatus(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
__init__(*args, **kwds)
failed = -1
failed_convergence = -2
failed_damping = -4
failed_factorization = -3
failed_too_long_simulation = -5
not_run = 0
success = 1
class amici.amici.SteadyStateStatusVector(*args)[source]
class amici.amici.StringDoubleMap(*args)[source]

Swig-Generated class templating Dict [str, float] to facilitate interfacing with C++ bindings.

class amici.amici.StringVector(*args)[source]

Swig-Generated class templating common python types including Iterable [str] and numpy.array [str] to facilitate interfacing with C++ bindings.

amici.amici.compiledWithOpenMP() bool[source]

AMICI extension was compiled with OpenMP?

Return type:

bool

amici.amici.getScaledParameter(unscaledParameter: float, scaling: int) float[source]

Apply parameter scaling according to scaling

Parameters:
  • unscaledParameter (float)

  • scaling (int) – parameter scaling

Return type:

float

Returns:

Scaled parameter

amici.amici.getUnscaledParameter(scaledParameter: float, scaling: int) float[source]

Remove parameter scaling according to scaling

Parameters:
  • scaledParameter (float) – scaled parameter

  • scaling (int) – parameter scaling

Return type:

float

Returns:

Unscaled parameter

amici.amici.parameterScalingFromIntVector(intVec: Sequence[int]) tuple[ParameterScaling][source]

Swig-Generated class, which, in contrast to other Vector classes, does not allow for simple interoperability with common Python types, but must be created using amici.amici.parameterScalingFromIntVector()

Return type:

tuple[amici.amici.ParameterScaling]

amici.amici.runAmiciSimulation(solver: Solver, edata: ExpData, model: Model, rethrow: bool = False) ReturnData[source]

Core integration routine. Initializes the solver and runs the forward and backward problem.

Parameters:
Return type:

amici.amici.ReturnData

Returns:

rdata pointer to return data object

amici.amici.runAmiciSimulations(solver: Solver, edatas: ExpDataPtrVector, model: Model, failfast: bool, num_threads: int) Iterable[ReturnData][source]

Same as runAmiciSimulation, but for multiple ExpData instances. When compiled with OpenMP support, this function runs multi-threaded.

Parameters:
Return type:

typing.Iterable[amici.amici.ReturnData]

Returns:

vector of pointers to return data objects

amici.amici.simulation_status_to_str(status: int) str[source]

Get the string representation of the given simulation status code (see ReturnData::status).

Parameters:

status (int) – Status code

Return type:

str

Returns:

Name of the variable representing this status code.