perCellQCMetrics: Per-cell quality control metrics

perCellQCMetricsR Documentation

Per-cell quality control metrics

Description

Compute per-cell quality control metrics for a count matrix or a SingleCellExperiment.

Usage

perCellQCMetrics(x, ...)

## S4 method for signature 'ANY'
perCellQCMetrics(
  x,
  subsets = NULL,
  percent.top = integer(0),
  threshold = 0,
  BPPARAM = SerialParam(),
  flatten = TRUE,
  percent_top = NULL,
  detection_limit = NULL
)

## S4 method for signature 'SummarizedExperiment'
perCellQCMetrics(x, ..., assay.type = "counts", exprs_values = NULL)

## S4 method for signature 'SingleCellExperiment'
perCellQCMetrics(
  x,
  subsets = NULL,
  percent.top = integer(0),
  ...,
  flatten = TRUE,
  assay.type = "counts",
  use.altexps = NULL,
  percent_top = NULL,
  exprs_values = NULL,
  use_altexps = NULL
)

Arguments

x

A numeric matrix of counts with cells in columns and features in rows.

Alternatively, a SummarizedExperiment or SingleCellExperiment object containing such a matrix.

...

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

For the SummarizedExperiment and SingleCellExperiment methods, further arguments to pass to the ANY method.

subsets

A named list containing one or more vectors (a character vector of feature names, a logical vector, or a numeric vector of indices), used to identify interesting feature subsets such as ERCC spike-in transcripts or mitochondrial genes.

percent.top

An integer vector specifying the size(s) of the top set of high-abundance genes. Used to compute the percentage of library size occupied by the most highly expressed genes in each cell.

threshold

A numeric scalar specifying the threshold above which a gene is considered to be detected.

BPPARAM

A BiocParallelParam object specifying how parallelization should be performed.

flatten

Logical scalar indicating whether the nested DataFrames in the output should be flattened.

percent_top, detection_limit, exprs_values, use_altexps

Soft deprecated equivalents to the arguments described above.

assay.type

A string or integer scalar indicating which assays in the x contains the count matrix.

use.altexps

Logical scalar indicating whether QC statistics should be computed for alternative Experiments in x. If TRUE, statistics are computed for all alternative experiments.

Alternatively, an integer or character vector specifying the alternative Experiments to use to compute QC statistics.

Alternatively NULL, in which case we only use alternative experiments that contain the specified assay.type.

Alternatively FALSE, in which case alternative experiments are not used.

Details

This function calculates useful QC metrics for identification and removal of potentially problematic cells. Obvious per-cell metrics are the sum of counts (i.e., the library size) and the number of detected features. The percentage of counts in the top features also provides a measure of library complexity.

If subsets is specified, these statistics are also computed for each subset of features. This is useful for investigating gene sets of interest, e.g., mitochondrial genes, Y chromosome genes. These statistics are stored as nested DataFrames in the subsets field of the output. For example, if the input subsets contained "Mito" and "Sex", the output would look like:

  output 
  |-- sum
  |-- detected
  |-- percent.top
  +-- subsets
      |-- Mito
      |   |-- sum
      |   |-- detected
      |   +-- percent
      +-- Sex 
          |-- sum
          |-- detected
          +-- percent

Here, the percent field contains the percentage of each cell's count sum assigned to each subset.

If use.altexps=TRUE, the same statistics are computed for each alternative experiment in x. This can also be an integer or character vector specifying the alternative Experiments to use. These statistics are also stored as nested DataFrames, this time in the altexps field of the output. For example, if x contained the alternative Experiments "Spike" and "Ab", the output would look like:

  output 
  |-- sum
  |-- detected
  |-- percent.top
  +-- altexps 
  |   |-- Spike
  |   |   |-- sum
  |   |   |-- detected
  |   |   +-- percent
  |   +-- Ab
  |       |-- sum
  |       |-- detected
  |       +-- percent
  +-- total 

The total field contains the total sum of counts for each cell across the main and alternative Experiments. The percent field contains the percentage of the total count in each alternative Experiment for each cell.

Note that the denominator for altexps$...$percent is not the same as the denominator for subset$...$percent. For example, if subsets contains a set of mitochondrial genes, the mitochondrial percentage would be computed as a fraction of the total endogenous coverage, while the altexps percentage would be computed as a fraction of the total coverage across all (endogenous and artificial) features.

If flatten=TRUE, the nested DataFrames are flattened by concatenating the column names with underscores. This means that, say, the subsets$Mito$sum nested field becomes the top-level subsets_Mito_sum field. A flattened structure is more convenient for end-users performing interactive analyses, but less convenient for programmatic access as artificial construction of strings is required.

Value

A DataFrame of QC statistics where each row corresponds to a column in x. This contains the following fields:

  • sum: numeric, the sum of counts for each cell.

  • detected: numeric, the number of observations above threshold.

If flatten=FALSE, the DataFrame will contain the additional columns:

  • percent.top: numeric matrix, the percentage of counts assigned to the top most highly expressed genes. Each column of the matrix corresponds to an entry of percent.top, sorted in increasing order.

  • subsets: A nested DataFrame containing statistics for each subset, see Details.

  • altexps: A nested DataFrame containing statistics for each alternative experiment, see Details. This is only returned for the SingleCellExperiment method.

  • total: numeric, the total sum of counts for each cell across main and alternative Experiments. This is only returned for the SingleCellExperiment method.

If flatten=TRUE, nested matrices and DataFrames are flattened to remove the hierarchical structure from the output DataFrame.

Author(s)

Aaron Lun

See Also

addPerCellQCMetrics, to add the QC metrics to the column metadata.

Examples

example_sce <- mockSCE()
stats <- perCellQCMetrics(example_sce)
stats

# With subsets.
stats2 <- perCellQCMetrics(example_sce, subsets=list(Mito=1:10), 
    flatten=FALSE)
stats2$subsets

# With alternative Experiments.
pretend.spike <- ifelse(seq_len(nrow(example_sce)) < 10, "Spike", "Gene")
alt_sce <- splitAltExps(example_sce, pretend.spike)
stats3 <- perCellQCMetrics(alt_sce, flatten=FALSE)
stats3$altexps



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