1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173
|
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/modelGeneVar.R
\name{modelGeneVar}
\alias{modelGeneVar}
\alias{modelGeneVar,ANY-method}
\alias{modelGeneVar,SingleCellExperiment-method}
\alias{modelGeneVar,SummarizedExperiment-method}
\title{Model the per-gene variance}
\usage{
\S4method{modelGeneVar}{ANY}(
x,
block = NULL,
design = NULL,
subset.row = NULL,
subset.fit = NULL,
...,
equiweight = TRUE,
method = "fisher",
BPPARAM = SerialParam()
)
\S4method{modelGeneVar}{SummarizedExperiment}(x, ..., assay.type = "logcounts")
}
\arguments{
\item{x}{A numeric matrix of log-normalized expression values where rows are genes and columns are cells.
Alternatively, a \linkS4class{SummarizedExperiment} containing such a matrix.}
\item{block}{A factor specifying the blocking levels for each cell in \code{x}.
If specified, variance modelling is performed separately within each block and statistics are combined across blocks.}
\item{design}{A numeric matrix containing blocking terms for uninteresting factors of variation.}
\item{subset.row}{See \code{?"\link{scran-gene-selection}"}, specifying the rows for which to model the variance.
Defaults to all genes in \code{x}.}
\item{subset.fit}{An argument similar to \code{subset.row}, specifying the rows to be used for trend fitting.
Defaults to \code{subset.row}.}
\item{...}{For the generic, further arguments to pass to each method.
For the ANY method, further arguments to pass to \code{\link{fitTrendVar}}.
For the \linkS4class{SummarizedExperiment} method, further arguments to pass to the ANY method.}
\item{equiweight}{A logical scalar indicating whether statistics from each block should be given equal weight.
Otherwise, each block is weighted according to its number of cells.
Only used if \code{block} is specified.}
\item{method}{String specifying how p-values should be combined when \code{block} is specified, see \code{\link{combineParallelPValues}}.}
\item{BPPARAM}{A \linkS4class{BiocParallelParam} object indicating whether parallelization should be performed across genes.}
\item{assay.type}{String or integer scalar specifying the assay containing the log-expression values.}
}
\value{
A \linkS4class{DataFrame} is returned where each row corresponds to a gene in \code{x} (or in \code{subset.row}, if specified).
This contains the numeric fields:
\describe{
\item{\code{mean}:}{Mean normalized log-expression per gene.}
\item{\code{total}:}{Variance of the normalized log-expression per gene.}
\item{\code{bio}:}{Biological component of the variance.}
\item{\code{tech}:}{Technical component of the variance.}
\item{\code{p.value, FDR}:}{Raw and adjusted p-values for the test against the null hypothesis that \code{bio<=0}.}
}
If \code{block} is not specified,
the \code{metadata} of the DataFrame contains the output of running \code{\link{fitTrendVar}} on the specified features,
along with the \code{mean} and \code{var} used to fit the trend.
If \code{block} is specified,
the output contains another \code{per.block} field.
This field is itself a DataFrame of DataFrames, where each internal DataFrame contains statistics for the variance modelling within each block and has the same format as described above.
Each internal DataFrame's \code{metadata} contains the output of \code{\link{fitTrendVar}} for the cells of that block.
}
\description{
Model the variance of the log-expression profiles for each gene,
decomposing it into technical and biological components based on a fitted mean-variance trend.
}
\details{
For each gene, we compute the variance and mean of the log-expression values.
A trend is fitted to the variance against the mean for all genes using \code{\link{fitTrendVar}}.
The fitted value for each gene is used as a proxy for the technical component of variation for each gene,
under the assumption that most genes exhibit a low baseline level of variation that is not biologically interesting.
The biological component of variation for each gene is defined as the the residual from the trend.
Ranking genes by the biological component enables identification of interesting genes for downstream analyses
in a manner that accounts for the mean-variance relationship.
We use log-transformed expression values to blunt the impact of large positive outliers and to ensure that large variances are driven by strong log-fold changes between cells rather than differences in counts.
Log-expression values are also used in downstream analyses like PCA, so modelling them here avoids inconsistencies with different quantifications of variation across analysis steps.
By default, the trend is fitted using all of the genes in \code{x}.
If \code{subset.fit} is specified, the trend is fitted using only the specified subset,
and the technical components for all other genes are determined by extrapolation or interpolation.
This could be used to perform the fit based on genes that are known to have low variance, thus weakening the assumption above.
Note that this does not refer to spike-in transcripts, which should be handled via \code{\link{modelGeneVarWithSpikes}}.
}
\section{Handling uninteresting factors}{
Setting \code{block} will estimate the mean and variance of each gene for cells in each level of \code{block} separately.
The trend is fitted separately for each level, and the variance decomposition is also performed separately.
Per-level statistics are then combined to obtain a single value per gene:
\itemize{
\item For means and variance components, this is done by averaging values across levels.
If \code{equiweight=FALSE}, a weighted average is used where the value for each level is weighted by the number of cells.
By default, all levels are equally weighted when combining statistics.
\item Per-level p-values are combined using \code{\link{combineParallelPValues}} according to \code{method}.
By default, Fisher's method is used to identify genes that are highly variable in any batch.
Whether or not this is responsive to \code{equiweight} depends on the chosen method.
\item Blocks with fewer than 2 cells are completely ignored and do not contribute to the combined mean, variance component or p-value.
}
Use of \code{block} is the recommended approach for accounting for any uninteresting categorical factor of variation.
In addition to accounting for systematic differences in expression between levels of the blocking factor,
it also accommodates differences in the mean-variance relationships.
Alternatively, uninteresting factors can be used to construct a design matrix to pass to the function via \code{design}.
In this case, a linear model is fitted to the expression profile for each gene and the residual variance is calculated.
This approach is useful for covariates or additive models that cannot be expressed as a one-way layout for use in \code{block}.
However, it assumes that the error is normally distributed with equal variance for all observations of a given gene.
Use of \code{block} and \code{design} together is currently not supported and will lead to an error.
}
\section{Computing p-values}{
The p-value for each gene is computed by assuming that the variance estimates are normally distributed around the trend, and that the standard deviation of the variance distribution is proportional to the value of the trend.
This is used to construct a one-sided test for each gene based on its \code{bio}, under the null hypothesis that the biological component is equal to zero.
The proportionality constant for the standard deviation is set to the \code{std.dev} returned by \code{\link{fitTrendVar}}.
This is estimated from the spread of per-gene variance estimates around the trend, so the null hypothesis effectively becomes \dQuote{is this gene \emph{more} variable than other genes of the same abundance?}
}
\examples{
library(scuttle)
sce <- mockSCE()
sce <- logNormCounts(sce)
# Fitting to all features.
allf <- modelGeneVar(sce)
allf
plot(allf$mean, allf$total)
curve(metadata(allf)$trend(x), add=TRUE, col="dodgerblue")
# Using a subset of features for fitting.
subf <- modelGeneVar(sce, subset.fit=1:100)
subf
plot(subf$mean, subf$total)
curve(metadata(subf)$trend(x), add=TRUE, col="dodgerblue")
points(metadata(subf)$mean, metadata(subf)$var, col="red", pch=16)
# With blocking.
block <- sample(LETTERS[1:2], ncol(sce), replace=TRUE)
blk <- modelGeneVar(sce, block=block)
blk
par(mfrow=c(1,2))
for (i in colnames(blk$per.block)) {
current <- blk$per.block[[i]]
plot(current$mean, current$total)
curve(metadata(current)$trend(x), add=TRUE, col="dodgerblue")
}
}
\seealso{
\code{\link{fitTrendVar}}, for the trend fitting options.
\code{\link{modelGeneVarWithSpikes}}, for modelling variance with spike-in controls.
}
\author{
Aaron Lun
}
|