Solvers

CPLEXParamSet(interface)[source]

This is a function which returns user specified CPLEX control parameters. It is not necessary to use a file like this if you want to use CPLEX default control parameters. It is intended to be a template for individual users to save with their own problem specific settings for CPLEX.

Usage

cpxControl = CPLEXParamSet(solver)

Input

  • interface – ‘tomlab_cplex’(default) or ‘ILOGcomplex’

Output

  • cpxControl – user specified CPLEX control parameters

Example

% (1) Paddy saves this file as CPLEXParamSetPaddyLPJob1
% (2) Paddy edits CPLEXParamSetPaddyLPJob1 in a problem specific way
% (3) Paddy then passes the name of this file to solveCobraLP_CPLEX using something like:

[solution, LPProblem] = solveCobraLP_CPLEX(LPProblem, [], [], [], 'CPLEXParamSetPaddyLPJob1');

Note

CPLEX consists of 4 different LP solvers which can be used to solve LP problems you can control which of the solvers, e.g. simplex or interior point solve using the CPLEX control parameter cpxControl.LPMETHOD

NLPobjPerFlux(fluxVector, Prob)[source]

Calculates the value of the objective - (Prob.osense * Prob.user.model.c)/sum(v.^2) based on a flux distribution This function is meant to be used with NLP solvers

Usage

value = NLPobjPerFlux(fluxVector, Prob)

Inputs

  • fluxVector – Flux vector
  • Prob – NLP problem structure

Output

  • value – -Objective flux / sum(v.^2)
buildCplexProblemFromCOBRAStruct(Problem)[source]

Build a cplex object from the given LP problem in COBRA Format Usage

cplexProblem = buildCplexProblemFromCOBRAStruct(LPproblem)

Input

  • LPproblem – A COBRA style Problem with the following fields:

    • .A - The equality and in equality matrix
    • .b - The right hand side values of the constraints
    • .ub - the upper bounds of the variables
    • .lb - the lower bounds of the variables
    • .osense - The objective sense (-1 for max, 1 for min)
    • .c - the objective coefficient vector for the linear part
    OPTIONAL:
    • .F - The objective coefficient matrix for the quadratic part.
    • .varType - The variable types for mixed integer problems (‘I’, integer, ‘C’, continous,’B’ binary)
    • .b_L - left hand sides of the constraints (will only be used if csense is empty)
    • .csense - The constraint senses, Default assumption is all ‘E’
    • .x0 - Basis to use
buildLPproblemFromModel(model, checked)[source]

Builds an COBRA Toolbox LP problem structure from a COBRA Toolbox model structure.

Usage

LPproblem = buildLPproblemFromModel(model)

Input

  • model – A COBRA model structure with at least the following fields
    • .S - The stoichiometric matrix
    • .c - Objective coeff vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector

Optional inputs

  • model – The model structure can also have these additional fields:
    • .b: accumulation/depletion vector (default 0 for each metabolite).
    • .osense: Objective sense (-1 means maximise (default), 1 means minimise)
    • .csense: a string with the constraint sense for each row in A (‘E’, equality(default), ‘G’ greater than, ‘L’ less than).
    • .C: the Constraint matrix;
    • .d: the right hand side vector for C;
    • .dsense: the constraint sense vector;
    • .E: the additional Variable Matrix
    • .evarub: the upper bounds of the variables from E;
    • .evarlb: the lower bounds of the variables from E;
    • .evarc: the objective coefficients of the variables from E;
    • .D: The matrix coupling additional Constraints (form C), with additional Variables (from E);
  • checked – Check the input (default: true);

Output

  • LPproblem – A COBRA LPproblem structure with the following fields:
    • .A: LHS matrix
    • .b: RHS vector
    • .c: Objective coeff vector
    • .lb: Lower bound vector
    • .ub: Upper bound vector
    • .osense: Objective sense (-1: maximise (default); 1: minimise)
    • .csense: string with the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).
changeCobraSolver(solverName, solverType, printLevel, validationLevel)[source]

Changes the Cobra Toolbox optimization solver(s)

Usage

solverOK = changeCobraSolver(solverName, solverType, printLevel, validationLevel)

Inputs

  • solverName – Solver name
  • solverType – Solver type, LP, MILP, QP, MIQP (opt, default LP, all). ‘all’ attempts to change all applicable solvers to solverName. This is purely a shorthand convenience.
  • printLevel – verbose level
    • if 0, warnings and errors are silenced
    • if > 0, warnings and errors are thrown. (default: 1)

Optional input

  • validationLevel – how much validation to use
    • -1: assign only the global variable. Do not assign any path.
    • 0: adjust solver paths but don’t validate the solver (default)
    • 1: validate but remove outputs
    • 2: validate and keep any outputs

Outputs

  • solverOKtrue if solver can be accessed, false if not
  • solverInstalledtrue if the solver is installed (not necessarily working)

Currently allowed LP solvers:

  • fully supported solvers

    cplex_direct CPLEX accessed directly through Tomlab cplex.m. This gives the user more control of solver parameters. e.g. minimising the Euclidean norm of the internal flux to get rid of net flux around loops
    dqqMinos DQQ solver
    glpk GLPK solver with Matlab mex interface (glpkmex)
    gurobi Gurobi solver
    ibm_cplex The IBM API for CPLEX using the CPLEX class
    matlab MATLAB’s linprog function
    mosek Mosek LP solver with Matlab API (using linprog.m from Mosek)
    pdco PDCO solver
    quadMinos quad solver
    tomlab_cplex CPLEX accessed through Tomlab environment (default)
  • legacy solvers:

    lindo_new Lindo API > v2.0
    lindo_legacy Lindo API < v2.0
    lp_solve lp_solve with Matlab API
    gurobi_mex Gurobi accessed through Matlab mex interface (Gurobi mex)
    opti CLP(recommended), CSDP, DSDP, OOQP and SCIP(recommended) solver installed and called with OPTI TB wrapper Lower level calls with installed mex files are possible but best avoided for all solvers

Currently allowed MILP solvers:

  • fully supported solvers:

    cplex_direct CPLEX accessed directly through Tomlab cplex.m. This gives the user more control of solver parameters. e.g. minimising the Euclidean norm of the internal flux to get rid of net flux around loops
    glpk glpk MILP solver with Matlab mex interface (glpkmex)
    gurobi Gurobi solver
    ibm_cplex The IBM API for CPLEX using the CPLEX class
    mosek Mosek LP solver with Matlab API (using linprog.m from Mosek)
    tomlab_cplex CPLEX MILP solver accessed through Tomlab environment
  • legacy solvers:

    gurobi_mex Gurobi accessed through Matlab mex interface (Gurobi mex)
    opti CLP(recommended), CSDP, DSDP, OOQP and SCIP(recommended) solver installed and called with OPTI TB wrapper Lower level calls with installed mex files are possible but best avoided for all solvers

