File: rlog.Rd

package info (click to toggle)
r-bioc-deseq2 1.14.1-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 2,816 kB
  • ctags: 9
  • sloc: cpp: 374; sh: 14; makefile: 2
file content (136 lines) | stat: -rw-r--r-- 6,259 bytes parent folder | download
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
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/rlog.R
\name{rlog}
\alias{rlog}
\alias{rlogTransformation}
\title{Apply a 'regularized log' transformation}
\usage{
rlog(object, blind = TRUE, intercept, betaPriorVar, fitType = "parametric")

rlogTransformation(object, blind = TRUE, intercept, betaPriorVar,
  fitType = "parametric")
}
\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 an 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{intercept}{by default, this is not provided and calculated automatically.
if provided, this should be a vector as long as the number of rows of object,
which is log2 of the mean normalized counts from a previous dataset.
this will enforce the intercept for the GLM, allowing for a "frozen" rlog
transformation based on a previous dataset.
You will also need to provide \code{mcols(object)$dispFit}.}

\item{betaPriorVar}{a single value, the variance of the prior on the sample
betas, which if missing is estimated from the data}

\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{
a \code{\link{DESeqTransform}} if a \code{DESeqDataSet} was provided,
or a matrix if a count matrix was provided as input.
Note that for \code{\link{DESeqTransform}} output, the matrix of
transformed values is stored in \code{assay(rld)}.
To avoid returning matrices with NA values, in the case of a row
of all zeros, the rlog transformation returns zeros
(essentially adding a pseudocount of 1 only to these rows).
}
\description{
This function transforms the count data to the log2 scale in a way 
which minimizes differences between samples for rows with small counts,
and which normalizes with respect to library size.
The rlog transformation produces a similar variance stabilizing effect as
\code{\link{varianceStabilizingTransformation}},
though \code{rlog} is more robust in the
case when the size factors vary widely.
The transformation is useful when checking for outliers
or as input for machine learning techniques
such as clustering or linear discriminant analysis.
\code{rlog} takes as input a \code{\link{DESeqDataSet}} and returns a
\code{\link{RangedSummarizedExperiment}} object.
}
\details{
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.

The regularization is on the log fold changes of the count for each sample
over an intercept, for each gene. As nearby count values for low counts genes
are almost as likely as the observed count, the rlog shrinkage is greater for low counts.
For high counts, the rlog shrinkage has a much weaker effect.
The fitted dispersions are used rather than the MAP dispersions
(so similar to the \code{\link{varianceStabilizingTransformation}}).

The prior variance for the shrinkag of log fold changes is calculated as follows: 
a matrix is constructed of the logarithm of the counts plus a pseudocount of 0.5,
the log of the row means is then subtracted, leaving an estimate of
the log fold changes per sample over the fitted value using only an intercept.
The prior variance is then calculated by matching the upper quantiles of the observed 
log fold change estimates with an upper quantile of the normal distribution.
A GLM fit is then calculated using this prior. It is also possible to supply the variance of the prior.
See the vignette for an example of the use and a comparison with \code{varianceStabilizingTransformation}.

The transformed values, rlog(K), are equal to
\eqn{rlog(K_{ij}) = \log_2(q_{ij}) = \beta_{i0} + \beta_{ij}}{rlog(K_ij) = log2(q_ij) = beta_i0 + beta_ij},
with formula terms defined in \code{\link{DESeq}}.

The parameters of the rlog transformation from a previous dataset
can be frozen and reapplied to new samples. See the 'Data quality assessment'
section of the vignette for strategies to see if new samples are
sufficiently similar to previous datasets. 
The frozen rlog is accomplished by saving the dispersion function,
beta prior variance and the intercept from a previous dataset,
and running \code{rlog} with 'blind' set to FALSE
(see example below).
}
\examples{

dds <- makeExampleDESeqDataSet(m=6,betaSD=1)
rld <- rlog(dds)
dists <- dist(t(assay(rld)))
plot(hclust(dists))

# run the rlog transformation on one dataset
design(dds) <- ~ 1
dds <- estimateSizeFactors(dds)
dds <- estimateDispersions(dds)
rld <- rlog(dds, blind=FALSE)

# apply the parameters to a new sample

ddsNew <- makeExampleDESeqDataSet(m=1)
mcols(ddsNew)$dispFit <- mcols(dds)$dispFit
betaPriorVar <- attr(rld,"betaPriorVar")
intercept <- mcols(rld)$rlogIntercept
rldNew <- rlog(ddsNew, blind=FALSE,
               intercept=intercept,
               betaPriorVar=betaPriorVar)
                           

}
\references{
Reference for regularized logarithm (rlog):

Michael I Love, Wolfgang Huber, Simon Anders: Moderated estimation of fold change and dispersion for RNA-seq data with DESeq2. Genome Biology 2014, 15:550. \url{http://dx.doi.org/10.1186/s13059-014-0550-8}
}
\seealso{
\code{\link{plotPCA}}, \code{\link{varianceStabilizingTransformation}}, \code{\link{normTransform}}
}