Welcome to framed!¶
Contents:
Installation¶
framed requires Python 2.7, which is available for all major operating systems. We recommend the Anaconda python distribution.
framed can be easily installed using pip:
pip install framed
Dependencies (automatically installed):
- libsbml (import/export of models in SBML format)
- sympy (parsing GPR associations)
- scipy (kinetic modeling)
- matplotlib (plotting methods)
Additional dependences (not automatically installed):
- seaborn (for some types of plots)
For constraint-based modeling you will need at least one solver installed. framed currently supports:
- Gurobi
- CPLEX
Getting started¶
Loading models¶
framed supports different kinds of metabolic models. In any case, loading a model is quite simple.
For constraint-based (cobra) models:
from framed import load_cbmodel
model = load_cbmodel('my_model.xml')
Note that different people have been using different modeling conventions to store constraint-based models in SBML format. framed handles this problem using flavors. We currently support two flavors:
- cobra: This is the format adopted by the cobra toolbox.
- fbc2: This is the new fbc2 extension for SBML.
If can optionally specify the flavor of the model you are loading and framed will automatically apply a set of rules to improve compatibility (for example, removing boundary metabolites):
model = load_cbmodel('my_model.xml', flavor='fbc2')
For kinetic models:
from framed import load_odemodel
model = load_odemodel('my_model.xml')
Saving models¶
If you apply any kind of changes to a model, you can easily save your model as follows:
from framed import save_sbml_model
save_sbml_model(model, 'my_output_file.xml')
For constraint-based models you can also specify a specific flavor (this is useful for converting between different flavors):
save_sbml_model(model, 'my_output_file.xml', flavor='cobra')
Human readable format¶
Additionally, framed introduces a new shorthand notation for import/export of models in an human-readable format
(You can see this format by typing print model
).
from framed import read_cbmodel_from_file, write_model_to_file
model = read_cbmodel_from_file('my_file.txt')
write_model_to_file(model, 'my_new_file.txt')
Note: This format is not yet available for kinetic models.
Model manipulation¶
Once you load a model you can easily access the model’s attributes perform different kinds of operations using the model object. Note that some operations might be specific to constraint-based or kinetic models.
You can quickly look at your model by printing it:
print model
R_ACALD: M_acald_c + M_coa_c + M_nad_c <-> M_accoa_c + M_h_c + M_nadh_c
R_ACALDt: M_acald_e <-> M_acald_c
R_ACKr: M_ac_c + M_atp_c <-> M_actp_c + M_adp_c
R_ACONTa: M_cit_c <-> M_acon_C_c + M_h2o_c
R_ACONTb: M_acon_C_c + M_h2o_c <-> M_icit_c
(...)
You can also inspect particular metabolites and reactions.
print model.reactions.R_PGI
R_PGI: M_g6p_c <-> M_f6p_c
print model.metabolites.M_g6p_c
D-Glucose-6-phosphate
print model.reactions.R_TKT1.gpr
(G_b2465 or G_b2935)
You can easily make changes to a model. For instance, let’s change the maximum glucose uptake rate:
model.reactions.R_EX_glc_e.lb = -10
But for programmatic access, you may want to use reactions as a dictionary instead:
for reaction in model.reactions.values():
reaction.lb = 0
For convenience, you can access some methods from the model class directly:
for r_id in model.reactions:
model.set_lower_bound(r_id, 0)
You can easily add a new reaction to a model (or replace an existing reaction with the same id):
model.add_reaction_from_str('R_new: A + 0.5 B --> 2 C')
For constraint-based models you can optionally define the flux bounds as well:
model.add_reaction_from_str('R_new: S + M_atp_c <-> P + M_adp_c [-10, 10]')
You can also remove reactions (as well as metabolites, genes or compartments):
model.remove_reaction('R_PGI')
model.remove_metabolites(['M_h2o_c', 'M_h_c'])
There is a lot more you can try! Just take a look into our API.
Model transformations¶
There are some utilities to make some overall model transformations. One obvious one is to remove blocked reactions and dead-end metabolites (only for constraint-based models):
from framed import simplify
simplify(model)
You can also create a fully irreversible model by splitting reversible reactions in two directions (useful for some applications):
from framed import make_irreversible
make_irreversible(model)
If you don’t want to change your original model, in both cases you can generate a new model:
new_model = simplify(model, inplace=False)
Constraint-based analysis¶
Simulation¶
framed implements several methods for simulating constraint-based models. Here are some examples:
from framed import FBA
solution = FBA(model)
print solution
Objective: 0.873921506968
Status: Optimal
All other simulation methods have a similar interface, just try:
from framed import pFBA, looplessFBA, MOMA, lMOMA, ROOM
If you have multiple solvers installed, you can tell framed to use a different solver:
from framed import set_default_solver
set_default_solver('cplex')
The simulation methods accept several additional arguments. Many of these arguments allow you to modify simulation parameters without changing the model.
For instance, you can easily change the model objective:
print FBA(model, objective={'R_ATPM': 1})
Objective: 175.0
Status: Optimal
You can also define additional constraints that will override the constraints in the model.
Note: The constraints can be specified as an interval or as a fixed value:
solution = pFBA(model, constraints={'R_EX_glc_e': (-5, 0), 'R_EX_o2_e': 0})
The solution object is interactive and you can use it for analysing different aspects of your results:
print solution.show_values(pattern='R_EX_')
R_EX_co2_e 22.8098
R_EX_glc_e -10
R_EX_h_e 17.5309
R_EX_h2o_e 29.1758
R_EX_nh4_e -4.76532
R_EX_o2_e -21.7995
R_EX_pi_e -3.2149
print solution.show_metabolite_balance('M_g6p_c', model)
[ --> o ] R_GLCpts 10
[ o --> ] R_Biomass_Ecoli_core_w_GAM -0.179154
[ o --> ] R_G6PDH2r -4.95998
[ o --> ] R_PGI -4.86086
print solution.show_metabolite_balance('M_g6p_c', model, percentage=True, sort=True)
[ --> o ] R_GLCpts 100.00%
[ o --> ] R_G6PDH2r -49.60%
[ o --> ] R_PGI -48.61%
[ o --> ] R_Biomass_Ecoli_core_w_GAM -1.79%
Flux variability analysis¶
To run flux variability analysis (FVA) just try:
from framed import FVA
result = FVA(model)
You can specify additional arguments such as required fraction of objective function
FVA(model, obj_percentage=0.9)
Since FVA is sensitive to the environmental conditions, you can easily override the model constraints:
FVA(model, constraints={'R_EX_o2': 0})
To speed-up computations, you can additionally specify a list of reactions that you are interested in analysing:
FVA(model, reactions=['R_PGI', 'R_PFK', 'R_PYK'])
If you just want to find and remove blocked reactions, then simply run:
from framed import blocked_reactions
blocked = blocked_reactions(model)
model.remove_reactions(blocked)
Better yet, if you can use the simplify method, which will also remove dead-end metabolites:
from framed import simplify
simplify(model)
Gene/reaction deletions¶
You can easily simulate gene and/or reaction deletions:
from framed import gene_deletion
genes = ['G_b2465', 'G_b2935']
solution = gene_deletion(model, genes)
from framed import reaction_deletion
reactions = ['R_PGI', 'R_PFK']
solution = reaction_deletion(model, reactions)
You can easily change the simulation method to any method available in framed:
gene_deletion(model, genes, method='pFBA')
gene_deletion(model, genes, method='MOMA')
gene_deletion(model, genes, method='lMOMA')
gene_deletion(model, genes, method='ROOM')
As always, you can easily override the model constraints without changing the model:
gene_deletion(model, genes, constraints={'R_EX_o2_e': 0})
Gene/reaction essentiality¶
You can calculate the set of essential genes (or reactions) as follows:
from framed import essential_genes
essential = essential_genes(model)
from framed import essential_reactions
essential = essential_reactions(model)
You can change the minimum growth threshold for which you consider a deletion to be lethal.
Let’s be more conservative and set it to 20% of the original growth rate:
essential_genes(model, min_growth=0.2)
Remember, gene (and reaction) essentiality is very much depedent on the environmental conditions.
Let’s try changing the carbon source:
essential_genes(model, constraints={'R_EX_glc_e': 0, 'R_EX_fru_e': (-10, 0)})
Omics-based methods¶
There are currently two simulation methods based on transcriptomics data implemented:
GIMME
from framed import GIMME
solution = GIMME(model, gene_expression)
Additional arguments (cutoff percentile, growth fraction) can be specified:
GIMME(model, gene_expression, cutoff=25, growth_frac=0.9)
E-Flux
from framed import eFlux:
solution = eFlux(model, gene_expression)
Additional arguments can be specified to scale the reaction rates to flux units using a measured reaction rate:
eFlux(model, gene_expression, scale_rxn='R_EX_glc_e', scale_value=11.5)
Strain design¶
framed doesn’t aim to be a strain design package. For that try Cameo instead.
Nonetheless, a few naive strain design methods are provided.
The first is a brute force approach that tries all possible combinations of gene/reaction deletions.
from framed import combinatorial_gene_deletion
objective = lambda x: x['R_EX_succ_e']
solutions = combinatorial_gene_deletion(model, objective, max_dels=3)
from framed import combinatorial_reaction_deletion
objective = lambda v: v['R_EX_succ_e']
solutions = combinatorial_reaction_deletion(model, objective, max_dels=3)
You can define more complex objective functions such as the BPCY, and you can easily change the simulation method:
biomass = model.detect_biomass_reaction()
BPCY = lambda v: v[biomass] * v['R_EX_succ_e'] / v['R_EX_glc_e']
solutions = combinatorial_gene_deletion(model, BPCY, method='MOMA', max_dels=3)
The other method is a heuristic hill-climbing approach:
from framed import greedy_gene_deletion
objective = lambda x: x['R_EX_succ_e']
solutions = greedy_gene_deletion(model, objective, max_dels=5)
from framed import greedy_reaction_deletion
objective = lambda x: x['R_EX_succ_e']
solutions = greedy_reaction_deletion(model, objective, max_dels=5)
You can fine-tune the balance between size of search space and speed by adjusting the population size:
solutions = greedy_gene_deletion(model, objective, max_dels=5, pop_size=20)
Plotting¶
framed provides some built-in plotting utilities. For instance, you can plot a flux envelope:
from framed import plot_flux_envelope
plot_flux_envelope(model, 'R_EX_o2_e', 'R_EX_glc_e')

