amici.amici.Model

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]

Initialize self. See help(type(self)) for accurate signature.

Methods Summary

__init__(*args, **kwargs)

Initialize self.

checkFinite(array, fun)

Check if the given array has only finite elements.

clone()

Clone this instance.

fJrz(nllh, iz, p, k, z, sigmaz)

Model specific implementation of fJrz

fJy(nllh, iy, p, k, y, sigmay, my)

Model specific implementation of fJy

fJz(nllh, iz, p, k, z, sigmaz, mz)

Model specific implementation of fJz

fdJrzdsigma(dJrzdsigma, iz, p, k, rz, sigmaz)

Model specific implementation of fdJrzdsigma

fdJrzdz(dJrzdz, iz, p, k, rz, sigmaz)

Model specific implementation of fdJrzdz

fdJydsigma(dJydsigma, iy, p, k, y, sigmay, my)

Model specific implementation of fdJydsigma

fdJzdsigma(dJzdsigma, iz, p, k, z, sigmaz, mz)

Model specific implementation of fdJzdsigma

fdJzdz(dJzdz, iz, p, k, z, sigmaz, mz)

Model specific implementation of fdJzdz

fdeltaqB(deltaqB, t, x, p, k, h, ip, ie, …)

Model specific implementation of fdeltaqB

fdeltasx(deltasx, t, x, p, k, h, w, ip, ie, …)

Model specific implementation of fdeltasx

fdeltax(deltax, t, x, p, k, h, ie, xdot, …)

Model specific implementation of fdeltax

fdeltaxB(deltaxB, t, x, p, k, h, ie, xdot, …)

Model specific implementation of fdeltaxB

fdrzdp(drzdp, ie, t, x, p, k, h, ip)

Model specific implementation of fdrzdp

fdrzdx(drzdx, ie, t, x, p, k, h)

Model specific implementation of fdrzdx

fdsigmaydp(dsigmaydp, t, p, k, ip)

Model specific implementation of fsigmay

fdsigmazdp(dsigmazdp, t, p, k, ip)

Model specific implementation of fsigmaz

fdydp(dydp, t, x, p, k, h, ip, w, dwdp)

Model specific implementation of fdydp

fdydx(dydx, t, x, p, k, h, w, dwdx)

Model specific implementation of fdydx

fdzdp(dzdp, ie, t, x, p, k, h, ip)

Model specific implementation of fdzdp

fdzdx(dzdx, ie, t, x, p, k, h)

Model specific implementation of fdzdx

frz(rz, ie, t, x, p, k, h)

Model specific implementation of frz

fsdx0()

Compute sensitivity of derivative initial states sensitivities sdx0.

fsigmay(sigmay, t, p, k)

Model specific implementation of fsigmay

fsigmaz(sigmaz, t, p, k)

Model specific implementation of fsigmaz

fsrz(srz, ie, t, x, p, k, h, sx, ip)

Model specific implementation of fsrz

fstau(stau, t, x, p, k, h, sx, ip, ie)

Model specific implementation of fstau

fsx0(sx0, t, x0, p, k, ip)

Model specific implementation of fsx0

fsx0_fixedParameters(sx0, t, x0, p, k, ip)

Model specific implementation of fsx0_fixedParameters

fsz(sz, ie, t, x, p, k, h, sx, ip)

Model specific implementation of fsz

fw(w, t, x, p, k, h, tcl)

Model specific implementation of fw

fx0(x0, t, p, k)

Model specific implementation of fx0

fx0_fixedParameters(x0, t, p, k)

Model specific implementation of fx0_fixedParameters

fy(y, t, x, p, k, h, w)

Model specific implementation of fy

fz(z, ie, t, x, p, k, h)

Model specific implementation of fz

getAlwaysCheckFinite()

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

getAmiciCommit()

Returns the amici commit that was used to generate the model

getAmiciVersion()

Returns the amici version that was used to generate the model

