Title: | Creates diffusion maps |
---|---|
Description: | Create and plot diffusion maps. |
Authors: | Philipp Angerer [cre, aut] , Laleh Haghverdi [ctb], Maren Büttner [ctb] , Fabian Theis [ctb] , Carsten Marr [ctb] , Florian Büttner [ctb] |
Maintainer: | Philipp Angerer <[email protected]> |
License: | GPL-3 |
Version: | 3.21.0 |
Built: | 2024-11-19 03:14:15 UTC |
Source: | https://github.com/bioc/destiny |
Convert a DiffusionMap
or DPT
object to other classes
## S4 method for signature 'DiffusionMap' as.data.frame(x, row.names = NULL, optional = FALSE, ...) fortify.DiffusionMap(model, data, ...) ## S4 method for signature 'DPT' as.data.frame(x, row.names = NULL, optional = FALSE, ...) fortify.DPT(model, data, ...) ## S4 method for signature 'DPT' as.matrix(x, ...)
## S4 method for signature 'DiffusionMap' as.data.frame(x, row.names = NULL, optional = FALSE, ...) fortify.DiffusionMap(model, data, ...) ## S4 method for signature 'DPT' as.data.frame(x, row.names = NULL, optional = FALSE, ...) fortify.DPT(model, data, ...) ## S4 method for signature 'DPT' as.matrix(x, ...)
x , model
|
A |
row.names |
NULL or a character vector giving the row names for the data frame. Missing values are not allowed. |
optional |
logical. If TRUE, setting row names and converting column names (to syntactic names: see make.names) is optional. |
... |
Passed to |
data |
ignored |
fortify is a ggplot2 generic allowing a diffusion map to be used as data
parameter in ggplot or qplot.
An object of the desired class
DiffusionMap accession methods, Extraction methods, DiffusionMap methods for more
library(Biobase) data(guo) dm <- DiffusionMap(guo) classes <- vapply(as.data.frame(dm), class, character(1L)) stopifnot(all(classes[paste0('DC', 1:20)] == 'numeric')) stopifnot(all(classes[featureNames(guo) ] == 'numeric')) stopifnot(all(classes[ varLabels(guo) ] == c('factor', 'integer')))
library(Biobase) data(guo) dm <- DiffusionMap(guo) classes <- vapply(as.data.frame(dm), class, character(1L)) stopifnot(all(classes[paste0('DC', 1:20)] == 'numeric')) stopifnot(all(classes[featureNames(guo) ] == 'numeric')) stopifnot(all(classes[ varLabels(guo) ] == c('factor', 'integer')))
Creates a color legend for a vector used to color a plot. It will use the current palette()
or the specified pal
as reference.
colorlegend( col, pal = palette(), log = FALSE, posx = c(0.9, 0.93), posy = c(0.05, 0.9), main = NULL, cex_main = par("cex.sub"), cex_axis = par("cex.axis"), col_main = par("col.sub"), col_lab = par("col.lab"), steps = 5, steps_color = 100, digit = 2, left = FALSE, ..., cex.main = NULL, cex.axis = NULL, col.main = NULL, col.lab = NULL )
colorlegend( col, pal = palette(), log = FALSE, posx = c(0.9, 0.93), posy = c(0.05, 0.9), main = NULL, cex_main = par("cex.sub"), cex_axis = par("cex.axis"), col_main = par("col.sub"), col_lab = par("col.lab"), steps = 5, steps_color = 100, digit = 2, left = FALSE, ..., cex.main = NULL, cex.axis = NULL, col.main = NULL, col.lab = NULL )
col |
Vector of factor, integer, or double used to determine the ticks. |
pal |
If |
log |
Use logarithmic scale? |
posx |
Left and right borders of the color bar relative to plot area (Vector of length 2; 0-1) |
posy |
Bottom and top borders of color bar relative to plot area (Vector of length 2; 0-1) |
main |
Legend title |
cex_main |
Size of legend title font (default: subtitle font size |
cex_axis |
Size of ticks/category labels (default: axis font size |
col_main |
Color of legend title (default: subtitle color |
col_lab |
Color of tick or category labels (default: axis color |
steps |
Number of labels in case of a continuous axis. If 0 or FALSE, draw no ticks |
steps_color |
Number of gradient samples in case of continuous axis |
digit |
Number of digits for continuous axis labels |
left |
logical. If TRUE, invert posx |
... |
Additional parameters for the text call used for labels |
cex.main , cex.axis , col.main , col.lab
|
For compatibility with |
When passed a factor or integer vector, it will create a discrete legend, whereas a double vector will result in a continuous bar.
This function is called for the side effect of adding a colorbar to a plot and returns nothing/NULL.
color_data <- 1:6 par(mar = par('mar') + c(0, 0, 0, 3)) plot(sample(6), col = color_data) colorlegend(color_data)
color_data <- 1:6 par(mar = par('mar') + c(0, 0, 0, 3)) plot(sample(6), col = color_data) colorlegend(color_data)
Creates a perceptually monotonously decreasing (or increasing) lightness color palette with different tones.
This was necessary in pre-viridis times, by now you can probably just use hcl.colors
cube_helix( n = 6, start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE ) scale_colour_cube_helix( ..., start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE, discrete = TRUE, guide = if (discrete) "legend" else "colourbar" ) scale_color_cube_helix( ..., start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE, discrete = TRUE, guide = if (discrete) "legend" else "colourbar" ) scale_fill_cube_helix( ..., start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE, discrete = TRUE, guide = if (discrete) "legend" else "colourbar" )
cube_helix( n = 6, start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE ) scale_colour_cube_helix( ..., start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE, discrete = TRUE, guide = if (discrete) "legend" else "colourbar" ) scale_color_cube_helix( ..., start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE, discrete = TRUE, guide = if (discrete) "legend" else "colourbar" ) scale_fill_cube_helix( ..., start = 0, r = 0.4, hue = 0.8, gamma = 1, light = 0.85, dark = 0.15, reverse = FALSE, discrete = TRUE, guide = if (discrete) "legend" else "colourbar" )
n |
Number of colors to return (default: 6) |
start |
Hue to start helix at ( |
r |
Number of rotations of the helix. Can be negative. (default: 0.4) |
hue |
Saturation. 0 means greyscale, 1 fully saturated colors (default: 0.8) |
gamma |
Emphasize darker (gamma < 1) or lighter (gamma > 1) colors (default: 1) |
light |
Lightest lightness (default: 0.85) |
dark |
Darkest lightness (default: 0.15) |
reverse |
logical. If TRUE, reverse lightness (default: FALSE) |
... |
parameters passed to |
discrete |
If TRUE, return a discrete scale, if FALSE a continuous one (default: TRUE) |
guide |
Type of scale guide to use. See |
A character
vector of hex colors with length n
palette(cube_helix()) image(matrix(1:6), col = 1:6, pch = 19, axes = FALSE) cr <- scales::colour_ramp(cube_helix(12, r = 3)) r <- runif(100) plot(1:100, r, col = cr(r), type = 'b', pch = 20)
palette(cube_helix()) image(matrix(1:6), col = 1:6, pch = 19, axes = FALSE) cr <- scales::colour_ramp(cube_helix(12, r = 3)) r <- runif(100) plot(1:100, r, col = cr(r), type = 'b', pch = 20)
The main function is DiffusionMap
, which returns an object you can plot
(plot.DiffusionMap
is then called).
Maintainer: Philipp Angerer [email protected] (ORCID)
Other contributors:
Laleh Haghverdi [email protected] [contributor]
Maren Büttner [email protected] (ORCID) [contributor]
Fabian Theis [email protected] (ORCID) [contributor]
Carsten Marr [email protected] (ORCID) [contributor]
Florian Büttner [email protected] (ORCID) [contributor]
Useful links:
Report bugs at https://github.com/theislab/destiny/issues
demo(destiny, ask = FALSE)
demo(destiny, ask = FALSE)
destiny provides several generic methods and implements them for the DiffusionMap
and Sigmas
classes.
eigenvalues(object) eigenvalues(object) <- value eigenvectors(object) eigenvectors(object) <- value sigmas(object) sigmas(object) <- value dataset(object) dataset(object) <- value distance(object) distance(object) <- value optimal_sigma(object)
eigenvalues(object) eigenvalues(object) <- value eigenvectors(object) eigenvectors(object) <- value sigmas(object) sigmas(object) <- value dataset(object) dataset(object) <- value distance(object) distance(object) <- value optimal_sigma(object)
object |
Object from which to extract or to which to assign a value |
value |
Value to assign within an object |
eigenvalues
retrieves the numeric eigenvalues
eigenvectors
retrieves the eigenvectors matrix
sigmas
retrieves the Sigmas
from an object utilizing it as kernel width
dataset
retrieves the data the object was created from
distance
retrieves the distance metric used to create the object, e.g. euclidean
optimal_sigma
retrieves the numeric value of the optimal sigma or local sigmas
DiffusionMap methods and Sigmas class for implementations
data(guo_norm) dm <- DiffusionMap(guo_norm) eigenvalues(dm) eigenvectors(dm) sigmas(dm) optimal_sigma(dm) dataset(dm) distance(dm)
data(guo_norm) dm <- DiffusionMap(guo_norm) eigenvalues(dm) eigenvectors(dm) sigmas(dm) optimal_sigma(dm) dataset(dm) distance(dm)
Get and set eigenvalues, eigenvectors, and sigma(s) of a DiffusionMap object.
## S4 method for signature 'DiffusionMap' eigenvalues(object) ## S4 replacement method for signature 'DiffusionMap' eigenvalues(object) <- value ## S4 method for signature 'DiffusionMap' eigenvectors(object) ## S4 replacement method for signature 'DiffusionMap' eigenvectors(object) <- value ## S4 method for signature 'DiffusionMap' sigmas(object) ## S4 replacement method for signature 'DiffusionMap' sigmas(object) <- value ## S4 method for signature 'DiffusionMap' dataset(object) ## S4 replacement method for signature 'DiffusionMap' dataset(object) <- value ## S4 method for signature 'DiffusionMap' distance(object) ## S4 replacement method for signature 'DiffusionMap' distance(object) <- value ## S4 method for signature 'DiffusionMap' optimal_sigma(object)
## S4 method for signature 'DiffusionMap' eigenvalues(object) ## S4 replacement method for signature 'DiffusionMap' eigenvalues(object) <- value ## S4 method for signature 'DiffusionMap' eigenvectors(object) ## S4 replacement method for signature 'DiffusionMap' eigenvectors(object) <- value ## S4 method for signature 'DiffusionMap' sigmas(object) ## S4 replacement method for signature 'DiffusionMap' sigmas(object) <- value ## S4 method for signature 'DiffusionMap' dataset(object) ## S4 replacement method for signature 'DiffusionMap' dataset(object) <- value ## S4 method for signature 'DiffusionMap' distance(object) ## S4 replacement method for signature 'DiffusionMap' distance(object) <- value ## S4 method for signature 'DiffusionMap' optimal_sigma(object)
object |
A DiffusionMap |
value |
Vector of eigenvalues or matrix of eigenvectors to get/set |
The assigned or retrieved value
Extraction methods, DiffusionMap methods, Coercion methods for more
data(guo) dm <- DiffusionMap(guo) eigenvalues(dm) eigenvectors(dm) sigmas(dm) dataset(dm) optimal_sigma(dm)
data(guo) dm <- DiffusionMap(guo) eigenvalues(dm) eigenvectors(dm) sigmas(dm) dataset(dm) optimal_sigma(dm)
Methods for external operations on diffusion maps
## S4 method for signature 'DiffusionMap' print(x) ## S4 method for signature 'DiffusionMap' show(object)
## S4 method for signature 'DiffusionMap' print(x) ## S4 method for signature 'DiffusionMap' show(object)
x , object
|
The DiffusionMap
object (print
), or NULL (show
), invisibly
DiffusionMap accession methods, Extraction methods, Coercion methods for more
data(guo) dm <- DiffusionMap(guo) print(dm) show(dm)
data(guo) dm <- DiffusionMap(guo) print(dm) show(dm)
The provided data can be a double matrix of expression data or a data.frame with all non-integer (double) columns being treated as expression data features (and the others ignored), an ExpressionSet, or a SingleCellExperiment.
DiffusionMap( data = stopifnot_distmatrix(distance), sigma = "local", k = find_dm_k(dataset_n_observations(data, distance) - 1L), n_eigs = min(20L, dataset_n_observations(data, distance) - 2L), density_norm = TRUE, ..., distance = c("euclidean", "cosine", "rankcor", "l2"), n_pcs = NULL, n_local = seq(to = min(k, 7L), length.out = min(k, 3L)), rotate = FALSE, censor_val = NULL, censor_range = NULL, missing_range = NULL, vars = NULL, knn_params = list(), verbose = !is.null(censor_range), suppress_dpt = FALSE )
DiffusionMap( data = stopifnot_distmatrix(distance), sigma = "local", k = find_dm_k(dataset_n_observations(data, distance) - 1L), n_eigs = min(20L, dataset_n_observations(data, distance) - 2L), density_norm = TRUE, ..., distance = c("euclidean", "cosine", "rankcor", "l2"), n_pcs = NULL, n_local = seq(to = min(k, 7L), length.out = min(k, 3L)), rotate = FALSE, censor_val = NULL, censor_range = NULL, missing_range = NULL, vars = NULL, knn_params = list(), verbose = !is.null(censor_range), suppress_dpt = FALSE )
data |
Expression data to be analyzed and covariates. Provide |
sigma |
Diffusion scale parameter of the Gaussian kernel. One of |
k |
Number of nearest neighbors to consider (default: a guess betweeen 100 and |
n_eigs |
Number of eigenvectors/values to return (default: 20) |
density_norm |
logical. If TRUE, use density normalisation |
... |
Unused. All parameters to the right of the |
distance |
Distance measurement method applied to |
n_pcs |
Number of principal components to compute to base calculations on. Using e.g. 50 DCs results in more regular looking diffusion maps.
The default NULL will not compute principal components, but use |
n_local |
If |
rotate |
logical. If TRUE, rotate the eigenvalues to get a slimmer diffusion map |
censor_val |
Value regarded as uncertain. Either a single value or one for every dimension (Optional, default: censor_val) |
censor_range |
Uncertainity range for censoring (Optional, default: none). A length-2-vector of certainty range start and end. TODO: also allow |
missing_range |
Whole data range for missing value model. Has to be specified if NAs are in the data |
vars |
Variables (columns) of the data to use. Specifying NULL will select all columns (default: All floating point value columns) |
knn_params |
Parameters passed to |
verbose |
Show a progressbar and other progress information (default: do it if censoring is enabled) |
suppress_dpt |
Specify TRUE to skip calculation of necessary (but spacious) information for |
A DiffusionMap object:
eigenvalues
Eigenvalues ranking the eigenvectors
eigenvectors
Eigenvectors mapping the datapoints to n_eigs
dimensions
sigmas
Sigmas object with either information about the find_sigmas heuristic run or just local or optimal_sigma.
data_env
Environment referencing the data used to create the diffusion map
eigenvec0
First (constant) eigenvector not included as diffusion component.
transitions
Transition probabilities. Can be NULL
d
Density vector of transition probability matrix
d_norm
Density vector of normalized transition probability matrix
k
The k parameter for kNN
n_pcs
Number of principal components used in kNN computation (NA if raw data was used)
n_local
The n_local
th nearest neighbor(s) is/are used to determine local kernel density
density_norm
Was density normalization used?
rotate
Were the eigenvectors rotated?
distance
Distance measurement method used
censor_val
Censoring value
censor_range
Censoring range
missing_range
Whole data range for missing value model
vars
Vars parameter used to extract the part of the data used for diffusion map creation
knn_params
Parameters passed to find_knn
DiffusionMap methods to get and set the slots. find_sigmas
to pre-calculate a fitting global sigma
parameter
data(guo) DiffusionMap(guo) DiffusionMap(guo, 13, censor_val = 15, censor_range = c(15, 40), verbose = TRUE) covars <- data.frame(covar1 = letters[1:100]) dists <- dist(matrix(rnorm(100*10), 100)) DiffusionMap(covars, distance = dists)
data(guo) DiffusionMap(guo) DiffusionMap(guo, 13, censor_val = 15, censor_range = c(15, 40), verbose = TRUE) covars <- data.frame(covar1 = letters[1:100]) dists <- dist(matrix(rnorm(100*10), 100)) DiffusionMap(covars, distance = dists)
Predict new data points using an existing DiffusionMap. The resulting matrix can be used in the plot method for the DiffusionMap
dm_predict(dm, new_data, ..., verbose = FALSE)
dm_predict(dm, new_data, ..., verbose = FALSE)
dm |
A |
new_data |
New data points to project into the diffusion map. Can be a matrix, data.frame, ExpressionSet, or SingleCellExperiment. |
... |
Passed to |
verbose |
Show progress messages? |
A matrix of projected diffusion components for the new data.
data(guo) g1 <- guo[, guo$num_cells != 32L] g2 <- guo[, guo$num_cells == 32L] dm <- DiffusionMap(g1) dc2 <- dm_predict(dm, g2) plot(dm, new_dcs = dc2)
data(guo) g1 <- guo[, guo$num_cells != 32L] g2 <- guo[, guo$num_cells == 32L] dm <- DiffusionMap(g1) dc2 <- dm_predict(dm, g2) plot(dm, new_dcs = dc2)
Treat DPT object as a matrix of cell-by-cell DPT distances.
## S4 method for signature 'DPT,index,index,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,index,missing,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,missing,index,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,missing,missing,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,index,index' x[[i, j, ...]] ## S4 method for signature 'DPT' nrow(x) ## S4 method for signature 'DPT' ncol(x) ## S4 method for signature 'DPT' dim(x)
## S4 method for signature 'DPT,index,index,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,index,missing,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,missing,index,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,missing,missing,logicalOrMissing' x[i, j, ..., drop = TRUE] ## S4 method for signature 'DPT,index,index' x[[i, j, ...]] ## S4 method for signature 'DPT' nrow(x) ## S4 method for signature 'DPT' ncol(x) ## S4 method for signature 'DPT' dim(x)
x |
|
i , j
|
|
... |
ignored |
drop |
If |
[
returns a dense matrix or (if applicable and isTRUE(drop)
) a vector.
[[
returns single distance value
nrow
and ncol
return the number of cells
dim
returns c(n_cells, n_cells)
data(guo_norm) dm <- DiffusionMap(guo_norm) dpt <- DPT(dm) set.seed(1) plot(dpt[random_root(dpt), ], Biobase::exprs(guo_norm)['DppaI', ])
data(guo_norm) dm <- DiffusionMap(guo_norm) dpt <- DPT(dm) set.seed(1) plot(dpt[random_root(dpt), ], Biobase::exprs(guo_norm)['DppaI', ])
Methods for the DPT class. branch_divide
subdivides branches for plotting (see the examples).
branch_divide(dpt, divide = integer(0L)) tips(dpt) ## S4 method for signature 'DPT' dataset(object) ## S4 replacement method for signature 'DPT' dataset(object) <- value
branch_divide(dpt, divide = integer(0L)) tips(dpt) ## S4 method for signature 'DPT' dataset(object) ## S4 replacement method for signature 'DPT' dataset(object) <- value
dpt , object
|
DPT object |
divide |
Vector of branch numbers to use for division |
value |
Value of slot to set |
branch_divide
and dataset<-
return the changed object, dataset
the extracted data, and tips
the tip indices.
plot.DPT uses branch_divide
for its divide
argument.
data(guo_norm) dpt <- DPT(DiffusionMap(guo_norm)) dpt_9_branches <- branch_divide(dpt, 1:3) plot(dpt_9_branches, col_by = 'branch')
data(guo_norm) dpt <- DPT(DiffusionMap(guo_norm)) dpt_9_branches <- branch_divide(dpt, 1:3) plot(dpt_9_branches, col_by = 'branch')
Create pseudotime ordering and assigns cell to one of three branches
DPT(dm, tips = random_root(dm), ..., w_width = 0.1)
DPT(dm, tips = random_root(dm), ..., w_width = 0.1)
dm |
A |
tips |
The cell index/indices from which to calculate the DPT(s) (integer of length 1-3) |
... |
Unused. All parameters to the right of the |
w_width |
Window width to use for deciding the branch cutoff |
Treat it as a matrix of pseudotime by subsetting ([ dim nrow ncol
as.matrix
), and as a list of pseudodime, and expression vectors ($ [[ names
as.data.frame
).
A DPT
object:
branch
matrix
(of integer
) recursive branch labels for each cell (row); NA
for undeceided. Use branch_divide
to modify this.
tips
matrix
(of logical
) indicating if a cell (row) is a tip of the corresponding banch level (col)
dm
DiffusionMap
used to create this DPT object
data(guo_norm) dm <- DiffusionMap(guo_norm) dpt <- DPT(dm) str(dpt)
data(guo_norm) dm <- DiffusionMap(guo_norm) dpt <- DPT(dm) str(dpt)
eigs
By default uses a random initialization vector that you can make deterministic using
set.seed
or override by specifying opts = list(initvec = ...)
.
eig_decomp(m, n_eigs, sym, ..., opts = list())
eig_decomp(m, n_eigs, sym, ..., opts = list())
m |
A matrix (e.g. from the Matrix package) or
a function (see |
n_eigs |
Number of eigenvectors to return. |
sym |
defunct and ignored. |
... |
Passed to |
opts |
Passed to |
see eigs
.
eig_decomp(cbind(c(1,0,-1), c(0,1,0), c(-1,0,1)), 2)
eig_decomp(cbind(c(1,0,-1), c(0,1,0), c(-1,0,1)), 2)
These functions present quick way to create ExpressionSet objects.
as.ExpressionSet(x, ...) ## S4 method for signature 'data.frame' as.ExpressionSet(x, annotation_cols = !sapply(x, is.double)) read.ExpressionSet(file, header = TRUE, ...)
as.ExpressionSet(x, ...) ## S4 method for signature 'data.frame' as.ExpressionSet(x, annotation_cols = !sapply(x, is.double)) read.ExpressionSet(file, header = TRUE, ...)
x |
data.frame to convert to an ExpressionSet. |
... |
Additional parameters to read.table |
annotation_cols |
The data.frame columns used as annotations. All others are used as expressions. (Logical, character or numerical index array) |
file |
File path to read ASCII data from |
header |
Specifies if the file has a header row. |
They work by using all continuous (double) columns as expression data, and all others as observation annotations.
an ExpressionSet object
read.table on which read.ExpressionSet
is based, and ExpressionSet.
library(Biobase) df <- data.frame(Time = seq_len(3), #integer column Actb = c(0.05, 0.3, 0.8), Gapdh = c(0.2, 0.03, 0.1)) set <- as.ExpressionSet(df) rownames(exprs(set)) == c('Actb', 'Gapdh') phenoData(set)$Time == 1:3
library(Biobase) df <- data.frame(Time = seq_len(3), #integer column Actb = c(0.05, 0.3, 0.8), Gapdh = c(0.2, 0.03, 0.1)) set <- as.ExpressionSet(df) rownames(exprs(set)) == c('Actb', 'Gapdh') phenoData(set)$Time == 1:3
Extract common information from objects.
Apart from the input data's branches,
you can extract diffusion components via $DCx
.
From DPT
objects, you can also extract the branch label via $Branch
,
or the diffusion pseudo time for a numbered cell via $DPTx
.
## S4 method for signature 'DiffusionMap' names(x) ## S4 method for signature 'DPT' names(x) ## S4 method for signature 'DiffusionMap,character,missing' x[[i, j, ...]] ## S4 method for signature 'DPT,character,missing' x[[i, j, ...]] ## S4 method for signature 'DiffusionMap' x$name ## S4 method for signature 'DPT' x$name
## S4 method for signature 'DiffusionMap' names(x) ## S4 method for signature 'DPT' names(x) ## S4 method for signature 'DiffusionMap,character,missing' x[[i, j, ...]] ## S4 method for signature 'DPT,character,missing' x[[i, j, ...]] ## S4 method for signature 'DiffusionMap' x$name ## S4 method for signature 'DPT' x$name
x |
|
i , name
|
Name of a diffusion component |
j |
N/A |
... |
ignored |
The names or data row, see respective generics.
Extract, names
for the generics. DiffusionMap accession methods, DiffusionMap methods, Coercion methods for more
data(guo) dm <- DiffusionMap(guo) dm$DC1 # A diffusion component dm$Actb # A gene expression vector dm$num_cells # Phenotype metadata dpt <- DPT(dm) dm$Branch dm$DPT1
data(guo) dm <- DiffusionMap(guo) dm$DC1 # A diffusion component dm$Actb # A gene expression vector dm$num_cells # Phenotype metadata dpt <- DPT(dm) dm$Branch dm$DPT1
The k
parameter for the k nearest neighbors used in DiffusionMap should be as big as possible while
still being computationally feasible. This function approximates it depending on the size of the dataset n
.
find_dm_k(n, min_k = 100L, small = 1000L, big = 10000L)
find_dm_k(n, min_k = 100L, small = 1000L, big = 10000L)
n |
Number of possible neighbors (nrow(dataset) - 1) |
min_k |
Minimum number of neighbors. Will be chosen for |
small |
Number of neighbors considered small. If/where |
big |
Number of neighbors considered big. If/where |
A vector of the same length as n
that contains suitable k
values for the respective n
curve(find_dm_k(n), 0, 13000, xname = 'n') curve(find_dm_k(n) / n, 0, 13000, xname = 'n')
curve(find_dm_k(n), 0, 13000, xname = 'n') curve(find_dm_k(n) / n, 0, 13000, xname = 'n')
Approximate k nearest neighbor search with flexible distance function.
find_knn( data, k, ..., query = NULL, distance = c("euclidean", "cosine", "rankcor", "l2"), method = c("covertree", "hnsw"), sym = TRUE, verbose = FALSE )
find_knn( data, k, ..., query = NULL, distance = c("euclidean", "cosine", "rankcor", "l2"), method = c("covertree", "hnsw"), sym = TRUE, verbose = FALSE )
data |
Data matrix |
k |
Number of nearest neighbors |
... |
Parameters passed to |
query |
Query matrix. Leave it out to use |
distance |
Distance metric to use. Allowed measures: Euclidean distance (default), cosine distance ( |
method |
Method to use. |
sym |
Return a symmetric matrix (as long as query is NULL)? |
verbose |
Show a progressbar? (default: FALSE) |
A list
with the entries:
index
A integer matrix containing the indices of the k nearest neighbors for each cell.
dist
A double matrix containing the distances to the k nearest neighbors for each cell.
dist_mat
A dgCMatrix
if sym == TRUE
,
else a dsCMatrix
().
Any zero in the matrix (except for the diagonal) indicates that the cells in the corresponding pair are close neighbors.
).The sigma with the maximum value in average dimensionality is close to the ideal one. Increasing step number gets this nearer to the ideal one.
find_sigmas( data, step_size = 0.1, steps = 10L, start = NULL, sample_rows = 500L, early_exit = FALSE, ..., censor_val = NULL, censor_range = NULL, missing_range = NULL, vars = NULL, verbose = TRUE )
find_sigmas( data, step_size = 0.1, steps = 10L, start = NULL, sample_rows = 500L, early_exit = FALSE, ..., censor_val = NULL, censor_range = NULL, missing_range = NULL, vars = NULL, verbose = TRUE )
data |
Data set with |
step_size |
Size of log-sigma steps |
steps |
Number of steps/calculations |
start |
Initial value to search from. (Optional. default: |
sample_rows |
Number of random rows to use for sigma estimation or vector of row indices/names to use. In the first case, only used if actually smaller than the number of available rows (Optional. default: 500) |
early_exit |
logical. If TRUE, return if the first local maximum is found, else keep running |
... |
Unused. All parameters to the right of the |
censor_val |
Value regarded as uncertain. Either a single value or one for every dimension |
censor_range |
Uncertainity range for censoring. A length-2-vector of certainty range start and end. TODO: also allow |
missing_range |
Whole data range for missing value model. Has to be specified if NAs are in the data |
vars |
Variables (columns) of the data to use. Specifying TRUE will select all columns (default: All floating point value columns) |
verbose |
logical. If TRUE, show a progress bar and plot the output |
Object of class Sigmas
Sigmas
, the class returned by this; DiffusionMap
, the class this is used for
data(guo) sigs <- find_sigmas(guo, verbose = TRUE) DiffusionMap(guo, sigs)
data(guo) sigs <- find_sigmas(guo, verbose = TRUE) DiffusionMap(guo, sigs)
Find tips in a DiffusionMap object
find_tips(dm_or_dpt, root = random_root(dm_or_dpt))
find_tips(dm_or_dpt, root = random_root(dm_or_dpt))
dm_or_dpt |
A |
root |
Root cell index from which to find tips. (default: random) |
An integer vector of length 3
data(guo) dm <- DiffusionMap(guo) is_tip <- l_which(find_tips(dm), len = ncol(guo)) plot(dm, col = factor(is_tip))
data(guo) dm <- DiffusionMap(guo) is_tip <- l_which(find_tips(dm), len = ncol(guo)) plot(dm, col = factor(is_tip))
featureNames <- ...
can be used to set the gene names used for plotting
(e.g. if the data contains hardly readably gene or transcript IDs).
dataset
gets the expressions used for the gene relevance calculations,
and distance
the distance measure.
## S4 method for signature 'GeneRelevance' print(x) ## S4 method for signature 'GeneRelevance' show(object) ## S4 method for signature 'GeneRelevance' featureNames(object) ## S4 replacement method for signature 'GeneRelevance,characterOrFactor' featureNames(object) <- value ## S4 method for signature 'GeneRelevance' dataset(object) ## S4 replacement method for signature 'GeneRelevance' dataset(object) <- value ## S4 method for signature 'GeneRelevance' distance(object) ## S4 replacement method for signature 'GeneRelevance' distance(object) <- value
## S4 method for signature 'GeneRelevance' print(x) ## S4 method for signature 'GeneRelevance' show(object) ## S4 method for signature 'GeneRelevance' featureNames(object) ## S4 replacement method for signature 'GeneRelevance,characterOrFactor' featureNames(object) <- value ## S4 method for signature 'GeneRelevance' dataset(object) ## S4 replacement method for signature 'GeneRelevance' dataset(object) <- value ## S4 method for signature 'GeneRelevance' distance(object) ## S4 replacement method for signature 'GeneRelevance' distance(object) <- value
x , object
|
|
value |
dataset
, distance
, and featureNames
return the stored properties.
The other methods return a GeneRelevance
object (print
, ... <- ...
),
or NULL (show
), invisibly
gene_relevance
, Gene Relevance plotting
data(guo_norm) dm <- DiffusionMap(guo_norm) gr <- gene_relevance(dm) stopifnot(distance(gr) == distance(dm)) featureNames(gr)[[37]] <- 'Id2 (suppresses differentiation)' # now plot it with the changed gene name(s)
data(guo_norm) dm <- DiffusionMap(guo_norm) gr <- gene_relevance(dm) stopifnot(distance(gr) == distance(dm)) featureNames(gr)[[37]] <- 'Id2 (suppresses differentiation)' # now plot it with the changed gene name(s)
The relevance map is cached insided of the DiffusionMap
.
gene_relevance( coords, exprs, ..., k = 20L, dims = 1:2, distance = NULL, smooth = TRUE, remove_outliers = FALSE, verbose = FALSE ) ## S4 method for signature 'DiffusionMap,missing' gene_relevance( coords, exprs, ..., k = 20L, dims = 1:2, distance = NULL, smooth = TRUE, remove_outliers = FALSE, verbose = FALSE ) ## S4 method for signature 'matrix,dMatrixOrMatrix' gene_relevance( coords, exprs, ..., pcs = NULL, knn_params = list(), weights = 1, k, dims, distance, smooth, remove_outliers, verbose )
gene_relevance( coords, exprs, ..., k = 20L, dims = 1:2, distance = NULL, smooth = TRUE, remove_outliers = FALSE, verbose = FALSE ) ## S4 method for signature 'DiffusionMap,missing' gene_relevance( coords, exprs, ..., k = 20L, dims = 1:2, distance = NULL, smooth = TRUE, remove_outliers = FALSE, verbose = FALSE ) ## S4 method for signature 'matrix,dMatrixOrMatrix' gene_relevance( coords, exprs, ..., pcs = NULL, knn_params = list(), weights = 1, k, dims, distance, smooth, remove_outliers, verbose )
coords |
A |
exprs |
An cells |
... |
Unused. All parameters to the right of the |
k |
Number of nearest neighbors to use |
dims |
Index into columns of |
distance |
Distance measure to use for the nearest neighbor search. |
smooth |
Smoothing parameters |
remove_outliers |
Remove cells that are only within one other cell's nearest neighbor, as they tend to get large norms. |
verbose |
If TRUE, log additional info to the console |
pcs |
A cell |
knn_params |
|
weights |
Weights for the partial derivatives. A vector of the same length as |
A GeneRelevance
object:
coords
A cells dims
matrix
or sparseMatrix
of coordinates (e.g. diffusion components), reduced to the dimensions passed as dims
exprs
A cells genes matrix of expressions
partials
Array of partial derivatives wrt to considered dimensions in reduced space
(genes cells
dimensions)
partials_norm
Matrix with norm of aforementioned derivatives. (n\_genes cells)
nn_index
Matrix of k nearest neighbor indices. (cells k)
dims
Column index for plotted dimensions. Can character
, numeric
or logical
distance
Distance measure used in the nearest neighbor search. See find_knn
smooth_window
Smoothing window used (see smth.gaussian
)
smooth_alpha
Smoothing kernel width used (see smth.gaussian
)
Gene Relevance methods, Gene Relevance plotting: plot_differential_map
/plot_gene_relevance
data(guo_norm) dm <- DiffusionMap(guo_norm) gr <- gene_relevance(dm) m <- t(Biobase::exprs(guo_norm)) gr_pca <- gene_relevance(prcomp(m)$x, m) # now plot them!
data(guo_norm) dm <- DiffusionMap(guo_norm) gr <- gene_relevance(dm) m <- t(Biobase::exprs(guo_norm)) gr_pca <- gene_relevance(prcomp(m)$x, m) # now plot them!
Gene expression data of 48 genes and an annotation column $num_cells
containing the cell stage at which the embryos were harvested.
data(guo) data(guo_norm)
data(guo) data(guo_norm)
An ExpressionSet with 48 features, 428 observations and 2 phenoData annotations.
The data is normalized using the mean of two housekeeping genes.
The difference between guo
and guo_norm
is the LoD being set to 10 in the former, making it usable with the censor_val
parameter of DiffusionMap.
an ExpressionSet with 48 features and 428 observations containing qPCR Ct values and a "num.cells" observation annotation.
Guoji Guo, Mikael Huss, Guo Qing Tong, Chaoyang Wang, Li Li Sun, Neil D. Clarke, Paul Robson [email protected]
http://www.sciencedirect.com/science/article/pii/S1534580710001103
Inverse of which. Converts an array of numeric or character indices to a logical index array. This function is useful if you need to perform logical operation on an index array but are only given numeric indices.
l_which(idx, nms = seq_len(len), len = length(nms), useNames = TRUE)
l_which(idx, nms = seq_len(len), len = length(nms), useNames = TRUE)
idx |
Numeric or character indices. |
nms |
Array of names or a sequence. Required if |
len |
Length of output array. Alternative to |
useNames |
Use the names of nms or idx |
Either nms
or len
has to be specified.
Logical vector of length len
or the same length as nms
all(l_which(2, len = 3L) == c(FALSE, TRUE, FALSE)) all(l_which(c('a', 'c'), letters[1:3]) == c(TRUE, FALSE, TRUE))
all(l_which(2, len = 3L) == c(FALSE, TRUE, FALSE)) all(l_which(c('a', 'c'), letters[1:3]) == c(TRUE, FALSE, TRUE))
plot(gene_relevance, 'Gene')
plots the differential map of this/these gene(s),
plot(gene_relevance)
a relevance map of a selection of genes.
Alternatively, you can use plot_differential_map
or plot_gene_relevance
on a GeneRelevance
or DiffusionMap
object, or with two matrices.
plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'matrix,matrix' plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'DiffusionMap,missing' plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'GeneRelevance,missing' plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) ## S4 method for signature 'matrix,matrix' plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) ## S4 method for signature 'DiffusionMap,missing' plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) ## S4 method for signature 'GeneRelevance,missing' plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'matrix,matrix' plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'DiffusionMap,missing' plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'GeneRelevance,missing' plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'GeneRelevance,character' plot(x, y, ...) ## S4 method for signature 'GeneRelevance,numeric' plot(x, y, ...) ## S4 method for signature 'GeneRelevance,missing' plot(x, y, ...)
plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'matrix,matrix' plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'DiffusionMap,missing' plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'GeneRelevance,missing' plot_differential_map( coords, exprs, ..., genes, dims = 1:2, pal = hcl.colors, faceter = facet_wrap(~Gene) ) plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) ## S4 method for signature 'matrix,matrix' plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) ## S4 method for signature 'DiffusionMap,missing' plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) ## S4 method for signature 'GeneRelevance,missing' plot_gene_relevance( coords, exprs, ..., iter_smooth = 2L, n_top = 10L, genes = NULL, dims = 1:2, pal = palette(), col_na = "grey", limit = TRUE ) plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'matrix,matrix' plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'DiffusionMap,missing' plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'GeneRelevance,missing' plot_gene_relevance_rank( coords, exprs, ..., genes, dims = 1:2, n_top = 10L, pal = c("#3B99B1", "#F5191C"), bins = 10L, faceter = facet_wrap(~Gene) ) ## S4 method for signature 'GeneRelevance,character' plot(x, y, ...) ## S4 method for signature 'GeneRelevance,numeric' plot(x, y, ...) ## S4 method for signature 'GeneRelevance,missing' plot(x, y, ...)
coords |
A |
exprs |
An cells |
... |
Passed to |
genes |
Genes to base relevance map on (vector of strings). You can also pass an index into the gene names (vector of numbers or logicals with length > 1). The default NULL means all genes. |
dims |
Names or indices of dimensions to plot. When not plotting a |
pal |
Palette. Either A colormap function or a list of colors. |
faceter |
A ggplot faceter like |
iter_smooth |
Number of label smoothing iterations to perform on relevance map. The higher the more homogenous and the less local structure. |
n_top |
Number the top n genes per cell count towards the score defining which genes to return and plot in the relevance map. |
col_na |
Color for cells that end up with no most relevant gene. |
limit |
Limit the amount of displayed gene labels to the amount of available colors in |
bins |
Number of hexagonal bins for |
x |
|
y |
Gene name(s) or index/indices to create differential map for. (integer or character) |
ggplot2 plot, when plotting a relevance map with a list member $ids
containing the gene IDs used.
gene_relevance
, Gene Relevance methods
data(guo_norm) dm <- DiffusionMap(guo_norm) gr <- gene_relevance(dm) plot(gr) # or plot_gene_relevance(dm) plot(gr, 'Fgf4') # or plot_differential_map(dm, 'Fgf4') guo_norm_mat <- t(Biobase::exprs(guo_norm)) pca <- prcomp(guo_norm_mat)$x plot_gene_relevance(pca, guo_norm_mat, dims = 2:3) plot_differential_map(pca, guo_norm_mat, genes = c('Fgf4', 'Nanog'))
data(guo_norm) dm <- DiffusionMap(guo_norm) gr <- gene_relevance(dm) plot(gr) # or plot_gene_relevance(dm) plot(gr, 'Fgf4') # or plot_differential_map(dm, 'Fgf4') guo_norm_mat <- t(Biobase::exprs(guo_norm)) pca <- prcomp(guo_norm_mat)$x plot_gene_relevance(pca, guo_norm_mat, dims = 2:3) plot_differential_map(pca, guo_norm_mat, genes = c('Fgf4', 'Nanog'))
If you want to plot the eigenvalues, simply plot(eigenvalues(dm)[start:end], ...)
plot.DiffusionMap( x, dims = 1:3, new_dcs = if (!is.null(new_data)) dm_predict(x, new_data), new_data = NULL, col = NULL, col_by = NULL, col_limits = NULL, col_new = "red", pal = NULL, pal_new = NULL, ..., ticks = FALSE, axes = TRUE, box = FALSE, legend_main = col_by, legend_opts = list(), interactive = FALSE, draw_legend = !is.null(col_by) || (length(col) > 1 && !is.character(col)), consec_col = TRUE, col_na = "grey", plot_more = function(p, ..., rescale = NULL) p ) ## S4 method for signature 'DiffusionMap,numeric' plot(x, y, ...) ## S4 method for signature 'DiffusionMap,missing' plot(x, y, ...)
plot.DiffusionMap( x, dims = 1:3, new_dcs = if (!is.null(new_data)) dm_predict(x, new_data), new_data = NULL, col = NULL, col_by = NULL, col_limits = NULL, col_new = "red", pal = NULL, pal_new = NULL, ..., ticks = FALSE, axes = TRUE, box = FALSE, legend_main = col_by, legend_opts = list(), interactive = FALSE, draw_legend = !is.null(col_by) || (length(col) > 1 && !is.character(col)), consec_col = TRUE, col_na = "grey", plot_more = function(p, ..., rescale = NULL) p ) ## S4 method for signature 'DiffusionMap,numeric' plot(x, y, ...) ## S4 method for signature 'DiffusionMap,missing' plot(x, y, ...)
x |
|
dims , y
|
Diffusion components (eigenvectors) to plot (default: first three components; 1:3) |
new_dcs |
An optional matrix also containing the rows specified with |
new_data |
A data set in the same format as |
col |
Single color string or vector of discrete or categoric values to be mapped to colors.
E.g. a column of the data matrix used for creation of the diffusion map. (default: |
col_by |
Specify a |
col_limits |
If |
col_new |
If |
pal |
Palette used to map the |
pal_new |
Palette used to map the |
... |
Parameters passed to plot, scatterplot3d, or plot3d (if |
ticks |
logical. If TRUE, show axis ticks (default: FALSE) |
axes |
logical. If TRUE, draw plot axes (default: Only if |
box |
logical. If TRUE, draw plot frame (default: TRUE or the same as |
legend_main |
Title of legend. (default: nothing unless |
legend_opts |
Other colorlegend options (default: empty list) |
interactive |
Use plot3d to plot instead of scatterplot3d? |
draw_legend |
logical. If TRUE, draw color legend (default: TRUE if |
consec_col |
If |
col_na |
Color for |
plot_more |
Function that will be called while the plot margins are temporarily changed
(its |
If you specify negative numbers as diffusion components (e.g. plot(dm, c(-1,2))
), then the corresponding components will be flipped.
The return value of the underlying call is returned, i.e. a scatterplot3d or rgl object.
data(guo) plot(DiffusionMap(guo))
data(guo) plot(DiffusionMap(guo))
Plots diffusion components from a Diffusion Map and the accompanying Diffusion Pseudo Time (DPT
)
plot.DPT( x, root = NULL, paths_to = integer(0L), dcs = 1:2, divide = integer(0L), w_width = 0.1, col_by = "dpt", col_path = rev(palette()), col_tip = "red", ..., col = NULL, legend_main = col_by ) ## S4 method for signature 'DPT,numeric' plot(x, y, ...) ## S4 method for signature 'DPT,missing' plot(x, y, ...)
plot.DPT( x, root = NULL, paths_to = integer(0L), dcs = 1:2, divide = integer(0L), w_width = 0.1, col_by = "dpt", col_path = rev(palette()), col_tip = "red", ..., col = NULL, legend_main = col_by ) ## S4 method for signature 'DPT,numeric' plot(x, y, ...) ## S4 method for signature 'DPT,missing' plot(x, y, ...)
x |
A |
paths_to |
Numeric Branch IDs. Are used as target(s) for the path(s) to draw. |
dcs |
The dimensions to use from the DiffusionMap |
divide |
If |
w_width |
Window width for smoothing the path (see |
col_by |
Color by 'dpt' (DPT starting at |
col_path |
Colors for the path or a function creating n colors |
col_tip |
Color for branch tips |
... |
Graphical parameters supplied to |
col |
See |
legend_main |
See |
y , root
|
Root branch ID. Will be used as the start of the DPT. (default: lowest branch ID)
(If longer than size 1, will be interpreted as |
The return value of the underlying call is returned, i.e. a scatterplot3d or rgl object for 3D plots.
data(guo_norm) dm <- DiffusionMap(guo_norm) dpt <- DPT(dm) plot(dpt) plot(dpt, 2L, col_by = 'branch') plot(dpt, 1L, 2:3, col_by = 'num_cells') plot(dpt, col_by = 'DPT3')
data(guo_norm) dm <- DiffusionMap(guo_norm) dpt <- DPT(dm) plot(dpt) plot(dpt, 2L, col_by = 'branch') plot(dpt, 1L, 2:3, col_by = 'num_cells') plot(dpt, col_by = 'DPT3')
Plot Sigmas object
## S4 method for signature 'Sigmas,missing' plot( x, col = par("fg"), col_highlight = "#E41A1C", col_line = "#999999", type = c("b", "b"), pch = c(par("pch"), 4L), only_dim = FALSE, ..., xlab = NULL, ylab = NULL, main = "" )
## S4 method for signature 'Sigmas,missing' plot( x, col = par("fg"), col_highlight = "#E41A1C", col_line = "#999999", type = c("b", "b"), pch = c(par("pch"), 4L), only_dim = FALSE, ..., xlab = NULL, ylab = NULL, main = "" )
x |
Sigmas object to plot |
col |
Vector of bar colors or single color for all bars |
col_highlight |
Color for highest bar. Overrides col |
col_line |
Color for the line and its axis |
type |
Plot type of both lines. Can be a vector of length 2 to specify both separately (default: 'b' aka “both lines and points”) |
pch |
Point identifier for both lines. Can be a vector of length 2 to specify both separately (default: |
only_dim |
logical. If TRUE, only plot the derivative line |
... |
Options passed to the call to plot |
xlab |
X label. NULL to use default |
ylab |
Either one y label or y labels for both plots. NULL to use both defauts, a NULL in a list of length 2 to use one default. |
main |
Title of the plot |
This method plots a Sigma object to the current device and returns nothing/NULL
data(guo) sigs <- find_sigmas(guo) plot(sigs)
data(guo) sigs <- find_sigmas(guo) plot(sigs)
Projection distance
projection_dist(dm, new_dcs = NULL, ..., new_data, verbose = FALSE)
projection_dist(dm, new_dcs = NULL, ..., new_data, verbose = FALSE)
dm |
A |
new_dcs |
Diffusion component matrix of which to calculate the distance to the data. |
... |
Passed to |
new_data |
New data points to project into the diffusion map. Can be a matrix, data.frame, ExpressionSet, or SingleCellExperiment. |
verbose |
If |
A vector of distances each new data point has to the existing data.
data(guo_norm) g2_32 <- guo_norm[, guo_norm$num_cells < 64] g64 <- guo_norm[, guo_norm$num_cells == 64] dm <- DiffusionMap(g2_32) d <- projection_dist(dm, new_data = g64)
data(guo_norm) g2_32 <- guo_norm[, guo_norm$num_cells < 64] g64 <- guo_norm[, guo_norm$num_cells == 64] dm <- DiffusionMap(g2_32) d <- projection_dist(dm, new_data = g64)
Finds a cell that has the maximum DPT distance from a randomly selected one.
random_root(dm_or_dpt)
random_root(dm_or_dpt)
dm_or_dpt |
A |
A cell index
data(guo) dm <- DiffusionMap(guo) random_root(dm)
data(guo) dm <- DiffusionMap(guo) random_root(dm)
Holds the information about how the sigma
parameter for a DiffusionMap was obtained,
and in this way provides a plotting function for the find_sigmas heuristic.
You should not need to create a Sigmas object yourself. Provide sigma
to DiffusionMap instead or use find_sigmas.
Sigmas(...) ## S4 method for signature 'Sigmas' optimal_sigma(object) ## S4 method for signature 'Sigmas' print(x) ## S4 method for signature 'Sigmas' show(object)
Sigmas(...) ## S4 method for signature 'Sigmas' optimal_sigma(object) ## S4 method for signature 'Sigmas' print(x) ## S4 method for signature 'Sigmas' show(object)
object , x
|
Sigmas object |
... |
See “Slots” below |
A Sigmas object is either created by find_sigmas or by specifying the sigma
parameter to DiffusionMap.
In the second case, if the sigma
parameter is just a number,
the resulting Sigmas
object has all slots except of optimal_sigma
set to NULL
.
Sigmas
creates an object of the same class
optimal_sigma
retrieves the numeric value of the optimal sigma or local sigmas
log_sigmas
Vector of length containing the
of the
s
dim_norms
Vector of length containing the average dimensionality
for the respective kernel widths
optimal_sigma
Multiple local sigmas or the mean of the two global s around the highest
(
c(optimal_idx, optimal_idx+1L)
)
optimal_idx
The index of the highest .
avrd_norms
Vector of length containing the average dimensionality for the corresponding sigma.
find_sigmas
, the function to determine a locally optimal sigma and returning this class
data(guo) sigs <- find_sigmas(guo, verbose = FALSE) optimal_sigma(sigs) print(sigs)
data(guo) sigs <- find_sigmas(guo, verbose = FALSE) optimal_sigma(sigs) print(sigs)
Handles DiffusionMap, Sigmas, and GeneRelevance.
## S4 method for signature 'DiffusionMap' updateObject(object, ..., verbose = FALSE) ## S4 method for signature 'Sigmas' updateObject(object, ..., verbose = FALSE) ## S4 method for signature 'GeneRelevance' updateObject(object, ..., verbose = FALSE)
## S4 method for signature 'DiffusionMap' updateObject(object, ..., verbose = FALSE) ## S4 method for signature 'Sigmas' updateObject(object, ..., verbose = FALSE) ## S4 method for signature 'GeneRelevance' updateObject(object, ..., verbose = FALSE)
object |
An object created with an older destiny release |
... |
ignored |
verbose |
tells what is being updated |
A DiffusionMap or Sigmas object that is valid when used with the current destiny release