sandbag: Cell cycle phase training

Description Usage Arguments Details Value Author(s) References See Also Examples

Description

Use gene expression data to train a classifier for cell cycle phase.

Usage

1
2
3
4
5
6
7
## S4 method for signature 'matrix'
sandbag(x, phases, gene.names=rownames(x), 
    fraction=0.5, subset.row=NULL)

## S4 method for signature 'SCESet'
sandbag(x, phases, subset.row=NULL, ..., 
    assay="counts", get.spikes=FALSE)

Arguments

x

A numeric matrix of gene expression values where rows are genes and columns are cells. Alternatively, a SCESet object containing such a matrix.

phases

A list of subsetting vectors specifying which cells are in each phase of the cell cycle. This should typically be of length 3, with elements named as "G1", "S" and "G2M".

gene.names

A character vector of gene names.

fraction

A numeric scalar specifying the minimum fraction to define a marker gene pair.

subset.row

A logical, integer or character scalar indicating the rows of x to use.

...

Additional arguments to pass to sandbag,matrix-method.

assay

A string specifying which assay values to use, e.g., counts or exprs.

get.spikes

A logical specifying whether spike-in transcripts should be used.

Details

This function implements the training step of the pair-based prediction method described by Scialdone et al. (2015). Pairs of genes (A, B) are identified from a training data set where in each pair, the fraction of cells in phase G1 with expression of A > B (based on expression values in training.data) and the fraction with B > A in each other phase exceeds fraction. These pairs are defined as the marker pairs for G1. This is repeated for each phase to obtain a separate marker pair set.

Pre-defined sets of marker pairs are provided for mouse and human (see Examples). The mouse set was generated as described by Scialdone et al. (2015), while the human training set was generated with data from Leng et al. (2015). Classification from test data can be performed using the cyclone function. For each cell, this involves comparing expression values between genes in each marker pair. The cell is then assigned to the phase that is consistent with the direction of the difference in expression in the majority of pairs.

For sandbag,SCESet-method, the matrix of counts is used but can be replaced with expression values by setting assays. By default, get.spikes=FALSE which means that any rows corresponding to spike-in transcripts will not be considered when picking markers. This is because the amount of spike-in RNA added will vary between experiments and will not be a robust predictor. Nonetheless, if all rows are required, users can set get.spikes=TRUE. Users can also manually select which rows to use via subset.row, which will override any setting of get.spikes.

While sandbag and its partner function cyclone were originally designed for cell cyclone phase classification, the same computational strategy can be used to classify cells into any mutually exclusive groupings. Any number and nature of groups can be specified in phases, e.g., differentiation lineages, activation states. Only the names of phases need to be modified to reflect the biology being studied.

Value

A named list of data.frames, where each data frame corresponds to a cell cycle phase and contains the names of the genes in each marker pair.

Author(s)

Antonio Scialdone, with modifications by Aaron Lun

References

Scialdone A, Natarajana KN, Saraiva LR et al. (2015). Computational assignment of cell-cycle stage from single-cell transcriptome data. Methods 85:54–61

Leng N, Chu LF, Barry C et al. (2015). Oscope identifies oscillatory genes in unsynchronized single-cell RNA-seq experiments. Nat. Methods 12:947–50

See Also

cyclone

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
ncells <- 50
ngenes <- 20
training <- matrix(rnorm(ncells*ngenes), ncol=ncells)
rownames(training) <- paste0("X", seq_len(ngenes))

is.G1 <- 1:20
is.S <- 21:30
is.G2M <- 31:50
out <- sandbag(training, list(G1=is.G1, S=is.S, G2M=is.G2M))
str(out)

# Getting pre-trained marker sets
mm.pairs <- readRDS(system.file("exdata", "mouse_cycle_markers.rds", package="scran"))
hs.pairs <- readRDS(system.file("exdata", "human_cycle_markers.rds", package="scran"))

Example output

Loading required package: BiocParallel
Loading required package: scater
Loading required package: Biobase
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, sd, var, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colMeans, colSums, colnames, do.call,
    duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
    lapply, lengths, mapply, match, mget, order, paste, pmax, pmax.int,
    pmin, pmin.int, rank, rbind, rowMeans, rowSums, rownames, sapply,
    setdiff, sort, table, tapply, union, unique, unsplit, which,
    which.max, which.min

Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

Loading required package: ggplot2

Attaching package: 'scater'

The following object is masked from 'package:stats':

    filter

List of 3
 $ G1 :'data.frame':	66 obs. of  2 variables:
  ..$ first : chr [1:66] "X1" "X3" "X5" "X6" ...
  ..$ second: chr [1:66] "X2" "X1" "X1" "X1" ...
 $ S  :'data.frame':	83 obs. of  2 variables:
  ..$ first : chr [1:83] "X1" "X1" "X1" "X1" ...
  ..$ second: chr [1:83] "X5" "X12" "X13" "X15" ...
 $ G2M:'data.frame':	86 obs. of  2 variables:
  ..$ first : chr [1:86] "X1" "X1" "X1" "X1" ...
  ..$ second: chr [1:86] "X3" "X6" "X7" "X8" ...
Warning message:
system call failed: Cannot allocate memory 

scran documentation built on July 19, 2017, 2:01 a.m.

Related to sandbag in scran...