optimizeProb-methods: Optimize Problem Object

Description Usage Arguments Details Value Methods Author(s) See Also Examples

Description

The generic optimizeProb performs the optimization of a mathematical programming object.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
## S4 method for signature 'modelorg'
optimizeProb(object,
             algorithm = SYBIL_SETTINGS("ALGORITHM"),
             gene = NULL,
             react = NULL,
             lb = NULL,
             ub = NULL,
             retOptSol = TRUE,
             obj_coef = NULL,
             lpdir = NULL,
             mtfobj = NULL,
             fldind = TRUE,
             prCmd = NA,
             poCmd = NA,
             prCil = NA,
             poCil = NA,
             ...)

## S4 method for signature 'sysBiolAlg'
optimizeProb(object,
             react = NULL,
             lb = NULL,
             ub = NULL,
             obj_coef = NULL,
             lpdir = NULL,
             fldind = TRUE,
             resetChanges = TRUE,
             prCmd = NA,
             poCmd = NA,
             prCil = NA,
             poCil = NA)

Arguments

object

An object of class modelorg or sysBiolAlg.

algorithm

A single character string giving the name of the algorithm to use. See parameter "ALGORITHM" in SYBIL_SETTINGS for possible values.
Default: SYBIL_SETTINGS("ALGORITHM").

gene

A character or integer vector containing gene id's or indices of gene id's in allGenes(model). If arguments lb and/or ub are additionally used (not NULL), upper and lower bounds will be applied to all fluxes on which the deletion of the genes given in gene have an effect. In this case, the first value in lb and ub is used. Default: NULL.

react

An object of class reactId, character or integer. Specifies the fluxes (variables) for which to change the upper and lower bound (see also arguments lb and ub) or objective coefficients (see also argument obj_coef). For class sysBiolAlg, it must be numeric. For class modelorg, setting react as no effect, if gene is also not NULL.
Default: NULL.

lb

Numeric vector, must have the same length as react. Contains the new values for the lower bounds of fluxes (variables) mentioned in react. If set to NULL, lower bounds for variables in react will be left unchanged. For class modelorg: if lb is of length one, lb is used for all elements in react.
Default: NULL.

ub

Same functionality as lb, but for upper bounds.
Default: NULL.

obj_coef

Numeric vector, must have the same length as react. Contains the new values for the objective coefficients of fluxes (variables) mentioned in react. All other objective coefficients stay untouched. If set to NULL, objective coefficients for variables in react will be left unchanged. For class modelorg: if obj_coef is of length one, obj_coef is used for all elements in react.
Default: NULL.

lpdir

Character value, direction of optimization. Can be set to "min" for minimization or "max" for maximization.
Default: SYBIL_SETTINGS("OPT_DIRECTION").

mtfobj

Only used, if argument algorithm is set to "mtf". A single numeric value giving a previously calculated optimized value of the objective function given in the model. The objective function of the model will be fixed to this value during optimization. If set to NULL, it will be computed by means of the "fba" algorithm. If additionally arguments solver and method are set, they will be used here too.
Default: NULL.

fldind

Boolean value. If set to TRUE, (default) indices in "react" are used only for reactions. If set to FALSE, indices in "react" are used for all variables during optimization, e.g. also for additional variables introduced by the mtf algorithm. Currently unused by class sysBiolAlg_room.
Default: TRUE.

resetChanges

Boolean value. If set to TRUE, (default) modifications of the problem object will be reset to their original values (e.g. changing upper and lower bounds for certain reactions). If set to FALSE, modifications will stay in the model.
Default: TRUE.

prCmd

A list of preprocessing commands. See Details below.
Default: NA.

poCmd

A list of postprocessing commands. See Details below.
Default: NA.

prCil

Can be used if optimizeProb is called several times (like in optimizer). The argument prCil gets the value of the loop variable and passes it to the preprocessing function. There, one can access it via the keyword “LOOP_VAR”. See also optimizer.
Default: NA.

poCil

Same as prCil, but for postprocessing.
Default: NA.

retOptSol

Boolean. Return an object of class optsol_optimizeProb or just a list containing the results.
Default: TRUE.

...

Only for the modelorg-method: further arguments passed to sysBiolAlg. See Details below.

Details

The arguments prCmd and poCmd can be used to execute R commands working on the problem object. All commands in prCmd are executed immediately before solving the problem; all commands in poCmd are executed after the problem has been solved. In all other aspects, the arguments work the same. The value of prCmd or poCmd are lists of character vectors (each list element is one command). Each command is a character vector and should be built as follows:

