Skip to content

Commit

Permalink
feat: assert solver status optimal
Browse files Browse the repository at this point in the history
The pattern of checking if solver.status != 'optimal' is not and then
raise appropriate exception started to appear in a number of places -
let's have a function for that.

Also tidy up the use of optlang string constants.
  • Loading branch information
hredestig committed Mar 31, 2017
1 parent 7c59314 commit 4863f78
Show file tree
Hide file tree
Showing 7 changed files with 46 additions and 29 deletions.
3 changes: 2 additions & 1 deletion cobra/core/solution.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,6 +9,7 @@
from warnings import warn

from numpy import zeros, asarray, nan
from optlang.interface import OPTIMAL
from pandas import Series

from cobra.util.solver import check_solver_status
Expand Down Expand Up @@ -96,7 +97,7 @@ def __init__(self, objective_value, status, reactions, fluxes,

def __repr__(self):
"""String representation of the solution instance."""
if self.status != "optimal":
if self.status != OPTIMAL:
return "<Solution {0:s} at 0x{1:x}>".format(self.status, id(self))
return "<Solution {0:.3f} at 0x{1:x}>".format(self.objective_value,
id(self))
Expand Down
19 changes: 7 additions & 12 deletions cobra/exceptions.py
Original file line number Diff line number Diff line change
Expand Up @@ -5,39 +5,34 @@
import optlang.interface


class SolveError(Exception):
class OptimizationError(Exception):
def __init__(self, message):
super(SolveError, self).__init__(message)
super(OptimizationError, self).__init__(message)


class Infeasible(SolveError):
class Infeasible(OptimizationError):
pass


class Unbounded(SolveError):
class Unbounded(OptimizationError):
pass


class FeasibleButNotOptimal(SolveError):
class FeasibleButNotOptimal(OptimizationError):
pass


class UndefinedSolution(SolveError):
class UndefinedSolution(OptimizationError):
pass


_OPTLANG_TO_EXCEPTIONS_DICT = dict((
OPTLANG_TO_EXCEPTIONS_DICT = dict((
(optlang.interface.INFEASIBLE, Infeasible),
(optlang.interface.UNBOUNDED, Unbounded),
(optlang.interface.FEASIBLE, FeasibleButNotOptimal),
(optlang.interface.UNDEFINED, UndefinedSolution)))


class OptimizationError(Exception):
def __init__(self, message):
super(OptimizationError, self).__init__(message)


class DefunctError(Exception):
"""Exception for retired functionality
Expand Down
3 changes: 2 additions & 1 deletion cobra/flux_analysis/phenotype_phase_plane.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@
from itertools import product
from multiprocessing import Pool
import pandas as pd
from optlang.interface import OPTIMAL

from numpy import (
nan, abs, arange, dtype, empty, int32, linspace, meshgrid, unravel_index,
Expand Down Expand Up @@ -409,7 +410,7 @@ def _envelope_for_points(model, reactions, grid, carbon_io):
for reaction, coordinate in zip(reactions, point):
reaction.bounds = (coordinate, coordinate)
model.solver.optimize()
if model.solver.status == 'optimal':
if model.solver.status == OPTIMAL:
for reaction, coordinate in zip(reactions, point):
results[reaction.id].append(coordinate)
results['direction'].append(direction)
Expand Down
7 changes: 4 additions & 3 deletions cobra/flux_analysis/single_deletion.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,7 @@

import pandas
from six import iteritems, string_types
from optlang.interface import OPTIMAL

import cobra.solvers as legacy_solvers
import cobra.util.solver as solvers
Expand Down Expand Up @@ -110,7 +111,7 @@ def single_reaction_deletion_fba(cobra_model, reaction_list, solver=None,
status = m.solver.status
status_dict[reaction.id] = status
growth_rate_dict[reaction.id] = m.solver.objective.value \
if status == "optimal" else 0.
if status == OPTIMAL else 0.
else:
# This entire block can be removed once the legacy solvers are
# deprecated
Expand Down Expand Up @@ -180,7 +181,7 @@ def single_reaction_deletion_moma(cobra_model, reaction_list, solver=None,
m.solver.optimize()
status = m.solver.status
status_dict[reaction.id] = status
if status == "optimal":
if status == OPTIMAL:
growth = m.variables.moma_old_objective.primal
else:
growth = 0.0
Expand Down Expand Up @@ -285,7 +286,7 @@ def single_gene_deletion_fba(cobra_model, gene_list, solver=None,
status = m.solver.status
status_dict[gene.id] = status
growth_rate_dict[gene.id] = m.solver.objective.value if \
status == "optimal" else 0.
status == OPTIMAL else 0.
else:
for gene in gene_list:
old_bounds = {}
Expand Down
5 changes: 2 additions & 3 deletions cobra/flux_analysis/variability.py
Original file line number Diff line number Diff line change
Expand Up @@ -158,9 +158,8 @@ def _fva_optlang(model, reaction_list, fraction, loopless):
prob = model.problem
with model as m:
m.solver.optimize()
if m.solver.status != "optimal":
raise ValueError("There is no optimal solution "
"for the chosen objective!")
sutil.assert_optimal(m, message="There is no optimal solution for "
"the chosen objective!")
# Add objective as a variable to the model than set to zero
# This also uses the fraction to create the lower bound for the
# old objective
Expand Down
10 changes: 6 additions & 4 deletions cobra/test/test_flux_analysis.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,6 +9,7 @@

import pytest
import numpy
from optlang.interface import OPTIMAL, INFEASIBLE
from six import StringIO, iteritems

import cobra.util.solver as sutil
Expand All @@ -18,6 +19,7 @@
from cobra.flux_analysis.sampling import ARCHSampler, OptGPSampler
from cobra.manipulation import convert_to_irreversible
from cobra.solvers import SolverNotFound, get_solver_name, solver_dict
from cobra.exceptions import Infeasible

try:
import scipy
Expand Down Expand Up @@ -144,7 +146,7 @@ def test_pfba(self, model, solver):
model.reactions.ATPM.lower_bound = 500
with warnings.catch_warnings():
warnings.simplefilter("error", UserWarning)
with pytest.raises((UserWarning, ValueError)):
with pytest.raises((UserWarning, Infeasible, ValueError)):
pfba(model, solver=solver)

@pytest.mark.parametrize("solver", all_solvers)
Expand Down Expand Up @@ -307,7 +309,7 @@ def test_fva_infeasible(self, model):
infeasible_model = model.copy()
infeasible_model.reactions.get_by_id("EX_glc__D_e").lower_bound = 0
# ensure that an infeasible model does not run FVA
with pytest.raises(ValueError):
with pytest.raises(Infeasible):
flux_variability_analysis(infeasible_model)

@pytest.mark.parametrize("solver", all_solvers)
Expand Down Expand Up @@ -385,8 +387,8 @@ def test_add_loopless(self, ll_test_model):
ll_test_model.reactions.v3.lower_bound = 1
ll_test_model.solver.optimize()
infeasible_status = ll_test_model.solver.status
assert feasible_status == "optimal"
assert infeasible_status == "infeasible"
assert feasible_status == OPTIMAL
assert infeasible_status == INFEASIBLE

def test_phenotype_phase_plane_benchmark(self, model, benchmark):
benchmark(calculate_phenotype_phase_plane,
Expand Down
28 changes: 23 additions & 5 deletions cobra/util/solver.py
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@
import optlang
import sympy

from cobra.exceptions import OptimizationError
from cobra.exceptions import OptimizationError, OPTLANG_TO_EXCEPTIONS_DICT
from cobra.util.context import get_context


Expand Down Expand Up @@ -378,8 +378,9 @@ def fix_objective_as_constraint(model, fraction=1, name='fixed_objective_{}'):
fix_objective_name = name.format(model.objective.name)
if fix_objective_name in model.constraints:
model.solver.remove(fix_objective_name)
solution = model.optimize()
objective_bound = solution.objective_value * fraction
model.solver.optimize()
assert_optimal(model)
objective_bound = model.solver.objective.value * fraction
if model.objective.direction == 'max':
ub, lb = None, objective_bound
else:
Expand All @@ -392,9 +393,9 @@ def fix_objective_as_constraint(model, fraction=1, name='fixed_objective_{}'):

def check_solver_status(status):
"""Perform standard checks on a solver's status."""
if status == "optimal":
if status == optlang.interface.OPTIMAL:
return
elif status == "infeasible":
elif status == optlang.interface.INFEASIBLE:
warn("solver status is '{}'".format(status), UserWarning)
elif status is None:
raise RuntimeError(
Expand All @@ -403,4 +404,21 @@ def check_solver_status(status):
raise OptimizationError("solver status is '{}'".format(status))


def assert_optimal(model, message='optimization failed'):
"""Assert model solver status is optimal.
Do nothing if model solver status is optimal, otherwise throw
appropriate exception depending on the status.
Parameters
----------
model : cobra.Model
The model to check the solver status for.
message : str (optional)
Message to for the exception if solver status was not optimal.
"""
if model.solver.status != optlang.interface.OPTIMAL:
raise OPTLANG_TO_EXCEPTIONS_DICT[model.solver.status](message)


import cobra.solvers as legacy_solvers # noqa

0 comments on commit 4863f78

Please sign in to comment.