cometspy package

Submodules

cometspy.comets module

The comets module runs COMETS simulations and stores output.

Generally, a comets object is created just before calling run(). Afterwards, any saved data (e.g. total_biomass) can be accessed from the object.

class cometspy.comets.comets(layout, parameters, relative_dir: str = '')

Bases: object

the main simulation object to run COMETS

The comets class stores a previously-made layout and params, and interacts with COMETS to run the simulation. It also stores any generated output, which could include total_biomass, biomass, media, or fluxes, as set in the params object. std_out from the COMETS simulation is saved in the attribute run_output and can be useful to examine to debug errors.

When creating a comets object, the optional relative_dir path is useful when one is to run multiple simulations simultaneously, otherwise temporary files may overwrite each other.

Parameters

layoutlayout

a cometspy.layout containing models and media information

parametersparams

a cometspy.params containing specified simulation parameters

relative_dirstr, optional

a directory to place temporary simulation files.

Attributes

layoutcometspy.layout

the layout containing cometspy.model[s] and media information

paramscometspy.params

the params containing simulation and default biological parameters

working_dirstr

the directory at which to save temporary sim files

GUROBI_HOMEstr

the directory where GUROBI exists on the system

COMETS_HOMEstr

the directory where COMETS exists on the system

VERSIONstr

the version of comets, read from the files at COMETS_HOME

classpath_piecesdict

classpath separated into library name (key) and location (value)

JAVA_CLASSPATHstr

a generated (overwritable) string containing the java classpath

run_outputstr

generated object containing text from COMETS sim’s std_out

run_errorsstr

generated object containing text from COMETS sim’s std_err

total_biomasspandas.DataFrame

generated object containing total biomass from simulation

biomasspandas.DataFrame

generated object containing spatially-explicit biomass from sim

specific_mediapandas.DataFrame

generated object containing information on selected media from sim

mediapandas.Dataframe

generated object containing spatially-explicit media from sim

fluxes_by_speciesdict{model_idpandas.DataFrame}

generated object containing each species’ spatial-explicit fluxes

genotypespandas.DataFrame

generated object containing genotypes if an evolution sim was run

Examples

>>> # assume layout and params have already been created
>>> # and cometspy was imported as c
>>> sim = c.comets(layout, params)
>>> sim.run() # this could take from seconds to hours
>>> print(sim.run_output) # to see the std_out
>>> sim.total_biomass.plot(x = "cycle")
>>> # assume params.all_params["writeBiomassLog"] == True, and that
>>> # params.all_params["BiomassLogRate"] = 1000, and that
>>> # params.all_params["maxCycles"] = 5000
>>> im = sim.get_biomass_image(layout.models[0].id, 4000)
>>> from matplotlib import pyplot as plt
>>> plt.imshow(im / np.max(im))
get_biomass_image(model_id: str, cycle: int) array

returns an image of biomass concentrations at a given cycle

This will only work if biomass was saved at the current cycle. This requires the following parameters to have been set: >>> params.set_param(“writeBiomassLog”,True) >>> params.set_param(“BiomassLogRate”, n) # n is in cycles

Parameters

model_idstr

the id of the model to get biomass data on

cycleint

the cycle to get the biomass data

Returns

A 2d numpy array which can be visualized like an image.

Examples

>>> sim.run()
>>> # assume sim.params.all_params["BiomassLogRate"] = 100 and that
>>> # sim.params.all_params["maxCycles"] = 1000
>>> im = sim.get_biomass_image("iJO1366", 500) # e.g. the ecoli model
>>> from matplotlib import pyplot as plt # may need to be installed
>>> plt.imshow(im)
get_flux_image(model_id: str, reaction_id: str, cycle: int) array

returns a 2d numpy array showing fluxes at a given cycle

This will only work if flux was saved at the current cycle. This requires the following parameters to have been set:

params.set_param(“writeFluxLog”,True) params.set_param(“FluxLogRate”, n) # n is in cycles

Parameters

model_idstr

the id of the model about which to get fluxes

reaction_idstr

the id of the reaction about which to get fluxes

cycleint

the cycle at which to get fluxes

Returns

a 2d numpy array which can be visualized like an image

Examples

>>> sim.run()
>>> # assume sim.params.all_params["FluxLogRate"] = 100 and that
>>> # sim.params.all_params["maxCycles"] = 1000
>>> # assume a model used was iJO1366
>>> im = sim.get_flux_image("iJO1366", "EX_ac_e", 500)
>>> from matplotlib import pyplot as plt # may need to be installed
>>> plt.imshow(im)
get_metabolite_image(met: str, cycle: int) array

