Exploration

biomassPrecursorCheck(model, checkCoupling, checkConservedQuantities, ATN)

Checks if biomass precursors are able to be synthesized.

USAGE:

[missingMets, presentMets, coupledMets, missingCofs, presentCofs] = biomassPrecursorCheck(model, checkCoupling, checkConservedQuantities)

INPUT:

model: COBRA model structure

OPTIONAL INPUT:
checkCoupling: Test, whether some compounds can only be produced

if there is a sink for other biomass precursors (Default: false)

checkConservedQuantities: true to check whether the cofactor pairs containing conserved moieties

(defined by the network structure) can be synthesized (e.g., ATP, NAD, NADPH, ACCOA, AA-charged tRNA, fatty acyl-ACP). They will otherwise be identified as missingMets (Default: false)

ATN: atom transition network outputed by buildAtomTransitionNetwork

If provided, true internal conserved moieties will be identified and used for checking conserved quantities (default [])

OUTPUTS:

missingMets: List of biomass precursors that are not able to be synthesized presentMets: List of biomass precursors that are able to be synthesized coupledMets: List of metabolites which need an exchange reaction for at least one other

biomass component because their production is coupled to it (returned only if checkCoupling = true)

missingCofs: List of cofactor pairs (defined by the network conserved moieties)

that are not able to be synthesized (returned only if checkConservedQuantities = true)

presentCofs: List of cofactor pairs that are able to be synthesized

(returned only if checkConservedQuantities = true)

essentialRxn4MultipleModels(allModels, objFun)

Performs single reactions deletions across multiple models and identifies which reactions have variable essentiality across all models for the chosen objective function (eg. it identifies reactions whos deletion would result in a zero flux through ATP consumption reaction - ATPM).

USAGE:

[essentialRxn4Models, dataStruct] = essentialRxn4MultipleModels(allModels, objFun)

INPUTS:

allModels: directory of the structure with multiple COBRA model structures objFun: string with objective function reaction (e.g. ‘ATPM’)

OUTPUT:
essentialRxn4Models: Table with reaction fluxes (within the objective function reaction)

after single deletion of model reaction (rows) across models (columns)

dataStruct: Structure with all models

Example

[essentialRxn4Models, dataStruct] = essentialRxn4MultipleModels(allModels, ‘ATPM’)

findBlockedIrrRxns(model, tol, varargin)
find all blocked irreversible reactions by solving one single LP problem:

min sum(z_pos + z_neg) s.t. Sv = 0

lb <= v <= ub v_j + tol*z_pos_j >= tol for each reaction j with lb_j >= 0 v_j - tol*z_neg_j <= -tol for each reaction j with ub_j <= 0

USAGE:

[blockedIrrRxns, fluxes] = findBlockedIrrRxns(model, tol, parameters)

INPUT:

model: COBRA model

OPTIONAL INPUTL

tol: tolerance for zeros (default feasTol*10, use default if the input is smaller) parameters: COBRA and solver-specific parameters, as a input structure or parameter/value inputs

OUTPUTS:

blockedIrrRxns: cell array of blocked irreversible reactions sol: solution structure from solveCobraLP tol: tolerance for zeros used (might be different from the input)

findBlockedReaction(model, method)

Determines those reactions which cannot carry any flux in the given simulation conditions.

USAGE:

BlockedReaction = findBlockedReaction(model)

INPUT:

model: COBRA model structure

OPTIONAL INPUT:
method: ‘FVA’ for flux variability analysis (default)

‘L2’ for 2-norm minimization

OUTPUT:

blockedReactions: List of blocked reactions

findCarbonRxns(model, nCarbonThr)

Returns the list of reactions that act of compounds which contain cabons greater than the threshold set.

USAGE:

[hiCarbonRxns, zeroCarbonRxns, nCarbon] = findCarbonRxns(model, nCarbonThr)

INPUTS:
model: Structure containing all necessary variables to described a

stoichiometric model

nCarbonThr: defines the min # of carbons that a metabolite, that is

acted on in a reaction, can have in the final list of reactions

OUTPUTS:
hiCarbonRxns: The list of reactions that act on metabolites with

greater than the threshold number of carbons

zeroCarbonRxns: Reactions with no carbon nCarbon: The number of carbons in each metabolite in the model

findEpistaticInteractions(model, doubleDeletionFitness, lethalFlag, minEffect)

Finds synthetic lethal and/or synthetic sick interactions based on double deletion analysis data

USAGE:

[interactions, epistaticEffect] = findEpistaticInteractions(model, doubleDeletionFitness, lethalFlag, minEffect)

INPUTS:

model: COBRA model structure doubleDeletionFitness: A matrix of fitness (or growth rate) values for

each of the double deletion strains. The diagonal of this matrix contains the single deletion fitness values.

OPTIONAL INPUTS:

lethalFlag: Only consider SL interactions (Default = false) minEffect: Minimum fitness effect considered to be significant

(Default = 1e-2)

OUTPUTS:
interactions: A sparse binary matrix indicating a SL or SS

interaction between two genes in the model

epistaticEffect: Magnitude of the epistatic interaction defined as

min(f1-f12, f2-f12) where f1 and f2 are the fitness values for the deletion strain of gene 1 and gene 2 respectively and f12 is the fitness value for the double deletion strain of genes 1 and 2

Note

The criteria for establishing a synthetic sick interaction are that the double deletion strain fitness must be at least minEffect lower than the fitness of either of the single deletion strains, i.e. f12 < f1 - minEffect and f12 < f2 - minEffect

The additional criterion for establishing a synthetic lethal interaction is that the double deletion fitness value is smaller than minEffect (i.e. essentially zero) f12 < minEffect

Note that the interactions matrix double counts all interactions

findExcRxns(model, inclObjFlag, irrevFlag)

Finds exchange and uptake rxns

USAGE:

[selExc, selUpt] = findExcRxns(model, inclObjFlag, irrevFlag)

INPUT:

model: COBRA model structure

OPTIONAL INPUTS:
inclObjFlag: Include objective rxns in the exchange rxn set (1) or not (0)

(Default = false)

irrevFlag: Model is in irreversible format (1) or not

(Default = false)

OUTPUTS:
selExc: Boolean vector indicating whether each reaction in

model is exchange or not

selUpt: Boolean vector indicating whether each reaction in

model is nutrient uptake or not

Note

Exchange reactions only have one non-zero (+1 / -1) element in the corresponding column of the stoichiometric matrix. Uptake reactions are exchange reactions are exchange reactions with negative lower bounds.

findGPRFromRxns(model, rxnIDs)

Get the Textual representations of the GPR rules for the indicated reactions. USAGE:

gprs = findGPRFromRxns(model,rxnIDs)

INPUTS:

model: The model to retrieve the GPR rules from rxnIDs: The reaction IDs to obtain the GPR rules for

OUTPUTS:

gprs: The textual GPR Rules

Authors:

Thomas Pfau Dec 2017

findGeneIDs(model, geneList)

Finds gene numbers in a model

USAGE:

geneID = findGeneIDs(model, geneList)

INPUTS:

model: COBRA model structure geneList: List of genes

OUTPUT:

geneID: List of gene IDs corresponding to geneList

findGenesFromMets(model, metList, param)

find the set of genes that correspond to a set of metabolites

USAGE:

[geneList] = findGenesFromMets(model, metabolites)

INPUTS:

model: COBRA model structure metList: List of metabolite identifiers in a model to find the corresponding genes for

OUTPUT:

geneList: List of genes corresponding to reactions

findGenesFromRxns(model, reactions)

Make a gene list of the genes that correspond to the selected reactions

USAGE:

[geneList] = findGenesFromRxns(model, reactions)

INPUTS:

model: COBRA model structure reactions: Reactions to find the corresponding genes for

OUTPUT:

geneList: List of genes corresponding to reactions

findHiCarbonRxns(model, nCarbonThr)

Returns the list of reactions that act on compounds which contain cabons greater than the threshold set.

USAGE:

[hiCarbonRxns, nCarbon] = findHiCarbonRxns(model, nCarbonThr)

INPUTS:
model: Structure containing all necessary variables to describe a

stoichiometric model

