Lengthy example of creating SBML models presented in the SBML specification.
54 ProgramName =
"createExampleModels" 55 ProgramVersion =
"1.0.0" 80 def createExampleEnzymaticReaction():
90 sbmlDoc = SBMLDocument(level, version)
98 model = sbmlDoc.createModel()
99 model.setId(
"EnzymaticReaction")
111 unitdef = model.createUnitDefinition()
112 unitdef.setId(
"per_second")
116 unit = unitdef.createUnit()
117 unit.setKind(UNIT_KIND_SECOND)
128 unitdef = model.createUnitDefinition()
129 unitdef.setId(
"litre_per_mole_per_second")
133 unit = unitdef.createUnit()
134 unit.setKind(UNIT_KIND_MOLE)
139 unit = unitdef.createUnit()
140 unit.setKind(UNIT_KIND_LITRE)
145 unit = unitdef.createUnit()
146 unit.setKind(UNIT_KIND_SECOND)
159 comp = model.createCompartment()
182 sp = model.createSpecies()
189 sp.setCompartment(compName)
204 sp.setInitialAmount(0)
210 sp = model.createSpecies()
211 sp.setCompartment(compName)
214 sp.setInitialAmount(0)
220 sp = model.createSpecies()
221 sp.setCompartment(compName)
224 sp.setInitialAmount(1e-20)
230 sp = model.createSpecies()
231 sp.setCompartment(compName)
234 sp.setInitialAmount(5e-21)
246 reaction = model.createReaction()
247 reaction.setId(
"veq")
253 spr = reaction.createReactant()
259 spr = reaction.createReactant()
267 spr = reaction.createProduct()
274 kl = reaction.createKineticLaw()
309 astCytosol = ASTNode(AST_NAME)
310 astCytosol.setName(
"cytosol")
312 astKon = ASTNode(AST_NAME)
313 astKon.setName(
"kon")
315 astKoff = ASTNode(AST_NAME)
316 astKoff.setName(
"koff")
318 astE = ASTNode(AST_NAME)
321 astS = ASTNode(AST_NAME)
324 astES = ASTNode(AST_NAME)
338 astTimes1 = ASTNode(AST_TIMES)
339 astTimes1.addChild(astKoff)
340 astTimes1.addChild(astES)
374 astTimes2 = ASTNode(AST_TIMES)
375 astTimes2.addChild(astE)
376 astTimes2.addChild(astS)
378 astTimes = ASTNode(AST_TIMES)
379 astTimes.addChild(astKon)
380 astTimes.addChild(astTimes2)
402 astMinus = ASTNode(AST_MINUS)
403 astMinus.addChild(astTimes)
404 astMinus.addChild(astTimes1)
430 astMath = ASTNode(AST_TIMES)
431 astMath.addChild(astCytosol)
432 astMath.addChild(astMinus)
448 para = kl.createParameter()
450 para.setValue(1000000)
451 para.setUnits(
"litre_per_mole_per_second")
455 para = kl.createParameter()
458 para.setUnits(
"per_second")
464 reaction = model.createReaction()
465 reaction.setId(
"vcat")
466 reaction.setReversible(
False)
475 spr = reaction.createReactant()
484 spr = reaction.createProduct()
489 spr = reaction.createProduct()
496 kl = reaction.createKineticLaw()
510 mathXMLString =
"""<math xmlns="http://www.w3.org/1998/Math/MathML"> 529 para = kl.createParameter()
532 para.setUnits(
"per_second")
548 def createExampleInvolvingUnits():
558 sbmlDoc = SBMLDocument(level, version)
564 sbmlDoc.getNamespaces().add(
"http://www.w3.org/1999/xhtml",
"xhtml")
572 model = sbmlDoc.createModel()
573 model.setId(
"unitsExample")
588 unitdef = model.createUnitDefinition()
589 unitdef.setId(
"substance")
593 unit = unitdef.createUnit()
594 unit.setKind(UNIT_KIND_MOLE)
605 unitdef = model.createUnitDefinition()
606 unitdef.setId(
"mmls")
610 unit = unitdef.createUnit()
611 unit.setKind(UNIT_KIND_MOLE)
616 unit = unitdef.createUnit()
617 unit.setKind(UNIT_KIND_LITRE)
622 unit = unitdef.createUnit()
623 unit.setKind(UNIT_KIND_SECOND)
630 unitdef = model.createUnitDefinition()
635 unit = unitdef.createUnit()
636 unit.setKind(UNIT_KIND_MOLE)
641 unit = unitdef.createUnit()
642 unit.setKind(UNIT_KIND_LITRE)
655 comp = model.createCompartment()
677 sp = model.createSpecies()
683 sp.setCompartment(compName)
695 sp.setInitialConcentration(1)
701 sp = model.createSpecies()
703 sp.setCompartment(compName)
704 sp.setInitialConcentration(1)
710 sp = model.createSpecies()
711 sp.setCompartment(compName)
713 sp.setInitialConcentration(1)
719 sp = model.createSpecies()
720 sp.setCompartment(compName)
722 sp.setInitialConcentration(1)
732 para = model.createParameter()
735 para.setUnits(
"mmls")
739 para = model.createParameter()
754 reaction = model.createReaction()
764 spr = reaction.createReactant()
773 spr = reaction.createProduct()
780 kl = reaction.createKineticLaw()
788 notesString =
"<xhtml:p> ((vm * s1)/(km + s1)) * cell </xhtml:p>" 789 kl.setNotes(notesString)
823 astMath = ASTNode(AST_TIMES)
825 astMath.addChild(ASTNode(AST_DIVIDE))
826 astDivide = astMath.getLeftChild()
828 astDivide.addChild(ASTNode(AST_TIMES))
829 astTimes = astDivide.getLeftChild()
831 astTimes.addChild(ASTNode(AST_NAME))
832 astTimes.getLeftChild().setName(
"vm")
834 astTimes.addChild(ASTNode(AST_NAME))
835 astTimes.getRightChild().setName(
"s1")
837 astDivide.addChild(ASTNode(AST_PLUS))
838 astPlus = astDivide.getRightChild()
840 astPlus.addChild(ASTNode(AST_NAME))
841 astPlus.getLeftChild().setName(
"km")
843 astPlus.addChild(ASTNode(AST_NAME))
844 astPlus.getRightChild().setName(
"s1")
846 astMath.addChild(ASTNode(AST_NAME))
847 astMath.getRightChild().setName(
"cell")
861 reaction = model.createReaction()
871 spr = reaction.createReactant()
880 spr = reaction.createProduct()
887 kl = reaction.createKineticLaw()
899 notesXMLNode = XMLNode(
900 XMLTriple(
"p",
"",
"xhtml"),
905 notesXMLNode.addChild(XMLNode(
" ((vm * s2)/(km + s2)) * cell "))
909 kl.setNotes(notesXMLNode)
923 mathXMLString =
"""<math xmlns="http://www.w3.org/1998/Math/MathML"> 951 reaction = model.createReaction()
961 spr = reaction.createReactant()
970 spr = reaction.createProduct()
977 kl = reaction.createKineticLaw()
981 notesString =
"<xhtml:p> ((vm * x1)/(km + x1)) * cell </xhtml:p>" 982 kl.setNotes(notesString)
988 mathXMLString =
"""<math xmlns="http://www.w3.org/1998/Math/MathML"> 1026 def createExampleInvolvingFunctionDefinitions():
1036 sbmlDoc = SBMLDocument(level, version)
1044 model = sbmlDoc.createModel()
1045 model.setId(
"functionExample")
1053 fdef = model.createFunctionDefinition()
1058 mathXMLString =
"""<math xmlns="http://www.w3.org/1998/Math/MathML"> 1073 fdef.setMath(astMath)
1081 compName =
"compartmentOne" 1085 comp = model.createCompartment()
1086 comp.setId(compName)
1106 sp = model.createSpecies()
1112 sp.setCompartment(compName)
1125 sp.setInitialConcentration(1)
1131 sp = model.createSpecies()
1133 sp.setCompartment(compName)
1134 sp.setInitialConcentration(0)
1144 para = model.createParameter()
1147 para.setUnits(
"second")
1159 reaction = model.createReaction()
1160 reaction.setId(
"reaction_1")
1161 reaction.setReversible(
False)
1170 spr = reaction.createReactant()
1171 spr.setSpecies(
"S1")
1179 spr = reaction.createProduct()
1180 spr.setSpecies(
"S2")
1186 kl = reaction.createKineticLaw()
1192 mathXMLString =
"""<math xmlns="http://www.w3.org/1998/Math/MathML"> 1201 <ci> compartmentOne </ci> 1235 def validateExampleSBML(sbmlDoc):
1237 print(
"validateExampleSBML: given a None SBML Document")
1240 consistencyMessages =
"" 1241 validationMessages =
"" 1243 numCheckFailures = 0
1244 numConsistencyErrors = 0
1245 numConsistencyWarnings = 0
1246 numValidationErrors = 0
1247 numValidationWarnings = 0
1254 numCheckFailures = sbmlDoc.checkInternalConsistency()
1255 if numCheckFailures > 0:
1257 for i
in range(0, numCheckFailures):
1258 sbmlErr = sbmlDoc.getError(i)
1259 if sbmlErr.isFatal()
or sbmlErr.isError():
1260 numConsistencyErrors += 1
1262 numConsistencyWarnings += 1
1263 sbmlDoc.printErrors()
1269 if numConsistencyErrors > 0:
1270 consistencyMessages +=
"Further validation aborted." 1272 numCheckFailures = sbmlDoc.checkConsistency()
1273 if numCheckFailures > 0:
1275 for i
in range(0, numCheckFailures):
1276 sbmlErr = sbmlDoc.getError(i)
1277 if sbmlErr.isFatal()
or sbmlErr.isError():
1278 numValidationErrors += 1
1280 numValidationWarnings += 1
1281 sbmlDoc.printErrors()
1286 if numConsistencyErrors > 0:
1288 if numConsistencyErrors > 1:
1290 print(
"ERROR: encountered " + str(
1291 numConsistencyErrors) +
" consistency error" + tmp +
" in model '" + sbmlDoc.getModel().getId() +
"'.")
1293 if numConsistencyWarnings > 0:
1295 if numConsistencyWarnings > 1:
1297 print(
"Notice: encountered " + str(numConsistencyWarnings)
1298 +
" consistency warning" + tmp
1299 +
" in model '" + sbmlDoc.getModel().getId() +
"'.")
1300 print(consistencyMessages)
1302 if numValidationErrors > 0:
1304 if numValidationErrors > 1:
1306 print(
"ERROR: encountered " + str(numValidationErrors) +
" validation error" + tmp
1307 +
" in model '" + sbmlDoc.getModel().getId() +
"'.")
1308 if numValidationWarnings > 0:
1310 if numValidationWarnings > 1:
1312 print(
"Notice: encountered " + str(numValidationWarnings)
1313 +
" validation warning" + tmp
1314 +
" in model '" + sbmlDoc.getModel().getId() +
"'.")
1315 print(validationMessages)
1317 return numConsistencyErrors == 0
and numValidationErrors == 0
1325 def writeExampleSBML(sbmlDoc, filename):
1328 print(
"Wrote file '" + filename +
"'")
1331 print(
"Failed to write '" + filename +
"'")
1356 sbmlDoc = createExampleEnzymaticReaction()
1357 SBMLok = validateExampleSBML(sbmlDoc)
1359 writeExampleSBML(sbmlDoc,
"enzymaticreaction.xml")
1367 sbmlDoc = createExampleInvolvingUnits()
1368 SBMLok = validateExampleSBML(sbmlDoc)
1370 writeExampleSBML(sbmlDoc,
"units.xml")
1378 sbmlDoc = createExampleInvolvingFunctionDefinitions()
1379 SBMLok = validateExampleSBML(sbmlDoc)
1381 writeExampleSBML(sbmlDoc,
"functiondef.xml")
1385 print(
"Unexpected exceptional condition encountered.")
1392 if __name__ ==
'__main__':