returns an image of metabolite concentrations at a given cycle

This will only work if media was saved at the current cycle. This requires the following parameters to have been set:

params.set_param(“writeMediaLog”,True) params.set_param(“MediaLogRate”, n) # n is in cycles

Notes

There may be a bug where cycles are + 1 what they should be. We will fix this soon but for now be aware you may need to +1 your desired cycle. The bug does not affect anything within the simulation.

Parameters

metstr

the name of the metabolite

cycleint

the cycle to get the metabolite data

Returns

A 2d numpy array which can be visualized like an image.

Examples

>>> sim.run()
>>> # assume sim.params.all_params["MediaLogRate"] = 100 and that
>>> # sim.params.all_params["maxCycles"] = 1000
>>> im = sim.get_metabolite_image("ac_e", 500)
>>> from matplotlib import pyplot as plt # may need to be installed
>>> plt.imshow(im)
get_metabolite_time_series(upper_threshold: float = 1000.0) DataFrame

returns a pandas DataFrame containing extracellular metabolite time series

Parameters

upper_thresholdfloat (optional)

metabolites ever above this are not returned

get_species_exchange_fluxes(model_id: str, threshold: float = 0.0) DataFrame

returns a pandas DataFrame containing one model’s exchange flux time series

Parameters

model_idstr

id of the model

thresholdfloat (optional)

abs(flux) must exceed this to be returned

run(delete_files: bool = True)

run a COMETS simulation

This runs the COMETS simulation specified by the layout and params objects supplied to this comets object. It creates the files needed for COMETS into the current directory or the optional one specified when this object was created.

Once complete (or if an error has occurred), this tries to read simulation logs including data, as well as the std_out in the run_output attribute.

If the optional delete_files is set to False, then temporary files and data log files are not deleted. They are deleted by default.

Parameters

delete_filesbool, optional

Whether to delete simulation and log files. The default is True.

Examples

>>> sim = c.comets(layout, params)
>>> sim.run(delete_files = True)
>>> print(sim.run_output)
>>> print(sim.total_biomass)
set_classpath(libraryname: str, path: str)

sets the full path to a required specified java library

This can be used to set non-default classpaths. Note that currently, it does not work for windows, because windows runs COMETS slightly differently than Unix.

Parameters

librarynamestr

the library for which the new path is being supplied.

pathstr

the full path, including file name, to the library

Examples

>>> sim = c.comets(layout, params)
>>> sim.set_classpath("jmatio", "/opt/jmatio/jmatio.jar")

cometspy.layout module

class cometspy.layout.layout(input_obj: list = None)

Bases: object

object containing model and environmental definitions for a COMETS sim.

The layout object is an essential object needed to start a comets simulation, along with a params object. For a simulation to run, the layout must contain model(s). Typically, the layout will be created after creating models, and these models will be given in a list during layout creation. Subsequently, methods like set_specific_metabolite() will be used to create environmental media definitions. Additionally, spatial information can be assigned, such as the spatial dimensions in lattice units (layout.grid), spatial barriers, etc.

Parameters

input_objlist(cometspy.model) or str, optional

either a list of cometspy models, or a path to a previously-written COMETS layout to load.

Examples

>>> import cometspy as c
>>> import cobra.test
>>> e_cobra = cobra.test.create_test_model('textbook')
>>> e = c.model(e_cobra)
>>> e.open_exchanges() 
>>> l = c.layout([e])
>>> # use the media from the textbook in the COMETS media
>>> for key in e_cobra.medium.keys():
>>>     l.set_specific_metabolite(key[3:], 10.)
>>> l.display_current_media()

Attributes

modelslist(cometspy.model)

a list of cometspy models. do not modify directly. use add_model()

gridlist of size two

the number of boxes in the x and y direction. default = [1, 1]

mediapandas.DataFrame

info on met initial values and diffusion. see set_specific_metabolite

refreshlist

info on met added per time. see set_specific_refresh

local_mediadict

info on spatial-explicit media. see set_specific_metabolite_at_location

local_refreshdict

info on met added per time to loc. see set_specific_refresh_at_location

local_staticdict

info on constant met value at loc. see set_specific_static_at_location

default_diff_cfloat

the baseline diffusion constant for all metabolites

barrierslist(tuple)

list of impenetrable locations. see add_barriers

region_mapnumpy.array or None

integer array specifying regions. see set_region_map

region_parametersdict

region-specific default diffusion params. see set_region_parameters

reactionslist

list of extracellular reactions. see add_external_reaction

periodic_medialist

list of media undergoing periodic cycles. see set_global_periodic_media

