Title: | Euclidean Distance Matrix Completion Tools |
---|---|
Description: | Implements various general algorithms to estimate missing elements of a Euclidean (squared) distance matrix. Includes optimization methods based on semi-definite programming found in Alfakih, Khadani, and Wolkowicz (1999)<doi:10.1023/A:1008655427845>, a non-convex position formulation by Fang and O'Leary (2012)<doi:10.1080/10556788.2011.643888>, and a dissimilarity parameterization formulation by Trosset (2000)<doi:10.1023/A:1008722907820>. When the only non-missing distances are those on the minimal spanning tree, the guided random search algorithm will complete the matrix while preserving the minimal spanning tree following Rahman and Oldford (2018)<doi:10.1137/16M1092350>. Point configurations in specified dimensions can be determined from the completions. Special problems such as the sensor localization problem, as for example in Krislock and Wolkowicz (2010)<doi:10.1137/090759392>, as well as reconstructing the geometry of a molecular structure, as for example in Hendrickson (1995)<doi:10.1137/0805040>, can also be solved. These and other methods are described in the thesis of Adam Rahman(2018)<https://hdl.handle.net/10012/13365>. |
Authors: | Adam Rahman [aut], R. Wayne Oldford [aut, cre, ths] |
Maintainer: | R. Wayne Oldford <[email protected]> |
License: | GPL-2 | GPL-3 |
Version: | 0.2.0 |
Built: | 2024-11-25 19:19:09 UTC |
Source: | https://github.com/cranhaven/cranhaven.r-universe.dev |
colamdR
returns the column approximate minimum degree permutation of a
sparse matrix S. The permutation of S, S[,p], will result in LU factors sparser
than S.
colamdR(M)
colamdR(M)
M |
A matrix to be permuted. |
This is an implementation of the colamd function available in SuiteSparse, and also implemented in Matlab.
A vector containing the column minimum degree permutation of the matrix M.
The authors of the code for "colamd" are Stefan I. Larimore and Timothy A. Davis ([email protected]), University of Florida.
M <- matrix(c(1,1,0,0,1,0,0,1,0,1,1,1,1,1,0,0,1,0,1,0), ncol=4) p <- colamdR(M) M[,p]
M <- matrix(c(1,1,0,0,1,0,0,1,0,1,1,1,1,1,0,0,1,0,1,0), ncol=4) p <- colamdR(M) M[,p]
dpf
returns a completed Euclidean Distance Matrix D, with dimension d,
from a partial Euclidean Distance Matrix using the methods of Trosset (2000)
dpf(D, d, toler = 1e-08, lower = NULL, upper = NULL, retainMST = FALSE)
dpf(D, d, toler = 1e-08, lower = NULL, upper = NULL, retainMST = FALSE)
D |
An nxn partial-distance matrix to be completed. D must satisfy a list of conditions (see details), with unkown entries set to NA |
d |
The dimension for the resulting completion. |
toler |
The convergence tolerance of the algorithm. Set to a default value of 1e-8 |
lower |
An nxn matrix containing the lower bounds for the unknown entries in D. If NULL, lower is set to be a matrix of 0s. |
upper |
An nxn matrix containing the upper bounds of the unknown entries in D. If NULL, upper[i,j] is set to be the shortest path between node i and node j. |
retainMST |
D logical input indicating if the current minimum spanning tree structure in D should be retained. If TRUE, a judicious choice of Lower is calculated internally such that the MST is retained. |
This is an implementation of the Dissimilarity Parameterization Formulation (DPF) for Euclidean Distance Matrix Completion, as proposed in 'Distance Matrix Completion by Numerical Optimization' (Trosset, 2000).
The method seeks to minimize the following:
where are the ordered eigenvalues of
. For details, see Trosset(2000)
The matrix D is a partial-distance matrix, meaning some of its entries are unknown. It must satisfy the following conditions in order to be completed:
diag(D) = 0
If is known,
If is unknown, so is
The graph of D must be connected. If D can be decomposed into two (or more) subgraphs, then the completion of D can be decomposed into two (or more) independent completion problems.
D |
The completed distance matrix with dimensionality d |
optval |
The minimum function value achieved during minimization (see details) |
Trosset, M.W. (2000). Distance Matrix Completion by Numerical Optimization.Computational Optimization and Applications, 17, 11–22, 2000.
set.seed(1337) D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) edmc(D, method="dpf", d=3, toler=1e-8)
set.seed(1337) D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) edmc(D, method="dpf", d=3, toler=1e-8)
edm2gram
Linear transformation of a Euclidean Distance Matrix to a Gram Matrix
edm2gram(D)
edm2gram(D)
D |
A Euclidean Distance Matrix |
While we specify that the input should be a Euclidean Distance Matrix (as this results in a Gram Matrix) the domain of edm2gram is the set of all real symmetric matrices. This function is particularly useful as it has the following property:
where is the space of symmetric, hollow matrices, negative definite on the space spanned by
and
is the space of centered positive definite matrices.
We can combine these two properties with a well known result: If D is a real symmetric matrix with 0 diagonal (call this matrix pre-EDM),
then D is a Euclidean Distance Matrix iff D is negative semi-definite on .
Using this result, combined with the properties of edm2gram we therefore have that
D is an EDM iff D is pre-EDM and is positive semi-definite.
G A Gram Matrix, where G = XX', and X is an nxp matrix containing the point configuration.
XY <- cbind(runif(100,0,1),runif(100,0,1)) D <- dist(XY) edm2gram(as.matrix(D))
XY <- cbind(runif(100,0,1),runif(100,0,1)) D <- dist(XY) edm2gram(as.matrix(D))
edm2psd
Convert an Euclidean Distance Matrix to a Positive Semi-definite Matrix
edm2psd(D, V = NULL)
edm2psd(D, V = NULL)
D |
A matrix in the set D_n^-. |
V |
A projection matrix satisfying V'1 = 0 and VV' = I |
For a matrix D in , edm2psd will be in the space of positive
semi-definite matrices. Therefore, if D also has zero diagonal, we have the following property:
D is a Euclidean Distance Matrix if and only if edm2psd is positive semi-definite.
This operator gives us another method to characterize the existence of a Euclidean distance matrix.
S A symmetric, positive semi-definite matrix
XY <- cbind(runif(100,0,1),runif(100,0,1)) D <- dist(XY) edm2psd(as.matrix(D))
XY <- cbind(runif(100,0,1),runif(100,0,1)) D <- dist(XY) edm2psd(as.matrix(D))
edmc
edmc(D, method = "dpf", ...)
edmc(D, method = "dpf", ...)
D |
An nxn partial-distance matrix to be completed, with unkown entries set to NA. |
method |
The algorithm to be used to complete the distance matrix D. One of sdp, npf, dpf, snl, or grs |
... |
The remaining input values required for the completion method specified in |
Depending on the method called, a number of input values are possible.
The return from edmc
depends on the method used. The help pages for each individual method
can be consulted for specific output.
set.seed(1337) D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) edmc(D,method = "dpf", d=3, toler=1e-8)
set.seed(1337) D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) edmc(D,method = "dpf", d=3, toler=1e-8)
getConfig
- given an nxn Euclidean distance matrix, produces a d-dimensional
point configuration of size n via eigendecomposition
getConfig(D, d)
getConfig(D, d)
D |
an nxn Euclidean distance matrix |
d |
the dimension for the configuration |
Given a distance matrix D, transform to a semi-definite matrix S using the linear transformation .
Using S, compute the eigen-decomposition
, where L is a diagonal matrix containing the singular-values of S,
and the columns of U contain the eigen-vectors. A point configuration X is then computed as:
To compute a configuration in d dimensions, the first d eigenvalues of S are used.
Y |
an nxd matrix containing the d-dimensional point configuration |
Accuracy |
the ratio of the sum of retained eigenvalues to the sum of all n eigenvalues obtained during decomposition |
set.seed(1337) D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) d <- 3 DStar <- dpf(D,d)$D getConfig(DStar,3)
set.seed(1337) D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) d <- 3 DStar <- dpf(D,d)$D getConfig(DStar,3)
gram2edm
Inverse Operator of edm2gram
gram2edm(B)
gram2edm(B)
B |
A centered, positive semi-definite matrix. |
The edm2gram function performs the following transformation:
where is the space of symmetric, hollow matrices, negative definite on the space spanned by
and
is the space of centered positive definite matrices.
The gram2edm function performs the inverse operation, taking a matrix in and transforming it to a matrix in
.
Therfore, gram2edm on is the inverse operator of edm2gram on
.
D A matrix in . If the input matrix B is a gram matrix, D is a Euclidean Distance Matrix.
X <- cbind(runif(100,0,1),runif(100,0,1)) G <- X %*% t(X) gram2edm(G)
X <- cbind(runif(100,0,1),runif(100,0,1)) G <- X %*% t(X) gram2edm(G)
grs
performs Euclidean Distance Matrix Completion using the guided random search algorithm
of Rahman & Oldford. Using this method will preserve the minimum spanning tree in the partial distance
matrix.
grs(D, d)
grs(D, d)
D |
An nxn partial-distance matrix to be completed. D must satisfy a list of conditions (see details), with unkown entries set to NA |
d |
The dimension for the resulting completion. |
The matrix D is a partial-distance matrix, meaning some of its entries are unknown. It must satisfy the following conditions in order to be completed:
diag(D) = 0
If is known,
If is unknown, so is
The graph of D must contain ONLY the minimum spanning tree distances
P |
The completed point configuration in dimension d |
D |
The completed Euclidean distance matrix |
Rahman, D., & Oldford, R.W. (2016). Euclidean Distance Matrix Completion and Point Configurations from the Minimal Spanning Tree.
#D matrix containing only the minimum spanning tree D <- matrix(c(0,3,NA,3,NA,NA, 3,0,1,NA,NA,NA, NA,1,0,NA,NA,NA, 3,NA,NA,0,1,NA, NA,NA,NA,1,0,1, NA,NA,NA,NA,1,0),byrow=TRUE, nrow=6) edmc(D, method="grs", d=3)
#D matrix containing only the minimum spanning tree D <- matrix(c(0,3,NA,3,NA,NA, 3,0,1,NA,NA,NA, NA,1,0,NA,NA,NA, 3,NA,NA,0,1,NA, NA,NA,NA,1,0,1, NA,NA,NA,NA,1,0),byrow=TRUE, nrow=6) edmc(D, method="grs", d=3)
mst
Compute a minimum spanning tree using Prim's algorithm
mst(D)
mst(D)
D |
A distance matrix |
MST a data frame object of 3 columns containing the parent nodes, child nodes, and corresponding weight of the MST edge
X <- runif(10,0,1) Y <- runif(10,0,1) D <- dist(cbind(X,Y)) mst(as.matrix(D))
X <- runif(10,0,1) Y <- runif(10,0,1) D <- dist(cbind(X,Y)) mst(as.matrix(D))
mstLB
Returns an nxn matrix containing the lower bounds for all unknown entries
in the partial distance matrix D such that the minimum spanning tree of the partial matrix D
is preserved upon completion.
mstLB(D)
mstLB(D)
D |
An nxn partial distance matrix to be completed |
The insight in constructing the lower bound is drawn from single-linkage clustering. Every edge in a spanning tree separates the vertices into two different groups, depending on which points remain connected to either one vertex or the other of that edge. Because the tree is a minimum spanning tree, if we select the largest edge, then the distance between any vertex of one group and any vertex of the other group must be at least as large as that of the the largest edge. This gives a lower bound for these distances that will preserve that edge in the minimum spanning tree. The same reasoning is applied recursively to each separate group, thus producing a lower bound on all edges.
The details of the algorithm can be found in Rahman & Oldford (2016).
Returns an nxn matrix containing the lower bound for the unknown entries in D
Rahman, D., & Oldford R.W. (2016). Euclidean Distance Matrix Completion and Point Configurations from the Minimal Spanning Tree.
D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) mstLB(D)
D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) mstLB(D)
mstUB
Compute the shortest path upper bound for all unknown entries in a partial distance matrix
mstUB(A)
mstUB(A)
A |
A (connected) partial distance matrix, with unknown entries set to Inf |
This function uses the shortest.paths() function, available in the igraph package.
UB A matrix containing the upper bounds for only the unknown entries. All other entries will be set to Inf.
A <- dist(cbind(rnorm(10,0,1),rnorm(10,0,1))) mstUB(as.matrix(A))
A <- dist(cbind(rnorm(10,0,1),rnorm(10,0,1))) mstUB(as.matrix(A))
npf
returns a completed Euclidean Distance Matrix D, with dimension d,
from a partial Euclidean Distance Matrix using the methods of Fang & O'Leary (2012)
npf( D, A = NA, d, dmax = (nrow(D) - 1), decreaseDim = 1, stretch = NULL, method = "Linear", toler = 1e-08 )
npf( D, A = NA, d, dmax = (nrow(D) - 1), decreaseDim = 1, stretch = NULL, method = "Linear", toler = 1e-08 )
D |
An nxn partial-distance matrix to be completed. D must satisfy a list of conditions (see details), with unkown entries set to NA. |
A |
a weight matrix, with |
d |
the dimension of the resulting completion |
dmax |
the maximum dimension to consider during dimension relaxation |
decreaseDim |
during dimension reduction, the number of dimensions to decrease each step |
stretch |
should the distance matrix be multiplied by a scalar constant? If no, stretch = NULL, otherwise stretch is a positive scalar |
method |
The method used for dimension reduction, one of "Linear" or "NLP". |
toler |
convergence tolerance for the algorithm |
This is an implementation of the Nonconvex Position Formulation (npf) for Euclidean Distance Matrix Completion, as proposed in 'Euclidean Distance Matrix Completion Problems' (Fang & O'Leary, 2012).
The method seeks to minimize the following:
where the function K() is that described in gram2edm, and the norm is Frobenius. Minimization is over X, the nxp matrix of node locations.
The matrix D is a partial-distance matrix, meaning some of its entries are unknown. It must satisfy the following conditions in order to be completed:
diag(D) = 0
If is known,
If is unknown, so is
The graph of D must be connected. If D can be decomposed into two (or more) subgraphs, then the completion of D can be decomposed into two (or more) independent completion problems.
D |
an nxn matrix of the completed Euclidean distances |
optval |
the minimum value achieved of the target function during minimization |
D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) A <- matrix(c(1,1,1,1,1,1, 1,1,1,0,1,0, 1,1,1,1,0,1, 1,0,1,1,1,0, 1,1,0,1,1,1, 1,0,1,0,1,1),byrow=TRUE, nrow=6) edmc(D, method="npf", d=3, dmax=5)
D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0),byrow=TRUE, nrow=6) A <- matrix(c(1,1,1,1,1,1, 1,1,1,0,1,0, 1,1,1,1,0,1, 1,0,1,1,1,0, 1,1,0,1,1,1, 1,0,1,0,1,1),byrow=TRUE, nrow=6) edmc(D, method="npf", d=3, dmax=5)
primPath
Given a starting node, creates the minimum spanning tree path through a point configuration.
primPath(A, start)
primPath(A, start)
A |
the distance matrix for which the minimum spanning tree path will be created |
start |
the starting node for the path |
Given a starting node, compute Prim's algorithm, resulting in the path taken to construct the minimum spanning tree.
return a 2x(n-1) matrix, where row 1 contains the parent nodes of the MST path, and row 2 contains the corresponding child nodes.
A <- dist(cbind(rnorm(100,0,1),rnorm(100,0,1))) primPath(as.matrix(A),1) primPath(as.matrix(A),2)
A <- dist(cbind(rnorm(100,0,1),rnorm(100,0,1))) primPath(as.matrix(A),1) primPath(as.matrix(A),2)
psd2edm
Transform a positive semi-definite matrix to a Euclidean Distance Matrix
psd2edm(S, V = NULL)
psd2edm(S, V = NULL)
S |
A symmetric, positive semi-definite matrix |
V |
A projection matrix satisfying V'1 = 0 and VV' = I |
The psd2edm function performs the inverse operation of the edm2psd function,
taking a matrix in and transforming it to a matrix in
.
Therefore, psd2edm on is the inverse operator of edm2psd on
.
For a symmetric positive semi-definite matrix S, psd2edm(S) will be in .
D A Euclidean Distance Matrix.
XY <- cbind(runif(100,0,1),runif(100,0,1)) S <- edm2psd(as.matrix(dist(XY))) D <- psd2edm(S)
XY <- cbind(runif(100,0,1),runif(100,0,1)) S <- edm2psd(as.matrix(dist(XY))) D <- psd2edm(S)
rgrs
Produce a point configuration given the edge lengths of the desired minimum spanning tree
rgrs( edges = NULL, d, n = NULL, theta = NULL, outlying = "N", skew = "N", stringy = "N" )
rgrs( edges = NULL, d, n = NULL, theta = NULL, outlying = "N", skew = "N", stringy = "N" )
edges |
A numeric vector containing the desired edge lengths of the minimum spanning tree. If n is specified, must be NULL. |
d |
the dimension of the resulting configuration. |
n |
the desired number of edge lengths to simulate. If edges is specified, must be set to NULL. |
theta |
Angle restriction during point proposal of the form (theta1,theta2,p), where p represents the probability of confining the proposal to [theta1,theta2]. Only used for d=2, otherwise NULL. See details for more in depth explanation. |
outlying |
One of "L", "M", or "H", specifying if the simulated edge lengths should have a Low, Medium, or High outlying scagnostic value. |
skew |
One of "L", "M", or "H", specifying if the simulated edge lengths should have a Low, Medium, or High skew scagnostic value. |
stringy |
One of "L", "M", or "H", specifying if the simulated edge lengths should have a Low, Medium, or High stringy scagnostic value. A numeric scalar specifying a value of stringy is also accepted. |
In 2-dimensions, when a new point is proposed, the position for the new point is determined by:
x <- x0 + r*sin(theta) y <- y0 + r*cos(theta)
where (x0,y0) is the base point, and r is the minimum spanning tree distance. theta is generated from a uniform distribution on (-pi,pi). By specifying the theta argument, the proposed theta is restricted, and is then generated from Uniform(theta1,theta2) or Uniform(-theta2,-theta1) with equal probability. This restriction allows the user to introduce striation into their point configuration.
An nxd matrix containing the d-dimensional locations of the points.
# An example where edge lengths are supplied EL <- runif(100,0,1) rgrs(edges = EL, d = 2) rgrs(edges = EL, d = 3) # An Example where edge lengths are simulated internally rgrs(d=2, n=100) rgrs(d=3, n=100) rgrs(d=2, n=100, outlying="H") rgrs(d=2, n=100, skew = "M") rgrs(d=2, n=100, stringy = "H") # An Example making use of theta rgrs(d=2, n=100, theta=c(pi/4,pi/3,.5))
# An example where edge lengths are supplied EL <- runif(100,0,1) rgrs(edges = EL, d = 2) rgrs(edges = EL, d = 3) # An Example where edge lengths are simulated internally rgrs(d=2, n=100) rgrs(d=3, n=100) rgrs(d=2, n=100, outlying="H") rgrs(d=2, n=100, skew = "M") rgrs(d=2, n=100, stringy = "H") # An Example making use of theta rgrs(d=2, n=100, theta=c(pi/4,pi/3,.5))
sdp
returns a completed Euclidean Distance Matrix D, with dimension d,
from a partial Euclidean Distance Matrix using the methods of Alfakih et. al. (1999)
sdp(D, A, toler = 1e-08)
sdp(D, A, toler = 1e-08)
D |
An nxn partial-distance matrix to be completed. D must satisfy a list of conditions (see details), with unkown entries set to NA. |
A |
a weight matrix, with |
toler |
convergence tolerance for the algorithm |
This is an implementation of the Semi-Definite Programming Algorithm (sdp) for Euclidean Distance Matrix Completion, as proposed in 'Solving Euclidean Distance Matrix Completion Problems via Semidefinite Programming' (Alfakih et. al., 1999).
The method seeks to minimize the following:
where the function psd2edm() is that described in psd2edm(), and the norm is Frobenius. Minimization is over S, a positive semidefinite matrix.
The matrix D is a partial-distance matrix, meaning some of its entries are unknown. It must satisfy the following conditions in order to be completed:
diag(D) = 0
If is known,
If is unknown, so is
The graph of D must be connected. If D can be decomposed into two (or more) subgraphs, then the completion of D can be decomposed into two (or more) independent completion problems.
D |
an nxn matrix of the completed Euclidean distances |
optval |
the minimum value achieved of the target function during minimization |
D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0), byrow=TRUE, nrow=6) A <- matrix(c(1,1,1,1,1,1, 1,1,1,0,1,0, 1,1,1,1,0,1, 1,0,1,1,1,0, 1,1,0,1,1,1, 1,0,1,0,1,1), byrow=TRUE, nrow=6) edmc(D, method="sdp", A=A, toler=1e-2)
D <- matrix(c(0,3,4,3,4,3, 3,0,1,NA,5,NA, 4,1,0,5,NA,5, 3,NA,5,0,1,NA, 4,5,NA,1,0,5, 3,NA,5,NA,5,0), byrow=TRUE, nrow=6) A <- matrix(c(1,1,1,1,1,1, 1,1,1,0,1,0, 1,1,1,1,0,1, 1,0,1,1,1,0, 1,1,0,1,1,1, 1,0,1,0,1,1), byrow=TRUE, nrow=6) edmc(D, method="sdp", A=A, toler=1e-2)
snl
solves the sensor network problem with
partial distance (squared) matrix D, and anchor positions anchors, in
dimension d.
snl(D, d, anchors = NULL)
snl(D, d, anchors = NULL)
D |
The partial distance matrix specifying the known distances between nodes. If anchors is specified (and is a pxr matrix), the p final columns and p final rows specify the distances between the anchors specified in anchors. |
d |
the dimension for the resulting completion |
anchors |
a pxr matrix specifying the d dimensional locations of the p anchors. If the anchorless problem is to be solved, anchors = NULL |
Set anchors=NULL to solve the anchorless (Euclidean distance matrix completion) problem in dimension d.
NOTE: When anchors is specified, the distances between the anchors must be in the bottom right corner of the matrix D, and anchors must have d columns.
X the d-dimensional positions of the localized sensors. Note that it may be the case that not all sensors could be localized, in which case X contains the positions of only the localized sensors.
Nathan Krislock and Henry Wolkowicz. Explicit sensor network localization using semidefinite representations and facial reductions. SIAM Journal on Optimization, 20(5):2679-2708, 2010.
D <- matrix(c(0,NA,.1987,NA,.0595,NA,.0159,.2251,.0036,.0875, NA,0,.0481,NA,NA,.0515,NA,.2079,.2230,NA, .1987,.0481,0,NA,NA,.1158,NA,NA,.1553,NA, NA,NA,NA,0,NA,NA,NA,.2319,NA,NA, .0595,NA,NA,NA,0,NA,.1087,.0894,.0589,.0159, NA,.0515,.1158,NA,NA,0,NA,NA,NA,NA, .0159,NA,NA,NA,.1087,NA,0,.3497,.0311,.1139, .2251,.2079,NA,.2319,.0894,NA,.3497,0,.1918,.1607, .0036,.2230,.1553,NA,.0589,NA,.0311,.1918,0,.1012, .0875,NA,NA,NA,.0159,NA,.1139,.1607,.1012,0),nrow=10, byrow=TRUE) anchors <- matrix(c(.5131,.9326, .3183,.3742, .5392,.7524, .2213,.7631), nrow=4,byrow=TRUE) d <- 2 #Anchorless Problem edmc(D, method="snl", d=2, anchors=NULL) #Anchored Problem edmc(D, method="snl", d=2, anchors=anchors)
D <- matrix(c(0,NA,.1987,NA,.0595,NA,.0159,.2251,.0036,.0875, NA,0,.0481,NA,NA,.0515,NA,.2079,.2230,NA, .1987,.0481,0,NA,NA,.1158,NA,NA,.1553,NA, NA,NA,NA,0,NA,NA,NA,.2319,NA,NA, .0595,NA,NA,NA,0,NA,.1087,.0894,.0589,.0159, NA,.0515,.1158,NA,NA,0,NA,NA,NA,NA, .0159,NA,NA,NA,.1087,NA,0,.3497,.0311,.1139, .2251,.2079,NA,.2319,.0894,NA,.3497,0,.1918,.1607, .0036,.2230,.1553,NA,.0589,NA,.0311,.1918,0,.1012, .0875,NA,NA,NA,.0159,NA,.1139,.1607,.1012,0),nrow=10, byrow=TRUE) anchors <- matrix(c(.5131,.9326, .3183,.3742, .5392,.7524, .2213,.7631), nrow=4,byrow=TRUE) d <- 2 #Anchorless Problem edmc(D, method="snl", d=2, anchors=NULL) #Anchored Problem edmc(D, method="snl", d=2, anchors=anchors)
sprosr
compute the three dimensional strucutre of a protein
molecule using its amino acid sequences using the semidefinite programming-based
protein structure determination (SPROS) method of Ramandi (2011)
sprosr( seq, aco, upl, hydrogen_omission = 1, f = c(10, 10, 10, 10, 10), in_max_res = NULL, in_min_res = NULL )
sprosr( seq, aco, upl, hydrogen_omission = 1, f = c(10, 10, 10, 10, 10), in_max_res = NULL, in_min_res = NULL )
seq |
A table containing the amino acid sequence of the protein in CYANA .seq format |
aco |
A table containing the angle constraint information in CYANA .aco format |
upl |
A table containing the distance constraint information in CYANA .upl format |
hydrogen_omission |
Should side-chain hydrogen atoms be omitted? TRUE/FALSE. Default is FALSE |
f |
Vector of length five detailing the multiplicative factors to be used. See details for more. |
in_max_res |
User overwrite of the maximum residue number. |
in_min_res |
User overwrite of the minimum residue number. |
The input files requires by sprosr follow the typical CYANA format. Each is a table with the following columns (no headers required).
Sequence File (seq)
column 1: amino acid residue name
column 2: residue number
Torsion Angle Restraint File (aco)
column 1: residue number (corresponding to seq file)
column 2: amino acid residue name
column 3: angle identifier, one of PHI or PSI
column 4: the lower limit of the angle specified in column 3
column 5: the upper limit of the angle specified in column 3
Distance Restraint File (upl)
column 1: residue number of the first atom (corresponding to seq file)
column 2: amino acid residue name of the first atom
column 3: atom name of the first atom
column 4: residue number of the second atom (corresponding to seq file)
column 5: amino acid residue name of the second atom
column 6: atom name of the second atom
column 7: upper distance limit (in Angstroms)
X |
Matrix containing the three dimensional point configuration of the protein structure. |
report |
A list containing the final violations of the protein |
Ramandi, Babak A., (2011). New Approaches to Protein NMR Automation. PhD Thesis. https://uwspace.uwaterloo.ca/bitstream/handle/10012/6389/Alipanahi_Ramandi_Babak.pdf;sequence=1
Demo Data - ACO
data(sprosr_aco)
data(sprosr_aco)
data.frame
Demo Data - SEQ
data(sprosr_seq)
data(sprosr_seq)
data.frame
Demo Data - UPL
data(sprosr_upl)
data(sprosr_upl)
data.frame