amici.importers.sbml

SBML Import

This module provides all necessary functionality to import a model specified in the Systems Biology Markup Language (SBML).

Functions

assignment_rules_to_observables(sbml_model)

Turn assignment rules into observables.

Classes

SbmlImporter(sbml_source, *[, ...])

Importer for models provided in the Systems Biology Markup Language (SBML).

Exceptions

SBMLException

Exception for SBML related errors.

exception amici.importers.sbml.SBMLException[source]

Exception for SBML related errors.

Raised, for example, for unsupported SBML features.

class amici.importers.sbml.SbmlImporter(sbml_source, *, show_sbml_warnings=False, from_file=True, discard_annotations=False, jax=False)[source]

Importer for models provided in the Systems Biology Markup Language (SBML).

__init__(sbml_source, *, show_sbml_warnings=False, from_file=True, discard_annotations=False, jax=False)[source]

Initialize.

Parameters:
  • sbml_source (str | pathlib.Path | libsbml.Model) – Either a path to SBML file where the model is specified, or a model string as created by sbml.sbmlWriter( ).writeSBMLToString() or an instance of libsbml.Model.

  • show_sbml_warnings (bool) – Indicates whether libSBML warnings should be displayed.

  • from_file (bool) – Whether sbml_source is a file name (True, default), or an SBML string

  • discard_annotations (bool) – discard information contained in AMICI SBML (spline) annotations (debug).

sbml2amici(model_name, output_dir=None, *, fixed_parameters=None, observation_model=None, verbose=40, assume_pow_positivity=False, compiler=None, allow_reinit_fixpar_initcond=True, compile=True, compute_conservation_laws=True, simplify=<function _default_simplify>, cache_simplify=False, generate_sensitivity_code=True, hardcode_symbols=None)[source]

Generate and compile AMICI C++ files for the model provided to the constructor.

The resulting model can be imported as a regular Python module (if compile=True), or used from C++ as described in the documentation of AMICI’s C++ interface.

Note that this generates model ODEs for changes in concentrations, not amounts unless the hasOnlySubstanceUnits attribute has been defined in the SBML model for a particular species.

Sensitivity analysis for local parameters is enabled by creating global parameters _{reactionId}_{localParameterName}.

Parameters:
  • model_name (str) – Name of the generated model package. Note that in a given Python session, only one model with a given name can be loaded at a time. The generated Python extensions cannot be unloaded. Therefore, make sure to choose a unique name for each model.

  • output_dir (str | pathlib.Path) – Directory where the generated model package will be stored. Defaults to amici.get_model_dir().

  • fixed_parameters (collections.abc.Iterable[str]) – SBML Ids to be excluded from sensitivity analysis

  • observation_model (list[amici.importers.utils.MeasurementChannel]) – The different measurement channels that make up the observation model, see MeasurementChannel. If None, default observables will be added (for all state variables of the model and all species, compartments, and assignment rule targets) and normally distributed measurement noise is assumed. If [], no observables will be added.

  • verbose (int | bool) – Verbosity level for logging, True/False defaults to logging.Error/logging.DEBUG.

  • assume_pow_positivity (bool) – if set to True, a special pow function is used to avoid problems with state variables that may become negative due to numerical errors

  • compiler (str) – Absolute path to the compiler executable to be used to build the Python extension, e.g. /usr/bin/clang.

  • allow_reinit_fixpar_initcond (bool) – Indicates whether reinitialization of initial states depending on fixed parameters is allowed for this model. Must be enabled if initial states are to be reset after pre-equilibration.

  • compile (bool) – If True, compile the generated Python package, if False, just generate code.

  • compute_conservation_laws (bool) – if set to True, conservation laws are automatically computed and applied such that the state-jacobian of the ODE right-hand-side has full rank. This option should be set to True when using the Newton algorithm to compute steadystate sensitivities. Conservation laws for constant species are enabled by default. Support for conservation laws for non-constant species is experimental and may be enabled by setting an environment variable AMICI_EXPERIMENTAL_SBML_NONCONST_CLS to either demartino to use the algorithm proposed by De Martino et al. (2014) https://doi.org/10.1371/journal.pone.0100750, or to any other value to use the deterministic algorithm implemented in conserved_moieties2.py. In some cases, the demartino may run for a very long time. This has been observed for example in the case of stoichiometric coefficients with many significant digits.

  • simplify (collections.abc.Callable | None) – If not None, this function will be used to simplify symbolic derivative expressions. Receives sympy expressions as only argument. To apply multiple simplifications, wrap them in a lambda expression.

  • cache_simplify (bool) – Whether to cache calls to the simplify method. Note that there are possible issues with PySB models: https://github.com/AMICI-dev/AMICI/pull/1672

  • generate_sensitivity_code (bool) – If False, the code required for sensitivity computation will not be generated.

  • hardcode_symbols (collections.abc.Sequence[str]) – List of SBML entity IDs that are to be hardcoded in the generated model. Their values cannot be changed anymore after model import. Currently, only parameters that are not targets of rules or initial assignments are supported.

