Skip to content

Commit c3acab3

Browse files
committed
Add made command
1 parent 16e1888 commit c3acab3

File tree

2 files changed

+351
-0
lines changed

2 files changed

+351
-0
lines changed

psamm/commands/made.py

Lines changed: 350 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,350 @@
1+
# This file is part of PSAMM.
2+
#
3+
# PSAMM is free software: you can redistribute it and/or modify
4+
# it under the terms of the GNU General Public License as published by
5+
# the Free Software Foundation, either version 3 of the License, or
6+
# (at your option) any later version.
7+
#
8+
# PSAMM is distributed in the hope that it will be useful,
9+
# but WITHOUT ANY WARRANTY; without even the implied warranty of
10+
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11+
# GNU General Public License for more details.
12+
#
13+
# You should have received a copy of the GNU General Public License
14+
# along with PSAMM. If not, see <http://www.gnu.org/licenses/>.
15+
#
16+
# Copyright 2014-2017 Jon Lund Steffensen <[email protected]>
17+
# Copyright 2016 Keith Dufault-Thompson <[email protected]>
18+
# Copyright 2016 Julie Cuddigan <[email protected]>
19+
20+
"""Metabolic Adjustment by Differential Expression (MADE) command."""
21+
22+
from __future__ import unicode_literals
23+
24+
import time
25+
import logging
26+
import math
27+
import csv
28+
29+
from six import iteritems
30+
31+
from ..command import SolverCommandMixin, MetabolicMixin, Command
32+
from ..util import MaybeRelative
33+
from ..lpsolver import lp
34+
from ..expression import boolean
35+
36+
# Module-level logging
37+
logger = logging.getLogger(__name__)
38+
39+
40+
class MadeFluxBalance(MetabolicMixin, SolverCommandMixin, Command):
41+
"""Run MADE flux balance analysis on the model.
42+
43+
Args:
44+
gene_var1 = Dictionary, key:value = gene expression objects:their new
45+
variable id, first set
46+
gene_var2 = Dictionary; key:value = gene expression objects:their new
47+
variable id, second set
48+
var_ineqvar1 = xi; Dictionary, key:value = new variable ids:their
49+
defined inequality variable, first set
50+
var_ineqvar2 = xi+1; Dictionary, key:value = new variable ids:their
51+
defined inequality variable, second set
52+
gene_pval = Dictionary, key:value = gene ID:gene fold change
53+
probability (pvalue)
54+
gene_diff = Dictionary, key:value = gene ID: binary up/down/constant
55+
regulation values
56+
gvdict = Dictionary, key:value = gene ID:defined variable ids from both
57+
sets (each key has 2 values)
58+
problem = Flux balance problem
59+
"""
60+
61+
@classmethod
62+
def init_parser(cls, parser):
63+
parser.add_argument(
64+
'--flux-threshold',
65+
help='Enter maximum objective flux as a decimal or percent',
66+
type=MaybeRelative, default=MaybeRelative('100%'))
67+
parser.add_argument(
68+
'--transc-file', help='Enter path to transcriptomic data file',
69+
metavar='FILE')
70+
parser.add_argument('--fva', help='Enable FVA', action='store_true')
71+
super(MadeFluxBalance, cls).init_parser(parser)
72+
73+
def run(self):
74+
"""Run MADE implementation."""
75+
gene_dict = self.get_gene_dict()
76+
77+
biomass_fun = self._model.biomass_reaction
78+
79+
# Create problem instance
80+
solver = self._get_solver(integer=True)
81+
prob = solver.create_problem()
82+
v_0 = prob.namespace()
83+
v_1 = prob.namespace()
84+
85+
# Define flux variables
86+
for reaction_id in self._mm.reactions:
87+
lower, upper = self._mm.limits[reaction_id]
88+
v_0.define([reaction_id], lower=lower, upper=upper)
89+
v_1.define([reaction_id], lower=lower, upper=upper)
90+
91+
# Create mass balance constraints for both conditions
92+
massbalance_0_lhs = {compound: 0 for compound in self._mm.compounds}
93+
massbalance_1_lhs = {compound: 0 for compound in self._mm.compounds}
94+
for (compound, reaction_id), value in iteritems(self._mm.matrix):
95+
massbalance_0_lhs[compound] += v_0(reaction_id) * value
96+
massbalance_1_lhs[compound] += v_1(reaction_id) * value
97+
for _, lhs in iteritems(massbalance_0_lhs):
98+
prob.add_linear_constraints(lhs == 0)
99+
for _, lhs in iteritems(massbalance_1_lhs):
100+
prob.add_linear_constraints(lhs == 0)
101+
102+
start_time = time.time()
103+
104+
# Set biomass flux threshold
105+
flux_threshold = self._args.flux_threshold
106+
if flux_threshold.relative:
107+
prob.set_objective(v_0(biomass_fun))
108+
result = prob.solve(lp.ObjectiveSense.Maximize)
109+
if not result:
110+
raise Exception('Failed to solve FBA')
111+
flux_threshold.reference = result.get_value(v_0(biomass_fun))
112+
113+
prob.add_linear_constraints(v_0(biomass_fun) >= float(flux_threshold))
114+
prob.add_linear_constraints(v_1(biomass_fun) >= float(flux_threshold))
115+
116+
gene_term_0 = prob.namespace()
117+
gene_term_1 = prob.namespace()
118+
119+
reaction_0 = prob.namespace(
120+
self._mm.reactions, types=lp.VariableType.Binary)
121+
reaction_1 = prob.namespace(
122+
self._mm.reactions, types=lp.VariableType.Binary)
123+
124+
for rxn_id, exp in sorted(iteritems(gene_dict)):
125+
create_gpr_constraints(
126+
prob, rxn_id, exp, reaction_0(rxn_id), gene_term_0)
127+
create_gpr_constraints(
128+
prob, rxn_id, exp, reaction_1(rxn_id), gene_term_1)
129+
130+
if self._args.transc_file is not None:
131+
con1, con2, gene_pval, gene_diff = idc(
132+
open_file(self._args.transc_file))
133+
134+
add_final_constraints(self._mm, prob, v_0, reaction_0)
135+
add_final_constraints(self._mm, prob, v_1, reaction_1)
136+
result = make_obj_fun(
137+
prob, gene_diff, gene_pval, gene_term_0, gene_term_1)
138+
139+
# Run FBA
140+
for reaction_id in sorted(self._mm.reactions):
141+
print('{}\t{}\t{}\t{}\t{}\t{}\t{}'.format(
142+
reaction_id, result.get_value(v_0(reaction_id)),
143+
result.get_value(v_1(reaction_id)),
144+
result.get_value(reaction_0(reaction_id)) > 0.5,
145+
result.get_value(reaction_1(reaction_id)) > 0.5,
146+
self._mm.get_reaction(reaction_id),
147+
gene_dict.get(reaction_id, '')))
148+
149+
logger.info('Solving took {:.2f} seconds'.format(
150+
time.time() - start_time))
151+
152+
def get_gene_dict(self):
153+
"""Using the reaction file called inside of the model file, it returns
154+
a dictionary with reaction IDs as keys and their associated
155+
gene-protein reaction (GPR) logic (i.e. (gene 1 and gene 2) or gene 3)
156+
as values of type str.
157+
"""
158+
gene_dict = {}
159+
for reaction in self._model.parse_reactions():
160+
if reaction.genes is not None:
161+
gene_dict[reaction.id] = boolean.Expression(reaction.genes)
162+
163+
return gene_dict
164+
165+
166+
def make_obj_fun(prob, gene_diff, gene_pval, gene_term_0, gene_term_1):
167+
"""Constructs the MADE objective funtion from dictionaries of LP variables.
168+
169+
Objective function consists of the summation of three functions dependent
170+
on the up/down regulation of gene expression between conditions. The
171+
functions contain a weighting function, and the difference between the
172+
binary representations of condition 1 and condition 2.
173+
"""
174+
i_vars = 0.0 # Increasing gene expression
175+
d_vars = 0.0 # Decreasing gene expression
176+
c_vars = 0.0 # Constant gene expression
177+
178+
def weight(p):
179+
return -math.log10(p)
180+
181+
for gene in gene_pval:
182+
# Comment by Julie/Matt?: Limitation of math.log()
183+
wp = max(2.2204460492e-16, gene_pval[gene])
184+
185+
x_0 = gene_term_0(boolean.Variable(gene))
186+
x_1 = gene_term_1(boolean.Variable(gene))
187+
# x_delta = x_0 XOR X_1
188+
prob.define(('xor', gene), types=lp.VariableType.Binary)
189+
x_delta = prob.var(('xor', gene))
190+
prob.add_linear_constraints(
191+
x_delta <= x_0 + x_1,
192+
x_delta >= x_0 - x_1, x_delta >= x_1 - x_0,
193+
x_delta <= 2 - x_0 - x_1)
194+
195+
if gene_diff[gene] == 1:
196+
i_vars += weight(wp) * (x_1 - x_0)
197+
elif gene_diff[gene] == -1:
198+
d_vars += weight(wp) * (x_0 - x_1)
199+
elif gene_diff[gene] == 0:
200+
c_vars += weight(wp) * x_delta
201+
202+
objective = i_vars + d_vars - c_vars
203+
204+
prob.set_objective(objective)
205+
result = prob.solve(lp.ObjectiveSense.Maximize)
206+
if not result:
207+
raise Exception('Unable to solve: ' + result.status)
208+
209+
obj_value = result.get_value(objective)
210+
logger.info('Objective: {}'.format(obj_value))
211+
212+
return result
213+
214+
215+
def create_gpr_constraints(prob, rxn_id, exp, reaction_var, gene_term):
216+
"""Opens all gene-logic containers, defines content, outputs the linear
217+
inequalities by calling bool_ineqs(). Sorts data into dictionaries that
218+
are used in other functions. Is recursive. No output.
219+
220+
Args:
221+
exp_obj: All of the expression objects (genes, AND, OR)
222+
var_gen: Counter used for relabeling the genes and arguments as
223+
variables
224+
new_var_id: Variable ID, also includes original reaction ID for first
225+
layer
226+
"""
227+
next_terms = iter([exp.root])
228+
current_type = None
229+
variable = reaction_var
230+
arguments = []
231+
stack = []
232+
while True:
233+
try:
234+
term = next(next_terms)
235+
except StopIteration:
236+
term = None
237+
238+
if gene_term.has_variable(term):
239+
arguments.append(gene_term(term))
240+
elif isinstance(term, boolean.Variable):
241+
gene_term.define([term], types=lp.VariableType.Binary)
242+
term_var = gene_term(term)
243+
arguments.append(term_var)
244+
elif isinstance(term, (boolean.And, boolean.Or)):
245+
stack.append((current_type, next_terms, variable, arguments))
246+
current_type = term.__class__
247+
next_terms = iter(term)
248+
arguments = []
249+
gene_term.define([term], types=lp.VariableType.Binary)
250+
variable = gene_term(term)
251+
else:
252+
# End of term
253+
if current_type is None:
254+
prob.add_linear_constraints(variable == arguments[0])
255+
break
256+
elif current_type == boolean.And:
257+
prob.add_linear_constraints(
258+
*and_constraints(variable, arguments))
259+
elif current_type == boolean.Or:
260+
prob.add_linear_constraints(
261+
*or_constraints(variable, arguments))
262+
263+
term_var = variable
264+
current_type, next_terms, variable, arguments = stack.pop()
265+
arguments.append(term_var)
266+
267+
268+
def and_constraints(var, arguments):
269+
"""Create constraints for boolean AND.
270+
271+
Creates constraints for: var = And(arguments) where var and each argument
272+
is a binary variable.
273+
"""
274+
var_sum = 0
275+
for arg in arguments:
276+
yield var <= arg
277+
var_sum += arg
278+
yield var >= var_sum - (len(arguments) - 1)
279+
280+
281+
def or_constraints(var, arguments):
282+
"""Create constraints for boolean OR.
283+
284+
Creates constraints for: var = Or(arguments) where var and each argument
285+
is a binary variable.
286+
"""
287+
var_sum = 0
288+
for arg in arguments:
289+
yield var >= arg
290+
var_sum += arg
291+
yield var <= var_sum
292+
293+
294+
def open_file(path):
295+
"""Returns the contents of model file in a tuple of dictionaries.
296+
File Form: tsv format, FOUR Columns: (1) Gene name, (2) Condition 1 Data,
297+
(3) Condition 2 Data, (4) P-value of the fold change for transition 1->2.
298+
"""
299+
file1 = open(path)
300+
con1_dict = {}
301+
con2_dict = {}
302+
pval_dict = {}
303+
304+
file1.readline()
305+
for row in csv.reader(file1, delimiter=str('\t')):
306+
con1_dict[row[0]] = float(row[1])
307+
con2_dict[row[0]] = float(row[2])
308+
if float(row[3]) == float(0.0):
309+
pval_dict[row[0]] = 1e-400
310+
else:
311+
pval_dict[row[0]] = float(row[3])
312+
313+
return con1_dict, con2_dict, pval_dict
314+
315+
316+
def idc(dicts):
317+
"""Used for accessing the list of dictionaries created in open_file()
318+
Creates a dictionary for the gene ID and a value = [-1, 0, +1]
319+
corresponding to decreasing, constant, and inreasing expression between the
320+
conditions.
321+
"""
322+
con1 = dicts[0]
323+
con2 = dicts[1]
324+
pval = dicts[2]
325+
diff = {}
326+
327+
for key in con1:
328+
if con2[key] == con1[key]:
329+
diff[key] = 0
330+
elif con2[key] > con1[key]:
331+
diff[key] = 1
332+
else:
333+
diff[key] = -1
334+
335+
return con1, con2, pval, diff
336+
337+
338+
def add_final_constraints(mm, prob, v, z):
339+
"""Takes the metabolic model, the LP Problem, and the binary
340+
dictionaries of each condition. Adds constraints connecting flux
341+
variables, reactions, and their flux bounds.
342+
"""
343+
for rxn in mm.reactions:
344+
vmin = mm.limits[rxn].lower
345+
vmax = mm.limits[rxn].upper
346+
flux_var = v(rxn)
347+
active = z(rxn)
348+
349+
prob.add_linear_constraints(
350+
active*vmax >= flux_var, flux_var >= active*vmin)

setup.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -63,6 +63,7 @@
6363
'gapcheck = psamm.commands.gapcheck:GapCheckCommand',
6464
'gapfill = psamm.commands.gapfill:GapFillCommand',
6565
'genedelete = psamm.commands.genedelete:GeneDeletionCommand',
66+
'made = psamm.commands.made:MadeFluxBalance',
6667
'masscheck = psamm.commands.masscheck:MassConsistencyCommand',
6768
'randomsparse = psamm.commands.randomsparse:RandomSparseNetworkCommand',
6869
'robustness = psamm.commands.robustness:RobustnessCommand',

0 commit comments

Comments
 (0)