Currently allowed QP solvers:

  • fully supported solvers:

    cplex_direct CPLEX accessed directly through Tomlab cplex.m. This gives the user more control of solver parameters. e.g. minimising the Euclidean norm of the internal flux to get rid of net flux around loops
    gurobi Gurobi solver
    ibm_cplex The IBM API for CPLEX using the CPLEX class
    mosek Mosek QP solver with Matlab API
    pdco PDCO solver
    tomlab_cplex CPLEX QP solver accessed through Tomlab environment
  • experimental support:

    qpng qpng QP solver with Matlab mex interface (in glpkmex package, only limited support for small problems)
  • legacy solvers:

    gurobi_mex Gurobi accessed through Matlab mex interface (Gurobi mex)
    opti CLP(recommended), CSDP, DSDP, OOQP and SCIP(recommended) solver installed and called with OPTI TB wrapper. Lower level calls with installed mex files are possible

Currently allowed MIQP solvers:

  • fully supported solvers:

    cplex_direct CPLEX accessed directly through Tomlab cplex.m. This gives the user more control of solver parameters. e.g. minimising the Euclidean norm of the internal flux to get rid of net flux around loops
    gurobi Gurobi solver
    ibm_cplex The IBM API for CPLEX using the CPLEX class
    tomlab_cplex CPLEX MIQP solver accessed through Tomlab environment
  • legacy solvers:

    gurobi_mex Gurobi accessed through Matlab mex interface (Gurobi mex)

Currently allowed NLP solvers:

  • fully supported solvers:

    matlab MATLAB’s fmincon.m
    quadMinos quad solver
  • experimental support:

    tomlab_snopt SNOPT solver accessed through Tomlab environment

Note

It is a good idea to put this function call into your startup.m file (usually matlabinstall/toolboxes/local/startup.m)

changeCobraSolverParams(solverType, paramName, paramValue)[source]

Changes parameters for the Cobra Toolbox optimization solver(s)

Usage

changeOK = changeCobraSolverParams(solverType, paramName, paramValue)

Inputs

  • solverType – Solver type, ‘LP’ or ‘MILP’ (opt, default, ‘LP’)
  • paramName – Parameter name
  • paramValue – Parameter value

Output

  • changeOK – Logical inicator that supplied parameter is allowed (= 1)

Explanation on parameters:

  • printLevel: Printing level
    • 0 - Silent
    • 1 - Warnings and Errors
    • 2 - Summary information (Default)
    • 3 - More detailed information
    • > 10 - Pause statements, and maximal printing (debug mode)
  • primalOnly: {(0), 1}; 1 = only return the primal vector (lindo solvers)
  • saveInput: Saves LPproblem to filename specified in field.
    i.e. parameters.saveInput = ‘LPproblem.mat’;
  • minNorm: {(0), scalar , n x 1 vector}, where [m, n] = size(S);
    If not zero then, minimise the Euclidean length of the solution to the LP problem. minNorm ~1e-6 should be high enough for regularisation yet maintain the same value for the linear part of the objective. However, this should be checked on a case by case basis, by optimization with and without regularisation.
  • optTol Optimality tolerance
  • feasTol Feasibility tolerance
  • timeLimit: Global solver time limit
  • intTol: Integrality tolerance
  • relMipGapTol: Relative MIP gap tolerance
  • logFile: Log file (for CPLEX)

Note

The available solver Parameters can be obtained by calling getSolverParamsOptionsForType(). If input argument minNorm is not zero, then minimise the Euclidean length of the solution to the LP problem. minNorm ~1e-6 should be high enough for regularisation yet maintain the same value for the linear part of the objective. However, this should be checked on a case by case basis, by optimization with and without regularisation.

checkGAMSSolvers(problemType)[source]

This function return the solvers that can be used in GAMS to solve the type of problem especified by the user

Usage

solvers = checkGAMSSolvers(problemType)

Input

  • problemType Type – string Description: string containing the problem type for which this function will search solvers. E.g.: problem type = ‘LP’ (Linear Programming)

Output

  • solvers Type – cell array for available GAMS solvers in your systems which allows the user to solve problems of type “problemType”

Example

solvers = checkGAMSSolvers('MIP')
% returns the GAMS solvers available to solve Mixed Integer Programming
% problems. You can see the entire list of problem types with the
% function getAvailableGAMSSolvers.m
checkSolFeas(LP, sol, maxInfeas, tol, internal)[source]

Returns the infeasibility of solutions given a COBRA model or LP structure, or a IBM-ILOG CPLEX class

Usage

[infeas, sol] = checkSolFeas(LP, sol, maxInfeas, tol)

Inputs

  • LP – COBRA model or LP structure, or a IBM-ILOG CPLEX class
  • sol – solution structure or columns of solution vectors. If LP is a CPLEX class with. Solution property sol can be omitted or empty.

Optional inputs

  • maxInfeas – if true (defaulted), infeas = maximum infeasiblity if false, infeas = struct of vectors of infeasibility with the following fields:

    • .con for infeasibility of constraints
    • .lb for infeasibility of lower bounds
    • .ub for infeasibility of upper bounds
    • .ind for infeasibility of indicator constraints. = -Inf if an indicator is not active. (CPLEX class only)
  • tol – feasibility tolerance (defaulted at the Cobra solver feasTol value). For determining if the input solution is indeed feasible. Used only if the input solution is a structure.

Outputs

  • infeas – maximum infeasibility (maxInfeas = true) or struct of vectors of infeasibility (maxInfeas = false)
  • sol – solution structure. Only available if the input solution is a structure. If \(infeas \leq tol\), sol.stat = 1. Otherwise no change.
configEnvVars(printLevel)[source]

Configures the global variables based on the system’s configuration First, all environment variables for each solver are defined together with all eventual solver paths. Then, there will be 4 methods marked that can be used to define the global variables:

  • 1: solver is on the path and at a standard location (*—)
  • 2: solver is on path but at a non-standard location (-*–)
  • 3: solver path is defined through environment variables (–*-)
  • 4: solver is not already on the path and the environment variable is not set, but the standard directory exists (—*)

If none of these 4 methods applies, the global solver path variable is not set and an appropriate message is returned

Usage

configEnvVars(printLevel)

Input

  • printLevel – default = 0, verbose level
getAvailableGAMSSolvers()[source]

This function return the GAMS solvers which are available in your system according to the license.

Usage

[summaryTable, booleanTable, problemTypes, solvers] = getAvailableGAMSSolvers

Outputs

  • summaryTable Type – cell array Description: matrix summarizing which problem can be solved by which solver. A “yes” in position (i,j) means that the GAMS solver “j” is available in your system to solve a problem of type “i”. Please note that headers are included in this summaryTable, so the first column and the first row have labels for solvers and problem types, respectively.
  • booleanTable Type – double matrix Description: matrix summarizing which problem can be solved by which solver. A “1” in position (i,j) means that the GAMS solver “j” is available in your system to solve a problem of type “i”
  • problemTypes Type – cell array of strings Description: list of problem types that can be solved in GAMS. E.g: LP, MIP, etc
  • solvers Type – cell array of strings Description: list of solvers available in GAMS
getAvailableSolversByType()[source]

Get the available solvers for the different solver types on the system.