add_barriers(barriers: list)

adds impenetrable barriers to the layout at specified locations

Barriers can be used in COMETS to create impenetrable spatial locations, for example to simulate rocks or to make a rounded simulation. Neither biomass nor metabolites can diffuse through barriers.

In COMETS, locations are zero-ordered.

Parameters

barrierslist(tuple)

A list containing tuples of integers of x,y locations.

Examples

>>> # make a diagonal line of barriers
>>> l = c.layout([model1, model2]) # assumes model1,2 are made
>>> l.grid = [5,5]
>>> l.add_barriers([(0,0),(1,1),(2,2),(3,3),(4,4)])
add_external_reaction(rxnName: str, metabolites: list, stoichiometry: list, **kwargs)

adds an extracellular reaction to the layout

External reactions are reactions not tied to a growing (living) model. They have rates defined in certain kwargs. There is no enyme concentration; it is assumed fixed.

Parameters

rxnNamestr

name of the external reaction

metaboliteslist(str)

list of metabolite names (strs)

stoichiometrylist(float)

stoichiometry of the metabolites. must be same order as metabolites

kwargs[‘Kcat’ or ‘Km’ or ‘K’]

rate parameters. Either Kcat and Km, or K must be set.

Examples

>>> # set an external reaction that converts lactose into glucose + gal
>>> # assumes import cometspy as c and a model called m has been made
>>> rxn_name = 'lactase'
>>> metabolites = ['lcts_e', 'glc__D_e', 'gal_e']
>>> stoichiometry = [-1., 1., 1.]
>>> K = 0.2 # mmol / hr
>>> l = c.layout([m]) # m is a cometspy.model that has exchanges for the mets above
>>> l.add_external_reaction(rxn_name, metabolites, stoichiometry,
                            K = K)
add_model(model)

add a cometspy model to this layout

This is the preferred way to add models to existing layout.

Parameters

modelcometspy.model

a cometspy.model with an initial_pop

Examples

>>> import cobra.test
>>> import cometspy as c
>>> layout = c.layout()
>>> ecoli = cobra.test.create_test_model("ecoli")
>>> salmonella = cobra.test.create_test_model("salmonella")
>>> ecoli = c.model(ecoli)
>>> salmonella = c.model(salmonella)
>>> ecoli.open_exchanges()
>>> salmonella.open_exchanges()
>>> ecoli.initial_pop = [0, 0, 1.e-7]
>>> salmonella.initial_pop = [0, 0, 1.e-7]
>>> layout.add_model(ecoli)
>>> layout.add_model(salmonella)
add_typical_trace_metabolites(amount: float = 1000.0, static: bool = True)

adds common BIGG-style metabolites to the environment. This only works if your models use metabolites which are formatted like ‘ca2_e’ and then should still only be used as a starting point.

Parameters

amountfloat, optional

the amount, in mmol, of these metabolites. Default is 1000.0.

staticbool, optional

The default is True. If false, these metabolites are depletable.

delete_model_files(working_dir)

deletes model files in specified directory

display_current_media()

a utility function to show current non-zero media amounts

get_model_ids() list

returns a list of the ids of the models in this layout.

read_comets_layout(input_obj: str)

load a comets layout from a file

If a COMETS layout file has been previously made, this function can load that file into the current cometspy layout object.

This is an unusual way of working with cometspy with rare uses.

Parameters

input_objstr

the path to the existing COMETS layout file to be loaded

set_global_periodic_media(metabolite: str, function: str, amplitude: float, period: float, phase: float, offset: float)

set cyclical changes in extracellular nutrient environment

Parameters

metabolitestr

name of the metabolite under a cycle

functionstr, one of [‘step’, ‘sin’, ‘cos’, ‘half_sin’, ‘half_cos’]

name of function that defines changing metabolite concentration

amplitudefloat

amplitude of the function in mmol

periodfloat

duration of the function period in hr

phasefloat

horizontal phase of the period (hr)

offsetfloat

vertical offset of the function (mmol)

set_metabolite_diffusion(diffusion_constant: float)

sets the default diffusion constant for all metabolites

Specific diffusion constants can be set subsequent to setting this.

Parameters

diffusion_constantfloat

the diffusion constant, in cm2 / s, that all metabolites use

set_models_pairs_frictions(pairs_frictions: list)

sets the values of the friction paremeters between the organisms and substrate

Parameters

modelint

the order of the model

locationtuple

the (x, y) location of the fixed metabolite amount. x, y are int

amountfloat

the amount, in mmol, that is in each box at each time step

set_models_substrate_frictions(friction: list)

