{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Simulating with FBA" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Scaling...\n", " A: min|aij| = 1.000e+00 max|aij| = 1.000e+00 ratio = 1.000e+00\n", "Problem data seem to be well scaled\n" ] } ], "source": [ "from cobra.io import load_model\n", "model = load_model(\"textbook\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Running FBA" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n" ] } ], "source": [ "solution = model.optimize()\n", "print(solution)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Model.optimize() function will return a Solution object. A solution object has several attributes:\n", "\n", " - `objective_value`: the objective value\n", " - `status`: the status from the linear programming solver\n", " - `fluxes`: 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.\n", " - `shadow_prices`: a pandas series with shadow price indexed by the metabolite identifier." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "0.8739215069684307" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "solution.objective_value" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CPU times: user 1.45 ms, sys: 8 µs, total: 1.46 ms\n", "Wall time: 1.47 ms\n" ] }, { "data": { "text/plain": [ "0.8739215069684307" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "%%time\n", "model.optimize().objective_value" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CPU times: user 86 µs, sys: 0 ns, total: 86 µs\n", "Wall time: 88.9 µs\n" ] }, { "data": { "text/plain": [ "0.8739215069684307" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "%%time\n", "model.slim_optimize()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Analyzing FBA solutions" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

Objective

1.0 Biomass_Ecoli_core = 0.8739215069684306

Uptake

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
MetaboliteReactionFluxC-NumberC-Flux
glc__D_eEX_glc__D_e106100.00%
nh4_eEX_nh4_e4.76500.00%
o2_eEX_o2_e21.800.00%
pi_eEX_pi_e3.21500.00%

Secretion

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
MetaboliteReactionFluxC-NumberC-Flux
co2_eEX_co2_e-22.811100.00%
h2o_eEX_h2o_e-29.1800.00%
h_eEX_h_e-17.5300.00%
" ], "text/plain": [ "" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model.summary()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

nadh_c

C21H27N7O14P2

Producing Reactions

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
PercentFluxReactionDefinition
13.14%5.064AKGDHakg_c + coa_c + nad_c --> co2_c + nadh_c + succoa_c
8.04%3.1Biomass_Ecoli_core1.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.02GAPDg3p_c + nad_c + pi_c <=> 13dpg_c + h_c + nadh_c
13.14%5.064MDHmal__L_c + nad_c <=> h_c + nadh_c + oaa_c
24.09%9.283PDHcoa_c + nad_c + pyr_c --> accoa_c + co2_c + nadh_c

Consuming Reactions

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
PercentFluxReactionDefinition
100.00%-38.53NADH164.0 h_c + nadh_c + q8_c --> 3.0 h_e + nad_c + q8h2_c
" ], "text/plain": [ "" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model.metabolites.nadh_c.summary()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Or to get a sense of the main energy production and consumption reactions" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

atp_c

C10H12N5O13P3

Producing Reactions

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
PercentFluxReactionDefinition
66.58%45.51ATPS4radp_c + 4.0 h_e + pi_c <=> atp_c + h2o_c + 3.0 h_c
23.44%16.02PGK3pg_c + atp_c <=> 13dpg_c + adp_c
2.57%1.758PYKadp_c + h_c + pep_c --> atp_c + pyr_c
7.41%5.064SUCOASatp_c + coa_c + succ_c <=> adp_c + pi_c + succoa_c

Consuming Reactions

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
PercentFluxReactionDefinition
12.27%-8.39ATPMatp_c + h2o_c --> adp_c + h_c + pi_c
76.46%-52.27Biomass_Ecoli_core1.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.2235GLNSatp_c + glu__L_c + nh4_c --> adp_c + gln__L_c + h_c + pi_c
10.94%-7.477PFKatp_c + f6p_c --> adp_c + fdp_c + h_c
" ], "text/plain": [ "" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model.metabolites.atp_c.summary()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Changing the Objectives" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [], "source": [ "biomass_rxn = model.reactions.get_by_id(\"Biomass_Ecoli_core\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Currently in the model, there is only one reaction in the objective (the biomass reaction), with an linear coefficient of 1." ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "{: 1.0}" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "from cobra.util.solver import linear_reaction_coefficients\n", "linear_reaction_coefficients(model)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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}`." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "{: 1.0}" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# change the objective to ATPM\n", "model.objective = \"ATPM\"\n", "\n", "# The upper bound should be 1000, so that we get\n", "# the actual optimal value\n", "model.reactions.get_by_id(\"ATPM\").upper_bound = 1000.\n", "linear_reaction_coefficients(model)" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "174.99999999999983" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model.optimize().objective_value" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can also have more complicated objectives including quadratic terms." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Running FVA" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [], "source": [ "from cobra.flux_analysis import flux_variability_analysis" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
minimummaximum
ACALD-2.623542e-140.000000e+00
ACALDt-2.623542e-140.000000e+00
ACKr-4.012477e-140.000000e+00
ACONTa2.000000e+012.000000e+01
ACONTb2.000000e+012.000000e+01
ACt2r-4.012477e-140.000000e+00
ADK10.000000e+001.705303e-13
AKGDH2.000000e+012.000000e+01
AKGt2r-1.451321e-140.000000e+00
ALCD2x-2.273737e-140.000000e+00
\n", "
" ], "text/plain": [ " minimum maximum\n", "ACALD -2.623542e-14 0.000000e+00\n", "ACALDt -2.623542e-14 0.000000e+00\n", "ACKr -4.012477e-14 0.000000e+00\n", "ACONTa 2.000000e+01 2.000000e+01\n", "ACONTb 2.000000e+01 2.000000e+01\n", "ACt2r -4.012477e-14 0.000000e+00\n", "ADK1 0.000000e+00 1.705303e-13\n", "AKGDH 2.000000e+01 2.000000e+01\n", "AKGt2r -1.451321e-14 0.000000e+00\n", "ALCD2x -2.273737e-14 0.000000e+00" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "flux_variability_analysis(model, model.reactions[:10])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Setting parameter `fraction_of_optimium=0.90` would give the flux ranges for reactions at 90% optimality." ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
minimummaximum
ACALD-2.6923080.0
ACALDt-2.6923080.0
ACKr-4.1176470.0
ACONTa8.46153820.0
ACONTb8.46153820.0
ACt2r-4.1176470.0
ADK10.00000017.5
AKGDH2.50000020.0
AKGt2r-1.4893620.0
ALCD2x-2.3333330.0
\n", "
" ], "text/plain": [ " minimum maximum\n", "ACALD -2.692308 0.0\n", "ACALDt -2.692308 0.0\n", "ACKr -4.117647 0.0\n", "ACONTa 8.461538 20.0\n", "ACONTb 8.461538 20.0\n", "ACt2r -4.117647 0.0\n", "ADK1 0.000000 17.5\n", "AKGDH 2.500000 20.0\n", "AKGt2r -1.489362 0.0\n", "ALCD2x -2.333333 0.0" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" } ], "source": [ "cobra.flux_analysis.flux_variability_analysis(\n", " model, model.reactions[:10], fraction_of_optimum=0.9)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
minimummaximum
FRD70.0980.0
SUCDi20.01000.0
\n", "
" ], "text/plain": [ " minimum maximum\n", "FRD7 0.0 980.0\n", "SUCDi 20.0 1000.0" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "loop_reactions = [model.reactions.FRD7, model.reactions.SUCDi]\n", "flux_variability_analysis(model, reaction_list=loop_reactions, loopless=False)" ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
minimummaximum
FRD70.00.0
SUCDi20.020.0
\n", "
" ], "text/plain": [ " minimum maximum\n", "FRD7 0.0 0.0\n", "SUCDi 20.0 20.0" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "flux_variability_analysis(model, reaction_list=loop_reactions, loopless=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Running FVA in summary methods" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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" ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

Objective

1.0 ATPM = 175.0

Uptake

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
MetaboliteReactionFluxRangeC-NumberC-Flux
glc__D_eEX_glc__D_e10[9.5; 10]6100.00%
o2_eEX_o2_e60[55.88; 60]00.00%

Secretion

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
MetaboliteReactionFluxRangeC-NumberC-Flux
ac_eEX_ac_e0[-2.059; 0]20.00%
acald_eEX_acald_e0[-1.346; 0]20.00%
akg_eEX_akg_e0[-0.7447; 0]50.00%
co2_eEX_co2_e-60[-60; -54.17]1100.00%
etoh_eEX_etoh_e0[-1.167; 0]20.00%
for_eEX_for_e0[-5.833; 0]10.00%
glu__L_eEX_glu__L_e0[-0.6731; 0]50.00%
h2o_eEX_h2o_e-60[-60; -54.17]00.00%
h_eEX_h_e0[-5.833; 0]00.00%
lac__D_eEX_lac__D_e0[-1.129; 0]30.00%
nh4_eEX_nh4_e0[0; 0.6731]00.00%
pi_eEX_pi_e0[0; 0.171]00.00%
pyr_eEX_pyr_e0[-1.346; 0]30.00%
succ_eEX_succ_e0[-0.875; 0]40.00%
" ], "text/plain": [ "" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model.optimize()\n", "model.summary(fva=0.95)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Similarly, variability in metabolite mass balances can also be checked with flux variability analysis." ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [ { "data": { "text/html": [ "

pyr_c

C3H3O3

Producing Reactions

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
PercentFluxRangeReactionDefinition
50.00%10[9.5; 10]GLCptsglc__D_e + pep_c --> g6p_c + pyr_c
0.00%0[-1.129; 0]LDH_Dlac__D_c + nad_c <=> h_c + nadh_c + pyr_c
0.00%0[0; 8.75]ME1mal__L_c + nad_c --> co2_c + nadh_c + pyr_c
0.00%0[0; 8.75]ME2mal__L_c + nadp_c --> co2_c + nadph_c + pyr_c
50.00%10[1.25; 18.75]PYKadp_c + h_c + pep_c --> atp_c + pyr_c
0.00%0[-1.346; 0]PYRt2h_e + pyr_e <=> h_c + pyr_c

Consuming Reactions

\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
PercentFluxRangeReactionDefinition
0.00%0[-0.1316; 0]Biomass_Ecoli_core1.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]PDHcoa_c + nad_c + pyr_c --> accoa_c + co2_c + nadh_c
0.00%0[-5.833; 0]PFLcoa_c + pyr_c --> accoa_c + for_c
0.00%0[-8.75; 0]PPSatp_c + h2o_c + pyr_c --> amp_c + 2.0 h_c + pep_c + pi_c
" ], "text/plain": [ "" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "model.metabolites.pyr_c.summary(fva=0.95)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Running pFBA" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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)](http://dx.doi.org/10.1038/msb.2010.47)." ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [], "source": [ "model.objective = 'Biomass_Ecoli_core'\n", "fba_solution = model.optimize()\n", "pfba_solution = cobra.flux_analysis.pfba(model)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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\"]`).\n", "\n", "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`:" ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "8.881784197001252e-16" ] }, "execution_count": 20, "metadata": {}, "output_type": "execute_result" } ], "source": [ "abs(fba_solution.fluxes[\"Biomass_Ecoli_core\"] - pfba_solution.fluxes[\n", " \"Biomass_Ecoli_core\"])" ] }, { "cell_type": "code", "execution_count": 21, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "True" ] }, "execution_count": 21, "metadata": {}, "output_type": "execute_result" } ], "source": [ "import numpy as np\n", "np.isclose(\n", " fba_solution.fluxes[\"Biomass_Ecoli_core\"], \n", " pfba_solution.fluxes[\"Biomass_Ecoli_core\"]\n", ")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Running geometric FBA" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "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)](http://dx.doi.org/10.1016/j.jtbi.2009.01.027)." ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [ { "data": { "text/html": [ "Optimal solution with objective value 0.000
\n", "\n", "\n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", " \n", "
fluxesreduced_costs
ACALD0.000000e+000.0
ACALDt0.000000e+000.0
ACKr7.454685e-150.0
ACONTa6.007250e+000.0
ACONTb6.007250e+000.0
.........
TALA1.496984e+000.0
THD20.000000e+000.0
TKT11.496984e+000.0
TKT21.181498e+000.0
TPI7.477382e+000.0
\n", "

95 rows × 2 columns

\n", "
" ], "text/plain": [ "" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" } ], "source": [ "geometric_fba_sol = cobra.flux_analysis.geometric_fba(model)\n", "geometric_fba_sol" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.7.1" } }, "nbformat": 4, "nbformat_minor": 1 }