nCarbonThr: defines the min # of carbons that a metabolite, that is

acted on in a reaction, can have in the final list of reactions

OUTPUTS:
hiCarbonRxns: The list of reactions that act on metabolites with

greater than the threshold number of carbons

zeroCarbonRxns: Reactions with no carbon nCarbon: The number of carbons in each metabolite in the model

findMatchingFieldEntries(field, searchString, isAnnotation, minSim)

Find elements of the field that are similar to the provided identifier

USAGE:

[matchingIDs, positons, similarities] = findMatchingFieldEntries(field, searchString, isAnnotation, minSim)

INPUTS:

field: A Cell Array of strings searchString: The item to search for isAnnotation: Whether the strings in the field are annotation

strings (multiple identifiers are concatenated with ; and will be separated)

minSim: minimum “similarity” between element and search

item.

OUTPUTS:

matchingIDs: All IDs which are matching to the search string. positions: The positions of these IDs in the supplied field. similarities: Similarities of the matchingIDs with the search

String

findMetFromCompartment(model, compartment)

Finds all the metabolites and their identifiers in a compartment of interest.

USAGE:

[compartmentMetabolites] = findMetFromCompartment(model,Compartment)

INPUTS:

model: COBRA model strcture compartment: compartment of interest (e.g.: ‘[m]’, ‘[n]’, ‘[e]’, etc.)

OUTPUT:

compartmentMetabolites: List of metabolites in the compartment of interest

findMetIDs(model, metList)

Finds metabolite numbers in a model

USAGE:

metID = findMetIds(model, metList)

INPUTS:

model: COBRA model structure metList: List of metabolites

OUTPUT:

metID: List of metabolite IDs corresponding to metList

findMetabolicJunctions(model, nRxnsJnc)

Finds metabolic branchpoints with different numbers of branches

USAGE:

validJunctionMets = findMetabolicJunctions(model, nRxnsJnc)

INPUTS:

model: COBRA model structure nRxnJnc: Number of reactions to be considered a junction

OUTPUT:

validJunctionMets: List of junction metabolites

findMetsFromRxns(model, reactions)

Get all Metabolites from a set of reactions.

USAGE:

[metList, stoichiometries] = findMetsFromRxns(model, reactions) [metList] = findMetsFromRxns(model, reactions)

INPUTS:

model: COBRA model structure reactions: Reaction IDs (Cell Array) or positions (Double array) to

find the corresponding metabolites for

OUTPUT:
metList: If only one output is requested, returns an ordered set of all

metabolites involved in the provided reactions. Otherwise, this is a Cell Array of cell arrays containing the metabolites involved in each of the provided reactions.

stoichiometries: this is a Cell array of double arrays of the

stoichiometric coefficients corresponding to the reactions in the order of provided reaction ids. If reactions not in the model are provided, those will be represented by empty arrays.

findNeighborRxns(model, rxns, asSingleArray, order, commonMets, withComp)

Identifies the reactions and the corresponding genes that are adjacent (having a common metabolite) to a reaction of interest. Useful for characterizing the network around an orphan reaction.

USAGE:

[neighborRxns, neighborGenes, mets] = findNeighborRxns(model, rxns, asSingleArray, order, commonMets, withComp)

INPUTS:

model: COBRA model structure rxns: the target reaction as a string or multiple reactions as cell array

OPTIONAL INPUTS:
asSingleArray: If false, then return cell array of cell arrays with neighbor reactions

for one particular connecting metabolite and input reaction combination. Else just return all neighbors in one cell array. (Default = false)

order: maximal order of neighbors to be returned (default = 1)

order >= 2 only works with asSingleArray = true Neighborhoods of order >= 2 will usually also return the input reactions.

commonMets: Cell array of common metabolites, that should not count as edges between reactions.

Use {‘’} if no such metabolite should be included (default = {‘atp’, ‘adp’, ‘h’, ‘h2o’, ‘pi’, ‘ppi’}).

withComp: if commonMets already have a compartment identifier, e.g. ‘atp[m]’, then true (default=false)

OUTPUTS:

neighborRxns: the neighboring rxns in the network, (having common metabolites) neighborGenes: the genes associated with the neighbor rxns mets: the metabolites in the target reaction

findOrphanRxns(model)

Finds all orphan reactions in model (reactions with no associated genes), not including exchange rxns

USAGE:

orphans = findOrphanRxns(model)

INPUT:

model: a COBRA model with GPRs

OUTPUT:

orphans: all orphan reactions in the model

findRootNPmets(model, findNCmets)

Finds the root no production (and no consumption) metabolites in a model, used by gapFind

USAGE:

gaps = findRootNPmets(model, findNCmets)

INPUT:

model a COBRA model

OPTIONAL INPUT:

findNCmets: find no consumption mets as well as no production (default = false)

OUTPUT:

gaps: all root no production metabolites

findRxnFromCompartment(model, compartment)

Finds all the reactions and their identifiers in a compartment of interest.

USAGE:

[compartmentReactions] = findRxnFromCompartment(model,Compartment)

INPUTS:

model: COBRA model strcture compartment: compartment of interest (e.g.: ‘[m]’, ‘[n]’, ‘[e]’, etc.)

OUTPUT:

compartmentMetabolites: List of reactions in the compartment of interest

findRxnIDs(model, rxnList)

Finds reaction indices in a model

USAGE:

rxnID = findRxnIDs(model, rxnList)

INPUTS:

model: COBRA model structure rxnList: cell array of reaction abbreviations

OUTPUT:

rxnID: indices for reactions corresponding to rxnList

findRxnsActiveWithGenes(model, genes)

Finds all reactions for which the provided genes show sufficient evidence (i.e. are usfficient according to the GPRs)

USAGE:

[Reaclist] = findRxnsActiveWithGenes(model, genes)

INPUT:

model: COBRA model structure genes: A list of gene identifiers

OUTPUTS:

Reaclist: Cell array of reactions which are supported by the provided genes

Note

Only reactions which do have a GPR are considered for this function. Reactions without GPRs are ignored, as we don’t have evidence.

findRxnsFromGenes(model, genes, numericFlag, ListResultsFlag)

Print every reaction associated with a gene of interest

USAGE:

[results, ListResults] = findRxnsFromGenes(model, genes, numericFlag, ListResultsFlag)

INPUTS:

model: COBRA model structure genes: string of single gene or cell array of multiple

genes for which rxns are desired.

OPTIONAL INPUTS:

numericFlag: 1 if using Human Recon (Default = 0) ListResultsFlag: 1 if you want to output ListResults (Default = 0)

OUTPUTS:
results: structure containing cell arrays for each gene.

Each cell array has one column of rxn abbreviations and one column containing the reaction formula

ListResults: same as above, but in a cell array

  • Ines Thiele Feb 2021 - Now also grRule is returned in the output

findRxnsFromMets(model, metList, varargin)

Returns a list of reactions in which at least one metabolite listed in metList participates.

USAGE:

[rxnList, rxnFormulaList] = findRxnsFromMets(model, metList, printFlag)

INPUTS:

model: COBRA model structure metList: Metabolite list

OPTIONAL INPUTS:

printFlag: Print reaction formulas to screen (Default = false) Property/Value: Allowed Properties are:

  • containsAll - If true only reactions containing all metabolites in metList are returned (Default: false)

  • printFlag - as above, will overwrite a printFlag set individually (Default: false)

  • producersOnly - Only return reactions which produce any of the given metabolites. (Default: false)

  • consumersOnly - Only return reactions which consume any of the given metabolites. (Default: false)

  • exclusive - Only return those reactions, which do not contain any metabolites but those given (Default: false)

OUTPUTS:

rxnList: List of reactions rxnFormulaList: Reaction formulas coresponding to rxnList

findRxnsFromSubSystem(model, subSystem)

Find all reactions which are part of a given subsystem.

USAGE:

[reactionNames,rxnPos] = findRxnsFromSubSystem(model,subSystem)

INPUT:
model: A COBRA model struct with at least rxns and

subSystems fields

subSystem: A String identifying a subsystem

OUTPUT:

reactionNames: A Cell array of reaction names rxnPos: A double array of positions of the reactions in