Usage

solvers = getAvailableSolversByType()

Output

  • solvers – struct containing one field per Problem type listing all solvers installed on the system for that problem type. Also contains an ALL field indicating all available solvers
getCobraSolverParams(solverType, paramNames, parameters)[source]

This function gets the specified parameters in paramNames from parameters, the global cobra paramters variable or default values set within this script. It will use values with the following priority

parameters > global parameters > default

The specified parameters will be delt to the specified output arguements. See examples below.

Usage

varargout = getCobraSolverParams(solverType, paramNames, parameters)

Inputs

  • solverType – Type of solver used: ‘LP’, ‘MILP’, ‘QP’, ‘MIQP’
  • paramNames – Cell array of strings containing parameter names OR one parameter name as string

Optional input

  • parameters – Structure with fields pertaining to parameter values that should be used in place of global or default parameters. parameters can be set to ‘default’ to use the default values set within this script.

Output

  • varargout – Variables which each value corresponding to paramNames is outputted to.

Example

parameters.saveInput = 'LPproblem.mat';
parameters.printLevel = 1;
[printLevel, saveInput] = getCobraSolverParams('LP', {'printLevel', 'saveInput'}, parameters);

%Example using default values
[printLevel, saveInput] = getCobraSolverParams('LP', {'printLevel','saveInput'}, 'default');
getCobraSolverParamsOptionsForType(solverType)[source]

This function returns the available optional parameters for the specified solver type.

Usage

paramnames = getCobraSolverParamsOptionsForType(solverType)

Input

  • solverType – One of the solver types available in the cobra Toolbox (‘LP’,’QP’,’MILP’,’MIQP’,’NLP’)