sets the values of the friction paremeters between the organisms and substrate

Parameters

modelint

the order of the model

locationtuple

the (x, y) location of the fixed metabolite amount. x, y are int

amountfloat

the amount, in mmol, that is in each box at each time step

set_region_map(region_map: array)

set a region_map dictating different regions in a spatial simulation

COMETS can have different regions with different substrate diffusivities and frictions. Here, you set the map defining the regions. Specifically, you provide either: 1) a numpy array whose shape == layout.grid, or 2) a list of lists whose first length is grid[0] and second len is grid[1]

Populating these objects should be integer values, beginning at 1 and incrementing only, that define the different grid areas. These are intimately connected to region_parameters, which are set with layout.set_region_parameters()

Parameters

region_mapnumpy.array(int) or list(list)

a 2d array of size layout.grid containing integers, or a list of lists which can be coerced to an integer array using np.array()

Examples

see set_region_parameters for an example

set_region_parameters(region: int, diffusion: list, friction: float)

set regions-specific diffusion and friction values

COMETS can have different regions with different substrate diffusivities and frictions. Here, you set those parameters.

This does not affect a simulation unless a region map is also set, using the layout.set_region_map() function.

Parameters

regionint

the region (set by set_region_map) that uses the other parameters

diffusionlist(float)

a list containing diffusion constants (cm2/s) for each metabolite

frictionfloat

a single value for the friction observed in this region

Examples

>>> # this example assumes you have already created a cometspy model
>>> # called "e" and have imported cometspy as c and numpy as np
>>> # Here we are making a 2x2 world with two regions in which region
>>> # 1 metabolites have diffusion constant 1.e-6 cm2/s and in region
>>> # 2 metabolites have diffusion constant 1.e-7 cm2/s, and in both
>>> # cases mets have default friction constant 1.
>>> l = c.layout([e])
>>> l.grid = [2,2]
>>> region_map = np.array([[1, 1], [2, 2]], dtype = int)
>>> l.set_region_map(region_map)
>>> n_mets = length(l.media) # how many metabolites there are
>>> l.set_region_parameters(1, [1.e-6] * n_mets, 1.)
>>> l.set_region_parameters(1, [1.e-7] * n_mets, 1.)
set_specific_metabolite(met: str, amount: float, static: bool = False)

sets the initial (or constant) value for a metabolite across space

This is the most common way to set a metabolite concentration. It sets the initial value of the named metabolite to ‘amount’ in every location. Optionally, if static = True, this concentration is fixed during the simulation.

Parameters

metstr

name of the metabolite

amountfloat

initial value for the metabolite in each box, in mmol

staticbool, optional

DEFAULT = False. If True, the ‘amount’ is fixed over time.

Examples

>>> l = c.layout([model])
>>> l.set_specific_metabolite("glc__D_e", 0.005)
>>> l.set_specific_metabolite("o2_e", 15, static = True)
>>> # a common thing to do is populate this from the cobra model, as 
>>> # shown here:
>>> import cobra.test
>>> import cometspy as c
>>> cobra_model = cobra.test.create_test_model("textbook")
>>> model = c.model(cobra_model)
>>> model.open_exchanges()
>>> l = c.layout([model])
>>> for key, value in cobra_model.medium.items():
>>>     l.set_specific_metabolite(key[3:], value)
set_specific_metabolite_at_location(met: str, location: tuple, amount: float)

set an initial value of a metabolite at a specific spatial location

In contrast to set_specific_metabolite, which sets initial values universally over space, this method sets a specific metabolite initial value at one location only.

Parameters

metstr

name of the metabolite

locationtuple

a tuple of integers specifying a specific location, e.g. (3,2)

amountfloat

the initial value for the metabolite at the location, in mmol

Examples

>>> l = c.layout([model])
>>> l.grid = [10,1]
>>> l.set_specific_metabolite_at_location("glc__D_e", (9,0), 0.005)

Notes

Remember that COMETS locations are zero-ordered.

set_specific_metabolite_diffusion(met: str, diffusion_constant: float)

sets the diffusion constant of a specific metabolite

Parameters

metstr

name of the metabolite

diffusion_constantfloat

the diffusion constant, in cm2/s, of the named metabolite

Examples

>>> l = c.layout([m])
>>> l.set_specific_metabolite_diffusion("ac_e", 5.2e-6)
>>> l.set_specific_metabolite_diffusion("gal_e", 4.e-6)
set_specific_refresh(met: str, amount: float)

set the amount to increment a metabolite continuously through time