getEventSigma(sigmaz, ie, nroots, t, edata)

Get event-resolved observable standard deviations.

getEventSigmaSensitivity(ssigmaz, ie, …)

Get sensitivities of event-resolved observable standard deviations.

getExpressionIds()

Get IDs of the expression.

getExpressionNames()

Get names of the expressions.

getFixedParameterById(par_id)

Get value of fixed parameter with the specified ID.

getFixedParameterByName(par_name)

Get value of fixed parameter with the specified name.

getFixedParameterIds()

Get IDs of the fixed model parameters.

getFixedParameterNames()

Get names of the fixed model parameters.

getFixedParameters()

Get values of fixed parameters.

getInitialStateSensitivities()

Get the initial states sensitivities.

getInitialStates()

Get the initial states.

getName()

Get the model name.

getObservableIds()

Get IDs of the observables.

getObservableNames()

Get names of the observables.

getObservableSigma(sigmay, it, edata)

Get time-resolved observable standard deviations

getObservableSigmaSensitivity(ssigmay, it, edata)

Sensitivity of time-resolved observable standard deviation.

getParameterById(par_id)

Get value of first model parameter with the specified ID.

getParameterByName(par_name)

Get value of first model parameter with the specified name.

getParameterIds()

Get IDs of the model parameters.

getParameterList()

Get the list of parameters for which sensitivities are computed.

getParameterNames()

Get names of the model parameters.

getParameterScale()

Get parameter scale for each parameter.

getParameters()

Get parameter vector.

getReinitializeFixedParameterInitialStates()

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

getSolver()

Retrieves the solver object

getStateIds()

Get IDs of the model states.

getStateIsNonNegative()

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

getStateNames()

Get names of the model states.

getSteadyStateSensitivityMode()

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

getTimepoint(it)

Get simulation timepoint for time index it.

getTimepoints()

Get the timepoint vector.

getUnobservedEventSensitivity(sz, ie)

Get sensitivity of z at final timepoint.

getUnscaledParameters()

Get parameters with transformation according to parameter scale applied.

hasCustomInitialStateSensitivities()

Return whether custom initial state sensitivities have been set.

hasCustomInitialStates()

Return whether custom initial states have been set.

hasExpressionIds()

Report whether the model has expression IDs set.

hasExpressionNames()

Report whether the model has expression names set.

hasFixedParameterIds()

Report whether the model has fixed parameter IDs set.

hasFixedParameterNames()

Report whether the model has fixed parameter names set.

hasObservableIds()

Report whether the model has observable IDs set.

hasObservableNames()

Report whether the model has observable names set.

hasParameterIds()

Report whether the model has parameter IDs set.

hasParameterNames()

Report whether the model has parameter names set.

hasQuadraticLLH()

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

hasStateIds()

Report whether the model has state IDs set.

hasStateNames()

Report whether the model has state names set.

isFixedParameterStateReinitializationAllowed()

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

k()

Get fixed parameters.

nMaxEvent()

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

ncl()

Get number of conservation laws.

nk()

Get number of constants

np()

Get total number of model parameters.

nplist()

Get number of parameters wrt to which sensitivities are computed.

nt()

Get number of timepoints.

nx_reinit()

Get number of solver states subject to reinitialization.

plist(pos)

Get entry in parameter list by index.

requireSensitivitiesForAllParameters()

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

setAllStatesNonNegative()

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

setAlwaysCheckFinite(alwaysCheck)

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

setFixedParameterById(par_id, value)

Set value of first fixed parameter with the specified ID.

setFixedParameterByName(par_name, value)

Set value of first fixed parameter with the specified name.

setFixedParameters(k)

Set values for constants.

setFixedParametersByIdRegex(par_id_regex, value)

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

setFixedParametersByNameRegex(…)

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

setInitialStateSensitivities(sx0)

Set the initial state sensitivities.

setInitialStates(x0)

