Last data update: 2014.03.03

R: Metric and non-metric Multidimensional Scaling for a distGPS...
mdsR Documentation

Metric and non-metric Multidimensional Scaling for a distGPS object.

Description

Generation of Multidimensional Scaling objects for the dissimilarities between elements given as an input in a distGPS object. Metric and non-metric algorithms are available, as well as an optimization algorithm for improving r-square correlation between observed and approximated distances. The MDS calculation for a given distance matrix can be splitted into smaller individual tasks and run in parallel, greatly improving CPU time and system memory usage. The S4 accessor functions getR2, getStress, getPoints retrieve R-square correlation, stress and points stored within a mds object respectively. The function is.adj is useful to know if a certain chroGPS MDS map has been adjusted by Procrustes or not (see help for procrustesAdj for details.)

Usage

mds(d, m = NULL, k = 2, type = "classic", add = FALSE, cor.method = "pearson", splitMDS = FALSE, split = 0.26, overlap = 0.025, stepSize=0.01, reshuffle = TRUE, set.seed = 149, mc.cores = 1, ...)
getR2(m)
getStress(m)
getPoints(m)

Arguments

d

Object of class distGPS with the pairwise observed dissimilarities between elements, a distance matrix.

m

(Optional). Object of class mds with a MDS object generated from the distances in d. Only MDS type "boostMDS" is available. The mds function performs an optimization of the approximated distances in m in order to improve r-square correlation between them and the observed dissimilarities en d, maximizing goodness of fit.

k

Dimensionality of the reconstructed space, typically set to 2 or 3.

type

Set to "classic" to perform classical MDS (uses function cmdscale from package stats). Set to "isoMDS" to use Kruskal's non-metric MDS (uses function isoMDS from package MASS) Set to "boostMDS" to perform r-square optimization of a pre-computed input MDS for that distance matrix.

add

Logical indicating if an additive constant c* should be computed, and added to the non-diagonal dissimilarities such that all n-1 eigenvalues are non-negative in cmdscale.

cor.method

A character string indicating which correlation coefficient (or covariance) is to be computed. One of "pearson" (default), "kendall", or "spearman", can be abbreviated.

splitMDS

Set to TRUE to perform computation of the MDS in parallel (see parameters below).

split

Proportion of elements to include in each (but last) distance submatrix.

overlap

Proportion of elements to be used as common anchor points between two adjacent distance submatrixes. These points will be used as spatial references to stitch each two adjacent MDS objects by Procrustes.

stepSize

Size for the quadratic search step to be used for R-square optimization if boostMDS is called, see specific help function for details.

reshuffle

Set to TRUE to perform random resampling of the input distance matrix before splitting it for parallel computation. This is often necessary to sufficiently capture the inherent variability of the data in each distance submatrix so that the stitching process can work properly, as the original data may present an arbitrary sorting of some kind. If a previous resampling of the data has been performed, this is not necessary.

set.seed

Random seed to perform the resampling.

mc.cores

Number of cores to be passed to the mclapply function from the parallel package, used to perform the parallel MDS computations.

...

Additional parameters passed to cmdscale, isoMDS or boostMDS, see each individual help file for details.

Value

The function returns a mds object. See help ("mds-Class") for details.

Methods

mds

signature(d = "distGPS", m = "missing"): Creates a mds object with points in a k-dimensional space approximating the pairwise distances in d.

mds

signature(d = "distGPS", m = "mds"): For the observed dissimilarities in d and a valid spatial representation of them in m, the function returns a mds object with an optimized representation of d in terms of R-square. The MDS stress measure is also returned. See help for boostMDS for details.

plot

signature(m = "mds"): S4 plot method for mds objects.

Author(s)

Oscar Reina

See Also

See functions cmdscale, isoMDS from package MASS.

Examples

x <- rbind(c(rep(0,15),rep(1,5)),c(rep(0,15),rep(1,5)),c(rep(0,19),1),c(rep(1,5),rep(0,15)))
rownames(x) <- letters[1:4]
d <- distGPS(x,metric='tanimoto',uniqueRows=TRUE)
mds1 <- mds(d)
mds1
plot(mds1)
#gps2xgmml(mds1, fname='chroGPS_factors.xgmml', fontSize=4,col=col2hex('red'), cex=8)

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(chroGPS)
Loading required package: IRanges
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: Biobase
Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

Loading required package: MASS
Loading required package: changepoint
Loading required package: zoo

Attaching package: 'zoo'

The following objects are masked from 'package:base':

    as.Date, as.Date.numeric

Successfully loaded changepoint package version 2.2.1
 NOTE: Predefined penalty values have changed.  Previous penalty values with a postfix 1 i.e. SIC1 are now without i.e. SIC and previous penalties without a postfix i.e. SIC are now with a postfix 0 i.e. SIC0. See NEWS and help files for further details.
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/chroGPS/mds.Rd_%03d_medium.png", width=480, height=480)
> ### Name: mds
> ### Title: Metric and non-metric Multidimensional Scaling for a distGPS
> ###   object.
> ### Aliases: mds mds-methods mds,distGPS,missing-method
> ###   mds,splitDistGPS,missing-method plot,mds,ANY-method getR2 getStress
> ###   getPoints is.adj getR2,mds-method getStress,mds-method
> ###   getPoints,mds-method is.adj,mds-method
> ### Keywords: graphs mds
> 
> ### ** Examples
> 
> x <- rbind(c(rep(0,15),rep(1,5)),c(rep(0,15),rep(1,5)),c(rep(0,19),1),c(rep(1,5),rep(0,15)))
> rownames(x) <- letters[1:4]
> d <- distGPS(x,metric='tanimoto',uniqueRows=TRUE)
> mds1 <- mds(d)
> mds1
Object of class MDS approximating distances between 3 objects 
R-squared= 1 Stress= 0 
> plot(mds1)
> #gps2xgmml(mds1, fname='chroGPS_factors.xgmml', fontSize=4,col=col2hex('red'), cex=8)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>