libSBML Python API  5.18.0
createSimpleModel.py

An example of creating a simple SBML Level 3 model.

1 #!/usr/bin/env python
2 ##
3 ## @file createSampleModel.py
4 ## @brief Creates a simple SBML model and prints it to stdout.
5 ## @author Michael Hucka
6 ##
7 ## <!--------------------------------------------------------------------------
8 ## This sample program is distributed under a different license than the rest
9 ## of libSBML. This program uses the open-source MIT license, as follows:
10 ##
11 ## Copyright (c) 2013-2018 by the California Institute of Technology
12 ## (California, USA), the European Bioinformatics Institute (EMBL-EBI, UK)
13 ## and the University of Heidelberg (Germany), with support from the National
14 ## Institutes of Health (USA) under grant R01GM070923. All rights reserved.
15 ##
16 ## Permission is hereby granted, free of charge, to any person obtaining a
17 ## copy of this software and associated documentation files (the "Software"),
18 ## to deal in the Software without restriction, including without limitation
19 ## the rights to use, copy, modify, merge, publish, distribute, sublicense,
20 ## and/or sell copies of the Software, and to permit persons to whom the
21 ## Software is furnished to do so, subject to the following conditions:
22 ##
23 ## The above copyright notice and this permission notice shall be included in
24 ## all copies or substantial portions of the Software.
25 ##
26 ## THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
27 ## IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
28 ## FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
29 ## THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
30 ## LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
31 ## FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
32 ## DEALINGS IN THE SOFTWARE.
33 ##
34 ## Neither the name of the California Institute of Technology (Caltech), nor
35 ## of the European Bioinformatics Institute (EMBL-EBI), nor of the University
36 ## of Heidelberg, nor the names of any contributors, may be used to endorse
37 ## or promote products derived from this software without specific prior
38 ## written permission.
39 ## ------------------------------------------------------------------------ -->
40 ##
41 
42 import sys
43 from libsbml import *
44 
45 
46 def check(value, message):
47  """If 'value' is None, prints an error message constructed using
48  'message' and then exits with status code 1. If 'value' is an integer,
49  it assumes it is a libSBML return status code. If the code value is
50  LIBSBML_OPERATION_SUCCESS, returns without further action; if it is not,
51  prints an error message constructed using 'message' along with text from
52  libSBML explaining the meaning of the code, and exits with status code 1.
53  """
54  if value is None:
55  raise SystemExit('LibSBML returned a null value trying to ' + message + '.')
56  elif type(value) is int:
57  if value == LIBSBML_OPERATION_SUCCESS:
58  return
59  else:
60  err_msg = 'Error encountered trying to ' + message + '.' \
61  + 'LibSBML returned error code ' + str(value) + ': "' \
62  + OperationReturnValue_toString(value).strip() + '"'
63  raise SystemExit(err_msg)
64  else:
65  return
66 
67 
68 def create_model():
69  """Returns a simple but complete SBML Level 3 model for illustration."""
70 
71  # Create an empty SBMLDocument object. It's a good idea to check for
72  # possible errors. Even when the parameter values are hardwired like
73  # this, it is still possible for a failure to occur (e.g., if the
74  # operating system runs out of memory).
75 
76  try:
77  document = SBMLDocument(3, 1)
78  except ValueError:
79  raise SystemExit('Could not create SBMLDocument object')
80 
81  # Create the basic Model object inside the SBMLDocument object. To
82  # produce a model with complete units for the reaction rates, we need
83  # to set the 'timeUnits' and 'extentUnits' attributes on Model. We
84  # set 'substanceUnits' too, for good measure, though it's not strictly
85  # necessary here because we also set the units for individual species
86  # in their definitions.
87 
88  model = document.createModel()
89  check(model, 'create model')
90  check(model.setTimeUnits("second"), 'set model-wide time units')
91  check(model.setExtentUnits("mole"), 'set model units of extent')
92  check(model.setSubstanceUnits('mole'), 'set model substance units')
93 
94  # Create a unit definition we will need later. Note that SBML Unit
95  # objects must have all four attributes 'kind', 'exponent', 'scale'
96  # and 'multiplier' defined.
97 
98  per_second = model.createUnitDefinition()
99  check(per_second, 'create unit definition')
100  check(per_second.setId('per_second'), 'set unit definition id')
101  unit = per_second.createUnit()
102  check(unit, 'create unit on per_second')
103  check(unit.setKind(UNIT_KIND_SECOND), 'set unit kind')
104  check(unit.setExponent(-1), 'set unit exponent')
105  check(unit.setScale(0), 'set unit scale')
106  check(unit.setMultiplier(1), 'set unit multiplier')
107 
108  # Create a compartment inside this model, and set the required
109  # attributes for an SBML compartment in SBML Level 3.
110 
111  c1 = model.createCompartment()
112  check(c1, 'create compartment')
113  check(c1.setId('c1'), 'set compartment id')
114  check(c1.setConstant(True), 'set compartment "constant"')
115  check(c1.setSize(1), 'set compartment "size"')
116  check(c1.setSpatialDimensions(3), 'set compartment dimensions')
117  check(c1.setUnits('litre'), 'set compartment size units')
118 
119  # Create two species inside this model, set the required attributes
120  # for each species in SBML Level 3 (which are the 'id', 'compartment',
121  # 'constant', 'hasOnlySubstanceUnits', and 'boundaryCondition'
122  # attributes), and initialize the amount of the species along with the
123  # units of the amount.
124 
125  s1 = model.createSpecies()
126  check(s1, 'create species s1')
127  check(s1.setId('s1'), 'set species s1 id')
128  check(s1.setCompartment('c1'), 'set species s1 compartment')
129  check(s1.setConstant(False), 'set "constant" attribute on s1')
130  check(s1.setInitialAmount(5), 'set initial amount for s1')
131  check(s1.setSubstanceUnits('mole'), 'set substance units for s1')
132  check(s1.setBoundaryCondition(False), 'set "boundaryCondition" on s1')
133  check(s1.setHasOnlySubstanceUnits(False), 'set "hasOnlySubstanceUnits" on s1')
134 
135  s2 = model.createSpecies()
136  check(s2, 'create species s2')
137  check(s2.setId('s2'), 'set species s2 id')
138  check(s2.setCompartment('c1'), 'set species s2 compartment')
139  check(s2.setConstant(False), 'set "constant" attribute on s2')
140  check(s2.setInitialAmount(0), 'set initial amount for s2')
141  check(s2.setSubstanceUnits('mole'), 'set substance units for s2')
142  check(s2.setBoundaryCondition(False), 'set "boundaryCondition" on s2')
143  check(s2.setHasOnlySubstanceUnits(False), 'set "hasOnlySubstanceUnits" on s2')
144 
145  # Create a parameter object inside this model, set the required
146  # attributes 'id' and 'constant' for a parameter in SBML Level 3, and
147  # initialize the parameter with a value along with its units.
148 
149  k = model.createParameter()
150  check(k, 'create parameter k')
151  check(k.setId('k'), 'set parameter k id')
152  check(k.setConstant(True), 'set parameter k "constant"')
153  check(k.setValue(1), 'set parameter k value')
154  check(k.setUnits('per_second'), 'set parameter k units')
155 
156  # Create a reaction inside this model, set the reactants and products,
157  # and set the reaction rate expression (the SBML "kinetic law"). We
158  # set the minimum required attributes for all of these objects. The
159  # units of the reaction rate are determined from the 'timeUnits' and
160  # 'extentUnits' attributes on the Model object.
161 
162  r1 = model.createReaction()
163  check(r1, 'create reaction')
164  check(r1.setId('r1'), 'set reaction id')
165  check(r1.setReversible(False), 'set reaction reversibility flag')
166  check(r1.setFast(False), 'set reaction "fast" attribute')
167 
168  species_ref1 = r1.createReactant()
169  check(species_ref1, 'create reactant')
170  check(species_ref1.setSpecies('s1'), 'assign reactant species')
171  check(species_ref1.setConstant(True), 'set "constant" on species ref 1')
172 
173  species_ref2 = r1.createProduct()
174  check(species_ref2, 'create product')
175  check(species_ref2.setSpecies('s2'), 'assign product species')
176  check(species_ref2.setConstant(True), 'set "constant" on species ref 2')
177 
178  math_ast = parseL3Formula('k * s1 * c1')
179  check(math_ast, 'create AST for rate expression')
180 
181  kinetic_law = r1.createKineticLaw()
182  check(kinetic_law, 'create kinetic law')
183  check(kinetic_law.setMath(math_ast), 'set math on kinetic law')
184 
185  # And we're done creating the basic model.
186  # Now return a text string containing the model in XML format.
187 
188  return writeSBMLToString(document)
189 
190 
191 if __name__ == '__main__':
192  print(create_model())