Set the initial states.

setNMaxEvent(nmaxevent)

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

setParameterById(*args)

Overload 1:

setParameterByName(*args)

Overload 1:

setParameterList(plist)

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

setParameterScale(*args)

Overload 1:

setParameters(p)

Set the parameter vector.

setParametersByIdRegex(par_id_regex, value)

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

setParametersByNameRegex(par_name_regex, value)

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

setReinitializeFixedParameterInitialStates(flag)

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

setStateIsNonNegative(stateIsNonNegative)

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

setSteadyStateSensitivityMode(mode)

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

setT0(t0)

Set simulation start time.

setTimepoints(ts)

Set the timepoint vector.

setUnscaledInitialStateSensitivities(sx0)

Set the initial state sensitivities.

t0()

Get simulation start time.

updateHeaviside(rootsfound)

Update the Heaviside variables h on event occurrences.

updateHeavisideB(rootsfound)

Updates the Heaviside variables h on event occurrences in the backward problem.

Attributes

app

AMICI application context

idlist

Flag array for DAE equations

lbw

Lower bandwidth of the Jacobian

nJ

Dimension of the augmented objective function for 2nd order ASA

ne

Number of events

nnz

Number of nonzero entries in Jacobian

nw

Number of common expressions

nx_rdata

Number of states

nx_solver

Number of states with conservation laws applied

nx_solver_reinit

Number of solver states subject to reinitialization

nxtrue_rdata

Number of states in the unaugmented system

nxtrue_solver

Number of states in the unaugmented system with conservation laws applied

ny

Number of observables

nytrue

Number of observables in the unaugmented system

nz

Number of event outputs

nztrue

Number of event outputs in the unaugmented system

o2mode

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

pythonGenerated

Flag indicating Matlab- or Python-based model generation

ubw

Upper bandwidth of the Jacobian

Methods

__init__(*args, **kwargs)[source]

Initialize self. See help(type(self)) for accurate signature.

checkFinite(array: Iterable[float], fun: str)int[source]

Check if the given array has only finite elements.

If not, try to give hints by which other fields this could be caused.

Parameters
  • array (typing.Iterable[float]) – Array to check

  • fun (str) – Name of the function that generated the values (for more informative messages).

Return type

int

Returns

amici::AMICI_RECOVERABLE_ERROR if a NaN/Inf value was found, amici::AMICI_SUCCESS otherwise

clone() → Iterable[amici.amici.Model][source]

Clone this instance.

Return type

Model

Returns

The clone

fJrz(nllh: Iterable[float], iz: int, p: Iterable[float], k: Iterable[float], z: Iterable[float], sigmaz: Iterable[float])None[source]

Model specific implementation of fJrz

Parameters
Return type

None

fJy(nllh: Iterable[float], iy: int, p: Iterable[float], k: Iterable[float], y: Iterable[float], sigmay: Iterable[float], my: Iterable[float])None[source]

Model specific implementation of fJy

Parameters
Return type

None

fJz(nllh: Iterable[float], iz: int, p: Iterable[float], k: Iterable[float], z: Iterable[float], sigmaz: Iterable[float], mz: Iterable[float])None[source]

Model specific implementation of fJz

Parameters
Return type

None

fdJrzdsigma(dJrzdsigma: Iterable[float], iz: int, p: Iterable[float], k: Iterable[float], rz: Iterable[float], sigmaz: Iterable[float])None[source]

Model specific implementation of fdJrzdsigma

Parameters
Return type

None

fdJrzdz(dJrzdz: Iterable[float], iz: int, p: Iterable[float], k: Iterable[float], rz: Iterable[float], sigmaz: Iterable[float])None[source]

Model specific implementation of fdJrzdz

Parameters
Return type

None

fdJydsigma(dJydsigma: Iterable[float], iy: int, p: Iterable[float], k: Iterable[float], y: Iterable[float], sigmay: Iterable[float], my: Iterable[float])None[source]

