Create & Annotate An SBML Model

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
# include <semsim/SemSim.h>

# include "sbml/SBMLTypes.h"

# include <iostream>

using namespace semsim;

int main() {
  // create SBML model
  LIBSBML_CPP_NAMESPACE_QUALIFIER SBMLDocument* d =
    new LIBSBML_CPP_NAMESPACE_QUALIFIER SBMLDocument(3, 2);
  LIBSBML_CPP_NAMESPACE_QUALIFIER Model* m = d->createModel();
  m->setId("beta_cell_model");

  // create a compartment to represent the cytosol
  LIBSBML_CPP_NAMESPACE_QUALIFIER Compartment* comp = m->createCompartment();
  comp->setId("cytosol");
  comp->setMetaId("cytosol");
  comp->setSize(1);
  comp->setConstant(true);

  // add a CV term for the cytosol
  LIBSBML_CPP_NAMESPACE_QUALIFIER CVTerm *cv_cytosol = new LIBSBML_CPP_NAMESPACE_QUALIFIER CVTerm();
  cv_cytosol->setQualifierType(LIBSBML_CPP_NAMESPACE_QUALIFIER BIOLOGICAL_QUALIFIER);
  cv_cytosol->setBiologicalQualifierType(LIBSBML_CPP_NAMESPACE_QUALIFIER BQB_IS);
  cv_cytosol->addResource(GO::get(5829).getURI().encode());
  comp->addCVTerm(cv_cytosol);

  LIBSBML_CPP_NAMESPACE_QUALIFIER UnitDefinition* unitdef = m->createUnitDefinition();
  unitdef->setId("molar");

  LIBSBML_CPP_NAMESPACE_QUALIFIER Unit* unit = unitdef->createUnit();
  unit->setKind(LIBSBML_CPP_NAMESPACE_QUALIFIER UNIT_KIND_MOLE);
  unit->setExponent(1);
  unit->setMultiplier(1);
  unit->setScale(1);

  unit = unitdef->createUnit();
  unit->setKind(LIBSBML_CPP_NAMESPACE_QUALIFIER UNIT_KIND_LITRE);
  unit->setExponent(-1);
  unit->setMultiplier(1);
  unit->setScale(1);

  // create the species for glucose
  LIBSBML_CPP_NAMESPACE_QUALIFIER Species* s = m->createSpecies();
  s->setCompartment("cytosol");
  s->setId("glucose");
  s->setMetaId("glucose");
  s->setInitialConcentration(0);
  s->setUnits("molar");
  s->setHasOnlySubstanceUnits(false);
  s->setConstant(false);
  s->setBoundaryCondition(false);

  // add a CV term to glucose
  LIBSBML_CPP_NAMESPACE_QUALIFIER CVTerm *cv_glucose = new LIBSBML_CPP_NAMESPACE_QUALIFIER CVTerm();
  cv_glucose->setQualifierType(LIBSBML_CPP_NAMESPACE_QUALIFIER BIOLOGICAL_QUALIFIER);
  cv_glucose->setBiologicalQualifierType(LIBSBML_CPP_NAMESPACE_QUALIFIER BQB_IS);
  cv_glucose->addResource(CHEBI::get(17234).getURI().encode());
  s->addCVTerm(cv_glucose);

  // create import reaction for glucose
  LIBSBML_CPP_NAMESPACE_QUALIFIER Reaction* reaction = m->createReaction();
  reaction->setId("glucose_import");
  reaction->setReversible(false);

  LIBSBML_CPP_NAMESPACE_QUALIFIER SpeciesReference* sr = reaction->createProduct();
  sr->setSpecies("glucose");
  sr->setConstant(false);
  LIBSBML_CPP_NAMESPACE_QUALIFIER KineticLaw* k = reaction->createKineticLaw();
  LIBSBML_CPP_NAMESPACE_QUALIFIER Parameter* p = k->createParameter();
  p->setId("glucose_import_rate");
  p->setValue(1);
  k->setMath(LIBSBML_CPP_NAMESPACE_QUALIFIER SBML_parseL3FormulaWithModel("glucose_import_rate", m));
  // alternate way to set formula
  // LIBSBML_CPP_NAMESPACE_QUALIFIER ASTNode* a =
  //   new LIBSBML_CPP_NAMESPACE_QUALIFIER ASTNode(LIBSBML_CPP_NAMESPACE_QUALIFIER AST_NAME);
  // a->setName("glucose_import_rate");
  // k->setMath(a);

  LIBSBML_CPP_NAMESPACE_QUALIFIER SBMLWriter sbml_writer;

  std::cerr << "**********************************************************\n";
  std::cerr << "SBML output:\n";
  std::cerr << "**********************************************************\n\n";

  std::cerr << sbml_writer.writeSBMLToString(d) << "\n";

  // importing the model into libSemSim will automatically
  // add metaids for any annotatable SBML elements that lack them
  SBMLImporter importer(d);
  SBMLModel& model = importer.getSBMLModel();

  std::cerr << "**********************************************************\n";

  std::cerr << "RDF serialization of annotations:\n";
  std::cerr << "**********************************************************\n\n";

  std::cerr << model.getRDF("./mymodel.xml", "turtle") << "\n";

  // this isn't bad - the RDF already contains a composite annotation
  // telling is that the species "glucose" is inside the "cytosol"
  // we currently have "cytosolic glucose"
  // however - what we want is "cytosolic glucose in a pancreatic beta cell"
  // to describe this, we need to annotate the compartment (not the species!)
  // and describe the cell type it is part of

  // the compartment currently has a *singular annotation*, which can't
  // represent membership in a "pancreatic beta cell", so we change this to
  // a composite annotation
  Component* c = model.getComponentForId("cytosol");
  c->setAnnotation(c->getAnnotation().makeComposite(OPB::get(523))); // 523 is spatial volume
  // now add a descriptor term specifying that the cytoplasm
  // is part of a pancreatic beta cell (CL:0000169)
  c->getCompositeAnnotation().addTerm(
    bqb::isPartOf, // the relation (the cytoplasm *is part of* ...)
    CL::get(169)   // the resource (pancreatic beta cell)
  );

  std::cerr << "**********************************************************\n";
  std::cerr << "RDF serialization with cell type annotation:\n";
  std::cerr << "**********************************************************\n\n";

  std::cerr << model.getRDF("./mymodel.xml", "turtle") << "\n";

  std::cerr << "**********************************************************\n";
  std::cerr << "Humanized version of the annotation information\n";
  std::cerr << "(ontology terms replaced by human-readable names):\n";
  std::cerr << "**********************************************************\n\n";

  std::cerr << model.humanize() << "\n";

  return 0;
}

