scaterPCA: Perform scater PCA on a SingleCellExperiment Object

View source: R/scater_PCA.R

scaterPCAR Documentation

Perform scater PCA on a SingleCellExperiment Object

Description

A wrapper to runPCA function to compute principal component analysis (PCA) from a given SingleCellExperiment object.

Usage

scaterPCA(
  inSCE,
  useAssay = "logcounts",
  useFeatureSubset = "hvg2000",
  scale = TRUE,
  reducedDimName = "PCA",
  nComponents = 50,
  ntop = 2000,
  useAltExp = NULL,
  seed = 12345,
  BPPARAM = BiocParallel::SerialParam()
)

Arguments

inSCE

Input SingleCellExperiment object.

useAssay

Assay to use for PCA computation. If useAltExp is specified, useAssay has to exist in assays(altExp(inSCE, useAltExp)). Default "logcounts"

useFeatureSubset

Subset of feature to use for dimension reduction. A character string indicating a rowData variable that stores the logical vector of HVG selection, or a vector that can subset the rows of inSCE. Default "hvg2000".

scale

Logical scalar, whether to standardize the expression values. Default TRUE.

reducedDimName

Name to use for the reduced output assay. Default "PCA".

nComponents

Number of principal components to obtain from the PCA computation. Default 50.

ntop

Automatically detect this number of variable features to use for dimension reduction. Ignored when using useReducedDim or using useFeatureSubset. Default 2000.

useAltExp

The subset to use for PCA computation, usually for the selected.variable features. Default NULL.

seed

Integer, random seed for reproducibility of PCA results. Default NULL.

BPPARAM

A BiocParallelParam object specifying whether the PCA should be parallelized.

Value

A SingleCellExperiment object with PCA computation updated in reducedDim(inSCE, reducedDimName).

Examples

data(scExample, package = "singleCellTK")
sce <- subsetSCECols(sce, colData = "type != 'EmptyDroplet'")
sce <- scaterlogNormCounts(sce, "logcounts")

# Example of ranking variable genes, selecting the top variable features,
# and running PCA. Make sure to increase the number of highly variable
# features (hvgNumber) and the number of principal components (nComponents)
# for real datasets
sce <- runModelGeneVar(sce, useAssay = "logcounts")
sce <- setTopHVG(sce, method = "modelGeneVar", hvgNumber = 100,
                 featureSubsetName = "hvf")
sce <- scaterPCA(sce, useAssay = "logcounts", scale = TRUE,
                 useFeatureSubset = "hvf", nComponents = 5)
                 
# Alternatively, let the scater PCA function select the top variable genes
sce <- scaterPCA(sce, useAssay = "logcounts", scale = TRUE,
                 useFeatureSubset = NULL, ntop = 100, nComponents = 5)

compbiomed/singleCellTK documentation built on Oct. 27, 2024, 3:26 a.m.