metabolites can be “refreshed” or “dripped in” during a COMETS simulation. This sets the amount added each hour, which is divided up by the number of time steps. It applies that refresh to all boxes in space.

Parameters

metstr

name of the metabolite refreshed

amountfloat

the amount added to (or subtracted from) each box, in mmol / hr.

set_specific_refresh_at_location(met: str, location: tuple, amount: float)

sets the amount to adjust a metabolite at a specific location in time

Like set_specific_refresh, but for a specific lacation. This Method is used to specify how much a metabolite is increased (or decreased) each hour in a specific location.

Parameters

metstr

the name of the metabolite

locationtuple

the (x, y) location to adjust, with x and y of type int

amountfloat

the amount added to (or subtracted from) the box, in mmol / hr

set_specific_static(met: str, amount: float)

sets a metabolite to a fixed value in every location

Parameters

metstr

name of the metabolite

amountfloat

the amount, in mmol, that is in each box at each time step

set_specific_static_at_location(met: str, location: tuple, amount: float)

sets a metabolite to a fixed value at a specific location

Parameters

metstr

name of the metabolite

locationtuple

the (x, y) location of the fixed metabolite amount. x, y are int

amountfloat

the amount, in mmol, that is in each box at each time step

update_models()

updates layout properties when models are added / removed

This is usually just called internally. However, if you manually change layout.models, then this method should be called to make other necessary changes.

write_layout(working_dir: str, to_append='')

writes just the COMETS layout file to the supplied path

Parameters

working_dirstr

the path to the directory where .current_layout will be written

write_model_files(working_dir='')

writes each model file

write_necessary_files(working_dir: str, to_append='')

writes the layout and the model files to file

This method is called by comets.run(). Alternatively it may be called directly if one wishes to inspect comets layout and model files.

A path must be specified if called directly.

Note: the layout file will be called “.current_layout” and therefore be hidden from file browsers by default.

Parameters

working_dirstr,

The directory where the files will be written.

to_appendstr,

String to append to written filenames

cometspy.model module

class cometspy.model.model(model: Model = None, randomtag: bool = False)

Bases: object

a COMETS metabolic model to use in a layout

A model contains information including the stoichiometric matrix, reaction boundaries and Michaelis-Menten parameters, founder populations, and various attributes related to motion in space.

A model is usually initiated by supplying a cobrapy model, and to be used in a COMETS simulation must minimally be given an initial population:

Parameters

modelcobra.Model or str, optional

Either a cobrapy Model, a path to a .cmd COMETS model, a path to a cobra sbml model which could be loaded with cobra.io.read_sbml_model() or None.

Attributes

initial_poplist

list of lists of spatial founder pops. e.g.[[x, y, grams], [x,y,grams]]

idstr

the name of the model. should be unique.

reactionspandas.DataFrame

DataFrame containing reaction information including boundaries.

smatpandas.DataFrame

DataFrame containing the stoichiometric matrix

metabolitespandas.DataFrame

single-column DataFrame containing the names of the metabolites

signalspandas.Dataframe

DataFrame on reactions whose bounds depend on metabolite concentrations

default_vmaxfloat

default uptake vmax in mmol / gDW /hr for Michaelis-Menten kinetics

default_kmfloat

default half-max metabolite concentration (M) for M-M kinetics

default_hillfloat

default hill coefficient for hill-like uptake kinetics

default_boundslist

list of two floats that are the default lower and upper reaction bounds

obj_stylestr

one of MAXIMIZE_OBJECTIVE_FLUX (fba) or MAX_OBJECTIVE_MIN_TOTAL (pfba)

optimizerstr

one of “GUROBI” or “GLPK”. not all functionality works with GLPK

Examples

>>> import cobra.test
>>> import cometspy as c
>>> ecoli = cobra.test.create_test_model("ecoli")
>>> model = c.model(ecoli)
>>> model.initial_pop = [0, 0, 1.e-12] # puts 1.e-12 gDW of biomass at 0,0
>>> print(model.id)
iJO1366
add_convection_parameters(packedDensity: float = 1.0, elasticModulus: float = 1.0, frictionConstant: float = 1.0, convDiffConstant: float = 1.0)

sets the parameters for biomass spread to use convection

This also requires one set the biomassMotionStyle to ‘Convection 2D’ in the comets.params object (see Example)

Parameters

packedDensity : float, optional elasticModulus : float, optional frictionConstant : float, optional convDiffConstant : float, optional

Example

import cobra.test import cometspy as c model = c.model(cobra.test.create_test_model(“ecoli”)) model.add_convection_parameters(1., 0.5, 2., 10.e-5) params = c.params() params.set_param(“biomassMotionStyle”, “Convection 2D”)