Expected Output

**********************************************************
SBML output:
**********************************************************

<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version2/core" level="3" version="2">
  <model id="beta_cell_model">
    <listOfUnitDefinitions>
      <unitDefinition id="molar">
        <listOfUnits>
          <unit kind="mole" exponent="1" scale="1" multiplier="1"/>
          <unit kind="litre" exponent="-1" scale="1" multiplier="1"/>
        </listOfUnits>
      </unitDefinition>
    </listOfUnitDefinitions>
    <listOfCompartments>
      <compartment metaid="cytosol" id="cytosol" size="1" constant="true">
        <annotation>
          <rdf:RDF xmlns:rdf="http://www.w3.org/1999/02/22-rdf-syntax-ns#" xmlns:dcterms="http://purl.org/dc/terms/" xmlns:vCard="http://www.w3.org/2001/vcard-rdf/3.0#" xmlns:vCard4="http://www.w3.org/2006/vcard/ns#" xmlns:bqbiol="http://biomodels.net/biology-qualifiers/" xmlns:bqmodel="http://biomodels.net/model-qualifiers/">
            <rdf:Description rdf:about="#cytosol">
              <bqbiol:is>
                <rdf:Bag>
                  <rdf:li rdf:resource="https://identifiers.org/GO:0005829"/>
                </rdf:Bag>
              </bqbiol:is>
            </rdf:Description>
          </rdf:RDF>
        </annotation>
      </compartment>
    </listOfCompartments>
    <listOfSpecies>
      <species metaid="glucose" id="glucose" compartment="cytosol" initialConcentration="0" substanceUnits="molar" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false">
        <annotation>
          <rdf:RDF xmlns:rdf="http://www.w3.org/1999/02/22-rdf-syntax-ns#" xmlns:dcterms="http://purl.org/dc/terms/" xmlns:vCard="http://www.w3.org/2001/vcard-rdf/3.0#" xmlns:vCard4="http://www.w3.org/2006/vcard/ns#" xmlns:bqbiol="http://biomodels.net/biology-qualifiers/" xmlns:bqmodel="http://biomodels.net/model-qualifiers/">
            <rdf:Description rdf:about="#glucose">
              <bqbiol:is>
                <rdf:Bag>
                  <rdf:li rdf:resource="http://identifiers.org/CHEBI:17234"/>
                </rdf:Bag>
              </bqbiol:is>
            </rdf:Description>
          </rdf:RDF>
        </annotation>
      </species>
    </listOfSpecies>
    <listOfReactions>
      <reaction id="glucose_import" reversible="false">
        <listOfProducts>
          <speciesReference species="glucose" constant="false"/>
        </listOfProducts>
        <kineticLaw>
          <math xmlns="http://www.w3.org/1998/Math/MathML">
            <ci> glucose_import_rate </ci>
          </math>
          <listOfLocalParameters>
            <localParameter id="glucose_import_rate" value="1"/>
          </listOfLocalParameters>
        </kineticLaw>
      </reaction>
    </listOfReactions>
  </model>
