Mechanistic Models

Mechanistic models in chi refer to any deterministic model that describes the evolution of quantities of interest in time. In systems biology such models are often inspired by biological mechanisms, which is why we go with the name chi.MechanisticModel. chi.MechanisticModel by no means have to be mechanism-based though, but may be any function of time that you may deem interesting.

Chi provides two ways to specify mechanistic models: 1. you can use the chi.MechanisticModel base class an implement its methods yourself; 2. you can specify the mechanistic model using the System Biology Markup Language (SBML) and instantiate the model using chi.SBMLModel. For detailed examples how either of those can be done, we refer to the Getting started.

Classes

Detailed API

class chi.MechanisticModel[source]

A base class for time series models of the form

\[\bar{y} = g(\bar{y}, \psi, t),\]

where \(\bar{y}\) are the model outputs, \(\psi\) are the model parameters and \(t\) is the time. \(g\) can be any (deterministic) function.

copy()[source]

Returns a deep copy of the mechanistic model.

enable_sensitivities(enabled, parameter_names=None)[source]

Enables the computation of the model output sensitivities to the model parameters if set to True.

The sensitivities of the model outputs are defined as the partial derviatives of the ouputs \(\bar{y}\) with respect to the model parameters \(\psi\)

Parameters:

enabled (bool) – A boolean flag which enables (True) / disables (False) the computation of sensitivities.

has_sensitivities()[source]

Returns a boolean indicating whether sensitivities have been enabled.

n_outputs()[source]

Returns the number of output dimensions.

By default this is the number of states.

n_parameters()[source]

Returns the number of parameters in the model.

Parameters of the model are initial state values and structural parameter values.

outputs()[source]

Returns the output names of the model.

parameters()[source]

Returns the parameter names of the model.

simulate(parameters, times)[source]

Returns the numerical solution of the model outputs (and optionally the sensitivites) for the specified parameters and times.

The model outputs are returned as a 2 dimensional NumPy array of shape (n_outputs, n_times). If sensitivities are enabled, a tuple is returned with the NumPy array of the model outputs and a NumPy array of the sensitivities of shape (n_times, n_outputs, n_parameters).

Parameters:
  • parameters (list, numpy.ndarray) – An array-like object with values for the model parameters.

  • times (list, numpy.ndarray) – An array-like object with time points at which the output values are returned.

Return type:

np.ndarray of shape (n_outputs, n_times) or (n_times, n_outputs, n_parameters)

supports_dosing()[source]

Returns a boolean whether dose administration with PKPDModel.set_dosing_regimen() is supported by the model.

class chi.SBMLModel(sbml_file)[source]

Instantiates a mechanistic model from a SBML specification.

Extends MechanisticModel.

Parameters:

sbml_file (str) – A path to the SBML model file that specifies the model.

copy()[source]

Returns a deep copy of the mechanistic model.

Note

Copying the model resets the sensitivity settings.

enable_sensitivities(enabled, parameter_names=None)[source]

Enables the computation of the model output sensitivities to the model parameters if set to True.

The sensitivities are computed using the forward sensitivities method, where an ODE for each sensitivity is derived. The sensitivities are returned together with the solution to the orginal system of ODEs when simulating the mechanistic model simulate().

The optional parameter names argument can be used to set which sensitivities are computed. By default the sensitivities to all parameters are computed.

Parameters:
  • enabled (bool) – A boolean flag which enables (True) / disables (False) the computation of sensitivities.

  • parameter_names (list[str], optional) – A list of parameter names of the model. If None sensitivities for all parameters are computed.

has_sensitivities()[source]

Returns a boolean indicating whether sensitivities have been enabled.

n_outputs()[source]

Returns the number of output dimensions.

By default this is the number of states.

n_parameters()[source]

Returns the number of parameters in the model.

Parameters of the model are initial state values and structural parameter values.

outputs()[source]

Returns the output names of the model.

parameters()[source]

Returns the parameter names of the model.

set_output_names(names)[source]

Assigns names to the model outputs. By default the myokit.Model names are assigned to the outputs.

Parameters:

names (dict[str, str]) – A dictionary that maps the current output names to new names.

set_outputs(outputs)[source]

Sets outputs of the model.

The outputs can be set to any quantifiable variable name of the myokit.Model, e.g. compartment.variable.

Note

Setting outputs resets the sensitivity settings (by default sensitivities are disabled.)

Parameters:

outputs (list[str]) – A list of output names.

set_parameter_names(names)[source]

Assigns names to the parameters. By default the myokit.Model names are assigned to the parameters.

Parameters:

names (dict[str, str]) – A dictionary that maps the current parameter names to new names.

simulate(parameters, times)[source]

Returns the numerical solution of the model outputs (and optionally the sensitivites) for the specified parameters and times.

The model outputs are returned as a 2 dimensional NumPy array of shape (n_outputs, n_times). If sensitivities are enabled, a tuple is returned with the NumPy array of the model outputs and a NumPy array of the sensitivities of shape (n_times, n_outputs, n_parameters).

Parameters:
  • parameters (list, numpy.ndarray) – An array-like object with values for the model parameters.

  • times (list, numpy.ndarray) – An array-like object with time points at which the output values are returned.