OUPTUT:
paramNames: The possible parameters that can be set for the
given solver Type (depends on the solver Type
getCobraSolverVersion(solverName, printLevel, rootPathSolver)[source]

detects the version of given COBRA solver

Usage

solverVersion = getCobraSolverVersion(solverName, printLevel, rootPathSolver)

Inputs

  • solverName – Name of the solver
  • printLevel – verbose level (default: 0)
  • rootPathSolver – Path to the solver installation

Output

  • solverVersion – string that contains the version number
getParamList(param, bottomFlag)[source]

Get all the end parameters and their paths for matching CPLEX parameters appropriately (e.g., if param.simplex.display is a parameter, then we will have ‘display’ in paramList and ‘param.simplex’ in paramPath)

Usage

[paramList, paramPath, addCur] = getParamList(param, bottomFlag)

Inputs

  • param – existing structure with parameters for Cplex
  • bottomFlag – boolean switch to extract parameters based on ‘Cur’ at the bottom level in the user-supplied parameter structure

Outputs

  • paramList – Cell with proper parameter names
  • paramPath – Cell with superseeding parameter structure
  • addCur – Structure with booleans whether or not .Cur has been added to the paramPath for user-supplied parameters
isCompatible(solverName, printLevel, specificSolverVersion)[source]

determine the compatibility status of a solver based on the file compatMatrix.rst

Usage

compatibleStatus = isCompatible(solverName, printLevel, specificSolverVersion)

Inputs

  • solverName – Name of the solver
  • printLevel – verbose level (default: 0)
  • specificSolverVersion – string with specific solver version (example: ‘12.7.1’ or ‘6.5.1’)

Output

  • compatibleStatus – compatibility status
    • 0: not compatible with the COBRA Toolbox (tested)
    • 1: compatible with the COBRA Toolbox (tested)
    • 2: unverified compatibility with the COBRA Toolbox (not tested)
liftModel(model, BIG, printLevel, fileName, directory)[source]

Lifts a COBRA model with badly-scaled stoichiometric and coupling constraints of the form: \(max c*v\) subject to: \(Sv = 0, x, Cv <= 0\) Converts it into a COBRA LPproblem structure, which can be used with solveCobraLP. Fluxes for the reactions should stay the same i.e. sol.full(1:nRxns) should yield an optimal flux vector.

Usage

LPproblem = liftModel(model, BIG, printLevel,fileName,directory)

Input

  • model – COBRA LPproblem Structure containing the original LP to be solved. The format of this struct is described in the documentation for solveCobraLP.m

Optional inputs

  • BIG – A parameter the controls the largest entries that appear in the reformulated problem (default = 1000).
  • printLevel – printLevel = 1 enables printing of problem statistics (default); printLevel = 0 silent
  • fileName – name of th file to load
  • directory – file directory (if model is empty, you can load it using fileName and directory)

Output

  • LPproblem – COBRA Structure contain the reformulated LP to be solved.
lp_solve(f, a, b, e, vlb, vub, xint, scalemode, keep)[source]

Solves mixed integer linear programming problems. Solves the MILP problem \(max v = f'*x\) \(a*x <> b\) \(vlb <= x <= vub\) x(int)` are integer

Usage

[obj, x, duals, stat] = lp_solve(f, a, b, e, vlb, vub, xint, scalemode, keep)

Inputs

  • fn vector of coefficients for a linear objective function.
  • am by n matrix representing linear constraints.
  • bm vector of right sides for the inequality constraints.
  • em vector that determines the sense of the inequalities:
    • e(i) = -1 ==> Less Than
    • e(i) = 0 ==> Equals
    • e(i) = 1 ==> Greater Than

Optional inputs

  • vlbn vector of lower bounds. If empty or omitted, then the lower bounds are set to zero.
  • vubn vector of upper bounds. May be omitted or empty.
  • xint – vector of integer variables. May be omitted or empty.
  • scalemode – scale flag. Off when 0 or omitted.
  • keep – Flag for keeping the lp problem after it’s been solved. If omitted, the lp will be deleted when solved.

Outputs

  • obj – Optimal value of the objective function.
  • x – Optimal value of the decision variables.
  • duals – solution of the dual problem.
  • stat – result of mxlpsolve function
optimizeCbModelNLP(model, varargin)[source]

Optimizes constraint-based model using a non-linear objective

Usage

[currentSol, allObjValues, allSolutions] = optimizeCbModelNLP(model, varargin)

Input

  • model – COBRA model structure

Optional inputs are parameter/value pairs

Optional inputs

  • objFunction – Name of the non-linear matlab function to be optimized (the corresponding m-file must be in the current matlab path) The function receives two arguments, the current flux vector, and the NLPProblem structure.
  • initFunction – Name of the matlab function used to generate random initial starting points. The function will be supplied with two arguments: the model and a cell array of input arguments (specified in the initArgs parameter)
  • osenseStr – Optimization direction (‘max’ or ‘min’), this will override any mention in the model.
  • nOpt – Number of independent optimization runs performed
  • objArgs – Cell array of arguments that are supplied to the objective function as objArguments in the NLPProblem structure (i.e. the second element, will have a field objArguments.)
  • initArgs – Cell array of arguments to the ‘initFunction’, will be provided as second input Argument to the initFunction
  • solveroptions – A Struct with options for the solver used. This is specific to the solver in question, but the fields should relate to options accepted by the solver.

Outputs

  • currentSol – Solution structure
  • allObjValues – Array of objective value of each iteration
  • allSolutions – Array of flux distribution of each iteration
optimizeTwoCbModels(model1, model2, osenseStr, minFluxFlag, verbFlag)[source]

Simultaneously solve two flux balance problems and minimize the difference between the two solutions

Usage

[solution1, solution2, totalFluxDiff] = optimizeTwoCbModels(model1, model2, osenseStr, minFluxFlag, verbFlag)

Inputs

  • model1 – The first COBRA model
  • model2 – The second COBRA model, where both models have mandatory fields:
    • S - Stoichiometric matrix
    • b - Right hand side = 0
    • c - Objective coefficients
    • lb - Lower bounds
    • ub - Upper bounds

Optional inputs

  • osenseStr – Maximize (‘max’)/minimize (‘min’) (Default = ‘max’)
  • minFluxFlag – Minimize the absolute value of fluxes in the optimal MOMA solution (Default = false)
  • verbFlag – Verbose output (Default = false)

Outputs

  • solution1 – Solution for the 1st model
  • solution2 – Solution for the 2nd model
  • totalFluxDiff – 1-norm of the difference between the flux vectors sum|v1-v2|

Example

solution
  f         Objective value
  x         Primal (flux vector)


First solves two separate FBA problems:
                             f1 = max/min c1'v1
                             subject to S1*v1 = b1
                                        lb1 <= v1 <= ub1
                             f2 = max/min c2'v2
                             subject to S2*v2 = b2
                                        lb2 <= v2 <= ub2

Then solves the following LP to obtain the two flux vectors with the
smallest possible 1-norm difference between them

                             min |v1-v2|
                               s.t. S1*v1 = b1
                                    c1'v1 = f1
                                    lb1 <= v1 <= ub1
                                    S2*v2 = b2
                                    c2'v2 = f2
                                    lb2 <= v2 <= ub2

Finally optionally minimizes the 1-norm of the flux vectors
parseSolverParameters(problemType, varargin)[source]

Parse the solver parameters for a specified problem

Usage

[cobraParams, solverParams] = parseSolverParameters(problemType,varargin)

Input

  • problemType – The type of the problem to get parameters for (‘LP’,’MILP’,’QP’,’MIQP’,’NLP’)

Optional input

  • varargin – Additional parameters either as parameter struct, or as parameter/value pairs. A combination is possible, if the parameter struct is either at the beginning or the end of the optional input. All fields of the struct which are not COBRA parameters (see getCobraSolverParamsOptionsForType) for this problem type will be passed on to the solver in a solver specific manner.

Outputs

  • cobraParams – The COBRA Toolbox specific parameters for this problem type given the provided parameters
  • solverParams – Additional parameters provided which are not part of the COBRA parameters and are assumed to be part of direct solver input structs.
reformulate(LPproblem, BIG, printLevel)[source]

Reformulates badly-scaled FBA program Transforms LPproblems with badly-scaled stoichiometric and coupling constraints of the form: \(max c*x\) subject to: math:Ax <= b

Eliminates the need for scaling and hence prevents infeasibilities after unscaling. After using PREFBA to transform a badly-scaled FBA program, please turn off scaling and reduce the aggressiveness of presolve.

Rransforms a badly-scaled LPproblem contained in the struct FBA and returns the transformed program in the structure FBA. reformulate assumes S and C do not contain very small entries and transforms constraints containing very large entries (entries larger than BIG). BIG should be set between 1000 and 10000 on double precision machines. printLevel = 1 or 0 enables/diables printing respectively.

Reformulation techniques are described in detail in: Y. Sun, R. M.T. Fleming, M. A. Saunders, I. Thiele, An Algorithm for Flux Balance Analysis of Multi-scale Biochemical Networks, submitted.

Usage

[LPproblem] = reformulate(LPproblem, BIG, printLevel)

Inputs

  • LPproblem – Structure contain the original LP to be solved. The format of this struct is described in the documentation for solveCobraLP.m
  • BIG – A parameter the controls the largest entries that appear in the reformulated problem.
  • printLevel – 1 enables printing of problem statistics; 0 = silent

Output

  • LPproblem – Structure contain the reformulated LP to be solved.
setCplexParam(LP, solverParams, verbFlag)[source]

Sets the parameters of the IBM ILOG CPLEX object according to the structure solverParams The solverParams structure has to contain the same structure as the Cplex.Param structue in a Cplex object. But values can be set by directly setting the respective parameter instead of the Cur value which would be necessary if directly modifying the Cplex.Param structure. Usage

LP = setCplexParam(LP, solverParams, verbFlag);

Inputs

  • LP – IBM-ILOG Cplex object

  • solverParams – parameter structure for Cplex. Check LP.Param. For example, [solverParams.simplex.display, solverParams.tune.display, solverParams.barrier.display, solverParams.sifting.display, solverParams.conflict.display] = deal(0); [solverParams.simplex.tolerances.optimality, solverParams.simplex.tolerances.feasibility] = deal(1e-9,1e-8);

    The full set of parameters can be obtained by calling ‘Cplex().Param’

  • verbFlag – true to show which parameter input is problematic if any (optional, default true)

setCplexParametersForProblem(cplexProblem, cobraParams, solverParams, problemType)[source]

Set the parameters for a specific problem from the COBRA Parameter structure and a solver specific parameter structre (latter has precedence). The cobra parameters structure contains fields as specified in getCobraSolverParamsOptionsForType, while solverParams needs to contain a structure compatible with setCplexParam. Usage

cplexProblem = setCplexParametersForProblem(cplexProblem, cobraParams, solverParams, ProblemType)

Inputs

  • cplexProblem – the Cplex() object to set the parameters
  • cobraParams – the COBRA parameter structure, with parameters as defined in getCobraSolverParamsOptionsForType
  • solverParams – the solver specific parameter structure has to be compatible with setCplexParam
  • problemType – The type of Problem (‘LP’,’MILP’,’QP’,’MIQP’).

set the printLevel to the cobra Parameters

solveCobraCPLEX(model, printLevel, basisReuse, conflictResolve, contFunctName, minNorm)[source]

Calls CPLEX to solve an LP or QP problem using the matlab API to cplex written by ILOG

Usage

[solution, model] = solveCobraCPLEX(model, printLevel, basisReuse, conflictResolve, contFunctName, minNorm)

Input

  • model – structure with mandatory and optional fields:
    • .A or .S: - m x n LHS matrix
    • .b - m x 1 RHS vector
    • .c - n x 1 Objective coeff vector
    • .lb - n x 1 Lower bound vector
    • .ub - n x 1 Upper bound vector
    • .osense - scalar Objective sense (-1 max, +1 min)
    • .rxns - (optional) cell array of reaction abbreviations (necessary for making a readable confilict resolution file).
    • .csense - (optional) Constraint senses, a string containting the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).
    • .LPBasis - (optional) Basis from previous solution of similar LP problem. See basisReuse

Optional inputs

  • printLevel – Printing level in the CPLEX m-file and CPLEX C-interface.

    • 0 - Silent
    • 1 - Warnings and Errors
    • 2 - Summary information (Default)
    • 3 - More detailed information
    • > 10 - Pause statements, and maximal printing (debug mode)
  • basisReuse – 0 - Use this for one of solution of an LP (Default); 1 - Returns a basis for reuse in the next LP i.e. outputs model.LPBasis

  • conflictResolve – 0 (Default); 1 If LP problem is proven to be infeasible by CPLEX, it will print out a ‘conflict resolution file’, which indicates the irreducible infeasible set of equaltiy & inequality constraints that together, combine to make the problem infeasible. This is useful for debugging an LP problem if you want to try to resolve a constraint conflict

  • contFunctName

    1. contFunctName = [] Use all default CLPEX control parameters, (Default);

    2. contFunctName = someString e.g. ‘someFunctionName’ uses the user specified control parameters defined in someFunctionName.m; (see template function CPLEXParamSet for details). 3. contFunctName = cpxControl structure (output from a file like CPLEXParamSet.m)

  • minNorm – {(0), 1 , n x 1 vector} If not zero then, minimise the Euclidean length of the solution to the LP problem. Gives the same objective, but minimises the square of flux. minNorm ~1e-6 should be high enough for regularisation yet keep the same objective

Output

  • solution – Structure containing the following fields describing a LP solution:
    • .full: Full LP solution vector
    • .obj: Objective value
    • .rcost: Lagrangian multipliers to the simple inequalties (Reduced costs)
    • .dual: Lagrangian multipliers to the equalities
    • .nInfeas: Number of infeasible constraints
    • .sumInfeas: Sum of constraint violation
    • .stat: COBRA Standardized solver status code:
      • 1 - Optimal solution
      • 2 - Unbounded solution
      • 0 - Infeasible
      • -1 - No solution reported (timelimit, numerical problem etc)
    • .origStat: CPLEX status code. Use cplexStatus(solution.origStat) for more information from the CPLEX solver
    • .solver solver used by cplex
    • .time time taken to solve the optimization problem

Optional output

  • model – with field:
    • .LPBasis - When input basisReuse=1, we return a basis for reuse in the next LP

CPLEX consists of 4 different LP solvers which can be used to solve sysbio optimization problems you can control which of the solvers, e.g. simplex vs interior point solver using the CPLEX control parameter cpxControl.LPMETHOD. At the moment, the solver is automatically chosen for you.

solveCobraLP(LPproblem, varargin)[source]

Solves constraint-based LP problems

Usage

solveCobraLP(LPproblem, varargin)

Input

  • LPproblem – Structure containing the following fields describing the LP problem to be solved
    • .A - LHS matrix
    • .b - RHS vector
    • .c - Objective coeff vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector
    • .osense - Objective sense (-1 means maximise (default), 1 means minimise)
    • .csense - Constraint senses, a string containting the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).

