project_cycle_space: Project data into the cell cycle pattern space

View source: R/project_cycle_space.R

project_cycle_spaceR Documentation

Project data into the cell cycle pattern space

Description

Project mouse and human single cell RNAseq data into a cell cycle embedding by a pre-learned reference projection matrix.

Usage

project_cycle_space(
  x,
  exprs_values = "logcounts",
  altexp = NULL,
  name = "tricycleEmbedding",
  ref.m = NULL,
  gname = NULL,
  gname.type = c("ENSEMBL", "SYMBOL"),
  species = c("mouse", "human"),
  AnnotationDb = NULL
)

Arguments

x

A numeric matrix of **log-expression** values where rows are features and columns are cells. Alternatively, a SummarizedExperiment or SingleCellExperiment containing such a matrix.

exprs_values

Integer scalar or string indicating which assay of x contains the **log-expression** values. Default: 'logcounts'

altexp

String or integer scalar specifying an alternative experiment containing the input data.

name

String specifying the name to be used to store the result in the reducedDims of the output. Default: 'tricycleEmbedding'

ref.m

A custom reference projection matrix to project the new data, where rows are features and columns are dimensions. Users need to use the same type of gname(or rownames of x) as for the ref.m. If no custom ref.m is given, the internal reference neuroRef will be used.

gname

Alternative rownames of x. If provided, this will be used to map genes within x with genes in ref.m. If not provided, the rownames of x will be used instead. Default: NULL

gname.type

The type of gene names as in gname or rownames of x. It can be either 'ENSEMBL' or 'SYMBOL'. If the user uses custom ref.m, this value will have no effect. Default: 'ENSEMBL'

species

The type of species in x. It can be either 'mouse' or 'human'. If the user uses custom ref.m, this value will have no effect. Default: 'mouse'

AnnotationDb

An AnnotationDb objects. If the user uses the internal reference to project human data, and provide rownames in the format of Ensembl IDs, this object will be used to map Ensembl IDs to gene SYMBOLs. If no AnnotationDb object being given, the function will use org.Hs.eg.db.

Details

The function will use pre-learned cell cycle pattern to project new data to show the cell cycle progression. If the user uses internal Neuropshere reference, the expression values must be **log-transformed**. Besides, we would assume the input data has been already preprocessed, library size normalized at least. The projection process is to take sum of weighted mean-centered expression of chosen genes, so the mean expression of a given gene could be affected without library size normalization.

Value

If the input is a numeric matrix or a SummarizedExperiment, a projection matrix with rows cells and column dimensions will be returned. The actual rotation matrix used to project the data is included in the attributes with name 'rotation'.

For SingleCellExperiment, an updated SingleCellExperiment is returned containing projection matrix in reducedDims(..., name).

Author(s)

Shijie C. Zheng

References

Zheng SC, et al. Universal prediction of cell cycle position using transfer learning. Genome Biology (2022) 23: 41 doi:10.1186/s13059-021-02581-y.

See Also

estimate_cycle_position, for inferring cell cycle position.

Examples

data(neurosphere_example, package = "tricycle")
neurosphere_example <- project_cycle_space(neurosphere_example)
reducedDimNames(neurosphere_example)
head(reducedDim(neurosphere_example, "tricycleEmbedding"))
plot(reducedDim(neurosphere_example, "tricycleEmbedding"))
names(attributes(reducedDim(neurosphere_example, "tricycleEmbedding")))

hansenlab/tricycle documentation built on March 19, 2022, 7:24 p.m.