</sbml>

**********************************************************
RDF serialization of annotations:
**********************************************************

@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> .
@prefix bqb: <http://biomodels.net/biology-qualifiers/> .

<#glucose_entity>
    bqb:is <http://identifiers.org/CHEBI:17234> ;
    bqb:occursIn <./mymodel#cytosol> .

<./mymodel#cytosol>
    bqb:is <https://identifiers.org/GO:0005829> .

<./mymodel#glucose>
    bqb:is <http://identifiers.org/opb/OPB_00425> ;
    bqb:isPropertyOf <#glucose_entity> .


**********************************************************
RDF serialization with cell type annotation:
**********************************************************

@prefix rdf: <http://www.w3.org/1999/02/22-rdf-syntax-ns#> .
@prefix bqb: <http://biomodels.net/biology-qualifiers/> .

<#cytosol_entity>
    bqb:is <https://identifiers.org/GO:0005829> ;
    bqb:isPartOf <#cytosol_entity_term1> .

<#cytosol_entity_term1>
    bqb:is <https://identifiers.org/CL:0000169> .

<#glucose_entity>
    bqb:is <http://identifiers.org/CHEBI:17234> ;
    bqb:occursIn <./mymodel#cytosol> .

<./mymodel#cytosol>
    bqb:is <http://identifiers.org/opb/OPB_00523> ;
    bqb:isPropertyOf <#cytosol_entity> .

<./mymodel#glucose>
    bqb:is <http://identifiers.org/opb/OPB_00425> ;
    bqb:isPropertyOf <#glucose_entity> .

**********************************************************
Humanized version of the annotation information
(ontology terms replaced by human-readable names):
**********************************************************

OPB:spatial volume -> (isPropertyOf) -> #cytosol -> (is) -> GO:cytoplasm -> (isPartOf) -> CL:type B pancreatic cell
OPB:chemical concentration -> (isPropertyOf) -> #glucose -> (is) -> CHEBI:glucose -> (occursIn) -> #cytosol