Return type:

amici._installation.amici.Model | None

Returns:

If compile is True and compilation was successful, an instance of the generated model class, otherwise None.

sbml2jax(model_name, output_dir=None, fixed_parameters=None, observation_model=None, verbose=40, compute_conservation_laws=True, simplify=<function _default_simplify>, cache_simplify=False, hybridization=None)[source]

Generate and compile AMICI jax files for the model provided to the constructor.

The resulting model can be imported as a regular Python module.

Note that this generates model ODEs for changes in concentrations, not amounts unless the hasOnlySubstanceUnits attribute has been defined for a particular species.

Parameters:
  • model_name (str) – Name of the generated model package. Note that in a given Python session, only one model with a given name can be loaded at a time. The generated Python extensions cannot be unloaded. Therefore, make sure to choose a unique name for each model.

  • output_dir (str | pathlib.Path) – Directory where the generated model package will be stored.

  • fixed_parameters (collections.abc.Iterable[str]) – SBML Ids to be excluded from sensitivity analysis

  • observation_model (list[amici.importers.utils.MeasurementChannel]) – The different measurement channels that make up the observation model, see amici.importers.utils.MeasurementChannel. Only time-resolved observables are supported here. If None, default observables will be added (for all state variables of the model and all species, compartments, and assignment rule targets) and normally distributed measurement noise is assumed. If [], no observables will be added.

  • verbose (int | bool) – verbosity level for logging, True/False default to logging.Error/logging.DEBUG

  • compute_conservation_laws (bool) – if set to True, conservation laws are automatically computed and applied such that the state-jacobian of the ODE right-hand-side has full rank. This option should be set to True when using the Newton algorithm to compute steadystate sensitivities. Conservation laws for constant species are enabled by default. Support for conservation laws for non-constant species is experimental and may be enabled by setting an environment variable AMICI_EXPERIMENTAL_SBML_NONCONST_CLS to either demartino to use the algorithm proposed by De Martino et al. (2014) https://doi.org/10.1371/journal.pone.0100750, or to any other value to use the deterministic algorithm implemented in conserved_moieties2.py. In some cases, the demartino may run for a very long time. This has been observed for example in the case of stoichiometric coefficients with many significant digits.

  • simplify (collections.abc.Callable | None) – If not None, this function will be used to simplify symbolic derivative expressions. Receives sympy expressions as only argument. To apply multiple simplifications, wrap them in a lambda expression.

  • cache_simplify (bool) – Whether to cache calls to the simplify method. Note that there are possible issues with PySB models: https://github.com/AMICI-dev/AMICI/pull/1672

  • hybridization (dict) – dict representation of the hybridization information in the PEtab YAML file, see https://petab-sciml.readthedocs.io/latest/format.html#problem-yaml-file

Return type:

None

sbml_model: libsbml.Model

the SBML model to import

amici.importers.sbml.assignment_rules_to_observables(sbml_model, filter_function=<function <lambda>>, as_dict=False)[source]

Turn assignment rules into observables.

The matching assignment rules and their targets are removed from the model.

Parameters:
  • sbml_model (libsbml.Model) – Model to operate on

  • filter_function (collections.abc.Callable) – Callback function taking assignment variable as input and returning True/False to indicate if the respective rule should be turned into an observable.

  • as_dict – If True, return a dict instead of a list, using the observable IDs as keys.

Return type:

list[amici.importers.utils.MeasurementChannel] | dict[str, amici.importers.utils.MeasurementChannel]

Returns:

A list or dict of measurement channels.