Optional inputs

  • varargin – Additional parameters either as parameter struct, or as parameter/value pairs. A combination is possible, if the parameter struct is either at the beginning or the end of the optional input. All fields of the struct which are not COBRA parameters (see getCobraSolverParamsOptionsForType) for this problem type will be passed on to the solver in a solver specific manner. Some optional parameters which can be passed to the function as parameter value pairs, or as part of the options struct are listed below:
  • printLevel – Printing level
    • 0 - Silent (Default)
    • 1 - Warnings and Errors
    • 2 - Summary information
    • 3 - More detailed information
    • > 10 - Pause statements, and maximal printing (debug mode)
  • saveInput – Saves LPproblem to filename specified in field. i.e. parameters.saveInput = ‘LPproblem.mat’;
  • minNorm – {(0), scalar , n x 1 vector}, where [m, n] = size(S); If not zero then, minimise the Euclidean length of the solution to the LP problem. minNorm ~1e-6 should be high enough for regularisation yet maintain the same value for the linear part of the objective. However, this should be checked on a case by case basis, by optimization with and without regularisation.
  • primalOnly – {(0), 1}; 1 = only return the primal vector (lindo solvers)
  • solverParams – solver-specific parameter structure. Formats supported are ILOG cplex and Tomlab parameter syntax. see example for details.

Output

  • solution – Structure containing the following fields describing a LP solution: * .full: Full LP solution vector * .obj: Objective value * .rcost: Reduced costs, dual solution to \(lb <= v <= ub\) * .dual: dual solution to A*v (‘E’ | ‘G’ | ‘L’) b * .solver: Solver used to solve LP problem * .algorithm: Algorithm used by solver to solve LP problem * .stat: Solver status in standardized form

    • 1 - Optimal solution
    • 2 - Unbounded solution
    • 3 - Partial success (OPTI-csdp) - will not give desired result from OptimizeCbModel
    • 0 - Infeasible
    • -1 - No solution reported (timelimit, numerical problem etc)
    • .origStat: Original status returned by the specific solver
    • .time: Solve time in seconds
    • .basis: (optional) LP basis corresponding to solution

Note

Optional parameters can also be set through the solver can be set through changeCobraSolver(‘LP’, value); changeCobraSolverParams(‘LP’, ‘parameter’, value) function. This includes the minNorm and the printLevel flags.

Example

%Optional parameters can be entered in three different ways {A,B,C}

%A) as a generic solver parameter followed by parameter value:
[solution] = solveCobraLP(LPCoupled, 'printLevel', 1);
[solution] = solveCobraLP(LPCoupled, 'printLevel', 1, 'feasTol', 1e-8);

%B) parameters structure with field names specific to a particular solvers
% internal parameter fields
[solution] = solveCobraLP(LPCoupled, parameters);

%C) as parameter followed by parameter value, with a parameter structure
%with field names specific to a particular solvers internal parameter,
%fields as the LAST argument
[solution] = solveCobraLP(LPCoupled, 'printLevel', 1, 'feasTol', 1e-6, parameters);
solveCobraLPCPLEX(LPProblem, printLevel, basisReuse, conflictResolve, contFunctName, minNorm, interface)[source]

Calls CPLEX to solve an LP problem By default, use the matlab interface to cplex written by TOMLAB, in preference to the one written by ILOG.

Usage

[solution, LPProblem] = solveCobraLPCPLEX(LPProblem, printLevel, basisReuse, conflictResolve, contFunctName, minNorm, interface)

Input

  • LPProblem – Structure containing the following fields describing the LP problem to be solved
    • .A - LHS matrix
    • .b - RHS vector
    • .c - Objective coeff vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector
    • .osense - Objective sense (-1 max, +1 min)
    • .rxns - (optional) cell array of reaction abbreviations (necessary for making a readable confilict resolution file).
    • .csense - (optional) Constraint senses, a string containting the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).
    • .LPBasis - (optional) Basis from previous solution of similar LP problem. See basisReuse

