5. Simulating with FBA¶
Simulations using flux balance analysis can be solved using Model.optimize()
. This will maximize or minimize (maximizing is the default) flux through the objective reactions.
[1]:
from cobra.io import load_model
model = load_model("textbook")
Scaling...
A: min|aij| = 1.000e+00 max|aij| = 1.000e+00 ratio = 1.000e+00
Problem data seem to be well scaled
5.1. Running FBA¶
[2]:
solution = model.optimize()
print(solution)
<Solution 0.874 at 0x7f993e396048>
The Model.optimize() function will return a Solution object. A solution object has several attributes:
objective_value
: the objective valuestatus
: the status from the linear programming solverfluxes
: a pandas series with flux indexed by reaction identifier. The flux for a reaction variable is the difference of the primal values for the forward and reverse reaction variables.shadow_prices
: a pandas series with shadow price indexed by the metabolite identifier.
For example, after the last call to model.optimize()
, if the optimization succeeds it’s status will be optimal. In case the model is infeasible an error is raised.
[3]:
solution.objective_value
[3]:
0.8739215069684307
The solvers that can be used with cobrapy are so fast that for many small to mid-size models computing the solution can be even faster than it takes to collect the values from the solver and convert to them python objects. With model.optimize
, we gather values for all reactions and metabolites and that can take a significant amount of time if done repeatedly. If we are only interested in the flux value of a single reaction or the objective, it is faster to instead use model.slim_optimize
which only does the optimization and returns the objective value leaving it up to you to fetch other values that you may need.
[4]:
%%time
model.optimize().objective_value
CPU times: user 1.45 ms, sys: 8 µs, total: 1.46 ms
Wall time: 1.47 ms
[4]:
0.8739215069684307
[5]:
%%time
model.slim_optimize()
CPU times: user 86 µs, sys: 0 ns, total: 86 µs
Wall time: 88.9 µs
[5]:
0.8739215069684307
5.1.1. Analyzing FBA solutions¶
Models solved using FBA can be further analyzed by using summary methods, which output printed text to give a quick representation of model behavior. Calling the summary method on the entire model displays information on the input and output behavior of the model, along with the optimized objective.
[6]:
model.summary()
[6]:
Objective
1.0 Biomass_Ecoli_core = 0.8739215069684306
Uptake
Metabolite | Reaction | Flux | C-Number | C-Flux |
---|---|---|---|---|
glc__D_e | EX_glc__D_e | 10 | 6 | 100.00% |
nh4_e | EX_nh4_e | 4.765 | 0 | 0.00% |
o2_e | EX_o2_e | 21.8 | 0 | 0.00% |
pi_e | EX_pi_e | 3.215 | 0 | 0.00% |
Secretion
Metabolite | Reaction | Flux | C-Number | C-Flux |
---|---|---|---|---|
co2_e | EX_co2_e | -22.81 | 1 | 100.00% |
h2o_e | EX_h2o_e | -29.18 | 0 | 0.00% |
h_e | EX_h_e | -17.53 | 0 | 0.00% |
In addition, the input-output behavior of individual metabolites can also be inspected using summary methods. For instance, the following commands can be used to examine the overall redox balance of the model
[7]:
model.metabolites.nadh_c.summary()
[7]:
nadh_c
C21H27N7O14P2
Producing Reactions
Percent | Flux | Reaction | Definition |
---|---|---|---|
13.14% | 5.064 | AKGDH | akg_c + coa_c + nad_c --> co2_c + nadh_c + succoa_c |
8.04% | 3.1 | Biomass_Ecoli_core | 1.496 3pg_c + 3.7478 accoa_c + 59.81 atp_c + 0.361 e4p_c + 0.0709 f6p_c + 0.129 g3p_c + 0.205 g6p_c + 0.2557 gln__L_c + 4.9414 glu__L_c + 59.81 h2o_c + 3.547 nad_c + 13.0279 nadph_c + 1.7867 oaa_c + 0.5191 pep_c + 2.8328 pyr_c + 0.8977 r5p_c --> 59.81 adp_c + 4.1182 akg_c + 3.7478 coa_c + 59.81 h_c + 3.547 nadh_c + 13.0279 nadp_c + 59.81 pi_c |
41.58% | 16.02 | GAPD | g3p_c + nad_c + pi_c <=> 13dpg_c + h_c + nadh_c |
13.14% | 5.064 | MDH | mal__L_c + nad_c <=> h_c + nadh_c + oaa_c |
24.09% | 9.283 | PDH | coa_c + nad_c + pyr_c --> accoa_c + co2_c + nadh_c |
Consuming Reactions
Percent | Flux | Reaction | Definition |
---|---|---|---|
100.00% | -38.53 | NADH16 | 4.0 h_c + nadh_c + q8_c --> 3.0 h_e + nad_c + q8h2_c |
Or to get a sense of the main energy production and consumption reactions
[8]:
model.metabolites.atp_c.summary()
[8]:
atp_c
C10H12N5O13P3
Producing Reactions
Percent | Flux | Reaction | Definition |
---|---|---|---|
66.58% | 45.51 | ATPS4r | adp_c + 4.0 h_e + pi_c <=> atp_c + h2o_c + 3.0 h_c |
23.44% | 16.02 | PGK | 3pg_c + atp_c <=> 13dpg_c + adp_c |
2.57% | 1.758 | PYK | adp_c + h_c + pep_c --> atp_c + pyr_c |
7.41% | 5.064 | SUCOAS | atp_c + coa_c + succ_c <=> adp_c + pi_c + succoa_c |
Consuming Reactions
Percent | Flux | Reaction | Definition |
---|---|---|---|
12.27% | -8.39 | ATPM | atp_c + h2o_c --> adp_c + h_c + pi_c |
76.46% | -52.27 | Biomass_Ecoli_core | 1.496 3pg_c + 3.7478 accoa_c + 59.81 atp_c + 0.361 e4p_c + 0.0709 f6p_c + 0.129 g3p_c + 0.205 g6p_c + 0.2557 gln__L_c + 4.9414 glu__L_c + 59.81 h2o_c + 3.547 nad_c + 13.0279 nadph_c + 1.7867 oaa_c + 0.5191 pep_c + 2.8328 pyr_c + 0.8977 r5p_c --> 59.81 adp_c + 4.1182 akg_c + 3.7478 coa_c + 59.81 h_c + 3.547 nadh_c + 13.0279 nadp_c + 59.81 pi_c |
0.33% | -0.2235 | GLNS | atp_c + glu__L_c + nh4_c --> adp_c + gln__L_c + h_c + pi_c |
10.94% | -7.477 | PFK | atp_c + f6p_c --> adp_c + fdp_c + h_c |
5.2. Changing the Objectives¶
The objective function is determined from the objective_coefficient attribute of the objective reaction(s). Generally, a “biomass” function which describes the composition of metabolites which make up a cell is used.
[9]:
biomass_rxn = model.reactions.get_by_id("Biomass_Ecoli_core")
Currently in the model, there is only one reaction in the objective (the biomass reaction), with an linear coefficient of 1.
[10]:
from cobra.util.solver import linear_reaction_coefficients
linear_reaction_coefficients(model)
[10]:
{<Reaction Biomass_Ecoli_core at 0x7f9941fce128>: 1.0}
The objective function can be changed by assigning Model.objective, which can be a reaction object (or just it’s name), or a dict
of {Reaction: objective_coefficient}
.
[11]:
# change the objective to ATPM
model.objective = "ATPM"
# The upper bound should be 1000, so that we get
# the actual optimal value
model.reactions.get_by_id("ATPM").upper_bound = 1000.
linear_reaction_coefficients(model)
[11]:
{<Reaction ATPM at 0x7f9941fce828>: 1.0}
[12]:
model.optimize().objective_value
[12]:
174.99999999999983
We can also have more complicated objectives including quadratic terms.
5.3. Running FVA¶
FBA will not give always give unique solution, because multiple flux states can achieve the same optimum. FVA (or flux variability analysis) finds the ranges of each metabolic flux at the optimum.
[13]:
from cobra.flux_analysis import flux_variability_analysis
[14]:
flux_variability_analysis(model, model.reactions[:10])
[14]:
minimum | maximum | |
---|---|---|
ACALD | -2.623542e-14 | 0.000000e+00 |
ACALDt | -2.623542e-14 | 0.000000e+00 |
ACKr | -4.012477e-14 | 0.000000e+00 |
ACONTa | 2.000000e+01 | 2.000000e+01 |
ACONTb | 2.000000e+01 | 2.000000e+01 |
ACt2r | -4.012477e-14 | 0.000000e+00 |
ADK1 | 0.000000e+00 | 1.705303e-13 |
AKGDH | 2.000000e+01 | 2.000000e+01 |
AKGt2r | -1.451321e-14 | 0.000000e+00 |
ALCD2x | -2.273737e-14 | 0.000000e+00 |
Setting parameter fraction_of_optimium=0.90
would give the flux ranges for reactions at 90% optimality.
[15]:
cobra.flux_analysis.flux_variability_analysis(
model, model.reactions[:10], fraction_of_optimum=0.9)
[15]:
minimum | maximum | |
---|---|---|
ACALD | -2.692308 | 0.0 |
ACALDt | -2.692308 | 0.0 |
ACKr | -4.117647 | 0.0 |
ACONTa | 8.461538 | 20.0 |
ACONTb | 8.461538 | 20.0 |
ACt2r | -4.117647 | 0.0 |
ADK1 | 0.000000 | 17.5 |
AKGDH | 2.500000 | 20.0 |
AKGt2r | -1.489362 | 0.0 |
ALCD2x | -2.333333 | 0.0 |
The standard FVA may contain loops, i.e. high absolute flux values that only can be high if they are allowed to participate in loops (a mathematical artifact that cannot happen in vivo). Use the loopless
argument to avoid such loops. Below, we can see that FRD7 and SUCDi reactions can participate in loops but that this is avoided when using the looplesss FVA.
[16]:
loop_reactions = [model.reactions.FRD7, model.reactions.SUCDi]
flux_variability_analysis(model, reaction_list=loop_reactions, loopless=False)
[16]:
minimum | maximum | |
---|---|---|
FRD7 | 0.0 | 980.0 |
SUCDi | 20.0 | 1000.0 |
[17]:
flux_variability_analysis(model, reaction_list=loop_reactions, loopless=True)
[17]:
minimum | maximum | |
---|---|---|
FRD7 | 0.0 | 0.0 |
SUCDi | 20.0 | 20.0 |
5.3.1. Running FVA in summary methods¶
Flux variability analysis can also be embedded in calls to summary methods. For instance, the expected variability in substrate consumption and product formation can be quickly found by
[17]:
model.optimize()
model.summary(fva=0.95)
[17]:
Objective
1.0 ATPM = 175.0
Uptake
Metabolite | Reaction | Flux | Range | C-Number | C-Flux |
---|---|---|---|---|---|
glc__D_e | EX_glc__D_e | 10 | [9.5; 10] | 6 | 100.00% |
o2_e | EX_o2_e | 60 | [55.88; 60] | 0 | 0.00% |
Secretion
Metabolite | Reaction | Flux | Range | C-Number | C-Flux |
---|---|---|---|---|---|
ac_e | EX_ac_e | 0 | [-2.059; 0] | 2 | 0.00% |
acald_e | EX_acald_e | 0 | [-1.346; 0] | 2 | 0.00% |
akg_e | EX_akg_e | 0 | [-0.7447; 0] | 5 | 0.00% |
co2_e | EX_co2_e | -60 | [-60; -54.17] | 1 | 100.00% |
etoh_e | EX_etoh_e | 0 | [-1.167; 0] | 2 | 0.00% |
for_e | EX_for_e | 0 | [-5.833; 0] | 1 | 0.00% |
glu__L_e | EX_glu__L_e | 0 | [-0.6731; 0] | 5 | 0.00% |
h2o_e | EX_h2o_e | -60 | [-60; -54.17] | 0 | 0.00% |
h_e | EX_h_e | 0 | [-5.833; 0] | 0 | 0.00% |
lac__D_e | EX_lac__D_e | 0 | [-1.129; 0] | 3 | 0.00% |
nh4_e | EX_nh4_e | 0 | [0; 0.6731] | 0 | 0.00% |
pi_e | EX_pi_e | 0 | [0; 0.171] | 0 | 0.00% |
pyr_e | EX_pyr_e | 0 | [-1.346; 0] | 3 | 0.00% |
succ_e | EX_succ_e | 0 | [-0.875; 0] | 4 | 0.00% |
Similarly, variability in metabolite mass balances can also be checked with flux variability analysis.
[18]:
model.metabolites.pyr_c.summary(fva=0.95)
[18]:
pyr_c
C3H3O3
Producing Reactions
Percent | Flux | Range | Reaction | Definition |
---|---|---|---|---|
50.00% | 10 | [9.5; 10] | GLCpts | glc__D_e + pep_c --> g6p_c + pyr_c |
0.00% | 0 | [-1.129; 0] | LDH_D | lac__D_c + nad_c <=> h_c + nadh_c + pyr_c |
0.00% | 0 | [0; 8.75] | ME1 | mal__L_c + nad_c --> co2_c + nadh_c + pyr_c |
0.00% | 0 | [0; 8.75] | ME2 | mal__L_c + nadp_c --> co2_c + nadph_c + pyr_c |
50.00% | 10 | [1.25; 18.75] | PYK | adp_c + h_c + pep_c --> atp_c + pyr_c |
0.00% | 0 | [-1.346; 0] | PYRt2 | h_e + pyr_e <=> h_c + pyr_c |
Consuming Reactions
Percent | Flux | Range | Reaction | Definition |
---|---|---|---|---|
0.00% | 0 | [-0.1316; 0] | Biomass_Ecoli_core | 1.496 3pg_c + 3.7478 accoa_c + 59.81 atp_c + 0.361 e4p_c + 0.0709 f6p_c + 0.129 g3p_c + 0.205 g6p_c + 0.2557 gln__L_c + 4.9414 glu__L_c + 59.81 h2o_c + 3.547 nad_c + 13.0279 nadph_c + 1.7867 oaa_c + 0.5191 pep_c + 2.8328 pyr_c + 0.8977 r5p_c --> 59.81 adp_c + 4.1182 akg_c + 3.7478 coa_c + 59.81 h_c + 3.547 nadh_c + 13.0279 nadp_c + 59.81 pi_c |
100.00% | -20 | [-28.75; -13] | PDH | coa_c + nad_c + pyr_c --> accoa_c + co2_c + nadh_c |
0.00% | 0 | [-5.833; 0] | PFL | coa_c + pyr_c --> accoa_c + for_c |
0.00% | 0 | [-8.75; 0] | PPS | atp_c + h2o_c + pyr_c --> amp_c + 2.0 h_c + pep_c + pi_c |
In these summary methods, the values are reported as a the center point +/- the range of the FVA solution, calculated from the maximum and minimum values.
5.4. Running pFBA¶
Parsimonious FBA (often written pFBA) finds a flux distribution which gives the optimal growth rate, but minimizes the total sum of flux. This involves solving two sequential linear programs, but is handled transparently by cobrapy. For more details on pFBA, please see Lewis et al. (2010).
[19]:
model.objective = 'Biomass_Ecoli_core'
fba_solution = model.optimize()
pfba_solution = cobra.flux_analysis.pfba(model)
These functions will likely give completely different objective values, because the objective value shown from pFBA is defined as sum(abs(pfba_solution.fluxes.values))
, while the objective value for standard FBA is defined as the weighted flux through the reactions being optimized for (e.g. fba_solution.fluxes["Biomass_Ecoli_core"]
).
Both pFBA and FBA should return identical results within solver tolerances for the objective being optimized. E.g. for a FBA problem which maximizes the reaction Biomass_Ecoli_core
:
[20]:
abs(fba_solution.fluxes["Biomass_Ecoli_core"] - pfba_solution.fluxes[
"Biomass_Ecoli_core"])
[20]:
8.881784197001252e-16
[21]:
import numpy as np
np.isclose(
fba_solution.fluxes["Biomass_Ecoli_core"],
pfba_solution.fluxes["Biomass_Ecoli_core"]
)
[21]:
True
5.5. Running geometric FBA¶
Geometric FBA finds a unique optimal flux distribution which is central to the range of possible fluxes. For more details on geometric FBA, please see K Smallbone, E Simeonidis (2009).
[22]:
geometric_fba_sol = cobra.flux_analysis.geometric_fba(model)
geometric_fba_sol
[22]:
fluxes | reduced_costs | |
---|---|---|
ACALD | 0.000000e+00 | 0.0 |
ACALDt | 0.000000e+00 | 0.0 |
ACKr | 7.454685e-15 | 0.0 |
ACONTa | 6.007250e+00 | 0.0 |
ACONTb | 6.007250e+00 | 0.0 |
... | ... | ... |
TALA | 1.496984e+00 | 0.0 |
THD2 | 0.000000e+00 | 0.0 |
TKT1 | 1.496984e+00 | 0.0 |
TKT2 | 1.181498e+00 | 0.0 |
TPI | 7.477382e+00 | 0.0 |
95 rows × 2 columns
[ ]: