| Title: | Landmark Multi-Dimensional Scaling |
|---|---|
| Description: | A fast dimensionality reduction method scaleable to large numbers of samples. Landmark Multi-Dimensional Scaling (LMDS) is an extension of classical Torgerson MDS, but rather than calculating a complete distance matrix between all pairs of samples, only the distances between a set of landmarks and the samples are calculated. |
| Authors: | Robrecht Cannoodt [aut, cre] (ORCID: <https://orcid.org/0000-0003-3641-729X>, github: rcannood), Wouter Saelens [aut] (ORCID: <https://orcid.org/0000-0002-7114-6248>, github: zouter) |
| Maintainer: | Robrecht Cannoodt <[email protected]> |
| License: | GPL-3 |
| Version: | 0.1.1 |
| Built: | 2026-05-29 10:33:29 UTC |
| Source: | https://github.com/dynverse/lmds |
A fast dimensionality reduction method scaleable to large numbers of samples. Landmark Multi-Dimensional Scaling (LMDS) is an extension of classical Torgerson MDS, but rather than calculating a complete distance matrix between all pairs of samples, only the distances between a set of landmarks and the samples are calculated.
Maintainer: Robrecht Cannoodt [email protected] (ORCID) (rcannood)
Authors:
Wouter Saelens [email protected] (ORCID) (zouter)
Useful links:
Report bugs at https://github.com/dynverse/lmds/issues
Perform MDS on landmarks and project other samples to the same space
cmdscale_landmarks(dist_2lm, ndim = 3, rescale = TRUE, ...)cmdscale_landmarks(dist_2lm, ndim = 3, rescale = TRUE, ...)
dist_2lm |
Distance matrix between the landmarks and all the samples in original dataset |
ndim |
The number of dimensions |
rescale |
Whether or not to rescale the final dimensionality reduction (recommended) |
... |
Extra params to pass to |
The dimensionality reduction in the form of a ncol(dist_2lm) by ndim matrix.
library(Matrix) x <- as.matrix(iris[,1:4]) dist_2lm <- select_landmarks(x) cmdscale_landmarks(dist_2lm)library(Matrix) x <- as.matrix(iris[,1:4]) dist_2lm <- select_landmarks(x) cmdscale_landmarks(dist_2lm)
A fast dimensionality reduction method scaleable to large numbers of samples. Landmark Multi-Dimensional Scaling (LMDS) is an extension of classical 'Torgerson MDS', but rather than calculating a complete distance matrix between all pairs of samples, only the distances between a set of landmarks and the samples are calculated.
lmds( x, ndim = 3, distance_method = c("euclidean", "pearson", "spearman", "cosine", "chisquared", "hamming", "kullback", "manhattan", "maximum", "canberra", "minkowski"), landmark_method = c("sample"), num_landmarks = 500 )lmds( x, ndim = 3, distance_method = c("euclidean", "pearson", "spearman", "cosine", "chisquared", "hamming", "kullback", "manhattan", "maximum", "canberra", "minkowski"), landmark_method = c("sample"), num_landmarks = 500 )
x |
A matrix, optionally sparse. |
ndim |
The number of dimensions |
distance_method |
The distance metric to use. Options are "euclidean" (default), "pearson", "spearman", "cosine", "manhattan". |
landmark_method |
The landmark selection method to use. Options are "sample" (default). |
num_landmarks |
The number of landmarks to use, |
The dimensionality reduction in the form of a nrow(x) by ndim matrix.
library(Matrix) x <- Matrix::rsparsematrix(1000, 1000, .01) lmds(x, ndim = 3)library(Matrix) x <- Matrix::rsparsematrix(1000, 1000, .01) lmds(x, ndim = 3)
In addition, the distances between the landmarks and all samples are calculated.
select_landmarks( x, distance_method = c("euclidean", "pearson", "spearman", "cosine", "chisquared", "hamming", "kullback", "manhattan", "maximum", "canberra", "minkowski"), landmark_method = c("sample"), num_landmarks = 500 )select_landmarks( x, distance_method = c("euclidean", "pearson", "spearman", "cosine", "chisquared", "hamming", "kullback", "manhattan", "maximum", "canberra", "minkowski"), landmark_method = c("sample"), num_landmarks = 500 )
x |
A matrix, optionally sparse. |
distance_method |
The distance metric to use. Options are "euclidean" (default), "pearson", "spearman", "cosine", "manhattan". |
landmark_method |
The landmark selection method to use. Options are "sample" (default). |
num_landmarks |
The number of landmarks to use, |
The distance matrix between the landmarks and all samples. In addition, an attribute "landmark_ix"
denotes the indices of landmarks that were sampled.
library(Matrix) x <- Matrix::rsparsematrix(1000, 1000, .01) select_landmarks(x)library(Matrix) x <- Matrix::rsparsematrix(1000, 1000, .01) select_landmarks(x)