public class MetabolicGraph extends org.jgrapht.graph.SimpleDirectedWeightedGraph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
Constructor and Description |
---|
MetabolicGraph(boolean reversible,
boolean compartments)
Creates an empty MetabolicGraph, a SimpleDirectedWeightedGraph.
|
MetabolicGraph(String version,
int numberOfCompounds,
boolean reversible,
boolean compartments)
Creates a MetabolicGraph using numberOfCompounds as initial compounds hash capacity.
|
Modifier and Type | Method and Description |
---|---|
org.jgrapht.graph.DefaultWeightedEdge |
addEdge(Vertex source,
Vertex target)
Adds an edge to connect source and target.
|
boolean |
addEdge(Vertex source,
Vertex target,
org.jgrapht.graph.DefaultWeightedEdge e)
Adds the edge to connect source and target.
|
boolean |
addVertex(Vertex vertex)
Adds the vertex to the graph and creates an entry in the compounds or
reactions HashMap.
|
void |
adjacencyMatrix(String file)
Deprecated.
|
double[] |
assortativities_old()
Deprecated.
|
double[] |
assortativities()
Calculates two assortativity coefficients on the directed MM-network: the pearson
correlation coefficient of (1) the compound in-degrees and the average in-degrees of
predecessors, (2) the compound out-degrees and the average out-degrees of successors.
|
BigDecimal |
averageReactionDegree()
Determines the average reaction degree, i.e., the average of the
in-degrees plus the out-degrees.
|
int[] |
balance(Vertex reaction)
Determines the mass balance of a reaction by summing over the substrate and product masses,
multiplied by the corresponding stoichiometric coefficients.
|
float |
characteristicPathLength()
Implementation of Floyd-Warshall's algorithm for computing all shortest paths in
a directed graph.
|
MetabolicGraph |
clone()
Creates an identical copy of the graph: reversibility is copied, vertices are
copied with names, mass and reversible reactions, edges are copied with weights.
|
double |
clusteringCoefficientUndirected()
Determines the undirected clustering coefficient in the corresponding
metabolite-metabolite network, i.e., the clustering coefficient in
the mm-network, where all edges are considered undirected.
|
boolean |
compare(MetabolicGraph graph)
Compares this MetabolicGraph to the given MetabolicGraph.
|
boolean |
compareEdges(MetabolicGraph graph)
Compares the edges of this MetabolicGraph to the edges of the given MetabolicGraph.
|
int |
compoundDeadEnds(String outputFile)
Prints the dead-end compounds to a file.
|
void |
compoundDegrees(String file,
boolean append)
Prints a tab-separated list of compound degrees.
|
void |
compoundInOutDegrees(String file,
boolean in,
boolean append)
Prints a tab-separated list of compound in- or out-degrees.
|
ArrayList<Integer> |
connectedComponents(boolean print)
Determines the sizes of all connected components in the graph by expanding
the transitive closure until every compound was visited.
|
double |
correlation(int[] x,
float meanX,
float[] y,
float meanY)
Calculates the empirical pearson correlation of x and y with given
means meanX and meanY.
|
int |
cycleCount(List<Set<Vertex>> components,
int n,
boolean printCycles,
boolean printComponents)
Given a set of strongly connected components, returns the number of n-cycles.
|
int |
distanceTo(MetabolicGraph graph)
Counts the differences in substrates and products between this graph
and the given graph.
|
double[] |
eigenvalueDecomposition(double[] matrix,
int m)
Calculates the full eigenvalue spectrum with double precision using the Fortran method dgeev.
|
float[] |
eigenvalueDecomposition(float[] matrix,
int m)
Calculates the full eigenvalue spectrum with float precision using the Fortran method sgeev.
|
double[] |
eigenvectorDouble(double[] A,
int m,
double d,
boolean forceDecomposition)
Calculates the eigenvector corresponding to the largest eigenvalue of A using the
Fortran routines dnaupd and dgemv.
|
float[] |
eigenvectorFloat(float[] A,
int m,
double d,
boolean forceDecomposition)
Calculates the eigenvector corresponding to the largest eigenvalue of A using the
Fortran routines snaupd and sgemv.
|
void |
equationSubstitutions(MetabolicGraph originalGraph,
HashMap<ArrayList<Integer>,EquivalenceClass> classes,
HashSet<Vertex> preserved,
String outputFile,
boolean doubles)
Writes each possible single and pair substitution of the reaction equations
of the original network to tab-separated files.
|
BigDecimal |
fixStoichiometry(Vertex reaction)
Fixes the stoichiometry of the reaction by finding the smallest factor
such that for all compounds i: factor*weight[i] is an integer.
|
HashMap<String,HashSet<Vertex>> |
getCompartments()
Returns a map of compartment strings and vertex sets.
|
HashMap<String,HashSet<Vertex>> |
getCompartments(Vertex reaction)
Returns the map of compartment strings and the substrates/products of the reaction which belong to
the compartment.
|
Vertex |
getCompound(String name,
String compartment)
Returns the compound vertex with the specified name and compartment.
|
Collection<Vertex> |
getCompounds()
Returns an immutable collection of the compound vertices.
|
Double |
getDeltaGDifference(Vertex reaction,
int side,
Vertex compound,
Vertex substitute,
int[] solve)
Calculates the difference in deltaG of the reaction after substituting compound
by substitute on the given side, without modifying the network.
|
Double |
getDeltaGDifference(Vertex reaction,
int side,
Vertex compound1,
Vertex compound2,
Vertex substitute1,
Vertex substitute2,
int[] solve)
Calculates the difference in deltaG of the reaction after substituting compound1
and compound1 by substitute1 and substitute2 on the given side.
|
Double |
getDeltaGn()
Calculates the sum of deltaGr over all reactions.
|
Double |
getDeltaGr(Vertex reaction,
boolean reversible)
Calculates the deltaGr of a reaction by summing over the negative deltaGf of
substrates and the positive deltaGf of products, multiplied by the stoichiometric
coefficients.
|
Vertex |
getReaction(String name)
Returns the reaction vertex with the specified name.
|
String |
getReactionEquation(MetabolicGraph graph,
Vertex reaction)
Returns the equation of the given reaction in the given graph as String.
|
Set<String> |
getReactionNames()
Returns the set of reactions names in the graph.
|
Collection<Vertex> |
getReactions()
Returns an immutable collection of the reaction vertices.
|
int |
getTransitionDegree(HashMap<ArrayList<Integer>,EquivalenceClass> classes,
HashSet<Vertex> preserved)
Calculates the number of possible substitutions for all reactions, i.e., the
degree of the node representing this graph in the transition graph Sigma_G.
|
boolean |
hasCompartments()
Indicates whether the graph has compartments, i.e., each compound vertex
is specified by a unique pair of a name and compartment.
|
boolean |
hasCycle(Vertex reaction)
Determines whether the given reaction has a cycle,
i.e., there is a substrates which is also a product.
|
boolean |
hasEdge(Vertex vertex1,
Vertex vertex2)
Tests whether the graph has an edge connecting vertex1 and vertex2 in either direction.
|
boolean |
isAdjacent(Vertex c1,
Vertex c2)
Tests if compound c1 is adjacent to compound c2, i.e., if there is a reaction which
has c1 as a substrate and c2 as a product.
|
boolean[] |
isConnected(ArrayList<String> compounds,
ArrayList<String> compartments,
ArrayList<Boolean> reversible)
Tests whether the given compounds are connected via paths with the given
reversibilities.
|
boolean |
isCyclic(Vertex reaction,
boolean print)
Determines whether the given reaction is a cyclic reaction,
i.e., the substrates and products are identical and have equal weights.
|
Vertex |
isDuplicate(Vertex reaction)
Returns the first found duplicate reaction.
|
boolean |
isNeighbour(Vertex v1,
Vertex v2)
Tests if the two compounds are neighbours, i.e., if they are connected
by a substrate-product or product-substrate relationship.
|
boolean |
isPath(ArrayList<String> compounds,
boolean bothDirections)
Checks if the given intermediary compounds are adjacent in any compartment.
|
boolean |
isReachable(Vertex c1,
Vertex c2,
HashSet<Vertex> visited)
Test if compound c2 is reachable from compound c1, i.e.,
if there exists a directed path from c1 to c2.
|
boolean |
isReversible()
Indicates whether the graph is reversible, i.e., for each reaction a reversed
reaction is added implicitly, which converts the products into the substrates.
|
HashMap<Vertex,Double[]> |
knockoutSet(Vertex knockout,
double threshold)
From a given knockout reaction, returns the set of reactions that are
affected by the knockout, and their knockout strength, up to the given threshold.
|
void |
leastIntegerMultiple(BigDecimal x)
Deprecated.
|
void |
leastIntegerMultiple(BigDecimal x,
BigDecimal y)
Deprecated.
|
double |
localCentrality(HashSet<Vertex> products,
HashSet<Vertex> substrates) |
double |
localCentrality(Vertex reaction1,
Vertex reaction2,
boolean reversible)
Calculates the local essentiality of reaction1 for reaction2, i.e., the
reciprocal of the smallest in-degree of a product of reaction1 which is
substrate of reaction2, or 0, if no such compound exists.
|
double |
localEssentiality(Vertex reaction,
boolean reversible)
Calculates the local reaction essentiality, i.e., the sum of essentialities of
the reaction for all its successor reactions.
|
void |
massMatrix(String outputFile,
boolean labelled)
Prints a m x n matrix containing the mass vectors of the compounds in the network,
where m is the number of compounds, and n the number of considered atomic elements
(specified as ELEMENTS in the configuration file).
|
void |
matrixSubstitutions(HashMap<ArrayList<Integer>,EquivalenceClass> classes,
HashSet<Vertex> preserved,
String outputFile,
boolean doubles,
boolean reversible)
Writes each possible single and pair substitution of the stoichiometric matrix
of the original network to tab-separated files.
|
int |
maxReactionDegree()
Returns the maximum undirected degree of a reaction in the network.
|
int |
maxReactionInDegree()
Returns the maximum in-degree of a reaction in the network.
|
int |
maxReactionOutDegree()
Returns the maximum out-degree of a reaction in the network.
|
HashSet<Vertex> |
neighborsOf(Vertex v)
Returns the set containing all second neighboring vertices.
|
HashSet<Vertex> |
predecessorsOf(Vertex v)
Returns the set containing all second predecessor vertices.
|
BigDecimal |
randomize(MetabolicGraph originalGraph,
HashMap<ArrayList<Integer>,EquivalenceClass> classes,
HashSet<Vertex> preserved,
float p,
boolean iterative,
float depth,
String deltaGFile,
String substitutabilityFile,
boolean append,
String logFile)
Performs mass balanced randomization on the graph by perturbing single compounds and pairs.
|
void |
randomScopeSizes(int[] seedSizes,
String scopeFile,
String scopeFile2,
boolean append)
Calculates random scopes for the given seed sizes and writes the scope sizes of each seed size
to a separate file.
|
TreeMap<Vertex,Number> |
reactionCentralities(MetabolicGraph originalGraph,
double d,
boolean doublePrecision,
boolean reversible)
Calculates the global propagation of the local reaction centralities.
|
void |
reactionDegrees(String file,
boolean append)
Prints a tab-separated list of reaction in-degrees and reaction out-degrees.
|
void |
readDeltaG(String file,
String logFile)
Reads the deltaGf of compounds from a tab-separated file and stores the values
with the compound vertices.
|
boolean |
removeEdge(org.jgrapht.graph.DefaultWeightedEdge e)
Removes the edge.
|
org.jgrapht.graph.DefaultWeightedEdge |
removeEdge(Vertex source,
Vertex target)
Removes the edge connecting source and target.
|
boolean |
removeVertex(Vertex vertex)
Removes the vertex from the graph and removes the entry from the
compounds or reactions HashMap.
|
HashSet<Vertex> |
scope_old(HashSet<Vertex> seed)
Deprecated.
|
Set<Vertex> |
scope(Set<Vertex> seed)
Returns the scope for the given set of seed compounds.
|
void |
setEdgeWeight(org.jgrapht.graph.DefaultWeightedEdge e,
double weight)
Sets the weight for the edge.
|
void |
stoichiometricMatrix(String file,
boolean reversible,
boolean labelled,
boolean sorted,
boolean log)
Exports the sparse stoichiometric matrix of the graph into a tab-separated file.
|
void |
substitutability(HashMap<ArrayList<Integer>,EquivalenceClass> classes,
String substitutabilityFile)
Calculate the sample space and substitutability class sizes for individual and pair substitutions.
|
HashSet<Vertex> |
successorsOf(Vertex v)
Returns the set containing all second successor vertices.
|
void |
switchRandomize(boolean strict)
Performs switch-randomization on the graph.
|
int |
validDeltaGr()
Counts the number of reactions containing only compounds
with a known deltaGf.
|
void |
weights(String file,
boolean append)
Prints the edge weights of the graph into a tab-separated file.
|
void |
write(String outputFile,
boolean printAllCoefficients)
Writes the graph to a tab- and delimiter-separated file readable by createGraph.
|
void |
write2(String outputFile,
boolean printAllCoefficients)
Writes the graph to an alternative tab- and delimiter-separated file format.
|
void |
writeBalance(String outputFile,
boolean append)
Writes the number of balanced reactions and the number of reactions which are fully annotated,
i.e., of which all substrates and products have an annotated mass.
|
containsEdge, containsVertex, degreeOf, edgeSet, edgesOf, getAllEdges, getEdge, getEdgeFactory, getEdgeSource, getEdgeTarget, getEdgeWeight, incomingEdgesOf, inDegreeOf, isAllowingLoops, isAllowingMultipleEdges, outDegreeOf, outgoingEdgesOf, setEdgeSetFactory, vertexSet
containsEdge, removeAllEdges, removeAllEdges, removeAllVertices, toString
equals, getClass, hashCode, notify, notifyAll, wait, wait, wait
public String version
public MetabolicGraph(boolean reversible, boolean compartments)
public MetabolicGraph(String version, int numberOfCompounds, boolean reversible, boolean compartments)
numberOfCompounds
- reversible
- public boolean isReversible()
public boolean hasCompartments()
public Vertex getCompound(String name, String compartment)
name
- compartment
- public Collection<Vertex> getCompounds()
public Vertex getReaction(String name)
name
- public Collection<Vertex> getReactions()
public Set<String> getReactionNames()
public String getReactionEquation(MetabolicGraph graph, Vertex reaction)
reaction
- public BigDecimal fixStoichiometry(Vertex reaction)
reaction
- public void leastIntegerMultiple(BigDecimal x)
x
- public void leastIntegerMultiple(BigDecimal x, BigDecimal y)
x
- y
- public int[] balance(Vertex reaction)
reaction
- @Deprecated public void adjacencyMatrix(String file) throws IOException
file
- IOException
public void stoichiometricMatrix(String file, boolean reversible, boolean labelled, boolean sorted, boolean log) throws IOException
file
- reversible
- labelled
- sorted
- log
- IOException
public void massMatrix(String outputFile, boolean labelled) throws IOException
outputFile
- The output file to write to.labelled
- If labelled==true, an additional header row and column are printed
containing the names of atomic elements, respectively compound names.IOException
public BigDecimal averageReactionDegree()
public int maxReactionInDegree()
public int maxReactionOutDegree()
public int maxReactionDegree()
public void compoundDegrees(String file, boolean append) throws IOException
file
- append
- IOException
public int compoundDeadEnds(String outputFile) throws IOException
outputFile
- IOException
public void reactionDegrees(String file, boolean append) throws IOException
file
- append
- IOException
public void compoundInOutDegrees(String file, boolean in, boolean append) throws IOException
file
- append
- IOException
public void weights(String file, boolean append) throws IOException
file
- IOException
public double clusteringCoefficientUndirected()
public boolean isNeighbour(Vertex v1, Vertex v2)
v1
- v2
- public boolean isAdjacent(Vertex c1, Vertex c2)
c1
- c2
- public boolean isReachable(Vertex c1, Vertex c2, HashSet<Vertex> visited)
c1
- c2
- visited
- public boolean[] isConnected(ArrayList<String> compounds, ArrayList<String> compartments, ArrayList<Boolean> reversible)
compounds
- compartments
- reversible
- public boolean isPath(ArrayList<String> compounds, boolean bothDirections)
compounds
- Intermediary compounds to check for connections.public boolean isCyclic(Vertex reaction, boolean print)
reaction
- print
- public boolean hasCycle(Vertex reaction)
reaction
- public HashMap<String,HashSet<Vertex>> getCompartments(Vertex reaction)
reaction
- public int cycleCount(List<Set<Vertex>> components, int n, boolean printCycles, boolean printComponents)
components
- n
- printCycles
- printComponents
- public void randomScopeSizes(int[] seedSizes, String scopeFile, String scopeFile2, boolean append) throws IOException
seedSizes
- scopeFile
- scopeFile2
- append
- IOException
public Set<Vertex> scope(Set<Vertex> seed)
public HashSet<Vertex> scope_old(HashSet<Vertex> seed)
seed
- public ArrayList<Integer> connectedComponents(boolean print)
public float characteristicPathLength()
Implementation of Floyd-Warshall's algorithm for computing all shortest paths in a directed graph. Adapted to calculate the average shortest path length. For reversible graphs the calculation could be faster by skipping the reverse path, but this implementation is already apx. 86 times faster than the naive algorithm. See http://algowiki.net/wiki/index.php/Floyd-Warshall%27s_algorithm..
NOTE: the result is (slightly) imprecise due to floating point representations.
public double[] assortativities()
public double[] assortativities_old()
public double correlation(int[] x, float meanX, float[] y, float meanY)
x
- meanX
- y
- meanY
- public double localCentrality(Vertex reaction1, Vertex reaction2, boolean reversible)
reaction1
- reaction2
- reversible
- public double localEssentiality(Vertex reaction, boolean reversible)
reaction
- public HashMap<Vertex,Double[]> knockoutSet(Vertex knockout, double threshold)
knockout
- threshold
- MetabolicGraph.localEssentiality(Vertex,boolean)
public TreeMap<Vertex,Number> reactionCentralities(MetabolicGraph originalGraph, double d, boolean doublePrecision, boolean reversible)
Calculates the global propagation of the local reaction centralities. These are given by the eigenvector corresponding to the largest eigenvalue of the reaction incidence matrix, normalized by the page rank transformation with damping factor d. The reaction incidence matrix is an m x m matrix containing the local essentiality scores between pairs of reactions.
If originalGraph==null, then the reaction centralities are calculated from this graph. Otherwise, the centralities are calculated repeatedly by preserving every reaction, one at a time. This allows calculating the centrality of a reaction randomized network while preserving the reaction itself, without additional randomization runs. originalGraph must be the original graph from which this network was randomized.
If reversible==true, then the local dependence between a source reaction and a target reaction is the same for both directions of a reversible source reaction (the knockout of the source affects both its directions equally). Otherwise, the local dependence considers both directions of a reversible source reaction independently. Note that both directions of a reversible target reaction are always affected independently (a knockout affects both directions independently).
originalGraph
- The originalGraph (not randomized) required for preserving
each reaction one at a time.d
- Damping factor for the PageRank transformation.doublePrecision
- If true, calculates the eigenvector in double precision,
otherwise float precision.reversible
- If true, both directions of a reversible source reaction are considered
as one, otherwise both directions are considered independently.localCentrality(Vertex,Vertex,boolean)
,
localCentrality(HashSet,HashSet)
,
eigenvectorFloat(float[],int,double,boolean)
,
eigenvectorDouble(double[],int,double,boolean)
public float[] eigenvectorFloat(float[] A, int m, double d, boolean forceDecomposition)
Calculates the eigenvector corresponding to the largest eigenvalue of A using the
Fortran routines snaupd and sgemv. If the fast calculation using the Arnoldi
method fails, then the eigenvector corresponding to the largest eigenvalue is
determined from the full eigenvalue spectrum using the LAPACK routine sgeev
(eigenvalueDecomposition(float[],int)
.
Before eigenvalue calculation, the matrix A is normalized according to page rank
with damping factor d:
A_i,j = e(r_i,r_j)*(A_i,j/colsum_j)*d+(1-d)/m)
where e(r_i,r_j) = max(1/inDegreeOf(c)) is the local essentiality for all c which are
a product of r_i and a substrate of r_j, or 0, if no such c exists;
colsum_j is the sum of all entries of column j;
m is the number of columns/rows.
The netlib-java Fortran interfaces to ARPACK were used for eigenvalue calculation, as they clearly outperformed five other java based libraries: Shared, JAMA, Colt, JBLAS, MTJ, as well as Fortran routines for calculating all eigenvalues, JLAPACK.
However, the ARPACK routines for directly calculating the largest eigenvalue
are unstable for some parameters. If an error occurs, the eigenvector corresponding
to the largest eigenvalue is determined from the full eigenvalue spectrum
using the LAPACK routines sgeev or dgeev (eigenvalueDecomposition(float[],int)
or eigenvalueDecomposition(float[],int)
).
Most critical parameters for the ARPACK routines are INFO, TOL, RESID and NCV. Certain values of TOL seem to cause numerical errors or a number overflow when using double precision and certain initial residual vectors RESID, either with INFO==0 or INFO==1. NCV causes error message if inadequate values are set.
The following parameters seem to be more or less stable for both float and double precision: d=0.85 or d=0.9 TOL = Float.MIN_VALUE, resp. Double.MIN_VALUE INFO = 0 RESID = new float[N], resp. new double[N]; (no initial residuals) NCV = Math.min(Math.max(2*NEV+1, 20), N);
The important parameter is d. The dgems/sgemv methods are highly unstable for other values of d, e.g. d=0.99, d=0.95, and d=0.995.
The following references were most informative for using ARPACK: - R. B. Lehoucq, D. C. Sorensen, C. Yang: ARPACK Users' Guide: Solution of Large Scale Eigenvalue Problems with Implicitly Restarted Arnoldi Methods. - Matlab implementation of eigs: eigs.m - SciPy Development Wiki, for snaupd.f and sneupd.f - BLAS (Basic Linear Algebra Subprograms) documentation, particularly for sgemv.f.
A
- Matrix for calculating the eigenvector in Fortran matrix style.m
- Number of rows/columns in A.d
- Damping factor.forceDecomposition
- If true and an error occurs, calculates the full eigenvalue spectrum.reactionCentralities(MetabolicGraph,double,boolean,boolean)
,
eigenvectorDouble(double[],int,double,boolean)
public double[] eigenvectorDouble(double[] A, int m, double d, boolean forceDecomposition)
Calculates the eigenvector corresponding to the largest eigenvalue of A using the
Fortran routines dnaupd and dgemv. If the fast calculation using the Arnoldi
method fails, then the eigenvector corresponding to the largest eigenvalue is
determined from the full eigenvalue spectrum using the LAPACK routine dgeev
(eigenvalueDecomposition(double[],int)
.
Before eigenvalue calculation, the matrix A is normalized according to page rank
with damping factor d:
A_i,j = e(r_i,r_j)*(A_i,j/colsum_j)*d+(1-d)/m)
where e(r_i,r_j) = max(1/inDegreeOf(c)) is the local essentiality for all c which are
a product of r_i and a substrate of r_j, or 0, if no such c exists;
colsum_j is the sum of all entries of column j;
m is the number of columns/rows.
See eigenvalueDecomposition(float[],int)
for additional
information.
A
- Matrix for calculating the eigenvector in Fortran matrix style.m
- Number of rows/columns in A.d
- Damping factor.forceDecomposition
- If true and an error occurs, calculates the full eigenvalue spectrum.reactionCentralities(MetabolicGraph,double,boolean,boolean)
,
eigenvectorFloat(float[],int,double,boolean)
public double[] eigenvalueDecomposition(double[] matrix, int m)
matrix
- The matrix in Fortran-style array format (column order).m
- Number of reactions.reactionCentralities(MetabolicGraph,double,boolean,boolean)
public float[] eigenvalueDecomposition(float[] matrix, int m)
matrix
- The matrix in Fortran-style array format (column order).m
- Number of reactions.reactionCentralities(MetabolicGraph,double,boolean,boolean)
public int getTransitionDegree(HashMap<ArrayList<Integer>,EquivalenceClass> classes, HashSet<Vertex> preserved)
public HashSet<Vertex> neighborsOf(Vertex v)
v
- public HashSet<Vertex> successorsOf(Vertex v)
v
- public HashSet<Vertex> predecessorsOf(Vertex v)
v
- public void readDeltaG(String file, String logFile) throws IOException
file
- IOException
public Double getDeltaGn()
public int validDeltaGr()
public Double getDeltaGr(Vertex reaction, boolean reversible)
reaction
- reversible
- public HashMap<String,HashSet<Vertex>> getCompartments()
public void substitutability(HashMap<ArrayList<Integer>,EquivalenceClass> classes, String substitutabilityFile) throws IOException
classes
- substitutabilityFile
- IOException
public BigDecimal randomize(MetabolicGraph originalGraph, HashMap<ArrayList<Integer>,EquivalenceClass> classes, HashSet<Vertex> preserved, float p, boolean iterative, float depth, String deltaGFile, String substitutabilityFile, boolean append, String logFile) throws IOException
Performs mass balanced randomization on the graph by perturbing single compounds and pairs. In a randomized reaction, substrates and products, and pairs thereof, are replaced by substitutes from the corresponding mass equivalence class.
The perturbations preserve the reaction degree with respect to compounds, i.e., the number of unique substrates/products is not changed. Mass balance is preserved: any reaction has the same numbers of chemical species on each side after perturbation. Generating duplicate reactions in the network is avoided.
If compartments==true, then compounds are only substituted within the same compartment. In addition, reactions which connect compartments and reactions with in- or out-degree of zero are skipped.
If preserveStrongComponents==true, then the sizes of strongly connected components are preserved in the following way: any substitution which changes the size of the component of the original compound is skipped. In that way, the creation of 'orphane' or 'dead end' metabolites is avoided.
classes
- Mass equivalence classes of the compounds.preserved
- A set of vertices (compounds or reactions) which are not randomized.p
- if p<1, then each reaction is only perturbed with probability p.depth
- factor indicating randomization strength. 1 indicates default values (see below).iterative
- if true, then each reaction is randomized (deg*deg+deg)/2 times, otherwise |V_r|*avgDeg^2 reaction sides are chosen at random and perturbed once.
reactions are randomized once.logFile
- IOException
public void matrixSubstitutions(HashMap<ArrayList<Integer>,EquivalenceClass> classes, HashSet<Vertex> preserved, String outputFile, boolean doubles, boolean reversible) throws IOException
stoichiometricMatrix(String,boolean,boolean,boolean,boolean)
with reversible==true). Otherwise, the indices correspond to the stoichiometric matrix
excluding reversed directions of reversible reactions (corresponding to a call of
stoichiometricMatrix(String,boolean,boolean,boolean,boolean)
with reversible==false).
Note: the order of substitutions is not the same as in equationSubstitutions(MetabolicGraph,HashMap,HashSet,String,boolean)
.
TODO include reaction equations from equationSubstitutions in output file.classes
- outputFile
- doubles
- reversible
- IOException
equationSubstitutions(MetabolicGraph,HashMap,HashSet,String,boolean)
public void equationSubstitutions(MetabolicGraph originalGraph, HashMap<ArrayList<Integer>,EquivalenceClass> classes, HashSet<Vertex> preserved, String outputFile, boolean doubles) throws IOException
matrixSubstitutions(HashMap,HashSet,String,boolean,boolean)
.originalGraph
- classes
- outputFile
- doubles
- IOException
matrixSubstitutions(HashMap,HashSet,String,boolean,boolean)
public Double getDeltaGDifference(Vertex reaction, int side, Vertex compound, Vertex substitute, int[] solve)
reaction
- side
- compound
- substitute
- public Double getDeltaGDifference(Vertex reaction, int side, Vertex compound1, Vertex compound2, Vertex substitute1, Vertex substitute2, int[] solve)
reaction
- side
- compound1
- compound2
- substitute1
- substitute2
- solve
- public int distanceTo(MetabolicGraph graph)
graph
- public Vertex isDuplicate(Vertex reaction)
public void switchRandomize(boolean strict)
Performs switch-randomization on the graph. Two (non-reversed) edges of the same type (substrate or product) are drawn at random and switched, i.e., the corresponding sources and targets of the edges are exchanged. The reaction degrees are preserved, because only vertices of the same type (substrate or product) are exchanged. If strict==true, then reversible reactions are never switched with irreversible reactions, such that the in- and out-degrees of compounds are preserved. If strict==false, then reversible reactions may be switched with irreversible reactions, in which case the in- and out-degrees of switched compounds increase respectively decrease by 1. In any case, the total number of chemical reactions any compound takes part in (reversible+irreversible) is preserved.
In case the drawn edges share a common source or target, the graph is not actually modified.
The number of perturbations is given by the number of distinct substrate edge pairs plus the number of distinct product edge pairs, i.e., s*(s-1)+p*(p-1), where s is the number of substrate edges and p the number of product edges.
Note that reversed reactions are not switched explicitly, as they are synchronized implicitly with
every modification (see addEdge(Vertex,Vertex)
).
strict
- public boolean hasEdge(Vertex vertex1, Vertex vertex2)
vertex1
- vertex2
- public boolean addVertex(Vertex vertex)
public boolean removeVertex(Vertex vertex)
public org.jgrapht.graph.DefaultWeightedEdge addEdge(Vertex source, Vertex target)
addEdge
in interface org.jgrapht.Graph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
addEdge
in class org.jgrapht.graph.AbstractBaseGraph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
public boolean addEdge(Vertex source, Vertex target, org.jgrapht.graph.DefaultWeightedEdge e)
addEdge
in interface org.jgrapht.Graph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
addEdge
in class org.jgrapht.graph.AbstractBaseGraph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
public org.jgrapht.graph.DefaultWeightedEdge removeEdge(Vertex source, Vertex target)
removeEdge
in interface org.jgrapht.Graph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
removeEdge
in class org.jgrapht.graph.AbstractBaseGraph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
public boolean removeEdge(org.jgrapht.graph.DefaultWeightedEdge e)
public void setEdgeWeight(org.jgrapht.graph.DefaultWeightedEdge e, double weight)
public MetabolicGraph clone()
clone
in class org.jgrapht.graph.AbstractBaseGraph<Vertex,org.jgrapht.graph.DefaultWeightedEdge>
public boolean compare(MetabolicGraph graph)
graph
- public boolean compareEdges(MetabolicGraph graph)
graph
- public void write(String outputFile, boolean printAllCoefficients) throws IOException
outputFile
- IOException
public void write2(String outputFile, boolean printAllCoefficients) throws IOException
outputFile
- IOException
public void writeBalance(String outputFile, boolean append) throws IOException
outputFile
- append
- IOException