Kinetic modeling¶
framed implements some basic (and experimental) support for working with kinetic models.
It now also supports models that contain assignment rules (see for example the Chassagnole 2002 E. coli model).
Simulation¶
Time-course¶
Running a simple time-course simulation (uses scipy odeint):
from framed import time_course
T, X = time_course(model, time=100)
You can change the number of integration steps:
time_course(model, time=100, steps=1000)
You can override model parameters without changing the model:
time_course(model, time=100, parameters={'Vmax1': 10.0, 'Vmax2': 20.0})
Steady-state¶
Find the steady-state metabolite concentrations and reaction rates:
from framed import find_steady_state
x_ss, v_ss = find_steady_state(model)
Again, you can easily override model parameters for simulation purposes:
find_steady_state(model, parameters={'Vmax1': 10.0})
Sampling¶
You can sample the steady-state solution space by manipulating the model parameters:
from framed import sample_kinetic_model
sample_kinetic_model(model, size=100)
You can define which parameters to sample:
sample_kinetic_model(model, size=100, parameters=['Vmax1', 'Vmax2'])
And you can define properties of the sampling distribution:
sample_kinetic_model(model, 100, distribution='normal', dist_args=(0, 1), log_scale=True)
Calibration¶
framed has basic support for calibrating kinetic models using time-course metabolomics data:
from framed import fit_from_metabolomics
fit_from_metabolomics(model, t_steps, data)
where data is a dictionary with time-course values for every measured metabolite.
You can specify which parameters you want to calibrate (by default it attempts to fit all parameters):
fit_from_metabolomics(model, t_steps, data, parameters=['Vmax1', 'Km1'])
and you can specify admissible bounds for the parameters:
fit_from_metabolomics(model, t_steps, data, bounds={'Vmax1': (0, 10), 'Km1': (0.1, 100)})
framed uses scipy’s minimize, function which implements several optimization methods. You can select the optimization method (see the scipy’s documentation for a list of available methods):
fit_from_metabolomics(model, t_steps, data, method='Nelder-Mead')
Plotting¶
framed implements some basic plotting utilities. For instance, you can plot a time-course simulation:
from framed import plot_timecourse
plot_timecourse(model, time=100)
You can define which metabolites you want to plot
plot_timecourse(model, time=100, metabolites=['S', 'P'])
You can overlay metabolomics data over your plots (for instance, after calibration):
plot_timecourse(model, time=100, data=my_data, data_steps=my_time_points)
As usual, you can override model parameters for a particular simulation without changing the model:
plot_timecourse(model, time=100, parameters={'Vmax': 10, 'Km': 0.1})
Here is a pratical example using the Chassagnole 2002 E. coli model:
plot_timecourse(model, 1e3, metabolites=['cglcex', 'cg6p', 'cpep', 'cpyr'],
xlabel='time', ylabel='concentration', parameters={'Dil': 0.2/3600})

Finally, you can plot flux sampling results for a chosen set of reactions.
Note: this can also be applied to flux sampling results obtained with constraint-based models.
from framed import plot_flux_sampling
plot_flux_sampling(model, sample, reactions=['vPGI', 'vPFK', 'vPK'])