Model specific implementation of fdJydsigma

Parameters
Return type

None

fdJzdsigma(dJzdsigma: Iterable[float], iz: int, p: Iterable[float], k: Iterable[float], z: Iterable[float], sigmaz: Iterable[float], mz: Iterable[float])None[source]

Model specific implementation of fdJzdsigma

Parameters
Return type

None

fdJzdz(dJzdz: Iterable[float], iz: int, p: Iterable[float], k: Iterable[float], z: Iterable[float], sigmaz: Iterable[float], mz: Iterable[float])None[source]

Model specific implementation of fdJzdz

Parameters
Return type

None

fdeltaqB(deltaqB: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], ip: int, ie: int, xdot: Iterable[float], xdot_old: Iterable[float], xB: Iterable[float])None[source]

Model specific implementation of fdeltaqB

Parameters
Return type

None

fdeltasx(deltasx: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], w: Iterable[float], ip: int, ie: int, xdot: Iterable[float], xdot_old: Iterable[float], sx: Iterable[float], stau: Iterable[float])None[source]

Model specific implementation of fdeltasx

Parameters
Return type

None

fdeltax(deltax: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], ie: int, xdot: Iterable[float], xdot_old: Iterable[float])None[source]

Model specific implementation of fdeltax

Parameters
Return type

None

fdeltaxB(deltaxB: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], ie: int, xdot: Iterable[float], xdot_old: Iterable[float], xB: Iterable[float])None[source]

Model specific implementation of fdeltaxB

Parameters
Return type

None

fdrzdp(drzdp: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], ip: int)None[source]

Model specific implementation of fdrzdp

Parameters
Return type

None

fdrzdx(drzdx: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float])None[source]

Model specific implementation of fdrzdx

Parameters
Return type

None

fdsigmaydp(dsigmaydp: Iterable[float], t: float, p: Iterable[float], k: Iterable[float], ip: int)None[source]

Model specific implementation of fsigmay

Parameters
Return type

None

fdsigmazdp(dsigmazdp: Iterable[float], t: float, p: Iterable[float], k: Iterable[float], ip: int)None[source]

Model specific implementation of fsigmaz

Parameters
Return type

None

fdydp(dydp: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], ip: int, w: Iterable[float], dwdp: Iterable[float])None[source]

Model specific implementation of fdydp

Parameters
Return type

None

fdydx(dydx: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], w: Iterable[float], dwdx: Iterable[float])None[source]

Model specific implementation of fdydx

Parameters
Return type

None

fdzdp(dzdp: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], ip: int)None[source]

Model specific implementation of fdzdp

Parameters
Return type

None

fdzdx(dzdx: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float])None[source]

Model specific implementation of fdzdx

Parameters
Return type

None

frz(rz: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float])None[source]

Model specific implementation of frz

Parameters
Return type

None

fsdx0()None[source]

Compute sensitivity of derivative initial states sensitivities sdx0.

Only necessary for DAEs.

Return type

None

fsigmay(sigmay: Iterable[float], t: float, p: Iterable[float], k: Iterable[float])None[source]

Model specific implementation of fsigmay

Parameters
Return type

None

fsigmaz(sigmaz: Iterable[float], t: float, p: Iterable[float], k: Iterable[float])None[source]

Model specific implementation of fsigmaz

Parameters
Return type

None

fsrz(srz: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], sx: Iterable[float], ip: int)None[source]

Model specific implementation of fsrz

Parameters
Return type

None

fstau(stau: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], sx: Iterable[float], ip: int, ie: int)None[source]

Model specific implementation of fstau

Parameters
Return type

None

fsx0(sx0: Iterable[float], t: float, x0: Iterable[float], p: Iterable[float], k: Iterable[float], ip: int)None

Model specific implementation of fsx0

Parameters
Return type

None