add_light(reaction: str, abs_coefficient: float, abs_base: float)

Causes a reaction to function in response to light

Parameters

reactionstr

the name of the reaction affected

abs_coefficientfloat

the absorption coefficient

abs_basefloat

absorption baseline

add_multitoxin(rxn_num, exch_ids, bound, vmax, kms, hills)

Adds a signaling relationship where >= 1 signal close a bound

This incorporates a hill-function reduction of (usually) an upper bound It is useful for incorporating multiple, additively-functioning signals.

bound = vmax * MULT((km_i^h_i) / (km_i^h_i + M_i^h_i)) where MULT is a multiplicative function, M_i is the molarity of the metabolite in the given exch_id, and km_i and h_i are the km and h for for that signal.

Parameters

rxn_numint

number of the affected reaction

exch_idslist[int]

exchange numbers of the signals

boundstr

“ub” or “lb”

vmaxfloat

maximum bound in absence of signals

kmslist[float]

kms for the signals

hillslist[float]

hill coefficients for the signals

Returns

None.

add_multspecies_convmodel_parameters(pressureKappa: float, pressureExponent: float, packBiomass: float, maxPressure: float)
add_neutral_drift_parameter(neutralDriftSigma)

sets the neutralDriftSigma parameter

add_noise_variance_parameter(noiseVariance: float)

sets the noise variance parameter

Parameters

noiseVariance : float

add_nonlinear_diffusion_parameters(zero: float = 1.0, n: float = 1.0, exponent: float = 1.0, hilln: float = 10.0, hillk: float = 0.9)

sets the model to use non-linear diffusion biomass spread.

This also requires one set the biomassMotionStyle to ‘ConvNonlin Diffusion 2D’ in the comets.params object (see Example)

Parameters

zero : float, optional n : float, optional exponent : float, optional hilln : float, optional hillk : float, optional

Example

>>> import cobra.test
>>> import cometspy as c
>>> model = c.model(cobra.test.create_test_model("ecoli"))
>>> model.add_nonlinear_diffusion_parameters(1., 1., 2., 5., 0.5)
>>> params = c.params()
>>> params.set_param("biomassMotionStyle", "ConvNonlin Diffusion 2D")
add_signal(rxn_num, exch_id, bound, function, parms)

adds a signal to the reaction rxn_num that changes bounds

Parameters

rxn_numint or ‘death’

reaction affected (int) or whether the signal causes death ‘death’

exch_idint

the number of the exchange reaction. see model.metabolites

boundstr ‘ub’ or ‘lb’

specifies whether the upper or lower bound is affected

functionstr ‘linear’ or ‘bounded_linear’ or ‘generalized_logistic’

specifies the function relating the metabolite conc. to the bound

parmslist (float)

a list of floats for the function

change_bounds(reaction: str, lower_bound: float, upper_bound: float)

changes the bounds of the specified reaction

Parameters

reactionstr

the name of the reaction

lower_boundfloat

the new value of the reaction’s lower bound in mmol / gDW / hr

upper_boundfloat

the new value of the reaction’s upper bound in mmol / gDW / hr

change_hill(reaction, hill)

changes the hill coefficient of the specified reaction.

Parameters

reactionstr

the name of the reaction

hillfloat

the new hill coef. for the reaction, for Hill-like kinetics

change_km(reaction, km)

changes the km of the specified reaction.

Parameters

reactionstr

the name of the reaction

kmfloat

the new km value for the reaction, for Monod (Michaelis-Menten) kinetics

change_optimizer(optimizer)

changes the optimizer to a specified one.

Parameters

optimizerstr

the name of the optimizer

change_vmax(reaction: str, vmax: float)

changes the vmax of the specified reaction.

Parameters

reactionstr

the name of the reaction

vmaxfloat

the new vmax value for the reaction, for Monod (Michaelis-Menten) kinetics

delete_comets_model(working_dir=None)

deletes a file version of this model if it exists.

Parameters

working_dirstr, optional

the directory where the comets model file exists

ensure_sinks_are_not_exchanges(suffix: str = '_c')

set exchange reactions ending in suffix (default = ‘_c’) to sink

Parameters

suffixstr, optional

the suffix to look for in exchange reactions. default = ‘_c’

get_bounds(reaction: str) tuple

returns a tuple (lb, ub) for the specified reaction

Parameters

reactionstr

the name of the reaction

get_exchange_metabolites() list

returns a list of the names of the exchange metabolites

get_reaction_names() list

returns a list of reaction names

