quickPerCellQC: Quick cell-level QC

quickPerCellQCR Documentation

Quick cell-level QC

Description

A convenient utility that identifies low-quality cells based on frequently used QC metrics.

Usage

quickPerCellQC(x, ...)

## S4 method for signature 'ANY'
quickPerCellQC(
  x,
  sum.field = "sum",
  detected.field = "detected",
  sub.fields = NULL,
  ...,
  lib_size = NULL,
  n_features = NULL,
  percent_subsets = NULL
)

## S4 method for signature 'SummarizedExperiment'
quickPerCellQC(
  x,
  ...,
  subsets = NULL,
  assay.type = "counts",
  other.args = list(),
  filter = TRUE
)

Arguments

x

A DataFrame containing per-cell QC statistics, as computed by perCellQCMetrics. Alternatively, a SummarizedExperiment object that can be used to create such a DataFrame via perCellQCMetrics.

...

For the generic, further arguments to pass to specific methods.

For the ANY method, further arguments to pass to isOutlier.

For the SummarizedExperiment method, further arguments to pass to the ANY method.

sum.field

String specifying the column of x containing the library size for each cell.

detected.field

String specifying the column of x containing the number of detected features per cell.

sub.fields

Character vector specifying the column(s) of x containing the percentage of counts in subsets of “control features”, usually mitochondrial genes or spike-in transcripts.

If set to TRUE, this will default to all columns in x with names following the patterns "subsets_.*_percent" and "altexps_.*_percent".

lib_size, n_features, percent_subsets

Soft-deprecated equivalents of the arguments above.

subsets, assay.type

Arguments to pass to the perCellQCMetrics function, exposed here for convenience.

other.args

A named list containing other arguments to pass to the perCellQCMetrics function.

filter

Logical scalar indicating whether to filter out low-quality cells from x.

Details

For DataFrame x, this function simply calls perCellQCFilters. The latter should be directly used in such cases; DataFrame inputs are soft-deprecated here.

For SummarizedExperiment x, this function is simply a convenient wrapper around perCellQCMetrics and perCellQCFilters.

Value

If filter=FALSE or x is a DataFrame, a DataFrame is returned with one row per cell and containing columns of logical vectors. Each column specifies a reason for why a cell was considered to be low quality, with the final discard column indicating whether the cell should be discarded.

If filter=TRUE, x is returned with the low-quality cells removed. QC statistics and filtering information for all remaining cells are stored in the colData.

Author(s)

Aaron Lun

See Also

perCellQCMetrics, for calculation of these metrics.

perCellQCFilters, to define filter thresholds based on those metrics.

Examples

example_sce <- mockSCE()

filtered_sce <- quickPerCellQC(example_sce, subsets=list(Mito=1:100),
    sub.fields=c("subsets_Mito_percent", "altexps_Spikes_percent"))
ncol(filtered_sce)

# Same result as the longer chain of expressions:
stats <- perCellQCMetrics(example_sce, subsets=list(Mito=1:100))
discard <- perCellQCFilters(stats, 
    sub.fields=c("subsets_Mito_percent", "altexps_Spikes_percent"))
filtered_sce2 <- example_sce[,!discard$discard]
ncol(filtered_sce2)


LTLA/scuttle documentation built on March 9, 2024, 11:16 a.m.