This Jupyter Notebook provides detailed, executable examples that demonstrate how to use pyomexmeta, the Python interface to the libOmexMeta library. This notebook refers to the OMEX Metadata Specification document which describes the complete specification for annotations of biosimulation models. Note that we anticipate that users of libOmexMeta will be programmers rather than biologists.
These examples demonstrate the construction of a variety of different types of annotations as defined in the OMEX Metadata Specification.
Section 1 demonstrates the creation of a singular annotations. Singular annotations are the lowest level of the pyomexmeta annotation API where the user wants to define the annotation triple directly.
Section 2 demonstrates the use of convenience methods for the creation of model-level annotations. Model-level annotations are knowledge about the whole model being annotated - for example, the author, related publication, or the relevante species of the model.
Section 3 demonstrates how to create the range of composite annotations defined in the OMEX Metadata Specification. Composite annotations are collections of singular annotations with a predefined structure and are used to construct precise physiological concepts that are not represented as singular terms in current knowledge resources. This section is divided into two parts, where we demonstrate these types of annotation first with models encoded in SBML and second with models encoded in CellML.
Note that in the following examples the SBML and CellML model strings are not aiming to be full models and may not be complete. They are purely defined to provide sufficient model snippets to enable the demonstration of pyomexmeta functionality.
In this section, we demonstrate how one can create a single annotation, which is an RDF triple object comprising a subject, a predicate and an object.
The SBML model that we use in this example is a basic toy model where A
transitions to B
using first order mass action kinetics.
# Example 1
from pyomexmeta import RDF, eUriType
# Define an SBML model to annotate
sbml1 = """<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version1/core" level="3" version="1">
<model metaid="ToyModel" id="ToyModel">
<listOfCompartments>
<compartment sboTerm="SBO:0000410" id="cell" metaid="cell" spatialDimensions="3" size="1" constant="true"/>
</listOfCompartments>
<listOfSpecies>
<species id="A" metaid="A" compartment="cell" initialConcentration="10" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="B" metaid="B" compartment="cell" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
</listOfSpecies>
<listOfParameters>
<parameter id="k1" value="0.1" constant="true"/>
</listOfParameters>
<listOfReactions>
<reaction id="R1" metaid="R1" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="A" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="B" stoichiometry="1" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> k1 </ci>
<ci> A </ci>
</apply>
</math>
</kineticLaw>
</reaction>
</listOfReactions>
</model>
</sbml>
"""
rdf_graph1 = RDF()
rdf_graph1.set_archive_uri("MyOmexArchive.omex")
rdf_graph1.set_model_uri("MyModel.sbml")
annot_editor = rdf_graph1.to_editor(sbml1, generate_new_metaids=False, sbml_semantic_extraction=False)
with annot_editor.new_singular_annotation() as singular_annotation:
# CHEBI:16236 = EtOH in chebi
singular_annotation.about("A") \
.predicate("bqbiol", "is") \
.resource_uri("CHEBI:16236")
print(rdf_graph1)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/MyOmexArchive.omex/MyModel.rdf#> . <http://omex-library.org/MyOmexArchive.omex/MyModel.sbml#A> bqbiol:is <https://identifiers.org/CHEBI:16236> .
This example shows the lowest level interface for creating an annotation. These three methods (about()
, predicate()
and resource_uri()
) demonstrate how to create the subject, predicate and resource (reference) parts of an RDF triple object respectively. Note, that there also exists a resource_literal()
and a resource_blank()
for creating RDF resources that are literal or blank nodes respectively.
This annotation states that the SBML species with metaid A
is the chemical entity that has the CHEBI identifier CHEBI:16236
, which is ethanol.
Here we demonstrate how to create annotations to indicate the creator, curator, taxon, publication ID, and creation date for an SBML model. We use the same SBML model from Section 1 above.
# example 2
from pyomexmeta import RDF, eUriType
# reuse SBML model from above
sbml2 = sbml1
rdf_graph2 = RDF()
rdf_graph2.set_archive_uri("MyOmexArchive.omex")
rdf_graph2.set_model_uri("MyModel.sbml")
annot_editor = rdf_graph2.to_editor(sbml2, generate_new_metaids=False, sbml_semantic_extraction=False)
annot_editor.add_creator("0000-0001-8254-4957") \
.add_curator("0000-0002-2390-6572") \
.add_taxon("9895") \
.add_pubmed("12334") \
.add_date_created("18-09-2020")
print(rdf_graph2)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix pubmed: <https://identifiers.org/pubmed:> . @prefix dc: <https://dublincore.org/specifications/dublin-core/dcmi-terms/> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix NCBI_Taxon: <https://identifiers.org/taxonomy:> . @prefix bqmodel: <http://biomodels.net/model-qualifiers/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/MyOmexArchive.omex/MyModel.rdf#> . <http://omex-library.org/MyOmexArchive.omex/MyModel.rdf#> dc:creator <https://orcid.org/0000-0002-2390-6572> . <http://omex-library.org/MyOmexArchive.omex/MyModel.sbml> bqbiol:hasTaxon <https://identifiers.org/taxonomy:9895> ; bqmodel:isDescribedBy <https://identifiers.org/pubmed:12334> ; dc:created [ dc:W3CDTF "18-09-2020" ] ; dc:creator <https://orcid.org/0000-0001-8254-4957> .
As described in the Omex Metadata v1.1 specification, it is also optionally possible to specify the people (curators, creators) with strings for names, email addresses, etc.
Composite annotations are collections of singular annotations with a predefined structure and are used to construct precise physiological concepts that are not represented as singular terms in current knowledge resources. Three types of composite annotations are defined in the OMEX Metadata Specification (version 1.1) and supported by libOmexMeta: 1) annotation for a physical entity; 2) a physical process; and 3) an energy differential. In this Section 3.1 we provide examples of each of these using models encoded in the SBML format and in Section 3.2 models encoded in the CellML format.
In many modeling situations, one may wish to annotate both the physical entity and the physical property that is being simulated, such as concentration, particle numbers, fluid volume or fluid pressure. This allows us to capture more complete biological semantics. For example, a cardiovascular fluid flow model might simulate the fluid volume or fluid pressure inside the left ventricle, or a chemical kinetic model might simulate the chemical concentration or the amount (in moles, perhaps) of a physical entity, such as glucose or ethanol.
In some modeling languages, only one of these aspects may be explicitly represented in the code. For example, in CellML, one has variables that represent chemical concentrations, or blood pressures, but the biological physical entity is implicit. Conversely, in SBML, the physical entity is specified as a “species”, but the physical property (usually chemical concentration or particle numbers) is left implicit. Thus, for annotation, we must allow for either situation.
A physical entity may be a complex structured object. For example, one may need to distinguish between nuclear and cytoplasmic
glucose. These are both the same entity (glucose), but they are part of different anatomic structures. In the following example we demonstrate how to annotate such a situation using another toy model depicting glucose in dynamic equilibrium between two
compartments, the nucleus
and the cytoplasm
.
# Example 3.1.1
from pyomexmeta import RDF, eUriType
sbml3_1_1 = """<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version1/core" level="3" version="1">
<model metaid="GlucoseTransport" id="GlucoseTransport">
<listOfCompartments>
<compartment id="cytoplasm" metaid="cytoplasm" spatialDimensions="3" size="1" constant="true"/>
<compartment id="nucleus" metaid="nucleus" spatialDimensions="3" size="1" constant="true"/>
</listOfCompartments>
<listOfSpecies>
<species id="glucose_c" metaid="glucose_c" compartment="cytoplasm" initialConcentration="10" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="glucose_n" metaid="glucose_n" compartment="nucleus" initialConcentration="100" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
</listOfSpecies>
<listOfParameters>
<parameter id="kimp" metaid="kimp" value="10" constant="true"/>
<parameter id="kexp" metaid="kexp" value="0.1" constant="true"/>
</listOfParameters>
<listOfReactions>
<reaction id="r1_imp" metaid="r1_imp" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="glucose_c" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="glucose_n" stoichiometry="1" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> cytoplasm </ci>
<ci> kimp </ci>
<ci> glucose_c </ci>
</apply>
</math>
</kineticLaw>
</reaction>
<reaction id="r2_exp" metaid="r2_exp" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="glucose_n" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="glucose_c" stoichiometry="1" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> nucleus </ci>
<ci> kexp </ci>
<ci> glucose_n </ci>
</apply>
</math>
</kineticLaw>
</reaction>
</listOfReactions>
</model>
</sbml>
"""
rdf_graph3_1_1 = RDF()
rdf_graph3_1_1.set_archive_uri("GlucoseTransport3_1_1.omex")
rdf_graph3_1_1.set_model_uri("glucose_transport3_1_1.sbml")
annot_editor = rdf_graph3_1_1.to_editor(sbml3_1_1, generate_new_metaids=False, sbml_semantic_extraction=False)
with annot_editor.new_physical_entity() as cytosolic_glucose:
# CHEBI:17234 = glucose
# GO:0005737 = GO cellular component term for cytoplasm.
# OPB:00340 = OPB term for chemical concentration
cytosolic_glucose.about("glucose_c", eUriType.MODEL_URI) \
.identity("CHEBI:17234")\
.is_part_of("GO:0005737") \
.has_property("OPB:OPB_00340")
with annot_editor.new_physical_entity() as nuclear_glucose:
# CHEBI:17234 = The GO cellular component term for nucleus.
nuclear_glucose.about("glucose_n", eUriType.MODEL_URI) \
.identity("CHEBI:17234")\
.is_part_of("GO:0005634") \
.has_property("OPB:OPB_00340")
print(rdf_graph3_1_1)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/GlucoseTransport3_1_1.omex/glucose_transport3_1_1.rdf#> . local:EntityProperty0000 bqbiol:isPropertyOf <http://omex-library.org/GlucoseTransport3_1_1.omex/glucose_transport3_1_1.sbml#glucose_c> ; bqbiol:isVersionOf <https://identifiers.org/OPB:OPB_00340> . local:EntityProperty0001 bqbiol:isPropertyOf <http://omex-library.org/GlucoseTransport3_1_1.omex/glucose_transport3_1_1.sbml#glucose_n> ; bqbiol:isVersionOf <https://identifiers.org/OPB:OPB_00340> . <http://omex-library.org/GlucoseTransport3_1_1.omex/glucose_transport3_1_1.sbml#glucose_c> bqbiol:is <https://identifiers.org/CHEBI:17234> ; bqbiol:isPartOf <https://identifiers.org/GO:0005737> . <http://omex-library.org/GlucoseTransport3_1_1.omex/glucose_transport3_1_1.sbml#glucose_n> bqbiol:is <https://identifiers.org/CHEBI:17234> ; bqbiol:isPartOf <https://identifiers.org/GO:0005634> .
In the above example we construct two physical entity annotations, one each for cytoplasmic and nuclear glucose. The argument to the about()
method is the metadata ID (metaid) for the SBML species element called glucose_c
. Internally, libOmexMeta creates a new RDF subject node local:EntityProperty0000
which is local to the annotation document we are creating.
The is_part_of
relation is clearly critical when an entity can exist in multiple compartments and the modeler needs to distinguish between them. However, this relation may be useful in other settings as well.
The has_property
relation may be used to indicate the property we are simulating. In this case, we are simulating the chemical concentration of the glucose entity. Often, this property is not made explicit; the use of the has_property
annotation is optional.
Perhaps the most prevalent type of biosimulation model is one that simulates biochemical reactions. Thus, we provide special attention to annotations for biological processes, such as reactions. Unlike physical entities, processes rarely have unique names. However, some processes can be categorized as a "version of" a Gene Ontology process term. Thus, in the main text, we used a "version of" alcohol dehydrogenase. Alcohol dehydrogenase ADH is the enzyme that catalyzes the conversion of ethanol EtOH into an aldehyde Aldehyde entity. Moreover, NAD+ gains an a H+ in the process as well as two electrons, thus neutralizing the charge on NAHD, the second product of this reaction. This reaction is depicted in the sbml model below with simple mass action kinetics.
# Example 3.1.2
from pyomexmeta import RDF, eUriType
sbml3_1_2_1 = """<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version1/core" level="3" version="1">
<model metaid="ADHModel" id="ADHModel">
<listOfCompartments>
<compartment id="cytosol" metaid="cytosol" spatialDimensions="3" size="1" constant="true"/>
</listOfCompartments>
<listOfSpecies>
<species id="NAD" metaid="NAD" compartment="cytosol" initialConcentration="10" hasOnlySubstanceUnits="false"
boundaryCondition="false" constant="false"/>
<species id="EtOH" metaid="EtOH" compartment="cytosol" initialConcentration="1"
hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="NADH" metaid="NADH" compartment="cytosol" initialConcentration="1"
hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="Aldehyde" metaid="Aldehyde" compartment="cytosol" initialConcentration="1"
hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="ADH" metaid="ADH" compartment="cytosol" initialConcentration="1"
hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
</listOfSpecies>
<listOfParameters>
<parameter id="vmax" metaid="k" value="0.5" constant="true"/>
</listOfParameters>
<listOfReactions>
<reaction id="ADHForwardReaction" metaid="ADHForwardReaction" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="EtOH" stoichiometry="1" constant="true"/>
<speciesReference species="NAD" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="Aldehyde" stoichiometry="1" constant="true"/>
<speciesReference species="NADH" stoichiometry="1" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci>k</ci>
<ci>ADH</ci>
<ci>EtOH</ci>
<ci>NAD</ci>
</apply>
</math>
</kineticLaw>
</reaction>
</listOfReactions>
</model>
</sbml>
"""
rdf_graph3_1_2_1 = RDF()
rdf_graph3_1_2_1.set_archive_uri("AHD.omex")
rdf_graph3_1_2_1.set_model_uri("AHDModel.sbml")
annot_editor = rdf_graph3_1_2_1.to_editor(sbml3_1_2_1, generate_new_metaids=False, sbml_semantic_extraction=False)
# OPB:OPB_00237 = Chemical flow rate
# GO:0004022 = GO term for Alcohol Dehydrogenase
with annot_editor.new_physical_process() as phy_process:
phy_process \
.about("ADHForwardReaction", eUriType.MODEL_URI) \
.is_version_of("GO:0004022") \
.add_source(physical_entity_reference="EtOH", uri_type=eUriType.MODEL_URI, multiplier=1.0) \
.add_source(physical_entity_reference="NAD", uri_type=eUriType.MODEL_URI, multiplier=1.0) \
.add_sink(physical_entity_reference="Aldehyde", uri_type=eUriType.MODEL_URI, multiplier=1.0) \
.add_sink(physical_entity_reference="NADH", uri_type=eUriType.MODEL_URI, multiplier=1.0) \
.add_mediator(physical_entity_reference="ADH", uri_type=eUriType.MODEL_URI) \
.has_property("OPB:OPB_00237") #Chemical flow rate
print(rdf_graph3_1_2_1)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix semsim: <http://bime.uw.edu/semsim/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/AHD.omex/AHDModel.rdf#> . local:MediatorParticipant0000 semsim:hasPhysicalEntityReference <http://omex-library.org/AHD.omex/AHDModel.sbml#ADH> . local:ProcessProperty0000 bqbiol:isPropertyOf <http://omex-library.org/AHD.omex/AHDModel.sbml#ADHForwardReaction> ; bqbiol:isVersionOf <https://identifiers.org/OPB:OPB_00237> . local:SinkParticipant0000 semsim:hasMultiplier "1"^^rdf:double ; semsim:hasPhysicalEntityReference <http://omex-library.org/AHD.omex/AHDModel.sbml#Aldehyde> . local:SinkParticipant0001 semsim:hasMultiplier "1"^^rdf:double ; semsim:hasPhysicalEntityReference <http://omex-library.org/AHD.omex/AHDModel.sbml#NADH> . local:SourceParticipant0000 semsim:hasMultiplier "1"^^rdf:double ; semsim:hasPhysicalEntityReference <http://omex-library.org/AHD.omex/AHDModel.sbml#EtOH> . local:SourceParticipant0001 semsim:hasMultiplier "1"^^rdf:double ; semsim:hasPhysicalEntityReference <http://omex-library.org/AHD.omex/AHDModel.sbml#NAD> . <http://omex-library.org/AHD.omex/AHDModel.sbml#ADHForwardReaction> semsim:hasMediatorParticipant local:MediatorParticipant0000 ; semsim:hasSinkParticipant local:SinkParticipant0000, local:SinkParticipant0001 ; semsim:hasSourceParticipant local:SourceParticipant0000, local:SourceParticipant0001 ; bqbiol:isVersionOf <https://identifiers.org/GO:0004022> .
In the above example, the reaction has two reactants, two products, and one mediating enzyme. Note that stoichiometry, where
appropriate (sources
and sinks
), can be specified via the multiplier
fields. As with the physical entity type composite annotations, the about method provides the link between the annotation document we are composing and the SBML source code. Similarly, the physical_entity_reference
arguments point to metaids on elements within the SBML document. The final optional clause has_property indicates that the reaction being simulated has the physical property of "chemical flow rate".
Finally, we provide an example for annotating the electrical potential across the cellular membrane, specifically resulting from calcium ions. In addition to source and sink, this is annotated with the OPB term for the Nernst, or reversal potential.
# Example 3_1_3
from pyomexmeta import RDF, eUriType
sbml3_1_3 = """<sbml xmlns="http://www.sbml.org/sbml/level3/version1/core" level="3" version="1">
<model metaid="NernstExample" id="NernstExample">
<listOfCompartments>
<compartment id="cytoplasm" metaid="cytoplasm" spatialDimensions="3" size="1" constant="true"/>
<compartment id="extracellular" metaid="extracellular" spatialDimensions="3" size="1" constant="true"/>
</listOfCompartments>
<listOfSpecies>
<species id="Ca_ex" metaid="Ca_ex" compartment="extracellular" initialConcentration="2" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="Ca_cyt" metaid="Ca_cyt" compartment="cytoplasm" initialConcentration="0.07" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
</listOfSpecies>
<listOfParameters>
<parameter id="NP" metaid="NernstPotential" value="137.04" constant="true"/>
</listOfParameters>
</model>
</sbml>"""
rdf_graph_3_1_3 = RDF()
rdf_graph_3_1_3.set_archive_uri("Example3_1_3.omex")
rdf_graph_3_1_3.set_model_uri("Example3_1_3.sbml")
annot_editor = rdf_graph_3_1_3.to_editor(sbml3_1_3, generate_new_metaids=False, sbml_semantic_extraction=False)
# Ca_cyt: Calcium Ions cytosol
# Ca_ex: Calcium Ions extracellular space
# NernstPotential: The metaID of the SBML reaction
# OPB/OPB_01581: Nernst reversal potential
with annot_editor.new_energy_diff() as energy_in:
energy_in \
.about("EnergyDiff000", eUriType.LOCAL_URI) \
.add_source(physical_entity_reference="Ca_ex", uri_type=eUriType.MODEL_URI) \
.add_sink(physical_entity_reference="Ca_cyt", uri_type=eUriType.MODEL_URI) \
.has_property(property_about="NernstPotential", about_uri_type=eUriType.MODEL_URI, is_version_of="OPB:OPB_01581")
print(rdf_graph_3_1_3)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix semsim: <http://bime.uw.edu/semsim/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/Example3_1_3.omex/Example3_1_3.rdf#> . local:EnergyDiff000 semsim:hasSinkParticipant local:SinkParticipant0000 ; semsim:hasSourceParticipant local:SourceParticipant0000 . local:SinkParticipant0000 semsim:hasPhysicalEntityReference <http://omex-library.org/Example3_1_3.omex/Example3_1_3.sbml#Ca_cyt> . local:SourceParticipant0000 semsim:hasPhysicalEntityReference <http://omex-library.org/Example3_1_3.omex/Example3_1_3.sbml#Ca_ex> . <http://omex-library.org/Example3_1_3.omex/Example3_1_3.sbml#NernstPotential> bqbiol:isPropertyOf local:EnergyDiff000 ; bqbiol:isVersionOf <https://identifiers.org/OPB:OPB_01581> .
As mentioned earlier, in CellML one has variables which represent various physical properties but the biological physical entity is implicit. In this example, we demonstrate how to annotate the Volume of blood in the Lumen of the left coronary artery using a snippet of CellML code.
# Example 3_2_1
from pyomexmeta import RDF, eUriType
cellml = """
<model xmlns="http://www.cellml.org/cellml/1.1#" xmlns:cmeta="http://www.cellml.org/metadata/1.0#"
name="annotation_examples" cmeta:id="annExamples">
<component name="main">
<variable cmeta:id="main.Volume" initial_value="100" name="Volume" units="dimensionless" />
<variable cmeta:id="main.MembraneVoltage" initial_value="-80" name="MembraneVoltage" units="dimensionless" />
<variable cmeta:id="main.ReactionRate" initial_value="1" name="ReactionRate" units="dimensionless" />
</component>
</model>
"""
rdf_graph_3_2_1 = RDF()
rdf_graph_3_2_1.set_archive_uri("physical_process.omex")
rdf_graph_3_2_1.set_model_uri("model.cellml")
annot_editor = rdf_graph_3_2_1.to_editor(cellml, generate_new_metaids=False, sbml_semantic_extraction=False)
# OPB:00154 = fluid volume
# FMA:9670 = portion of blood
# FMA:18228 = Lumen of the left coronary artery
with annot_editor.new_physical_entity() as arterial_volume:
arterial_volume.about("local-entity-0", eUriType.LOCAL_URI) \
.identity("FMA:9670") \
.is_part_of("FMA:18228") \
.has_property("main.Volume", eUriType.MODEL_URI, "opb:OPB_00154")
print(rdf_graph_3_2_1)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/physical_process.omex/model.rdf#> . <http://omex-library.org/physical_process.omex/model.cellml#main.Volume> bqbiol:isPropertyOf local:local-entity-0 ; bqbiol:isVersionOf <https://identifiers.org/opb:OPB_00154> . local:local-entity-0 bqbiol:is <https://identifiers.org/FMA:9670> ; bqbiol:isPartOf <https://identifiers.org/FMA:18228> .
In this example we construct a physical entity annotation for the arterial blood volume. The argument to the about method is the cmetaid for the CellML variable element called main.Volume
.
<variable cmeta:id="main.Volume" initial_value="100" name="Volume" units="dimensionless" />
Internally, libOmexMeta creates a new entity (local:local-entity-0
) which is local to the annotation document we are creating that links the variable with cmetaid main.Volume
to the information that annotates this physical entity.
Here we describe the annotation of the CellML version of the conversion of ethanol
into aldehyde
as described in 3.1.2. In CellML this reaction is represented by the molar flow rate from the source
ethanol
to the sink
acetaldehyde
via the mediator alcohol dehydrogenase 1A
within the cytosol of a hepatocyte
.
# Example 3_2_2
from pyomexmeta import RDF, eUriType
cellml = """
<model xmlns="http://www.cellml.org/cellml/1.1#" xmlns:cmeta="http://www.cellml.org/metadata/1.0#"
name="annotation_examples" cmeta:id="annExamples">
<component name="main">
<variable cmeta:id="main.Volume" initial_value="100" name="Volume" units="dimensionless" />
<variable cmeta:id="main.MembraneVoltage" initial_value="-80" name="MembraneVoltage" units="dimensionless" />
<variable cmeta:id="main.ReactionRate" initial_value="1" name="ReactionRate" units="dimensionless" />
</component>
</model>
"""
rdf_graph_3_2_2 = RDF()
rdf_graph_3_2_2.set_archive_uri("physical_process.omex")
rdf_graph_3_2_2.set_model_uri("model.cellml")
annot_editor = rdf_graph_3_2_2.to_editor(cellml, generate_new_metaids=False, sbml_semantic_extraction=False)
# fma:14515 = hepatocyte
# fma:66836 = portion of cytosol
with annot_editor.new_physical_entity() as cytosol:
cytosol \
.about("cytosol", eUriType.LOCAL_URI) \
.identity("FMA:66836") \
.is_part_of("FMA:14515")
# http://purl.obolibrary.org/obo/PR_000003767 = alcohol dehydrogenase 1A
with annot_editor.new_physical_entity() as mediator:
mediator \
.about("mediator", eUriType.LOCAL_URI) \
.identity("http://purl.obolibrary.org/obo/PR_000003767") \
.is_part_of("cytosol", eUriType.LOCAL_URI)
# CHEBI:16236 = ethanol
with annot_editor.new_physical_entity() as source:
source \
.about("source", eUriType.LOCAL_URI) \
.identity("CHEBI:16236") \
.is_part_of("cytosol", eUriType.LOCAL_URI)
# CHEBI:15343 = acetaldehyde
with annot_editor.new_physical_entity() as sink:
sink \
.about("sink", eUriType.LOCAL_URI) \
.identity("CHEBI:15343") \
.is_part_of("cytosol", eUriType.LOCAL_URI)
# opb:OPB_00592 = chemical molar flow rate
# GO:0004022 = alcohol dehydrogenase (NAD+) activity
with annot_editor.new_physical_process() as reaction_rate:
reaction_rate \
.about("process", eUriType.LOCAL_URI) \
.is_version_of("GO:0004022") \
.add_source("source", eUriType.LOCAL_URI, multiplier=1) \
.add_sink("sink", eUriType.LOCAL_URI, multiplier=1) \
.add_mediator("mediator", eUriType.LOCAL_URI) \
.has_property(property_about="main.ReactionRate", about_uri_type=eUriType.MODEL_URI, is_version_of="opb:OPB_00592")
print(rdf_graph_3_2_2)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix semsim: <http://bime.uw.edu/semsim/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/physical_process.omex/model.rdf#> . <http://omex-library.org/physical_process.omex/model.cellml#main.ReactionRate> bqbiol:isPropertyOf local:process ; bqbiol:isVersionOf <https://identifiers.org/opb:OPB_00592> . local:MediatorParticipant0000 semsim:hasPhysicalEntityReference local:mediator . local:SinkParticipant0000 semsim:hasMultiplier "1"^^rdf:double ; semsim:hasPhysicalEntityReference local:sink . local:SourceParticipant0000 semsim:hasMultiplier "1"^^rdf:double ; semsim:hasPhysicalEntityReference local:source . local:cytosol bqbiol:is <https://identifiers.org/FMA:66836> ; bqbiol:isPartOf <https://identifiers.org/FMA:14515> . local:mediator bqbiol:is <http://purl.obolibrary.org/obo/PR_000003767> ; bqbiol:isPartOf local:cytosol . local:process semsim:hasMediatorParticipant local:MediatorParticipant0000 ; semsim:hasSinkParticipant local:SinkParticipant0000 ; semsim:hasSourceParticipant local:SourceParticipant0000 ; bqbiol:isVersionOf <https://identifiers.org/GO:0004022> . local:sink bqbiol:is <https://identifiers.org/CHEBI:15343> ; bqbiol:isPartOf local:cytosol . local:source bqbiol:is <https://identifiers.org/CHEBI:16236> ; bqbiol:isPartOf local:cytosol .
We show here an example of annotating an energy differential, specifically the electrical potential across a membrane of a pancreatic islet.
# Example 3_2_3
from pyomexmeta import RDF, eUriType
cellml = """
<model xmlns="http://www.cellml.org/cellml/1.1#" xmlns:cmeta="http://www.cellml.org/metadata/1.0#"
name="annotation_examples" cmeta:id="annExamples">
<component name="main">
<variable cmeta:id="main.Volume" initial_value="100" name="Volume" units="dimensionless" />
<variable cmeta:id="main.MembraneVoltage" initial_value="-80" name="MembraneVoltage" units="dimensionless" />
<variable cmeta:id="main.ReactionRate" initial_value="1" name="ReactionRate" units="dimensionless" />
</component>
</model>
"""
rdf_graph_3_2_3 = RDF()
rdf_graph_3_2_3.set_archive_uri("physical_process.omex")
rdf_graph_3_2_3.set_model_uri("model.cellml")
annot_editor = rdf_graph_3_2_3.to_editor(cellml, generate_new_metaids=False, sbml_semantic_extraction=False)
# FMA:16016 = Pancreatic islet
# FMA:9672 = Intercellular matrix
with annot_editor.new_physical_entity() as matrix:
matrix \
.about("intercellular-matrix", eUriType.LOCAL_URI) \
.identity("FMA:9672") \
.is_part_of("FMA:16016")
# CHEBI:24870 - ion
with annot_editor.new_physical_entity() as sink:
sink \
.about("sink", eUriType.LOCAL_URI) \
.identity("CHEBI:24870") \
.is_part_of("intercellular-matrix", eUriType.LOCAL_URI)
# FMA:70586 = Type B cell of pancreatic islet
# FMA:66836 = Portion of cytosol
with annot_editor.new_physical_entity() as cytosol:
cytosol \
.about("cytosol", eUriType.LOCAL_URI) \
.identity("FMA:66836") \
.is_part_of("FMA:70586")
# CHEBI:24870 - ion
with annot_editor.new_physical_entity() as source:
source \
.about("source", eUriType.LOCAL_URI) \
.identity("CHEBI:24870") \
.is_part_of("cytosol", eUriType.LOCAL_URI)
# opb:OPB_00506 = Electrical potential
with annot_editor.new_energy_diff() as energy_diff:
energy_diff \
.about("energy-diff", eUriType.LOCAL_URI) \
.add_source("source", eUriType.LOCAL_URI) \
.add_sink("sink", eUriType.LOCAL_URI) \
.has_property("main.MembraneVoltage", eUriType.MODEL_URI, "opb:OPB_00506")
print(rdf_graph_3_2_3)
@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> . @prefix bqbiol: <http://biomodels.net/biology-qualifiers/> . @prefix semsim: <http://bime.uw.edu/semsim/> . @prefix OMEXlib: <http://omex-library.org/> . @prefix local: <http://omex-library.org/physical_process.omex/model.rdf#> . <http://omex-library.org/physical_process.omex/model.cellml#main.MembraneVoltage> bqbiol:isPropertyOf local:energy-diff ; bqbiol:isVersionOf <https://identifiers.org/opb:OPB_00506> . local:SinkParticipant0000 semsim:hasPhysicalEntityReference local:sink . local:SourceParticipant0000 semsim:hasPhysicalEntityReference local:source . local:cytosol bqbiol:is <https://identifiers.org/FMA:66836> ; bqbiol:isPartOf <https://identifiers.org/FMA:70586> . local:energy-diff semsim:hasSinkParticipant local:SinkParticipant0000 ; semsim:hasSourceParticipant local:SourceParticipant0000 . local:intercellular-matrix bqbiol:is <https://identifiers.org/FMA:9672> ; bqbiol:isPartOf <https://identifiers.org/FMA:16016> . local:sink bqbiol:is <https://identifiers.org/CHEBI:24870> ; bqbiol:isPartOf local:intercellular-matrix . local:source bqbiol:is <https://identifiers.org/CHEBI:24870> ; bqbiol:isPartOf local:cytosol .