load_cobra_model(curr_m: Model, randomtag: bool = False)

creates the COMETS model from the supplied cobra model

This is usually used internally when creating a COMETS model.

Parameters

curr_mcobra.Model

the cobra model object to be converted to the COMETS model object

Example

>>> import cometspy as c
>>> import cobra.test
>>> ecoli = cobra.test.create_test_model('ecoli')
>>> model = c.model()
>>> model.load_cobra_model(ecoli)
open_exchanges(lower_bound: float = -1000.0, upper_bound: float = 1000.0)

sets all exchange boundaries to the specifed values

This method is useful when first making a COMETS model so that medium definitions are not carried over from cobra into COMETS.

Parameters

lower_boundfloat, optional

default = -1000. in units of mmol / gDW / hr

upper_boundfloat, optional

default = 1000. in units of mmol / gDW / hr

read_cobra_model(path: str, randomtag: bool = False)

reads a cobra model from a file and loads it into this model

This is an alternative way to initialize a COMETS model, by supplying it with the path to a cobra model that can be read using cobra.io.read_sbml_model.

Parameters

pathstr

path to a cobra model in sbml format

Example

>>> import cometspy as c
>>> path_to_file = "./my_new_model.xml"
>>> model = c.model()
>>> model.read_cobra_model(path_to_file)
read_comets_model(path: str, randomtag: bool = False)

an alternative way to create a COMETS model by reading from a file

This allows a user to load previously-saved COMETS models. The contents populate this object’s attributes.

Parameters

pathstr

the path to the COMETS model file

Example

>>> import cometspy as c
>>> path_to_model = "./iJO1366.cmd" # this must actually exist
>>> model = c.model()
>>> model.read_comets_model(path_to_model)
write_comets_model(working_dir: str = None)

writes the COMETS model object to a file

This writes the current model object in COMETS format to file, either in the current working directory or in the directory specified by the optional parameter working_dir. It is mostly used by the comets object to run simulations, though a user may with to write these to examine directly or to save.

Parameters

working_dirstr, optional

a directory path to put COMETS files. for example “./data_files/”

cometspy.params module

class cometspy.params.params(global_params: str = None, package_params: str = None)

Bases: object

object storing simulation-related and some default biological parameters

The params object is an essential object needed to start a comets simulation, along with a layout object containing models. The params object stores simulation information such as timeStep and maxCycles, as well as many default biological parameters such as defaultKm and defaultDiffC (the default metabolite diffusion constant).

All of the parameters are stored in a dictionary called ‘all_params’ and can either be adjusted directly or using set_param().

When params are written to a COMETS object for COMETS to process, two files are generated (‘.current_global’ and ‘current_package’)

Parameters

global_paramsstr, optional

optional path to an existing ‘global_params’ file to load

package_paramsstr, optional

optional path to an existing ‘package_params’ file to load

Examples

>>> # make a params object and change two params
>>> import cometspy as c
>>> params = c.params()
>>> params.set_param("timeStep", 0.01) # hours
>>> params.set_param("spaceWidth", 0.5) # cm
>>> # view all the params
>>> params.show_params()
>>> # access the params directly
>>> params.all_params # this is a dict
>>> # use a params object in a comets simulation
>>> import cobra.test
>>> model = c.model(cobra.test.create_test_model("textbook"))
>>> model.initial_pop = [0, 0, 1.e-7]
>>> model.open_exchanges()
>>> layout = c.layout([model])
>>> # normally we'd now alter the media in the layout
>>> layout.add_typical_trace_metabolites()
>>> layout.set_specific_metabolite('lcts_e', 0.05) # mmol
>>> sim = c.comets(layout, params)
>>> sim.run()

Attributes

all_paramsdict

contains every possible param as keys, with their value as the value

params_unitsdict

dictionary containing the units for each possible parameter

get_param(name: str) object

returns the value associated with the given parameter name

Parameters

namestr

the name of the parameter whose value is desired

Returns

object

the type of object depends on the parameter requested

set_param(name: str, value)

alters a specific parameter value

See show_params() for a list of possible parameters to adjust.

Parameters

namestr

the name of a specific parameter to change

valuevariable

the type of the value depends on the parameter.

Examples

>>> p = cometspy.params()
>>> p.set_param("writeBiomassLog", True)
>>> p.set_param("maxCycles", 100)
>>> p.set_param("BiomassLogRate",5)
>>> p.set_param("defaultVmax", 10.24)
show_params()

utility function to show all possible parameters and their units

write_params(out_glb: str, out_pkg: str)

writes params data to the specified files

Usually this is only used internally, though a user could use it to pre-generate params files.

