corral: corral: Correspondence analysis on a single matrix

View source: R/corral.R

corral_matR Documentation

corral: Correspondence analysis on a single matrix

Description

corral can be used for dimension reduction to find a set of low-dimensional embeddings for a count matrix.

corral is a wrapper for corral_mat and corral_sce, and can be called on any of the acceptable input types.

Usage

corral_mat(
  inp,
  method = c("irl", "svd"),
  ncomp = 30,
  row.w = NULL,
  col.w = NULL,
  rtype = c("standardized", "indexed", "hellinger", "freemantukey", "pearson"),
  vst_mth = c("none", "sqrt", "freemantukey", "anscombe"),
  ...
)

corral_sce(
  inp,
  method = c("irl", "svd"),
  ncomp = 30,
  whichmat = "counts",
  fullout = FALSE,
  subset_row = NULL,
  ...
)

corral(inp, ...)

## S3 method for class 'corral'
print(x, ...)

Arguments

inp

matrix (any type), SingleCellExperiment, or SummarizedExperiment. If using SingleCellExperiment or SummarizedExperiment, then include the whichmat argument to specify which slot to use (defaults to counts).

method

character, the algorithm to be used for svd. Default is irl. Currently supports 'irl' for irlba::irlba or 'svd' for stats::svd

ncomp

numeric, number of components; Default is 30

row.w

numeric vector; the row weights to use in chi-squared scaling. Defaults to 'NULL', in which case row weights are computed from the input matrix.

col.w

numeric vector; the column weights to use in chi-squared scaling. For instance, size factors could be given here. Defaults to 'NULL', in which case column weights are computed from the input matrix.

rtype

character indicating what type of residual should be computed; options are '"indexed"', '"standardized"' (or '"pearson"' is equivalent), '"freemantukey"', and '"hellinger"'; defaults to '"standardized"' for corral and '"indexed"' for corralm. '"indexed"', '"standardized"', and '"freemantukey"' compute the respective chi-squared residuals and are appropriate for count data. The '"hellinger"' option is appropriate for continuous data.

vst_mth

character indicating whether a variance-stabilizing transform should be applied prior to calculating chi-squared residuals; defaults to '"none"'

...

(additional arguments for methods)

whichmat

character; defaults to counts, can also use logcounts or normcounts if stored in the sce object

fullout

boolean; whether the function will return the full corral output as a list, or a SingleCellExperiment; defaults to SingleCellExperiment (FALSE). To get back the corral_mat-style output, set this to TRUE.

subset_row

numeric, character, or boolean vector; the rows to include in corral, as indices (numeric), rownames (character), or with booleans (same length as the number of rows in the matrix). If this parameter is NULL, then all rows will be used.

x

(print method) corral object; the list output from corral_mat

Value

When run on a matrix, a list with the correspondence analysis matrix decomposition result:

d

a vector of the diagonal singular values of the input mat (from SVD output)

u

a matrix of with the left singular vectors of mat in the columns (from SVD output)

v

a matrix of with the right singular vectors of mat in the columns. When cells are in the columns, these are the cell embeddings. (from SVD output)

eigsum

sum of the eigenvalues for calculating percent variance explained

SCu and SCv

standard coordinates, left and right, respectively

PCu and PCv

principal coordinates, left and right, respectively

When run on a SingleCellExperiment, returns a SCE with the embeddings (PCv from the full corral output) in the reducedDim slot corral (default). Also can return the same output as corral_mat when fullout is set to TRUE.

For matrix and SummarizedExperiment input, returns list with the correspondence analysis matrix decomposition result (u,v,d are the raw svd output; SCu and SCv are the standard coordinates; PCu and PCv are the principal coordinates)

For SummarizedExperiment input, returns the same as for a matrix.

.

Examples

mat <- matrix(sample(0:10, 5000, replace=TRUE), ncol=50)
result <- corral_mat(mat)
result <- corral_mat(mat, method = 'irl', ncomp = 5)

library(DuoClustering2018)
sce <- sce_full_Zhengmix4eq()[1:100,1:100]
result_1 <- corral_sce(sce)
result_2 <- corral_sce(sce, method = 'svd')
result_3 <- corral_sce(sce, method = 'irl', ncomp = 30, whichmat = 'logcounts')


library(DuoClustering2018)
sce <- sce_full_Zhengmix4eq()[1:100,1:100]
corral_sce <- corral(sce,whichmat = 'counts')

mat <- matrix(sample(0:10, 500, replace=TRUE), ncol=25)
corral_mat <- corral(mat, ncomp=5)

mat <- matrix(sample(1:100, 10000, replace = TRUE), ncol = 100)
corral(mat)

laurenhsu1/corral documentation built on Feb. 19, 2023, 10:37 p.m.