reactionNames in the model (same order).

USAGE:
%Obtain all reactions with Glycolysis in their respective subSystems

field.

[reactionNames,rxnPos] = findRxnsFromSubSystem(model,’Glycolysis’)

findRxnsInActiveWithGenes(model, genes)

Finds all reactions for which the provided genes show sufficient evidence of their absence. (i.e. make the corresponding GPRs always be zero)

USAGE:

[Reaclist] = findRxnsActiveWithGenes(model, genes)

INPUTS:

model: COBRA model structure genes: A list of gene identifiers

OUTPUT:

Reaclist: Cell array of reactions which are supported by the provided genes

Note

Only reactions which do have a GPR are considered for this function. Reactions without GPRs are ignored, as we don’t have evidence.

findSubSysGen(model)

Lists the subsystem that a reaction occurs in encoded by a gene. Returns list of subsystems. If multiple reactions are associated with gene, subsystem of first occurance will be listed.

USAGE:

[GenSubSystem] = findSubSysGen(model)

INPUT:

model: COBRA model structure

OUTPUT:

GenSubSystem: array listing genes and subsystmes

findSubSystemsFromGenes(model, genes, varargin)

Returns the subsystems associated with the provided genes.

USAGE:

[geneSubSystems,singleList] = findSubSystemsFromGenes(model,genes,…)

INPUT:

model: COBRA model structure

OPTIONAL INPUTS:
genes: The genes to find subSystems for

(Default: model.genes)

varargin: Optional arguments as parameter/value pairs, or parameter struct.
  • structResult - return individual Gene Associations as a struct instead of a Cell Array with each gene represented by a field. Gene names might be adapted to fit to field names and all genes will be prefixed with gene_.

  • onlyOneSub - Only return one subSystem for each gene.

OUTPUT:

geneSubSystems: All subsystems associated with the provided genes.

OPTIONAL OUTPUT:
singleList: A Cell array of size(nGenes,2), with the first column

representing the gene and the second column containing the subSystems associated with the gene.

findTrspRxnFromMet(model, metList, compFlag)

Find transport reactions for defined metabolites. Option available to define the compartment of the transport

USAGE:

[TrspRxns] = findTrspRxnFromMet(model, metList, compFlag)

INPUTS:

model: COBRA model structure metList: metabolites list

OPTIONAL INPUT:

compFlag: compartment of the transport (e.g. ‘c’, cytosol)

OUTPUT:

TrspRxns: List of transporters reactions

getMIRIAMAnnotations(model, type, varargin)

Get the MIRIAM Annotations a struct array of annotation structs for the given type

USAGE:

annotations = getMIRIAMAnnotations(model,field,…)

INPUTS:

