Adding More Complex Constraints On A Model
Adding More Complex Constraints On A Model
Reviewers:
INTRODUCTION
The COBRA Toolbox offers the possibility to add additional constraints to a model that are not direct flux constraints. One
example is already indicated in the constraining models tutorial which introduces coupling constraints[1], but there are plenty of
other possibilities that can be achieved using constraints.
The tools introduced in this tutorial are mainly aimed at developers who want to implement complex algorithms or formalism
within the COBRA Toolbox environment, but the basics are also useful to integrate specific properties of a model or specific
literature data. The section "The COBRA Model structure" is aimed at developers and can be skipped by users who don't to go
into the details.
C - The matrix of additional constraints on the reactions, which are not derived from the steady state condition and
Overall A COBRA model will be translated into a Linear problem by combining these matrices in the following way:
The Cobra toolbox allows one sided inequality constraints and equality constraints on a model. i.e. there is currently no
mechanism to add a constraint like:
and
which is commonly what solvers will translate these constraints into anyways.
It is best practice to not manually alter the sizes of these matrices, but to use the appropriate model manipulation functions
(addReaction/addMultipleReactions, addMetabolite/addMultipleMetabolites, addCOBRAVariables, addCOBRAConstraints). These
functions will ensure, that the model structure stays in sync, and all necessary fields are updated appropriately. There is a
convenience function to create the corresponding LPproblem struct to the given model: buildLPproblemFromModel(). This
function will create and initialize all necessary fields and builds the corresponding LP problem. Any additional modifications of
the LP problem imposed by an algorithm should be done on this LP and not on the original model structure.
PROCEDURE
Initially, we will load a model that we want to modify by adding a few additional constraints. The model used will be the simple
E.Coli core model:
initCobraToolbox(false)
_____ _____ _____ _____ _____ |
/ ___| / _ \ | _ \ | _ \ / ___ \ | COnstraint-Based Reconstruction and Analysis
| | | | | | | |_| | | |_| | | |___| | | The COBRA Toolbox - 2018
| | | | | | | _ { | _ / | ___ | |
| |___ | |_| | | |_| | | | \ \ | | | | | Documentation:
\_____| \_____/ |_____/ |_| \_\ |_| |_| | http://opencobra.github.io/cobratoolbox
|
Now, there are two genes which code for aconitase in the ecoli core model: b0118 (aconB) and b1276 (aconA). Both the citrate
hydratase and the aconitatde dehydratase have the same GPR rule: (b0118 or b1276), so they can both use either enzyme.
aconAgene = 'b1276';
aconBgene = 'b0118';
We would like to add a constraint that not only restricts the activity of these two reactions, but also ensures that the turnover
rates are considered. One mg of aconA can support a total flux of ~0.367 mmol/h of ACONTa
printRxnFormula(model,'rxnAbbrList',{'ACONTa'},'gprFlag', true);
However, it will not be able to do both at the same time. Therefore, we can assume that in addition to its normal metabolites,
the reaction also consumes some of the enzyme (for this time step). However, both reactions can also be catalysed by aconB,
so they might actually not use any of aconA. Therefore, we need additional constraints that represent the activity through these
two reactions. In addition, we need variables that represent the efficiency of one mg of protein to catalyse the reactions.
The function takes an existing model, and a set of variable ids (idList) and generates all required fields for these variables. You
can modify the properties of the variables by the following parameters which can either be supplied as parameter/value pairs:
addCOBRAVariables(model,{'newVar'},'lb', -5,'ub',3);
or as a parameter struct:
params = struct();
params.lb = -5;
params.ub = 3;
addCOBRAVariables(model,{'newVar'},params);
We will now add the variables 'aconA' and 'aconB' with lower bounds 0 and upper bounds according to the values determined
above.
aconVars = {'aconA','aconB'};
model = addCOBRAVariables(model,aconVars,'lb',[0;0],'ub',[aconAAmount;aconBAmount]);
We further need a conversion between the used amount of aconA and the potential flux through ACONTa. We also need this for
ACONTb and the same for aconB.
linkedReactions = {'ACONTa','ACONTb'};
for enzyme = 1:numel(aconVars)
for linkedReaction = 1:numel(linkedReactions)
model = addCOBRAVariables(model,{strcat(aconVars{enzyme},'to',linkedReactions{linkedReaction})},'lb',0);
end
end
This function will add constraints to the given model using a list of reactions (idList) involved in the constraint and right hand
side values for each constraint (d). By default, the constraint is assumed to be a limiting constraint, and the coefficients of the
reactions is assumed to be one.
params = struct();
params.c = [0.2,0.4];
params.dsense = 'L';
constMod = addCOBRAConstraints(model,{'FBP','FBA'}, 5, params);
'c' - The coefficient matrix with one column per reaction id and one row per added constraint (default: 1 for each element
in idList)
'dsense' - the sense vector with one element per added constraint ('E' for equality, 'L' for lower than, 'G' for greater
than), or one element which is used for all constraints (default: 'L').
'ConstraintID' - a cell array of strings with one element for each added constraint (default: ConstraintXYZ, with XYZ
being the position in the ctrs vector)
checkDuplicates - Whether to check for duplicate Constraints (they don't hurt, but they don't help). Note that duplicate
IDs are still not allowed.
Coming back to our example, the first constraint we add is for the amount of aconA, which should be steady (i.e. not more than
made available by the aconA variable). More precisely, the amount enzyme made available for ACONTa and ACONTb, should be
balanced with the amount of enzyme made available by the aconA variable.
Finally, we have a a system, in which the two aconitase reactions are competing for the available enzymes.
orig_sol = optimizeCbModel(model_orig)
Modifying constraints allows the adjustment of the coefficients of the constraint (parameter 'c'), the directionality (parameter
'dsense') and right hand side of the constraint ( parameter 'd') as well as the name (same as for variables). There are several
ways to update the coefficients: Either a full row (including both values for the C and the E matrix - if any) needs to be provided,
or the IDs along with the coefficients (similar to the situation when generating the constraints) need to be provided. It is
important to note that ALL coefficients will be reset if a constraint is modified using the changeCOBRAConstraint function.
model = changeCOBRAVariable(model,'aconA','ub',aconAAmount*2);
less_restricted_sol = optimizeCbModel(model)
model = changeCOBRAConstraints(model,'ACONTaFlux','idList',{'aconAtoACONTa','aconBtoACONTa','ACONTa'},...
'c',[aconACit,aconBCit*0.5,-1]);
less_efficient_aconA = optimizeCbModel(model)
less_efficient_aconA = struct with fields:
full: [101×1 double]
obj: 0.0159
rcost: [101×1 double]
dual: [76×1 double]
slack: [76×1 double]
solver: 'gurobi'
algorithm: 'default'
stat: 1
origStat: 'OPTIMAL'
time: 0.0030
basis: [1×1 struct]
x: [95×1 double]
f: 0.0159
y: [72×1 double]
w: [95×1 double]
v: [95×1 double]
We can again see the objective drop.
While we only used a very simple example for this tutorial, this interplay can improve predictive qualities substantially (for more
have a look at e.g. [4])
References
[1] Thiele I, Fleming RM, Bordbar A, Schellenberger J, Palsson BØ. Functional characterization of alternate optimal solutions of
Escherichia coli's transcriptional and translational machinery. Biophys J. 98(10):2072-81 (2010).
[2] The UniProt Consortium, UniProt: the universal protein knowledgebase, Nucleic Acids Res. 45: D158-D169 (2017)
[3]Jacek R. Wiśniewski, Dariusz Rakus, Quantitative analysis of the Escherichia coli proteome, Data in Brief 1, 7-11, (2014)
[4] Sánchez et al, Improving the phenotype predictions of a yeast genome-scale metabolic model by incorporating enzymatic
constraints, Mol Sys Biol, 13:935 (2017)