fsx0_fixedParameters(sx0: Iterable[float], t: float, x0: Iterable[float], p: Iterable[float], k: Iterable[float], ip: int)None

Model specific implementation of fsx0_fixedParameters

Parameters
Return type

None

fsz(sz: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], sx: Iterable[float], ip: int)None[source]

Model specific implementation of fsz

Parameters
Return type

None

fw(w: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], tcl: Iterable[float])None[source]

Model specific implementation of fw

Parameters
Return type

None

fx0(x0: Iterable[float], t: float, p: Iterable[float], k: Iterable[float])None

Model specific implementation of fx0

Parameters
Return type

None

fx0_fixedParameters(x0: Iterable[float], t: float, p: Iterable[float], k: Iterable[float])None

Model specific implementation of fx0_fixedParameters

Parameters
Return type

None

fy(y: Iterable[float], t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float], w: Iterable[float])None[source]

Model specific implementation of fy

Parameters
Return type

None

fz(z: Iterable[float], ie: int, t: float, x: Iterable[float], p: Iterable[float], k: Iterable[float], h: Iterable[float])None[source]

Model specific implementation of fz

Parameters
Return type

None

getAlwaysCheckFinite()bool[source]

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

Return type

boolean

Returns

that

getAmiciCommit()str

Returns the amici commit that was used to generate the model

Return type

string

Returns

ver amici commit string

getAmiciVersion()str

Returns the amici version that was used to generate the model

Return type

string

Returns

ver amici version string

getEventSigma(sigmaz: Iterable[float], ie: int, nroots: int, t: float, edata: amici.amici.ExpData)None[source]

Get event-resolved observable standard deviations.

Parameters
Return type

None

getEventSigmaSensitivity(ssigmaz: Iterable[float], ie: int, nroots: int, t: float, edata: amici.amici.ExpData)None[source]

Get sensitivities of event-resolved observable standard deviations.

Total derivative (only forward sensitivities).

Parameters
  • ssigmaz (typing.Iterable[float]) – Output buffer (dimension: nz x nplist, row-major)

  • ie (int) – Event index

  • nroots (int) – Event occurrence

  • t (float) – Timepoint

  • edata (amici.amici.ExpData) – Pointer to experimental data (optional, pass nullptr to ignore)

Return type

None

getExpressionIds()amici.amici.StringVector[source]

Get IDs of the expression.

Return type

StringVector

Returns

Expression IDs

getExpressionNames()amici.amici.StringVector[source]

Get names of the expressions.

Return type

StringVector

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()amici.amici.StringVector[source]

Get IDs of the fixed model parameters.

Return type

StringVector

Returns

Fixed parameter IDs

getFixedParameterNames()amici.amici.StringVector[source]

Get names of the fixed model parameters.

Return type

StringVector

Returns

Fixed parameter names

getFixedParameters()amici.amici.DoubleVector[source]

Get values of fixed parameters.

Return type

DoubleVector

Returns

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

getInitialStateSensitivities()amici.amici.DoubleVector[source]

Get the initial states sensitivities.

Return type

DoubleVector

Returns

vector of initial state sensitivities

getInitialStates()amici.amici.DoubleVector[source]

Get the initial states.

Return type

DoubleVector

Returns

Initial state vector

getName()str[source]

Get the model name.

Return type

string

Returns

Model name

getObservableIds()amici.amici.StringVector[source]

Get IDs of the observables.

Return type

StringVector

Returns

Observable IDs

getObservableNames()amici.amici.StringVector[source]

Get names of the observables.

Return type

StringVector

Returns

Observable names

getObservableSigma(sigmay: Iterable[float], it: int, edata: amici.amici.ExpData)None[source]

Get time-resolved observable standard deviations

Parameters
Return type

None

getObservableSigmaSensitivity(ssigmay: Iterable[float], it: int, edata: amici.amici.ExpData)None[source]

Sensitivity of time-resolved observable standard deviation.