The result will be an object of class ppProc. A few examples for arguments prCmd or poCmd (all arguments must be lists, see examples section below):

1
2

will be translated to the command

1
2

with LP_PROB being the placeholder for the variable name of the problem object. The vector

1
2
    c("writeProb", "LP_PROB", "'Ec_core.lp'", "'lp'")
  

will be translated to the command

1
2
    writeProb(LP_PROB, 'Ec_core.lp', 'lp')
  

The first element will be the function name and the others the arguments to that function. The list of commands

1
2
3
4
5
6
    list("sensitivityAnalysis",
         c("getDjCPLEX", "LP_PROB@oobj@env",
           "LP_PROB@oobj@lp", "0", "react_num(Ec_core)-1"
         )
    )
  

will be translated to the commands

1
2
3
4
    sensitivityAnalysis(LP_PROB)
    getDjCPLEX(LP_PROB@oobj@env, LP_PROB@oobj@lp,
               0, react_num(Ec_core)-1)
  

For more information on the usage of prCmd and poCmd, see the examples section below.

The method optimizeProb for class modelorg generates a subclass of class sysBiolAlg and calls optimizeProb for that object again. Argument MoreArgs is used to transport arguments to the second optimizeProb call. Argument ... instead is used to transport arguments to the constructor function sysBiolAlg, for example algorithm, solver, method and solverParm. See SYBIL_SETTINGS for possible values.

Arguments gene, react, lb, ub and react cause changes in the problem object (object of class optObj, slot problem of class sysBiolAlg). These changes will be reset immediately after optimization if argument resetChanges is set to TRUE, otherwise changes will persist.

Value

Calls to optimizeProb returns either an object of class optsol_optimizeProb of length one if argument retOptSol is set to TRUE and object is of class modelorg, or a list containing the results of the optimization:

ok

Return value of the optimizer (e.g. “solution process was successful” or “time limit exceeded”).

obj

Value of the objective function after optimization.

stat

Status value of the optimization (e.g. “solution is optimal” or “no feasible solution exists”).

fluxes

The resulting flux distribution.

fldind

Pointers to columns (variables) representing a flux (reaction) in the original network. The variable fldind[i] in the solution object represents reaction i in the original network.

preP

An object of class ppProc if a preprocessing command was given.

postP

An object of class ppProc if a postprocessing command was given.

Methods

signature(object = "modelorg")

Translates the object of class modelorg into an object of class sysBiolAlg and calls optimizeProb again.

signature(object = "sysBiolAlg")

Run optimization with the given problem object.

Author(s)

Gabriel Gelius-Dietrich <geliudie@uni-duesseldorf.de>

Maintainer: Mayo Roettger <mayo.roettger@hhu.de>

See Also

modelorg, applyChanges and sysBiolAlg.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
## Not run: 
## The examples here require the package glpkAPI to be
## installed. If that package is not available, you have to set
## the argument 'solver' (the default is: solver = SYBIL_SETTINGS("SOLVER")).

## load the example data set
data(Ec_core)

## run optimizeProb(), Ec_sf will be an object of
## class optsol_optimizeProb
Ec_sf <- optimizeProb(Ec_core)

## run optimizeProb(), Ec_sf will be a list
Ec_sf <- optimizeProb(Ec_core, retOptSol = FALSE)

## do FBA, change the upper and lower bounds for the reactions
## "ATPM" and "PFK".
optimizeProb(Ec_core, react = c("ATPM", "PFK"),
             lb = c(3, -3), ub = c(5, 6))
          
## do FBA, perform sensitivity analysis after optimization
optimizeProb(Ec_core, poCmd = list("sensitivityAnalysis"))

## do FBA, write the problem object to file in lp-format
optimizeProb(Ec_core,
             poCmd = list(c("writeProb", "LP_PROB",
                            "'Ec_core.lp'", "'lp'")))

## do FBA, use "cplexAPI" as lp solver. Get all lower bounds before
## solving the problem. After solving, perform a sensitivity
## analysis and retrieve the reduced costs
opt <- optimizeProb(Ec_core, solver = "cplexAPI",
                    prCmd = list(c("getColsLowBnds", "LP_PROB", "1:77")),
                    poCmd = list("sensitivityAnalysis",
                                 c("getDjCPLEX",
                                 "LP_PROB@oobj@env",
                                 "LP_PROB@oobj@lp",
                                 "0", "react_num(Ec_core)-1")))
## get lower bounds
preProc(opt)
## get results of sensitivity analysis
postProc(opt)

## End(Not run)

sybil documentation built on May 31, 2021, 5:08 p.m.