trainSingleR | R Documentation |
Train the SingleR classifier on one or more reference datasets with known labels.
trainSingleR(
ref,
labels,
test.genes = NULL,
genes = "de",
sd.thresh = NULL,
de.method = c("classic", "wilcox", "t"),
de.n = NULL,
de.args = list(),
aggr.ref = FALSE,
aggr.args = list(),
recompute = TRUE,
restrict = NULL,
assay.type = "logcounts",
check.missing = TRUE,
approximate = FALSE,
num.threads = bpnworkers(BPPARAM),
BNPARAM = NULL,
BPPARAM = SerialParam()
)
ref |
A numeric matrix of expression values where rows are genes and columns are reference samples (individual cells or bulk samples). Each row should be named with the gene name. In general, the expression values are expected to be log-transformed, see Details. Alternatively, a SummarizedExperiment object containing such a matrix. Alternatively, a list or List of SummarizedExperiment objects or numeric matrices containing multiple references. |
labels |
A character vector or factor of known labels for all samples in Alternatively, if |
test.genes |
Character vector of the names of the genes in the test dataset, i.e., the row names of |
genes |
A string containing Alternatively, if
If |
sd.thresh |
Deprecated and ignored. |
de.method |
String specifying how DE genes should be detected between pairs of labels.
Defaults to |
de.n |
An integer scalar specifying the number of DE genes to use when |
de.args |
Named list of additional arguments to pass to |
aggr.ref |
Logical scalar indicating whether references should be aggregated to pseudo-bulk samples for speed, see |
aggr.args |
Further arguments to pass to |
recompute |
Deprecated and ignored. |
restrict |
A character vector of gene names to use for marker selection.
By default, all genes in |
assay.type |
An integer scalar or string specifying the assay of |
check.missing |
Logical scalar indicating whether rows should be checked for missing values. If true and any missing values are found, the rows containing these values are silently removed. |
approximate |
Deprecated, use |
num.threads |
Integer scalar specifying the number of threads to use for index building. |
BNPARAM |
A BiocNeighborParam object specifying how the neighbor search index should be constructed. |
BPPARAM |
A BiocParallelParam object specifying how parallelization should be performed.
Relevant for marker detection if |
This function uses a training data set to select interesting features and construct nearest neighbor indices in rank space.
The resulting objects can be re-used across multiple classification steps with different test data sets via classifySingleR
.
This improves efficiency by avoiding unnecessary repetition of steps during the downstream analysis.
The automatic marker detection (genes="de"
) identifies genes that are differentially expressed between labels.
This is done by identifying the median expression within each label, and computing differences between medians for each pair of labels.
For each label, the top de.n
genes with the largest differences compared to another label are chosen as markers to distinguish the two labels.
The expression values are expected to be log-transformed and normalized.
Classification with classifySingleR
assumes that the test dataset contains all marker genes that were detected from the reference.
If the test and reference datasets do not have the same genes in the same order, we can set test.genes
to the row names of the test dataset.
This will instruct trainSingleR
to only consider markers that are present in the test dataset.
Any subsequent call to classifySingleR
will also check that test.genes
is consistent with rownames(test)
.
On a similar note, if restrict
is specified, marker selection will only be performed using the specified subset of genes.
This can be convenient for ignoring inappropriate genes like pseudogenes or predicted genes.
It has the same effect as filtering out undesirable rows from ref
prior to calling trainSingleR
.
Unlike test.genes
, setting restrict
does not introduce further checks on rownames(test)
in classifySingleR
.
For a single reference, a List is returned containing:
built
:An external pointer to various indices in C++ space.
Note that this cannot be serialized and should be removed prior to any saveRDS
step.
ref
:The reference expression matrix.
This may have fewer columns than the input ref
if aggr.ref = TRUE
.
markers
:A list containing unique
, a character vector of all marker genes used in training;
and full
, a list of list of character vectors containing the markers for each pairwise comparison between labels.
labels
:A list containing unique
, a character vector of all unique reference labels;
and full
, a character vector containing the assigned label for each column in ref
.
For multiple references, a List of Lists is returned where each internal List corresponds to a reference in ref
and has the same structure as described above.
Rather than relying on the in-built feature selection, users can pass in their own features of interest to genes
.
The function expects a named list of named lists of character vectors, with each vector containing the DE genes between a pair of labels.
For example:
genes <- list( A = list(A = character(0), B = "GENE_1", C = c("GENE_2", "GENE_3")), B = list(A = "GENE_100", B = character(0), C = "GENE_200"), C = list(A = c("GENE_4", "GENE_5"), B = "GENE_5", C = character(0)) )
If we consider the entry genes$A$B
, this contains marker genes for label "A"
against label "B"
.
That is, these genes are upregulated in "A"
compared to "B"
.
The outer list should have one list per label, and each inner list should have one character vector per label.
(Obviously, a label cannot have markers against itself, so this is just set to character(0)
.)
Alternatively, genes
can be a named list of character vectors containing per-label markers.
For example:
genes <- list( A = c("GENE_1", "GENE_2", "GENE_3"), B = c("GENE_100", "GENE_200"), C = c("GENE_4", "GENE_5") )
The entry genes$A
represent the genes that are upregulated in A
compared to some or all other labels.
This allows the function to handle pre-defined marker lists for specific cell populations.
However, it obviously captures less information than marker sets for the pairwise comparisons.
If genes
is manually passed, ref
can be the raw counts or any monotonic transformation thereof.
There is no need to supply (log-)normalized expression values for the benefit of the automatic marker detection.
Similarly, for manual genes
, the values of de.method
, de.n
and sd.thresh
have no effect.
The default SingleR policy for dealing with multiple references is to perform the classification for each reference separately and combine the results
(see ?combineRecomputedResults
for an explanation).
To this end, if ref
is a list with multiple references, marker genes are identified separately within each reference if genes = NULL
.
Rank calculation and index construction is then performed within each reference separately.
The result is identical to lapply
ing over a list of references and runing trainSingleR
on each reference.
Alternatively, genes
can still be used to explicitly specify marker genes for each label in each of multiple references.
This is achieved by passing a list of lists to genes
,
where each inner list corresponds to a reference in ref
and can be of any format described in “Custom feature specification”.
Thus, it is possible for genes
to be - wait for it - a list (per reference) of lists (per label) of lists (per label) of character vectors.
The default marker selection is based on log-fold changes between the per-label medians and is very much designed with bulk references in mind.
It may not be effective for single-cell reference data where it is not uncommon to have more than 50% zero counts for a given gene such that the median is also zero for each group.
Users are recommended to either set de.method
to another DE ranking method, or detect markers externally and pass a list of markers to genes
(see Examples).
In addition, it is generally unnecessary to have single-cell resolution on the reference profiles.
We can instead set aggr.ref=TRUE
to aggregate per-cell references into a set of pseudo-bulk profiles using aggregateReference
.
This improves classification speed while using vector quantization to preserve within-label heterogeneity and mitigate the loss of information.
Note that any aggregation is done after marker gene detection; this ensures that the relevant tests can appropriately penalize within-label variation.
Users should also be sure to set the seed as the aggregation involves randomization.
Aaron Lun, based on the original SingleR
code by Dvir Aran.
classifySingleR
, where the output of this function gets used.
combineRecomputedResults
, to combine results from multiple references.
rebuildIndex
, to rebuild the index after external memory is invalidated.
# Making up some data for a quick demonstration.
ref <- .mockRefData()
# Normalizing and log-transforming for automated marker detection.
ref <- scuttle::logNormCounts(ref)
trained <- trainSingleR(ref, ref$label)
trained
length(trained$markers$unique)
# Alternatively, computing and supplying a set of label-specific markers.
by.t <- scran::pairwiseTTests(assay(ref, 2), ref$label, direction="up")
markers <- scran::getTopMarkers(by.t[[1]], by.t[[2]], n=10)
trained <- trainSingleR(ref, ref$label, genes=markers)
length(trained$markers$unique)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.