Total derivative (can be used with both adjoint and forward sensitivity).

Parameters
  • ssigmay (typing.Iterable[float]) – Buffer (dimension: ny x nplist, row-major)

  • it (int) – Timepoint index

  • edata (amici.amici.ExpData) – Pointer to experimental data instance (optional, pass nullptr to ignore)

Return type

None

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()amici.amici.StringVector[source]

Get IDs of the model parameters.

Return type

StringVector

Returns

Parameter IDs

getParameterList()amici.amici.IntVector[source]

Get the list of parameters for which sensitivities are computed.

Return type

IntVector

Returns

List of parameter indices

getParameterNames()amici.amici.StringVector[source]

Get names of the model parameters.

Return type

StringVector

Returns

The parameter names

getParameterScale()amici.amici.ParameterScalingVector[source]

Get parameter scale for each parameter.

Return type

ParameterScalingVector >

Returns

Vector of parameter scales

getParameters()amici.amici.DoubleVector[source]

Get parameter vector.

Return type

DoubleVector

Returns

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

getReinitializeFixedParameterInitialStates()bool[source]

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

Return type

boolean

Returns

flag true / false

getSolver()amici.amici.Solver

Retrieves the solver object

Return type

Solver

Returns

The Solver instance

getStateIds()amici.amici.StringVector[source]

Get IDs of the model states.

Return type

StringVector

Returns

Sate IDs

getStateIsNonNegative()amici.amici.BoolVector[source]

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

Return type

BoolVector

Returns

Vector of flags

getStateNames()amici.amici.StringVector[source]

Get names of the model states.

Return type

StringVector

Returns

State names

getSteadyStateSensitivityMode()amici.amici.SteadyStateSensitivityMode[source]

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

Return type

int

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()amici.amici.DoubleVector[source]

Get the timepoint vector.

Return type

DoubleVector

Returns

Timepoint vector

getUnobservedEventSensitivity(sz: Iterable[float], ie: int)None[source]

Get sensitivity of z at final timepoint.

Ignores sensitivity of timepoint. Total derivative.

Parameters
Return type

None

getUnscaledParameters()amici.amici.DoubleVector[source]

Get parameters with transformation according to parameter scale applied.

Return type

DoubleVector

Returns

Unscaled parameters

hasCustomInitialStateSensitivities()bool[source]

Return whether custom initial state sensitivities have been set.

Return type

boolean

Returns

true if has custom initial state sensitivities, otherwise false.

hasCustomInitialStates()bool[source]

Return whether custom initial states have been set.

Return type

boolean

Returns

true if has custom initial states, otherwise false

hasExpressionIds()bool[source]

Report whether the model has expression IDs set.

Return type

boolean

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

boolean

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

boolean

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

boolean

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

boolean

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

boolean

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

boolean

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

boolean

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

boolean

Returns

boolean flag

hasStateIds()bool[source]

Report whether the model has state IDs set.

Return type

boolean

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

boolean

Returns

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

isFixedParameterStateReinitializationAllowed()bool

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

Return type

boolean

Returns

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

k() → Iterable[float][source]

Get fixed parameters.

Return type

float

Returns

Pointer to constants array

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

nk()int[source]

Get number of constants

Return type

int

Returns

Length of constant vector

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

nt()int[source]

Get number of timepoints.

Return type

int

Returns

Number of timepoints

nx_reinit()int[source]

Get number of solver states subject to reinitialization.

Return type

int

Returns

Model member nx_solver_reinit

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

requireSensitivitiesForAllParameters()None[source]

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

NOTE: Resets initial state sensitivities.

Return type

None

setAllStatesNonNegative()None[source]

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

Return type

None

setAlwaysCheckFinite(alwaysCheck: bool)None[source]

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

Parameters

alwaysCheck (bool) –

Return type

None

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

Set value of first fixed parameter with the specified ID.