Optional inputs

  • printLevel – Printing level in the CPLEX m-file and CPLEX C-interface.

    • 0 - Silent
    • 1 - Warnings and Errors
    • 2 - Summary information (Default)
    • 3 - More detailed information
    • > 10 - Pause statements, and maximal printing (debug mode)
  • basisReuse – 0 - Use this for one of solution of an LP (Default); 1 - Returns a basis for reuse in the next LP i.e. outputs model.LPBasis

  • conflictResolve – 0 (Default); 1 If LP problem is proven to be infeasible by CPLEX, it will print out a ‘conflict resolution file’, which indicates the irreducible infeasible set of equaltiy & inequality constraints that together, combine to make the problem infeasible. This is useful for debugging an LP problem if you want to try to resolve a constraint conflict

  • contFunctName – structure or function with parameters (only for tomlab_cplex or ILOGcomplex)

    • when using the tomlab_cplex interface

      1. contFunctName = [] Use all default CLPEX control parameters, (Default);
      2. contFunctName = someString e.g. ‘someFunctionName’ uses the user specified control parameters defined in someFunctionName.m (see template function CPLEXParamSet for details).
      3. contFunctName = cpxControl structure (output from a file like CPLEXParamSet.m)
    • when using the ILOGcomplex interface (parameter structure for Cplex). The full set of parameters can be obtained by calling Cplex().Param. For example:

      • [solverParams.simplex.display, solverParams.tune.display, solverParams.barrier.display, solverParams.sifting.display, solverParams.conflict.display] = deal(0);
      • [solverParams.simplex.tolerances.optimality, solverParams.simplex.tolerances.feasibility] = deal(1e-9, 1e-8);
  • minNorm – {(0), 1 , n x 1 vector} If not zero then, minimise the Euclidean length of the solution to the LP problem. Gives the same objective, but minimises the square of flux. minNorm ~1e-6 should be high enough for regularisation yet keep the same objective

  • interface – {‘ILOGcomplex’, ‘ILOGsimple’, ‘tomlab_cplex’} Default is the tomlab_cplex interface

Output

  • solution – Structure containing the following fields describing a LP solution:
    • .full: Full LP solution vector
    • .obj: Objective value
    • .rcost: Lagrangian multipliers to the simple inequalties (Reduced costs)
    • .dual: Lagrangian multipliers to the equalities
    • .nInfeas: Number of infeasible constraints
    • .sumInfeas: Sum of constraint violation
    • .stat: COBRA Standardized solver status code:
      • 1 - Optimal solution
      • 2 - Unbounded solution
      • 0 - Infeasible
      • -1 - No solution reported (timelimit, numerical problem etc)
    • .origStat: CPLEX status code. Use cplexStatus(solution.origStat) for more information from the CPLEX solver
    • .solver solver used by cplex
    • .time time taken to solve the optimization problemtime taken to solve the optimization problem

Optional output

  • LPProblem – with field:
    • .LPBasis: When input basisReuse = 1, we return a basis for reuse in the next LP

CPLEX consists of 4 different LP solvers which can be used to solve sysbio optimization problems you can control which of the solvers, e.g. simplex vs interior point solver using the CPLEX control parameter cpxControl.LPMETHOD. At the moment, the solver is automatically chosen for you

solveCobraLPCPLEXcard(LPProblem, printLevel, basisReuse, conflictResolve, contFunctName, minNorm, theNorm)[source]

Calls CPLEX to solve an LP problem By default, use the matlab interface to cplex written by TOMLAB, in preference to the one written by ILOG.

Usage

solveCobraLPCPLEXcard(LPProblem, printLevel, basisReuse, conflictResolve, contFunctName, minNorm, theNorm)

Input

  • model – structure with mandatory and optional fields:
    • .A or .S: - m x n LHS matrix
    • .b - m x 1 RHS vector
    • .c - n x 1 Objective coeff vector
    • .lb - n x 1 Lower bound vector
    • .ub - n x 1 Upper bound vector
    • .osense - scalar Objective sense (-1 max, +1 min)
    • .rxns - (optional) cell array of reaction abbreviations (necessary for making a readable confilict resolution file).
    • .csense - (optional) Constraint senses, a string containting the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).
    • .LPBasis - (optional) Basis from previous solution of similar LP problem. See basisReuse

Optional inputs

  • printLevel – Printing level in the CPLEX m-file and CPLEX C-interface.
    • 0 - Silent
    • 1 - Warnings and Errors
    • 2 - Summary information (Default)
    • 3 - More detailed information
    • > 10 - Pause statements, and maximal printing (debug mode)
  • basisReuse – 0 - Use this for one of solution of an LP (Default); 1 - Returns a basis for reuse in the next LP i.e. outputs model.LPBasis
  • conflictResolve – 0 (Default); 1 If LP problem is proven to be infeasible by CPLEX, it will print out a ‘conflict resolution file’, which indicates the irreducible infeasible set of equaltiy & inequality constraints that together, combine to make the problem infeasible. This is useful for debugging an LP problem if you want to try to resolve a constraint conflict
  • contFunctName – = [] Use all default CLPEX control parameters, (Default); = someString e.g. ‘someFunctionName’ uses the user specified control parameters defined in someFunctionName.m (see template function CPLEXParamSet for details). = cpxControl structure (output from a file like CPLEXParamSet.m)
  • minNorm – {(0), 1 , n x 1 vector} If not zero then, minimise the Euclidean length of the solution to the LP problem. Gives the same objective, but minimises the square of flux. minNorm ~1e-6 should be high enough for regularisation yet keep the same objective
  • theNorm – {‘zero’, ‘one’, (‘two’)} Controls which norm is minimized. ‘zero’ minimizes cardinality for nonzero entries in minNorm, ‘one’ minimizes taxicab norm for nonzero entries in minNorm (not implemented), ‘two’ minimizes Euclidean norm for nonzero entries in minNorm (default)

Output

  • solution – Structure containing the following fields describing a LP solution:
    • .full: Full LP solution vector
    • .obj: Objective value
    • .rcost: Lagrangian multipliers to the simple inequalties (Reduced costs)
    • .dual: Lagrangian multipliers to the equalities
    • .nInfeas: Number of infeasible constraints
    • .sumInfeas: Sum of constraint violation
    • .stat: COBRA Standardized solver status code:
      • 1 - Optimal solution
      • 2 - Unbounded solution
      • 0 - Infeasible
      • -1 - No solution reported (timelimit, numerical problem etc)
    • .origStat: CPLEX status code. Use cplexStatus(solution.origStat) for more information from the CPLEX solver
    • .solver solver used by cplex
    • .time time taken to solve the optimization problem

Optional output

  • LPProblem – with field:
    • .LPBasis: When input basisReuse = 1, we return a basis for reuse in the next LP

CPLEX consists of 4 different LP solvers which can be used to solve sysbio optimization problems you can control which of the solvers, e.g. simplex vs interior point solver using the CPLEX control parameter cpxControl.LPMETHOD. At the moment, the solver is automatically chosen for you.

solveCobraLPLindo(A, b, c, csense, lb, ub, osense, primalOnlyFlag, oldAPIFlag, verbLevel, method)[source]

Solves a LP problem using Lindo

Usage

[obj, x, y, w, s, solStatus] = solveCobraLPLindo(A, b, c, csense, lb, ub, osense, primalOnlyFlag, oldAPIFlag, verbLevel, method)

