amici.sbml_import

SBML Import

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

Functions

assignmentRules2observables(sbml_model[, ...])

Turn assignment rules into observables.

get_species_initial(species)

Extract the initial concentration from a given species

replace_logx(math_str)

Replace logX(.) by log(., X) since sympy cannot parse the former

Classes

SbmlImporter(sbml_source[, ...])

Class to generate AMICI C++ files for a model provided in the Systems Biology Markup Language (SBML).

class amici.sbml_import.SbmlImporter(sbml_source, show_sbml_warnings=False, from_file=True, discard_annotations=False)[source]

Class to generate AMICI C++ files for a model provided in the Systems Biology Markup Language (SBML).

Variables:
  • show_sbml_warnings – indicates whether libSBML warnings should be displayed

  • symbols – dict carrying symbolic definitions

  • sbml_reader

    The libSBML sbml reader

    Warning

    Not storing this may result in a segfault.

  • sbml_doc

    document carrying the sbml definition

    Warning

    Not storing this may result in a segfault.

  • sbml – SBML model to import

  • compartments – dict of compartment ids and compartment volumes

  • stoichiometric_matrix – stoichiometric matrix of the model

  • flux_vector – reaction kinetic laws

  • flux_ids – identifiers for elements of flux_vector

  • _local_symbols – model symbols for sympy to consider during sympification see locals`argument in `sympy.sympify

  • species_assignment_rules – Assignment rules for species. Key is symbolic identifier and value is assignment value

  • compartment_assignment_rules – Assignment rules for compartments. Key is symbolic identifier and value is assignment value

  • parameter_assignment_rules – assignment rules for parameters, these parameters are not permissible for sensitivity analysis

  • initial_assignments – initial assignments for parameters, these parameters are not permissible for sensitivity analysis

  • sbml_parser_settings – sets behaviour of SBML Formula parsing

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

Create a new Model instance.

Parameters:
  • sbml_source (typing.Union[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 annotations (debug).

add_d_dt(d_dt, variable, variable0, name)[source]

Creates or modifies species, to implement rate rules for compartments and species, respectively.

Parameters:
Return type:

None

add_local_symbol(key, value)[source]

Add local symbols with some sanity checking for duplication which would indicate redefinition of internals, which SBML permits, but we don’t.

Parameters:
check_event_support()[source]

Check possible events in the model, as AMICI does currently not support :rtype: None

  • delays in events

  • priorities of events

  • events fired at initial time

Furthermore, event triggers are optional (e.g., if an event is fired at initial time, no trigger function is necessary). In this case, warn that this event will have no effect.

check_support()[source]

Check whether all required SBML features are supported. Also ensures that the SBML contains at least one reaction, or rate rule, or assignment rule, to produce change in the system over time.

Return type:

None

is_assignment_rule_target(element)[source]

Checks if an element has a valid assignment rule in the specified model.

Parameters:

element (libsbml.SBase) – SBML variable

Return type:

bool

Returns:

boolean indicating truth of function name

is_rate_rule_target(element)[source]

Checks if an element has a valid assignment rule in the specified model.

Parameters:

element (libsbml.SBase) – SBML variable

Return type:

bool

Returns:

boolean indicating truth of function name

sbml2amici(model_name, output_dir=None, observables=None, event_observables=None, constant_parameters=None, sigmas=None, event_sigmas=None, noise_distributions=None, event_noise_distributions=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, log_as_log10=True, 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 Matlab or C++ as described in the documentation of the respective AMICI interface.

Note that this generates model ODEs for changes in concentrations, not amounts unless the hasOnlySubstanceUnits attribute has been defined 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 (typing.Union[str, pathlib.Path]) – Directory where the generated model package will be stored.

  • observables (dict[str, dict[str, str]]) – Observables to be added to the model: dictionary( observableId:{'name':observableName (optional), 'formula':formulaString)}).

  • event_observables (dict[str, dict[str, str]]) – Event observables to be added to the model: dictionary( eventObservableId:{'name':eventObservableName (optional), 'event':eventId, 'formula':formulaString)})

  • constant_parameters (collections.abc.Iterable[str]) – list of SBML Ids identifying constant parameters

  • sigmas (dict[str, typing.Union[str, float]]) – dictionary(observableId: sigma value or (existing) parameter name)

  • event_sigmas (dict[str, typing.Union[str, float]]) – dictionary(eventObservableId: sigma value or (existing) parameter name)

  • noise_distributions (dict[str, typing.Union[str, typing.Callable]]) – dictionary(observableId: noise type). If nothing is passed for some observable id, a normal model is assumed as default. Either pass a noise type identifier, or a callable generating a custom noise string. For noise identifiers, see amici.import_utils.noise_distribution_to_cost_function().

  • event_noise_distributions (dict[str, typing.Union[str, typing.Callable]]) – dictionary(eventObservableId: noise type). If nothing is passed for some observable id, a normal model is assumed as default. Either pass a noise type identifier, or a callable generating a custom noise string. For noise identifiers, see amici.import_utils.noise_distribution_to_cost_function().

  • verbose (typing.Union[int, bool]) – verbosity level for logging, True/False default 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) – see amici.de_export.ODEExporter

  • 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 (typing.Optional[typing.Callable]) – see amici.ODEModel._simplify

  • cache_simplify (bool) – see amici.ODEModel.__init__()

  • log_as_log10 (bool) – If True, log in the SBML model will be parsed as log10 (default), if False, log will be parsed as natural logarithm ln.

  • 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:

None

amici.sbml_import.assignmentRules2observables(sbml_model, filter_function=<function <lambda>>)[source]

Turn assignment rules into observables.

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

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

Returns:

A dictionary(observableId:{ ‘name’: observableName, ‘formula’: formulaString })

amici.sbml_import.get_species_initial(species)[source]

Extract the initial concentration from a given species

Parameters:

species (libsbml.Species) – species index

Return type:

sympy.core.expr.Expr

Returns:

initial species concentration

amici.sbml_import.replace_logx(math_str)[source]

Replace logX(.) by log(., X) since sympy cannot parse the former

Parameters:

math_str (typing.Union[str, float, None]) – string for sympification

Return type:

typing.Union[str, float, None]

Returns:

sympifiable string