Parameters
  • par_id (str) – Fixed parameter id

  • value (float) – Fixed parameter value

Return type

None

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

Set value of first fixed parameter with the specified name.

Parameters
  • par_name (str) – Fixed parameter ID

  • value (float) – Fixed parameter value

Return type

None

setFixedParameters(k: amici.amici.DoubleVector)None[source]

Set values for constants.

Parameters

k (amici.amici.DoubleVector) – Vector of fixed parameters

Return type

None

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: amici.amici.DoubleVector)None[source]

Set the initial state sensitivities.

Parameters

sx0 (amici.amici.DoubleVector) – vector of initial state sensitivities with chainrule applied. This could be a slice of ReturnData::sx or ReturnData::sx0

Return type

None

setInitialStates(x0: amici.amici.DoubleVector)None[source]

Set the initial states.

Parameters

x0 (amici.amici.DoubleVector) – Initial state vector

Return type

None

setNMaxEvent(nmaxevent: int)None[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

Return type

None

setParameterById(*args)None[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) – 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 (string) – Parameter ID

  • value (float) – Parameter value

Return type

None

setParameterByName(*args)None[source]

Overload 1:

Set value of first model parameter with the specified name.

Parameters
  • par_name (string) – 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) – 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

Return type

None

setParameterList(plist: amici.amici.IntVector)None[source]

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

NOTE: Resets initial state sensitivities.

Parameters

plist (amici.amici.IntVector) – List of parameter indices

Return type

None

setParameterScale(*args)None[source]

Overload 1:

Set parameter scale for each parameter.

NOTE: Resets initial state sensitivities.

Parameters

pscale (int) – Scalar parameter scale to be set for all parameters


Overload 2:

Set parameter scale for each parameter.

NOTE: Resets initial state sensitivities.

Parameters

pscaleVec (ParameterScalingVector >) – Vector of parameter scales

Return type

None

setParameters(p: amici.amici.DoubleVector)None[source]

Set the parameter vector.

Parameters

p (amici.amici.DoubleVector) – Vector of parameters

Return type

None

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

setReinitializeFixedParameterInitialStates(flag: bool)None[source]

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

Parameters

flag (bool) – Fixed parameters reinitialized?

Return type

None

setStateIsNonNegative(stateIsNonNegative: amici.amici.BoolVector)None[source]

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

Parameters

stateIsNonNegative (amici.amici.BoolVector) – Vector of flags

Return type

None

setSteadyStateSensitivityMode(mode: amici.amici.SteadyStateSensitivityMode)None[source]

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

Parameters

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

Return type

None

setT0(t0: float)None[source]

Set simulation start time.

Parameters

t0 (float) – Simulation start time

Return type

None

setTimepoints(ts: amici.amici.DoubleVector)None[source]

Set the timepoint vector.

Parameters

ts (amici.amici.DoubleVector) – New timepoint vector

Return type

None

setUnscaledInitialStateSensitivities(sx0: amici.amici.DoubleVector)None[source]

Set the initial state sensitivities.

Parameters

sx0 (amici.amici.DoubleVector) – Vector of initial state sensitivities without chainrule applied. This could be the readin from a model.sx0data saved to HDF5.

Return type

None

t0()float[source]

Get simulation start time.

Return type

float

Returns

Simulation start time

updateHeaviside(rootsfound: amici.amici.IntVector)None[source]

Update the Heaviside variables h on event occurrences.

Parameters

rootsfound (amici.amici.IntVector) – Provides the direction of the zero-crossing, so adding it will give the right update to the Heaviside variables (zero if no root was found)

Return type

None

updateHeavisideB(rootsfound: Iterable[int])None[source]

Updates the Heaviside variables h on event occurrences in the backward problem.

Parameters

rootsfound (typing.Iterable[int]) – Provides the direction of the zero-crossing, so adding it will give the right update to the Heaviside variables (zero if no root was found)

Return type

None