Inputs

  • A – LHS matrix
  • b – RHS vector
  • c – Objective coeff vector

Optional inputs

  • csense – Constraint senses
  • lb – Lower bound vector
  • ub – Upper bound vector
  • osense – Objective sense
  • primalOnlyFlag – Get the primal soln only
  • oldAPIFLag – should be true if Lindo API <2.0 is used and false for newer versions of the API
  • verbLevel – verbose level
  • method – default = 0, chooses the algorithm

Outputs

  • obj – Objective value
  • xn vector: value of each variable in x
  • ym vector: dual variables
  • w – vector of reduced costs
  • sm vector: value of each slack in s
  • solStatus – Original status returned by the specific solver
solveCobraMILP(MILPproblem, varargin)[source]

Solves constraint-based MILP problems The solver is defined in the CBT_MILP_SOLVER global variable (set using changeCobraSolver). Solvers currently available are ‘tomlab_cplex’ and ‘glpk’

Usage

solution = solveCobraMILP(MILPproblem, parameters)

Input

  • MILPproblem – Structure containing the following fields describing the LP problem to be solved
    • .A - LHS matrix
    • .b - RHS vector
    • .c - Objective coeff vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector
    • .osense - Objective sense (-1 max, +1 min)
    • .csense - Constraint senses, a string containting the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).
    • .vartype - Variable types (‘C’ continuous, ‘I’ integer, ‘B’ binary)
    • .x0 - Initial solution

Optional parameters can be entered using parameters structure or as parameter followed by parameter value: i.e. ,’printLevel’, 3) Setting parameters = ‘default’ uses default setting set in getCobraSolverParameters.

Optional inputs

  • varargin – Additional parameters either as parameter struct, or as parameter/value pairs. A combination is possible, if the parameter struct is either at the beginning or the end of the optional input. All fields of the struct which are not COBRA parameters (see getCobraSolverParamsOptionsForType) for this problem type will be passed on to the solver in a solver specific manner. Some optional parameters which can be passed to the function as parameter value pairs, or as part of the options struct are listed below:
  • timeLimit – Global solver time limit
  • intTol – Integrality tolerance
  • relMipGapTol – Relative MIP gap tolerance
  • logFile – Log file (for CPLEX)
  • printLevel – Printing level
    • 0 - Silent (Default)
    • 1 - Warnings and Errors
    • 2 - Summary information
    • 3 - More detailed information
    • > 10 - Pause statements, and maximal printing (debug mode)
  • saveInput – Saves LPproblem to filename specified in field. i.e. parameters.saveInput = ‘LPproblem.mat’;

Output

  • solution – Structure containing the following fields describing a MILP solution
    • .cont: Continuous solution
    • .int: Integer solution
    • .full: Full MILP solution vector
    • .obj: Objective value
    • .solver: Solver used to solve MILP problem
    • .stat: Solver status in standardized form (see below)
      • 1 - Optimal solution found
      • 2 - Unbounded solution
      • 0 - Infeasible MILP
      • -1 - No integer solution exists
      • 3 - Other problem (time limit etc, but integer solution exists)
    • .origStat: Original status returned by the specific solver
    • .time: Solve time in seconds
solveCobraMIQP(MIQPproblem, varargin)[source]

Solves constraint-based QP problems The solver defined in the CBT_MIQP_SOLVER global variable (set using changeCobraSolver). Solvers currently available are ‘tomlab_cplex’

Usage

solution = solveCobraMIQP(MIQPproblem, varargin)

Solves problems of the type \(min osense * 0.5 x' * F * x + osense * c' * x\) s/t. \(lb <= x <= ub\) \(A * x <=/=/>= b\) xi = integer

Input

  • MIQPproblem – Structure containing the following fields describing the MIQP problem to be solved
    • .A - LHS matrix
    • .b - RHS vector
    • .F - F matrix for quadratic objective (see above)
    • .c - Objective coeff vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector
    • .osense - Objective sense (-1 max, +1 min)
    • .csense - Constraint senses, a string containting the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).

Optional parameters can be entered using parameters structure or as parameter followed by parameter value: i.e. ,’printLevel’, 3) Setting parameters = ‘default’ uses default setting set in getCobraSolverParameters.

Optional inputs

  • varargin – Additional parameters either as parameter struct, or as parameter/value pairs. A combination is possible, if the parameter struct is either at the beginning or the end of the optional input. All fields of the struct which are not COBRA parameters (see getCobraSolverParamsOptionsForType) for this problem type will be passed on to the solver in a solver specific manner. Some optional parameters which can be passed to the function as parameter value pairs, or as part of the options struct are listed below:
  • printLevel – Print level for solver
  • saveInput – Saves LPproblem to filename specified in field.

Output

  • solution – Structure containing the following fields describing a MIQP solution
    • .full: Full MIQP solution vector
    • .obj: Objective value
    • .solver: Solver used to solve MIQP problem
    • .stat: Solver status in standardized form (see below)
      • 1 - Optimal solution found
      • 2 - Unbounded solution
      • 0 - Infeasible MIQP
      • -1 - No optimal solution found (time limit etc)
      • 3 - Solution exists but with problems
    • .origStat: Original status returned by the specific solver
    • .time: Solve time in seconds
solveCobraNLP(NLPproblem, varargin)[source]

Solves a COBRA non-linear (objective and/or constraints) problem. Solves a problem of the following form: optimize objFunction(x) or c’*x st. \(A*x <=> b or b_L < A*x < b_U\) and \(d_L < d(x) < d_U\) where A is a matrix, d(x) is an optional function and the objective is either a general function or a linear function.

Usage

solution = solveCobraNLP(NLPproblem, varargin)

Input

  • NLPproblem – Non-linear optimization problem (fields up to ‘c’ are mandatory, below c are optional)
    • .A - LHS matrix
    • .b - RHS vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector
    • .osense - Objective sense (-1 for maximisation, 1 for minimisation)
    • .csense - Constraint senses (‘L’,’E’,’G’)
    • .objFunction - Function to evaluate as the objective (The function will receive two inputs, First the flux vector to evaluate and second the NLPproblem struct. The function should be provided as a string (or c)
    • .c - linear objective such that c*x is optimized.
    • .x0 - Initial solution
    • .gradFunction - Name of the function that computes the n x 1 gradient vector (ignored if ‘d’ is set).
    • .H - Name of the function that computes the n x n Hessian matrix
    • .fLowBnd - A lower bound on the function value at optimum.
    • .d - Name of function that computes the mN nonlinear constraints
    • .dd - Name of function that computes the constraint Jacobian mN x n
    • .d2d - Name of function that computes the second part of the Lagrangian function (only needed for some solvers)
    • .d_L - Lower bound vector in nonlinear constraints
    • .d_U - Upper bound vector in nonlinear constraints
    • .user - Solver specific user parameters structure

Note that ‘b_L’ and ‘b_U’ can be used in place of ‘b’ and ‘csense’

Optional parameters can be entered using parameters structure or as parameter followed by parameter value: i.e. ,’printLevel’, 3) Setting parameters = ‘default’ uses default setting set in getCobraSolverParameters.