model: The COBRA model structure. type: the basic field to get annotations for (e.g. rxn, met, or

model

OPTIONAL INPUTS:
varagin: Additional arguments as parameter/value pairs.
  • bioQualifiers - A Cell array of BioQualifiers to look for if not provided, or empty, all bioQualifiers defined in getBioQualifiers() will be used

  • ids - A specific ID or IDs to get the annotation data for. Cannot be combined with type model.(Default: model.([type ‘s’])).

OUTPUT:
annotations: A struct array with the following structure:
  • .annotation -> A Struct array with one element per element of the given type.

  • .annotation.cvterms -> A Struct array of controlled vocabulary terms with one element per qualifier used

  • .annotation.cvterms.qualifier -> the bioQualifier for this all ressourcesof this cvterm

  • .annotation.cvterms.qualifierType -> the Qualifier type (modelQualifier or bioQualifier) for all ressources of this cvterm

  • .annotation.cvterms.ressources -> struct with the following fields:

  • .annotation.cvterms.ressources.database -> the database for this ressource.

  • .annotation.cvterms.ressources.id-> The ID in the database for this ressource.

getModelSizes(model)

Get the sizes of the basic fields of the model structure

USAGE:

[nMets, nRxns, nCtrs, nVars, nGenes, nComps] = getModelSizes(model)

INPUT:

model: A COBRA model structure

OUTPUTS:

nMets: The number of metabolites in the model nRxns: The number of reactions in the model nCtrs: The number of constraints in the model nVars: The number of variables in the model nGenes: The number of genes in the model nComps: The number of compartments in the model

getModelSubSystems(model)

Get unique set of subSystems present in a model

USAGE:

[subSystems] = getModelSubSystems(model)

INPUT:
model: A COBRA model struct with at least the

subSystems fields

OUTPUT:
subSystems: A Cell Array of strings containing all

subSystems in the model

USAGE:

%Get all subSystems present in the model. [subSystems] = getModelSubSystems(model)

getObjectiveSense(model)

Get the objective sense of the model (both the osenseStr (‘max’ or ‘min’) and the correspdoning osense value (-1 or 1)

USAGE:

[osenseStr,osense] = getObjectiveSense(model)

INPUTS:
model: The model to obtain the sense for. If the model has a

osenseStr field, it has to be either ‘min’ or ‘max’ (case insensitive) otherwise this function will error.

OUTPUS:

osenseStr: ‘min’ or ‘max’ osense: 1 or -1;

identifyFastBlockedRxns(model, rxnList, printLevel, sTol)

This function evaluates the presence of blocked reactions in a metabolic model

USAGE:

[BlockedRxns] = identifyFastBlockedRxns(model,rxnList, printLevel,sTol)

INPUTS:

organisms: model in COBRA model structure format rxnList: nx1 cell array with reactions to test printLevel: Verbose level (default: printLevel = 1) sTol: Solver tolerance for flux (default: 1e-6)

OUTPUT:

BlockedRxns: nx1 cell array containing blocked reactions

isReactionInSubSystem(model, reactions, subSystem)

Determine whether a reaction is in a given subSystem.

USAGE:

[subSystems] = getModelSubSystems(model)

INPUT:
model: A COBRA model struct with at least rxns and

subSystems fields

reactions: Either a string identifying a reaction, or a

cell array of strings identifying multiple reactions, or a double vector identifying the positions.

subSystem: A String identifying a subsystem.

OUTPUT:

present: a boolean vector for each provided reaction.

USAGE:

%Get all subSystems present in the model. [subSystems] = getModelSubSystems(model)

plotEssentialRxns(essentialRxn4Models, essentialityRange, numModelsPresent)

Identifies and plots (heatmap) a set of essential reactions (reactions of interest) based on conditional inputs.

USAGE:

rxnInterest4Models = plotEssentialRxns(essentialRxn4Models, essentialityRange, numModelsPresent)

INPUTS:
essentialRxn4Models: Table with reaction fluxes (within the objective function reaction)

after single deletion of model reaction (rows) across models (columns) This input can be obtained from essentialRxn4MultipleModels.m

essentialityRange: Range of fluxes (e.g. [-100,100]) numModelsPresent: Minimum number of models where a reaction is essential in order to be ploted.

OUTPUT:
rxnInterest4Models: Table with reaction fluxes (within the objective function reaction)

for reactions of interest (rows) across models (columns)

Example

rxnInterest4Models = plotEssentialRxns(essentialRxn4Models, essentialityRange, numModelsPresent)

printConstraints(model, minInf, maxInf, rxnSelection, modelAfter, printLevel)

Print all network constraints that are between -Inf (minInf) or +Inf (maxInf) inclusive

USAGE:

printConstraints(model, minInf, maxInf)

INPUTS:

model: COBRA model structure

OPTIONAL INPUTS:

minInf: value that is considered as -Inf (or desired minimum cutoff value) maxInf: value that is considered as +Inf (or desired maximum cutoff value) rxnSelection: boolean vector or cell array of reaction abbreviations for the reactions to be printed modelAfter: model after some perturbation to the bounds printLevel:

printCouplingConstraints(model, printLevel)

Print all reaction flux coupling constraints

USAGE:

printConstraints(model)

INPUTS:

model: COBRA model structure printLevel: Printing the constraints: rxn1 + rxn2 + rxn3 - rxn4 >= d

if printLevel>0, and the formula for the reactions if printLevel>1 (Default = 0).

printFluxBounds(model, rxns)

Prints the reactionID and upper/lower flux bounds.

USAGE:

printFluxBounds(model, rxns)

INPUTS:

model: The model to print

OPTIONAL INPUTS:
rxns: a string array of reaction ids for which flux bounds need to

be printed

printFluxVector(model, fluxData, nonZeroFlag, excFlag, sortCol, fileName, headerRow, formulaFlag, gprFlag)

Prints a flux vector with reaction labels

USAGE:

printFluxVector(model, fluxData, nonZeroFlag, excFlag, sortCol, fileName, headerRow, formulaFlag, gprFlag)

INPUTS:

model: COBRA model structure fluxData: n x k Data matrix/vector (for example, solution.v)

OPTIONAL INPUTS:

nonZeroFlag: Only print nonzero rows (Default = false) excFlag: Only print exchange fluxes (Default = false) sortCol: Column used for sorting (-1, none; 0, labels; >0, data

columns;) (Default = -1)

fileName: Name of output file (Default = []) headerRow: k x 1 cell array Header (Default = []) formulaFlag: Print reaction formulas (Default = false) gprFlag: Print reaction GPR (Default = false)

printGPRForRxns(model, rxnIDs, printLevel)

Print the GPRs (in textual format) for the given RxnIDs USAGE:

gprs = printGPRForRxns(model,rxnID)

INPUTS:

model: The model to retrieve the GPR rules from rxnIDs: The reaction IDs to obtain the GPR rules for

OPTIONAL INPUTS:

printLevel: Whether to print out the GPRs. If printLevel is 0,

the function will only return the Strings in a cell array.

OUTPUTS:

gprs: The textual GPR Rules

Authors:

Thomas Pfau Dec 2017

printLabeledData(labels, data, nonzeroFlag, sortCol, fileName, headerRow, sortMode)

Print a matrix of data with labels

USAGE:

printLabeledData(labels, data, nonzeroFlag, sortCol, fileName, headerRow, sortMode)

INPUTS:

labels: Row labels data: Data matrix/vector nonzeroFlag: Only print nonzero rows (opt) sortCol: Column used for sorting (-1, none; 0, labels; >0, data columns; opt) fileName: Name of output file (opt) headerRow: Header (opt) sortMode: Sort mode, ‘ascend’ or ‘descend’ (opt, default ‘ascend’)

printMatrix(A, format, file)

Prints matrix into a file or screen

USAGE:

printMatrix(A, format, file)

INPUTS:

A: Matrix format: Format string (opt, default ‘%6.4ft’) file: File name (opt)

printRxnFormula(model, varargin)

Prints the reaction formulas for a list of reactions

Reactions that have an upperbound <= 0 and lowerbound < 0 will have its directionality reversed unless directionFlag = false.

USAGE:

formulas = printRxnFormula(model, varargin)

INPUTS:

model: COBRA model structure

OPTIONAL INPUTS:

varargin: Optional Inputs provided as ‘ParameterName’, Value

pairs. the following parameternames are available:

  • rxnAbbrList: Cell array of rxnIDs to be printed (Default = print all reactions)

  • printFlag: Print formulas or just return them (Default = true)

  • lineChangeFlag: Append a line change at the end of each line

    (Default = true)

  • metNameFlag: Print full met names instead of abbreviations

    (Default = false)

  • fid: Optional file identifier for printing in files

    (default 1, i.e. stdout)

  • directionFlag: Checks directionality of reaction. See Note.

    (Default = false)

  • gprFlag: Print gene protein reaction association

    (Default = false)

  • proteinFlag: Print the protein names associated with the genes in the

    GPRs associated with the reactions. (Default = false)

  • printBounds: Print the upper and lower Bounds of the reaction (Default = false)

OUTPUT:

formulas: Cell array containing formulas of specified reactions

Examples

  1. print only ATPM and TKT1: printRxnFormula(model,’rxnAbbrList’,{‘ATPM’,’TKT1’});

  2. print the reactions of the model with the metabolite Names instead of ids. printRxnFormula(model, ‘metNameFlag’,true);

  3. print all reactions with the metabolite names to a given fileID printRxnFormula(model, ‘metNameFlag’,true, ‘fid’, fileID);

printUptakeBound(model)

Returns the indices to the substrate uptake reactions and prints substrate uptake bounds

USAGE:

uptakeInd = printUptakeBound(model)

INPUTS:

model: CORBRA model structure

OUTPUTS:

upInd: Vector containing indecies of uptake reactions

searchModel(model, searchTerm, varargin)

Search for the specified term in the given model. The performed search is fuzzy if similarity is lower than 1.

USAGE:

results = searchModel(model,searchTerm,…)

INPUTS:

model: The model to search in searchTerm: The term to search for varargin: Additional parameters as parameter/value pairs or value struct.

Available parameters are:
  • printLevel - Printout switch, 0 silent, 1 print results. (Default = 1) .

  • similarity - Minimum similarity (as provided in calcSim) for matches.

OUTPUTS:

results: The results struct build as follows:

  • .field - a struct array containing of the basic fields (e.g. rxns, mets etc) of the model in which matches were found with each field indicating the matches as detailed below.

  • .field.id - The id of the matching element

  • .field.matches - information about the matches of the respective ID.

  • .field.matches.value - the value that matched.

  • .field.matches.source - the field that contained the matching value.

surfNet(varargin)

A simple and convenient tool to navigate through the metabolic network interactively, possibly with given flux vectors, in the Matlab command window using mouse once the starting point is called by command.

USAGE:

surfNet(model, object, metNameFlag, flux, nonzeroFluxFlag, showMets, printFields, charPerLine, similarity) surfNet(model, object, … , ‘parameter’, ‘value’, …)

INPUT:

model: COBRA model

OPTIONAL INPUTS:
object: met, rxn or gene names, can be a cell array of multiple rxns, mets and genes, or simply a string

If the input is none of the above, will search related mets, rxns and genes using searchModel.m (default: objective reactions or a random reaction)

(the arguments below can also be inputted as parameter/value pairs and partial matching is supported)

metNameFlag: print model.metNames in reaction formulas. (default: false) flux: flux vectors for the model, a #reactions x #columns matrix. (default: [])

If given, the producing and consuming reactions displayed are in accordance with the flux direction.

nonzeroFluxFlag: show only reactions with nonzero fluxes if flux is given (default: true) showMets: show metabolites in a list when printing reactions (default: true) printFields: a cell array of field names of the fields to be printed for each met and rxn.

Default: {‘metNames’,’metFormulas’,’rxnNames’,’lb’,’ub’,’grRules’} Use, e.g., {‘default’, ‘rxnNotes’}, to print model.rxnNotes in addition to the defaulted fields.

charPerLine: max. no. of characters per line for reaction formulas (default 0, the current width of the command window) thresholdSim: threshold for the similarity of the results to be printed from searching the model for object

if object is not a met, rxn or gene. Default 0.8. See the parameter similarity in searchModel.m

Examples

surfNet(model) % for starting the navigation at the objective reaction surfNet(model, ‘glc-D[c]’) % for starting the navigation at D-glucose surfNet(model, ‘glc-D[c]’, ‘metNameFlag’, 1) % for printing model.metNames in reaction formulas surfNet(model, ‘glc-D[c]’, ‘m’, 1) % support partial mapping of the parameter name surfNet(model, ‘glc-D[c]’, ‘f’, flux) % for starting at glucose and show reactions with non-zero flux in ‘flux’ only surfNet(model, ‘EX_glc-D(e)’, ‘s’, 0) % to show the reaction formulas only but not the details of metabolites during navigation surfNet(model, ‘EX_glc-D(e)’, ‘p’, {‘default’, ‘metKEGGID’}) % to print also the charge for each metabolite printed surfNet(model, ‘EX_glc-D(e)’, ‘p’, {‘d’, ‘metK’}) % support unambiguous partial matching surfNet(model, {‘glc-D[c]’; ‘fru[c]’; ‘b1779’}) % to view several mets, rxns, genes surfNet(model, ‘glucose’) % search the model for mets, rxns, or genes similar to the term ‘glucose’ surfNet(model, ‘glucose’, ‘t’, 0.6) % search with a lower similarity threshold