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
|
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/vst.R
\name{varianceStabilizingTransformation}
\alias{varianceStabilizingTransformation}
\alias{getVarianceStabilizedData}
\title{Apply a variance stabilizing transformation (VST) to the count data}
\usage{
varianceStabilizingTransformation(object, blind = TRUE, fitType = "parametric")
getVarianceStabilizedData(object)
}
\arguments{
\item{object}{a DESeqDataSet or matrix of counts}
\item{blind}{logical, whether to blind the transformation to the experimental
design. blind=TRUE should be used for comparing samples in a manner unbiased by
prior information on samples, for example to perform sample QA (quality assurance).
blind=FALSE should be used for transforming data for downstream analysis,
where the full use of the design information should be made.
blind=FALSE will skip re-estimation of the dispersion trend, if this has already been calculated.
If many of genes have large differences in counts due to
the experimental design, it is important to set blind=FALSE for downstream
analysis.}
\item{fitType}{in case dispersions have not yet been estimated for \code{object},
this parameter is passed on to \code{\link{estimateDispersions}} (options described there).}
}
\value{
\code{varianceStabilizingTransformation} returns a
\code{\link{DESeqTransform}} if a \code{DESeqDataSet} was provided,
or returns a a matrix if a count matrix was provided.
Note that for \code{\link{DESeqTransform}} output, the matrix of
transformed values is stored in \code{assay(vsd)}.
\code{getVarianceStabilizedData} also returns a matrix.
}
\description{
This function calculates a variance stabilizing transformation (VST) from the
fitted dispersion-mean relation(s) and then transforms the count data (normalized
by division by the size factors or normalization factors), yielding a matrix
of values which are now approximately homoskedastic (having constant variance along the range
of mean values). The transformation also normalizes with respect to library size.
The \code{\link{rlog}} is less sensitive
to size factors, which can be an issue when size factors vary widely.
These transformations are useful when checking for outliers or as input for
machine learning techniques such as clustering or linear discriminant analysis.
}
\details{
For each sample (i.e., column of \code{counts(dds)}), the full variance function
is calculated from the raw variance (by scaling according to the size factor and adding
the shot noise). We recommend a blind estimation of the variance function, i.e.,
one ignoring conditions. This is performed by default, and can be modified using the
'blind' argument.
Note that neither rlog transformation nor the VST are used by the
differential expression estimation in \code{\link{DESeq}}, which always
occurs on the raw count data, through generalized linear modeling which
incorporates knowledge of the variance-mean dependence. The rlog transformation
and VST are offered as separate functionality which can be used for visualization,
clustering or other machine learning tasks. See the transformation section of the
vignette for more details.
The transformation does not require that one has already estimated size factors
and dispersions.
A typical workflow is shown in Section \emph{Variance stabilizing transformation}
in the package vignette.
If \code{\link{estimateDispersions}} was called with:
\code{fitType="parametric"},
a closed-form expression for the variance stabilizing
transformation is used on the normalized
count data. The expression can be found in the file \file{vst.pdf}
which is distributed with the vignette.
\code{fitType="local"},
the reciprocal of the square root of the variance of the normalized counts, as derived
from the dispersion fit, is then numerically
integrated, and the integral (approximated by a spline function) is evaluated for each
count value in the column, yielding a transformed value.
\code{fitType="mean"}, a VST is applied for Negative Binomial distributed counts, 'k',
with a fixed dispersion, 'a': ( 2 asinh(sqrt(a k)) - log(a) - log(4) )/log(2).
In all cases, the transformation is scaled such that for large
counts, it becomes asymptotically (for large values) equal to the
logarithm to base 2 of normalized counts.
The variance stabilizing transformation from a previous dataset
can be "frozen" and reapplied to new samples.
The frozen VST is accomplished by saving the dispersion function
accessible with \code{\link{dispersionFunction}}, assigning this
to the \code{DESeqDataSet} with the new samples, and running
varianceStabilizingTransformation with 'blind' set to FALSE.
Then the dispersion function from the previous dataset will be used
to transform the new sample(s).
Limitations: In order to preserve normalization, the same
transformation has to be used for all samples. This results in the
variance stabilizition to be only approximate. The more the size
factors differ, the more residual dependence of the variance on the
mean will be found in the transformed data. \code{\link{rlog}} is a
transformation which can perform better in these cases.
As shown in the vignette, the function \code{meanSdPlot}
from the package \pkg{vsn} can be used to see whether this is a problem.
}
\examples{
dds <- makeExampleDESeqDataSet(m=6)
vsd <- varianceStabilizingTransformation(dds)
dists <- dist(t(assay(vsd)))
# plot(hclust(dists))
}
\references{
Reference for the variance stabilizing transformation for counts with a dispersion trend:
Simon Anders, Wolfgang Huber: Differential expression analysis for sequence count data. Genome Biology 2010, 11:106. \url{http://dx.doi.org/10.1186/gb-2010-11-10-r106}
}
\seealso{
\code{\link{plotPCA}}, \code{\link{rlog}}, \code{\link{normTransform}}
}
\author{
Simon Anders
}
|