Return type:

np.ndarray of shape (n_outputs, n_times) or a Tuple of two np.ndarray, the first of shape (n_outputs, n_times) and the second of shape (n_times, n_outputs, n_parameters)

supports_dosing()

Returns a boolean whether dose administration with PKPDModel.set_dosing_regimen() is supported by the model.

time_unit()[source]

Returns the model’s unit of time.

class chi.PKPDModel(sbml_file)[source]

Instantiates a PKPD model from a SBML specification.

Extends SBMLModel.

Parameters:

sbml_file (str) – A path to the SBML model file that specifies the PKPD model.

administration()[source]

Returns the mode of administration in form of a dictionary.

The dictionary has the keys ‘compartment’ and ‘direct’. The former provides information about which compartment is dosed, and the latter whether the dose is administered directly ot indirectly to the compartment.

copy()[source]

Returns a deep copy of the mechanistic model.

Note

Copying the model resets the sensitivity settings.

dosing_regimen()[source]

Returns the dosing regimen of the compound in form of a myokit.Protocol. If the protocol has not been set, None is returned.

enable_sensitivities(enabled, parameter_names=None)[source]

Enables the computation of the model output sensitivities to the model parameters if set to True.

The sensitivities are computed using the forward sensitivities method, where an ODE for each sensitivity is derived. The sensitivities are returned together with the solution to the orginal system of ODEs when simulating the mechanistic model simulate().

The optional parameter names argument can be used to set which sensitivities are computed. By default the sensitivities to all parameters are computed.

Parameters:
  • enabled (bool) – A boolean flag which enables (True) / disables (False) the computation of sensitivities.

  • parameter_names (list[str], optional) – A list of parameter names of the model. If None sensitivities for all parameters are computed.

has_sensitivities()

Returns a boolean indicating whether sensitivities have been enabled.

n_outputs()

Returns the number of output dimensions.

By default this is the number of states.

n_parameters()

Returns the number of parameters in the model.

Parameters of the model are initial state values and structural parameter values.

outputs()

Returns the output names of the model.

parameters()

Returns the parameter names of the model.

set_administration(compartment, amount_var='drug_amount', direct=True)[source]

Sets the route of administration of the compound.

The compound is administered to the selected compartment either directly or indirectly. If it is administered directly, a dose rate variable is added to the drug amount’s rate of change expression

\[\frac{\text{d}A}{\text{d}t} = \text{RHS} + r_d,\]

where \(A\) is the drug amount in the selected compartment, RHS is the rate of change of \(A\) prior to adding the dose rate, and \(r_d\) is the dose rate.

The dose rate can be set by set_dosing_regimen().

If the route of administration is indirect, a dosing compartment is added to the model, which is connected to the selected compartment. The dose rate variable is then added to the rate of change expression of the dose amount variable in the dosing compartment. The drug amount in the dosing compartment flows at a linear absorption rate into the selected compartment

\[\begin{split}\frac{\text{d}A_d}{\text{d}t} = -k_aA_d + r_d \\ \frac{\text{d}A}{\text{d}t} = \text{RHS} + k_aA_d,\end{split}\]

where \(A_d\) is the amount of drug in the dose compartment and \(k_a\) is the absorption rate.

Setting an indirect administration route changes the number of parameters of the model, because an initial dose compartment drug amount and a absorption rate parameter are added.

Parameters:
  • compartment (str) – Compartment to which doses are either directly or indirectly administered.

  • amount_var (str, optional) – Drug amount variable in the compartment. By default the drug amount variable is assumed to be ‘drug_amount’.

  • direct (bool, optional) – A boolean flag that indicates whether the dose is administered directly or indirectly to the compartment.

set_dosing_regimen(dose, start=0, duration=0.01, period=None, num=None)[source]

Sets the dosing regimen with which the compound is administered.

The route of administration can be set with set_administration(). However, the type of administration, e.g. bolus injection or infusion, may be controlled with the duration input.

By default the dose is administered as a bolus injection (duration on a time scale that is 100 fold smaller than the basic time unit). To model an infusion of the dose over a longer time period, the duration can be adjusted to the appropriate time scale.

By default the dose is administered once. To apply multiple doses provide a dose administration period.

Parameters:
  • dose (float or myokit.Protocol) – The amount of the compound that is injected at each administration, or a myokit.Protocol instance that defines the dosing regimen.

  • start (float, optional) – Start time of the treatment. By default the administration starts at t=0.

  • duration (float, optional) – Duration of dose administration. By default the duration is set to 0.01 of the time unit (bolus).

  • period (float, optional) – Periodicity at which doses are administered. If None the dose is administered only once.

  • num (int, optional) – Number of administered doses. If None and the periodicity of the administration is not None, doses are administered indefinitely.

set_output_names(names)

Assigns names to the model outputs. By default the myokit.Model names are assigned to the outputs.

Parameters:

names (dict[str, str]) – A dictionary that maps the current output names to new names.

set_outputs(outputs)

Sets outputs of the model.

The outputs can be set to any quantifiable variable name of the myokit.Model, e.g. compartment.variable.

