amici.sbml_import
SBML Import
This module provides all necessary functionality to import a model specified in the Systems Biology Markup Language (SBML).
Functions
|
Turn assignment rules into observables. |
|
Extract the initial concentration from a given species |
|
Replace logX(.) by log(., X) since sympy cannot parse the former |
Classes
|
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 (
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 stringdiscard_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:
d_dt (
sympy.core.expr.Expr
) – The rate rule (or, right-hand side of an ODE).variable (
sympy.core.symbol.Symbol
) – The subject of the rate rule.variable0 (
float
|sympy.core.expr.Expr
) – The initial value of the variable.name (
str
) – Species name, only applicable if this function generates a new species
- Return type:
- 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:
key (
str
) – local symbol keyvalue (
sympy.core.expr.Expr
) – local symbol value
- 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:
- 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:
- 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:
- 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 in the SBML model for a particular species.
Sensitivity analysis for local parameters is enabled by creating global parameters
_{reactionId}_{localParameterName}
.Note
When providing expressions for (event) observables and their sigmas as strings (see below), those will be passed to
sympy.sympify()
. The supported grammar is not well defined. Note there can be issues with, for example,==
or n-ary (n>2) comparison operators. Also note that operator precedence and function names may differ from SBML L3 formulas or PEtab math expressions. Passing a sympy expression directly will be the safer option for more complex expressions.Note
When passing sympy expressions, all Symbols therein must have the
real=True
assumption.Note
In any math expressions passed to this function,
time
will be interpreted as the time symbol.- 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.observables (
dict
[str
,dict
[str
,str
|sympy.core.expr.Expr
]]) –Observables to be added to the model:
dict( observableId: { 'name': observableName, # optional 'formula': formulaString or sympy expression, } )
If the observation function is passed as a string, it will be passed to
sympy.sympify()
(see note above).event_observables (
dict
[str
,dict
[str
,str
|sympy.core.expr.Expr
]]) –Event observables to be added to the model:
dict( eventObservableId: { 'name': eventObservableName, # optional 'event':eventId, 'formula': formulaString or sympy expression, } )
If the formula is passed as a string, it will be passed to
sympy.sympify()
(see note above).constant_parameters (
collections.abc.Iterable
[str
]) – list of SBML Ids identifying constant parameterssigmas (
dict
[str
,str
|float
|sympy.core.expr.Expr
]) –Expression for the scale parameter of the noise distribution for each observable. This can be a numeric value, a sympy expression, or an expression string that will be passed to
sympy.sympify()
.{observableId: sigma}
event_sigmas (
dict
[str
,str
|float
|sympy.core.expr.Expr
]) –Expression for the scale parameter of the noise distribution for each observable. This can be a numeric value, a sympy expression, or an expression string that will be passed to
sympy.sympify()
.{eventObservableId: sigma}
noise_distributions (
dict
[str
,str
|collections.abc.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, seeamici.import_utils.noise_distribution_to_cost_function()
.event_noise_distributions (
dict
[str
,str
|collections.abc.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, seeamici.import_utils.noise_distribution_to_cost_function()
.verbose (
int
|bool
) – Verbosity level for logging,True
/False
defaults tologging.Error
/logging.DEBUG
.assume_pow_positivity (
bool
) – if set toTrue
, a special pow function is used to avoid problems with state variables that may become negative due to numerical errorscompiler (
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
) – Seeamici.de_export.DEExporter
.compile (
bool
) – IfTrue
, compile the generated Python package, ifFalse
, just generate code.compute_conservation_laws (
bool
) – if set toTrue
, 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 toTrue
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 variableAMICI_EXPERIMENTAL_SBML_NONCONST_CLS
to eitherdemartino
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 inconserved_moieties2.py
. In some cases, thedemartino
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
) – Seeamici.DEModel._simplify
.cache_simplify (
bool
) – Seeamici.DEModel.__init__()
.log_as_log10 (
bool
) – IfTrue
, log in the SBML model will be parsed aslog10
(default), ifFalse
, log will be parsed as natural logarithmln
.generate_sensitivity_code (
bool
) – IfFalse
, 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:
- sbml2jax(model_name, output_dir=None, observables=None, sigmas=None, noise_distributions=None, verbose=40, compute_conservation_laws=True, simplify=<function _default_simplify>, cache_simplify=False, log_as_log10=True)[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.observables (
dict
[str
,dict
[str
,str
]]) – Observables to be added to the model:dictionary( observableId:{'name':observableName (optional), 'formula':formulaString)})
.sigmas (
dict
[str
,str
|float
]) – dictionary(observableId: sigma value or (existing) parameter name)noise_distributions (
dict
[str
,str
|collections.abc.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, seeamici.import_utils.noise_distribution_to_cost_function()
.verbose (
int
|bool
) – verbosity level for logging,True
/False
default tologging.Error
/logging.DEBUG
compute_conservation_laws (
bool
) – if set toTrue
, 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 toTrue
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 variableAMICI_EXPERIMENTAL_SBML_NONCONST_CLS
to eitherdemartino
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 inconserved_moieties2.py
. In some cases, thedemartino
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
) – seeamici.DEModel._simplify
cache_simplify (
bool
) – seeamici.DEModel.__init__()
log_as_log10 (
bool
) – IfTrue
, log in the SBML model will be parsed aslog10
(default), ifFalse
, log will be parsed as natural logarithmln
.
- Return type:
- amici.sbml_import.assignmentRules2observables(sbml_model, filter_function=<function <lambda>>)[source]
Turn assignment rules into observables.
- Parameters:
sbml_model (
libsbml.Model
) – Model to operate onfilter_function (
collections.abc.Callable
) – Callback function taking assignment variable as input and returningTrue
/False
to indicate if the respective rule should be turned into an observable.
- Returns:
A dictionary(observableId:{ ‘name’: observableName, ‘formula’: formulaString })