modelGeneCV2 | R Documentation |
Model the squared coefficient of variation (CV2) of the normalized expression profiles for each gene, fitting a trend to account for the mean-variance relationship across genes.
modelGeneCV2(x, ...)
## S4 method for signature 'ANY'
modelGeneCV2(
x,
size.factors = NULL,
block = NULL,
subset.row = NULL,
subset.fit = NULL,
...,
equiweight = TRUE,
method = "fisher",
BPPARAM = SerialParam()
)
## S4 method for signature 'SummarizedExperiment'
modelGeneCV2(x, ..., assay.type = "counts")
## S4 method for signature 'SingleCellExperiment'
modelGeneCV2(x, size.factors = NULL, ...)
x |
A numeric matrix of counts where rows are genes and columns are cells. Alternatively, a SummarizedExperiment or SingleCellExperiment containing such a matrix. |
... |
For the generic, further arguments to pass to each method. For the ANY method, further arguments to pass to For the SummarizedExperiment method, further arguments to pass to the ANY method. For the SingleCellExperiment method, further arguments to pass to the SummarizedExperiment method. |
size.factors |
A numeric vector of size factors for each cell in |
block |
A factor specifying the blocking levels for each cell in |
subset.row |
See |
subset.fit |
An argument similar to |
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 |
method |
String specifying how p-values should be combined when |
BPPARAM |
A BiocParallelParam object indicating whether parallelization should be performed across genes. |
assay.type |
String or integer scalar specifying the assay containing the counts. |
For each gene, we compute the CV2 and mean of the counts after scaling them by size.factors
.
A trend is fitted to the CV2 against the mean for all genes using fitTrendCV2
.
The fitted value for each gene is used as a proxy for the technical noise, assuming that most genes exhibit a low baseline level of variation that is not biologically interesting.
The ratio of the total CV2 to the trend is used as a metric to rank interesting genes, with larger ratios being indicative of strong biological heterogeneity.
By default, the trend is fitted using all of the genes in x
.
If subset.fit
is specified, the trend is fitted using only the specified subset,
and the values of trend
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 modelGeneCV2WithSpikes
.
If no size factors are supplied, they are automatically computed depending on the input type:
If size.factors=NULL
for the ANY method, the sum of counts for each cell in x
is used to compute a size factor via the librarySizeFactors
function.
If size.factors=NULL
for the SingleCellExperiment method, sizeFactors(x)
is used if available.
Otherwise, it defaults to library size-derived size factors.
If size.factors
are supplied, they will override any size factors present in x
.
A DataFrame is returned where each row corresponds to a gene in x
(or in subset.row
, if specified).
This contains the numeric fields:
mean
:Mean normalized expression per gene.
total
:Squared coefficient of variation of the normalized expression per gene.
trend
:Fitted value of the trend.
ratio
:Ratio of total
to trend
.
p.value, FDR
:Raw and adjusted p-values for the test against the null hypothesis that ratio<=1
.
If block
is not specified,
the metadata
of the DataFrame contains the output of running fitTrendCV2
on the specified features,
along with the mean
and cv2
used to fit the trend.
If block
is specified,
the output contains another 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 metadata
contains the output of fitTrendCV2
for the cells of that block.
The p-value for each gene is computed by assuming that the CV2 estimates are normally distributed around the trend, and that the standard deviation of the CV2 distribution is proportional to the value of the trend.
This is used to construct a one-sided test for each gene based on its ratio
, under the null hypothesis that the ratio is equal to or less than 1.
The proportionality constant for the standard deviation is set to the std.dev
returned by fitTrendCV2
.
This is estimated from the spread of per-gene CV2 values around the trend, so the null hypothesis effectively becomes “is this gene more variable than other genes of the same abundance?”
Setting block
will estimate the mean and variance of each gene for cells in each level of 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:
For means and CV2 values, this is done by taking the geometric mean across blocking levels.
If 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.
Per-level p-values are combined using combineParallelPValues
according to 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 equiweight
depends on the chosen method.
Blocks with fewer than 2 cells are completely ignored and do not contribute to the combined mean, variance component or p-value.
Aaron Lun
library(scuttle)
sce <- mockSCE()
# Simple case:
spk <- modelGeneCV2(sce)
spk
plot(spk$mean, spk$total, pch=16, log="xy")
curve(metadata(spk)$trend(x), add=TRUE, col="dodgerblue")
# With blocking:
block <- sample(LETTERS[1:2], ncol(sce), replace=TRUE)
blk <- modelGeneCV2(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, pch=16, log="xy")
curve(metadata(current)$trend(x), add=TRUE, col="dodgerblue")
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.