Note

Setting outputs resets the sensitivity settings (by default sensitivities are disabled.)

Parameters:

outputs (list[str]) – A list of output names.

set_parameter_names(names)

Assigns names to the parameters. By default the myokit.Model names are assigned to the parameters.

Parameters:

names (dict[str, str]) – A dictionary that maps the current parameter names to new names.

simulate(parameters, times)

Returns the numerical solution of the model outputs (and optionally the sensitivites) for the specified parameters and times.

The model outputs are returned as a 2 dimensional NumPy array of shape (n_outputs, n_times). If sensitivities are enabled, a tuple is returned with the NumPy array of the model outputs and a NumPy array of the sensitivities of shape (n_times, n_outputs, n_parameters).

Parameters:
  • parameters (list, numpy.ndarray) – An array-like object with values for the model parameters.

  • times (list, numpy.ndarray) – An array-like object with time points at which the output values are returned.

Return type:

np.ndarray of shape (n_outputs, n_times) or a Tuple of two np.ndarray, the first of shape (n_outputs, n_times) and the second of shape (n_times, n_outputs, n_parameters)

supports_dosing()[source]

Returns a boolean whether dose administration with PKPDModel.set_dosing_regimen() is supported by the model.

time_unit()

Returns the model’s unit of time.

class chi.ReducedMechanisticModel(mechanistic_model)[source]

A wrapper class for a MechanisticModel instance that can be used to fix model parameters to fixed values.

Extends MechanisticModel.

Parameters:

mechanistic_model (chi.MechanisticModel) – A mechanistic model.

copy()[source]

Returns a deep copy of the reduced model.

Note

Copying the model resets the sensitivity settings.

dosing_regimen()[source]

Returns the dosing regimen of the compound in form of a myokit.Protocol. If the protocol has not been set, None is returned.

If the model does not support dose administration, None is returned.

enable_sensitivities(enabled)[source]

Enables the computation of the output sensitivities with respect to the free model parameters.

fix_parameters(name_value_dict)[source]

Fixes the value of model parameters, and effectively removes them as a parameter from the model. Fixing the value of a parameter at None, sets the parameter free again.

Parameters:

name_value_dict – A dictionary with model parameter names as keys, and parameter values as values.

has_sensitivities()[source]

Returns a boolean indicating whether sensitivities have been enabled.

mechanistic_model()[source]

Returns the original mechanistic model.

n_fixed_parameters()[source]

Returns the number of fixed model parameters.

n_outputs()[source]

Returns the number of output dimensions.

By default this is the number of states.

n_parameters()[source]

Returns the number of parameters in the model.

Parameters of the model are initial state values and structural parameter values.

outputs()[source]

Returns the output names of the model.

parameters()[source]

Returns the parameter names of the model.

set_dosing_regimen(dose, start=0, duration=0.01, period=None, num=None)[source]

Sets the dosing regimen with which the compound is administered.

The route of administration can be set with set_administration(). However, the type of administration, e.g. bolus injection or infusion, may be controlled with the duration input.

By default the dose is administered as a bolus injection (duration on a time scale that is 100 fold smaller than the basic time unit). To model an infusion of the dose over a longer time period, the duration can be adjusted to the appropriate time scale.

By default the dose is administered once. To apply multiple doses provide a dose administration period.

Parameters:
  • dose – The amount of the compound that is injected at each administration.

  • start – Start time of the treatment.

  • duration – Duration of dose administration. For a bolus injection, a dose duration of 1% of the time unit should suffice. By default the duration is set to 0.01 (bolus).

  • period – Periodicity at which doses are administered. If None the dose is administered only once.

  • num – Number of administered doses. If None and the periodicity of the administration is not None, doses are administered indefinitely.

set_output_names(names)[source]

Assigns names to the outputs. By default the myokit.Model names are assigned to the outputs.

Parameters:

names – A dictionary that maps the current output names to new names.

set_outputs(outputs)[source]

Sets outputs of the model.

Parameters:

outputs – A list of quantifiable variable names of the myokit.Model, e.g. compartment.variable.

set_parameter_names(names)[source]

Assigns names to the parameters. By default the myokit.Model names are assigned to the parameters.

Parameters:

names – A dictionary that maps the current parameter names to new names.

simulate(parameters, times)[source]

Returns the numerical solution of the model outputs (and optionally the sensitivites) for the specified parameters and times.

The model outputs are returned as a 2 dimensional NumPy array of shape (n_outputs, n_times). If sensitivities are enabled, a tuple is returned with the NumPy array of the model outputs and a NumPy array of the sensitivities of shape (n_times, n_outputs, n_parameters).

Parameters:
  • parameters (list, numpy.ndarray) – An array-like object with values for the model parameters.

  • times (list, numpy.ndarray) – An array-like object with time points at which the output values are returned.

Return type:

np.ndarray of shape (n_outputs, n_times) or (n_times, n_outputs, n_parameters)

supports_dosing()[source]

Returns a boolean whether dose administration with PKPDModel.set_dosing_regimen() is supported by the model.

time_unit()[source]

Returns the model’s unit of time.