Parameters

out_glbstr

the path and name of the ‘global’ parameters file to be written

out_pkgstr

the path and name of the ‘package’ parameters file to be written

cometspy.utils module

The utils module contains helper functions generating spatial patterns.

cometspy.utils.chemostat(models: list, reservoir_media: dict, dilution_rate: float) tuple

helper function to let a user skip some steps when generating a chemostat

This sets relevant simulation parameters (e.g. deathRate, metaboliteDilutionRate) and layout values (e.g. refresh media) based upon a “reservoir” definition and a dilution rate.

It generates a layout that has the reservoir_media as the initial values, as well as set it to drip in / out based upon the dilution rate.

The returned layout and params can be further modified before supplying to a comets object if desired.

Parameters

modelslist(cometspy.model)

list of cometspy.model(s) with initial_pop set to use in the sim

reservoir_mediadict

media definition with metabolite names as keys as mmol amounts as values

dilution_ratefloat

the dilution rate of the chemostat, in 1/hr

Returns

tuple (layout, params)

a cometspy.layout object and a cometspy.params object

Examples

>>> import cobra.test
>>> import cometspy as c
>>> from cometspy.utils import chemostat
>>> # make a model from a cobra model, open exchange reactions, and give a pop
>>> tb = cobra.test.create_test_model("textbook")
>>> m = c.model(tb)
>>> m.initial_pop = [0, 0, 1.e-4]
>>> m.open_exchanges()
>>> reservoir = {'glc__D_e' : 0.01, 'nh4_e' : 1000., 'pi_e' : 1000.}
>>> layout, params = chemostat([m], reservoir, 0.1)
>>> params.set_param("maxCycles", 100)
>>> sim = c.comets(layout, params)
>>> sim.run()
>>> print(sim.total_biomass)
cometspy.utils.grow_rocks(n: int, xrange: tuple, yrange: tuple, mean_size: int) list

grows simple simulated rocks by adding random adjacent points from seeds

n number of seed points are generated first with pick_random_locations. Then, mean_size * n - n additional points are added to these seed locations. For each new point, one random location out of the set of all possible unoccupied locations next to occupied locations are chosen. Only lattice points directly to the NSEW are considered. This process is repeated until all new points are assigned. This usually results in ‘rocks’ of different shapes and sizes.

This function can be very slow (> 1 min) when n * mean_size > 2000

Parameters

nint

number of seed points to generate rocks from

xrangetuple

x range possible, e.g. (0, 5)

yrangetuple

y range possible, e.g. (0, 10)

mean_sizeint

average size in lattice units of a generated rock

Returns

list

list of all points generated including seed points

cometspy.utils.pick_random_locations(n: int, xrange: tuple, yrange: tuple, forbidden_locs: set = {}) list

returns a list of n x,y tuples corresponding to locations in the range given

Parameters

nint

number of locations desired

xrangetuple

the x-range (min, max) of the x range possible

yrangetuple

the y-range (min, max) of the y range possible

forbidden_locsset, optional

A list of tuples that cannot be chosen.

Returns

list

a list of (x,y) values.

Examples

>>> from cometspy.utils import pick_random_locations
>>> locs = pick_random_locations(3, (0, 10), (0,10))
>>> locs

Module contents

a python package for making and running COMETS simulations.

cometspy is a straight-forward python interface to the COMETS program. Use for dynamic flux-balance analysis, with multiple models, heterogenous spatio-temporal environments, and evolution.

To use this package, you must also have the actual COMETS program installed.

cometspy workflow:

models -> layout, layout + params -> comets, comets.run()

cometspy hello.world

>>> import cobra.test
>>> import cometspy as c
>>> # make a model from a cobra model, open exchange reactions, and give a pop
>>> tb = cobra.test.create_test_model("textbook")
>>> m = c.model(tb)
>>> m.initial_pop = [0, 0, 1.e-4]
>>> m.open_exchanges()
>>> # make a layout with the model, and give it three nutrients
>>> l = c.layout([m])
>>> l.set_specific_metabolite("glc__D_e", 0.01)
>>> l.set_specific_metabolite("nh4_e", 1000, static = True)
>>> l.set_specific_metabolite("pi_e", 1000, static = True)
>>> # make params and change one default
>>> p = c.params()
>>> p.set_param("maxCycles", 100)
>>> # make a simulation and run it!
>>> sim = c.comets(l, p)
>>> sim.run()
>>> print(sim.total_biomass)
>>> # saved data objects are pandas.DataFrames, and therefore can be plotted
>>> # sim.total_biomass.plot(x = "cycle")