Optional input

  • varargin – Additional parameters either as parameter struct, or as parameter/value pairs. A combination is possible, if the parameter struct is either at the beginning or the end of the optional input. All fields of the struct which are not COBRA parameters (see getCobraSolverParamsOptionsForType) for this problem type will be passed on to the solver in a solver specific manner.

Output

  • solution – Structure containing the following fields describing a NLP solution:
    • .full: Full NLP solution vector
    • .obj: Objective value
    • .rcost: Reduced costs
    • .dual: Dual solution
    • .solver: Solver used to solve NLP problem
    • .stat: Solver status in standardized form
      • 1 - Optimal solution
      • 2 - Unbounded solution
      • 0 - Infeasible
      • -1 - No solution reported (timelimit, numerical problem etc)
    • .origStat: Original status returned by the specific solver
    • .time: Solve time in seconds
    • .origSolStruct Original solution structure%
solveCobraQP(QPproblem, varargin)[source]

Solves constraint-based QP problems

The solver is defined in the CBT_MILP_SOLVER global variable (set using changeCobraSolver). Solvers currently available are ‘tomlab_cplex’, ‘mosek’ and ‘qpng’ (limited support for small problems)

Solves problems of the type \(min 0.5 x' * F * x + osense * c' * x\) s/t \(lb <= x <= ub\) \(A * x <=/=/>= b\)

Usage

solution = solveCobraQP(QPproblem, varargin)

Input

  • QPproblem – Structure containing the following fields describing the QP
    • .A - LHS matrix
    • .b - RHS vector
    • .F - F matrix for quadratic objective (see above)
    • .c - Objective coeff vector
    • .lb - Lower bound vector
    • .ub - Upper bound vector
    • .osense - Objective sense for the linear part (-1 max, +1 min)
    • .csense - Constraint senses, a string containing the constraint sense for each row in A (‘E’, equality, ‘G’ greater than, ‘L’ less than).

Optional parameters can be entered using parameters structure or as parameter followed by parameter value: i.e. ,’printLevel’, 3) Setting parameters = ‘default’ uses default setting set in getCobraSolverParameters.

Optional inputs

  • varargin – Additional parameters either as parameter struct, or as parameter/value pairs. A combination is possible, if the parameter struct is either at the beginning or the end of the optional input. All fields of the struct which are not COBRA parameters (see getCobraSolverParamsOptionsForType) for this problem type will be passed on to the solver in a solver specific manner. Some optional parameters which can be passed to the function as parameter value pairs, or as part of the options struct are listed below:
  • printLevel – Print level for solver
  • saveInput – Saves LPproblem to filename specified in field.

Output

  • solution – Structure containing the following fields describing a QP solution
    • .full: Full QP solution vector
    • .rcost: Reduced costs, dual solution to \(lb <= x <= ub\)
    • .dual: dual solution to \(A*x <=/=/>= b\)
    • .slack: slack variable such that \(A*x + s = b\)
    • .obj: Objective value
    • .solver: Solver used to solve QP problem
    • .origStat: Original status returned by the specific solver
    • .time: Solve time in seconds
    • .stat: Solver status in standardized form (see below)
      • 1 - Optimal solution
      • 2 - Unbounded solution
      • 0 - Infeasible
      • -1 - No solution reported (timelimit, numerical problem etc)
tuneParam(LPProblem, contFunctName, timelimit, nrepeat, printLevel)[source]

Optimizes cplex parameters to make model resolution faster. Particularly interetsing for large-scale MILP models and repeated runs of optimisation. While, it won’t optimize memory space nor model constraints for numerical infeasibilities, tuneParam will provide the optimal set of solver parameters for feasbile models. It requires IBM ILOG cplex (for now).

Usage

optimalParameters = tuneParam(LPProblem,contFunctName,timelimit,nrepeat,printLevel);

Inputs

  • LPProblem – MILP as COBRA LP problem structure
  • contFunctName – Parameters structure containing the name and value. A set of routine parameters will be added by the solver but won’t be reported.
  • timelimit – default is 10000 second
  • nrepeat – number of row/column permutation of the original problem, reports robust results. sets the CPX_PARAM_TUNINGREPEAT parameter High values of nrepeat would require consequent memory and swap.
  • printLevel – 0/(1)/2/3

Output

  • optimParam – structure of optimal parameter values directly usable as contFunctName argument in solveCobraLP function
tuneParamForModel(model, varargin)[source]

Optimizes cplex parameters to make model resolution faster. Particularly interetsing for large-scale MILP models and repeated runs of optimisation. While, it won’t optimize memory space nor model constraints for numerical infeasibilities, tuneParam will provide the optimal set of solver parameters for feasbile models. It requires IBM ILOG cplex (for now).

Usage

optimalParameters = tuneParam(model,contFunctName,1000,1000,0);

Inputs

  • model – A COBRA model struct.
  • contFunctName – Parameters structure containing the name and value. A set of routine parameters will be added by the solver but won’t be reported.
  • timelimit – default is 10000 second
  • nrepeat – number of row/column permutation of the original problem, reports robust results. sets the CPX_PARAM_TUNINGREPEAT parameter High values of nrepeat would require consequent memory and swap.
  • printLevel – 0/1/2/3

Output

  • optimParam – structure of optimal parameter values directly usable as contFunctName argument in solveCobraLP function

Note

This is just a wrapper function that calls the tuneParam function using a Cobra model structure converted to a LP problem.

verifyCobraProblem(XPproblem, x, tol, verbose)[source]

Verifies dimensions of fields in XPproblem and determines if they are valid LP, QP, MILP, MIQP problems. Also checks inputs for NaN. If x is provided, it will see if x is a valid solution to tolerance (tol).

Usage

[statusOK, invalidConstraints, invalidVars, objective] = verifyCobraProblem(XPproblem, x, tol, verbose)

Input

  • XPproblem – struct containing:
    • .A - Constraints matrix
    • .b - rhs
    • .csense - vector of ‘E’, ‘L’, ‘G’ for equality, Less than and Greater than constraint
    • .lb, .ub - lower and upper bound on variables
    • .c - objective coefficients
    • .F - quadratic objective (optional, only used for QP, MIQP problems)
    • .vartype - vector of ‘C’, ‘I’, ‘B’ for ‘continuous’, ‘integer’, ‘binary’ variables (optional, only used for MILP, MIQP problems).

Optional inputs

  • x – Vector. Function will determine if x satisfies XPproblem
  • tol – numerical tolerance to which all constraints should be verified to. (default = 1e-8)
  • verbose – Controls whether results are printed to screen.(Default = true)

Outputs

  • statusOK – Returns -1 if any field in XPproblem has an error, returns 0 if the x vector is not valid for XPproblem and returns 1 if at least one problem type is satisfied
  • invalidConstraints – Vector which lists a 1 for any constaint that is invalid
  • invalidVars – Vector which lists a 1 for any variable that is invalid
  • objective – Objective of XPproblem