Title: | Wrench normalization for sparse count data |
---|---|
Description: | Wrench is a package for normalization sparse genomic count data, like that arising from 16s metagenomic surveys. |
Authors: | Senthil Kumar Muthiah [aut], Hector Corrada Bravo [aut, cre] |
Maintainer: | Hector Corrada Bravo <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.25.0 |
Built: | 2024-11-19 04:42:11 UTC |
Source: | https://github.com/bioc/Wrench |
Obtain robust means. .
.estimSummary(res, estim.type = "s2.w.mean", ...)
.estimSummary(res, estim.type = "s2.w.mean", ...)
res |
result structure of |
estim.type |
estimator type |
... |
other parameters |
a chosen summary statistic
Log Postive-conditional weight computations for wrench estimators.
.getCondLogWeights(res)
.getCondLogWeights(res)
res |
result structure of |
inverse variance weights when using positive conditional models.
Postive-conditional weight computations for wrench estimators.
.getCondWeights(res)
.getCondWeights(res)
res |
result structure of |
positive conditional weights for each sample
This function is used to derive weights for feature-wise compositional estimates. Our (default) intention is to derive these based on average occurrences across the dataset, as just a function of sample depth, and not with particular relevance to groups.
.getHurdle(mat, hdesign = model.matrix(~-1 + log(colSums(mat))), pres.abs.mod = TRUE, thresh = FALSE, thresh.val = 1e-08, ...)
.getHurdle(mat, hdesign = model.matrix(~-1 + log(colSums(mat))), pres.abs.mod = TRUE, thresh = FALSE, thresh.val = 1e-08, ...)
mat |
count matrix |
hdesign |
design matrix for the logistic; the default is usually sufficient. |
pres.abs.mod |
TRUE if glm regression is for presence or absence. FALSE if glm regression is for counts. |
thresh |
TRUE if numerically one/zero probability occurrences must be thresholded |
thresh.val |
if thresh is true, the numerically one/zero probability occurrences is thresholded to this value |
... |
other parameters |
A list with components:
pi0.fit - list with feature-wise glm.fit objects
pi0 - matrix with fitted probabilities
Marginal weight computations for wrench estimators.
.getMargWeights(res, z.adj, ...)
.getMargWeights(res, z.adj, ...)
res |
result structure of |
z.adj |
TRUE if the result structure was generated with |
... |
other parameters |
inverse marginal variances for robust mean computing
This function generates the reference.
.getReference(mat, ref.est = "sw.means", ...)
.getReference(mat, ref.est = "sw.means", ...)
mat |
count matrix; rows are features and columns are samples |
ref.est |
reference estimate method |
... |
other parameters |
the reference to be used for normalization
Obtain variances of logged counts.
.gets2(mat, design = model.matrix(mat[1, ] ~ 1), plot = FALSE, ebs2 = TRUE, smoothed = FALSE, ...)
.gets2(mat, design = model.matrix(mat[1, ] ~ 1), plot = FALSE, ebs2 = TRUE, smoothed = FALSE, ...)
mat |
count matrix; rows are features and columns are samples. |
design |
model matrix for the count matrix |
plot |
if the mean-variance trend function (the same as that of voom) needs to be plot. |
ebs2 |
if regularization of variances needs to be performed. |
smoothed |
TRUE if all the variance estimates must be based on the mean-variance trend function. |
... |
other parameters |
a vector with variance estimates for logged feature-wise counts.
Get weighted means for matrix
.getWeightedMean(mat, w = rep(1, nrow(mat)))
.getWeightedMean(mat, w = rep(1, nrow(mat)))
mat |
input matrix |
w |
weights |
column-wise weighted means.
Get weighted median for matrix
.getWeightedMedian(mat, w = rep(1, nrow(mat)))
.getWeightedMedian(mat, w = rep(1, nrow(mat)))
mat |
input matrix |
w |
weights |
column-wise weighted means.
Obtain normalization factors for sparse, under-sampled count data that often arise with metagenomic count data.
wrench(mat, condition, etype = "w.marg.mean", ebcf = TRUE, z.adj = FALSE, phi.adj = TRUE, detrend = FALSE, ...)
wrench(mat, condition, etype = "w.marg.mean", ebcf = TRUE, z.adj = FALSE, phi.adj = TRUE, detrend = FALSE, ...)
mat |
count matrix; rows are features and columns are samples |
condition |
a vector with group information on the samples |
etype |
weighting strategy with the following options:
|
ebcf |
TRUE if empirical bayes regularization of ratios needs to be performed. Default recommended. |
z.adj |
TRUE if the feature-wise ratios need to be adjusted by hurdle probabilities (arises when taking marginal expectation). Default recommended. |
phi.adj |
TRUE if estimates need to be adjusted for variance terms (arises when considering positive-part expectations). Default recommended. |
detrend |
FALSE if any linear dependence between sample-depth and compositional factors needs to be removed. (setting this to TRUE reduces variation in compositional factors and can improve accuracy, but requires an extra assumption that no linear dependence between compositional factors and sample depth is present in samples). |
... |
other parameters |
a list
with components:
nf
, normalization factors for samples passed.
Samples with zero total counts are removed from output.
ccf
, compositional correction factors.
Samples with zero total counts are removed from output.
others
, a list
with results from intermediate computations.
qref
, reference chosen.
design
, design matrix used for computation of positive-part parameters.
s2
, feature-wise variances of logged count data.
r
, (regularized) ratios of feature-wise proportions.
radj
, adjustments made to the regularized ratios based
on z.adj and phi.adj settings.
M. Senthil Kumar
#Obtain counts matrix and some group information require(metagenomeSeq) data(mouseData) cntsMatrix <- MRcounts(mouseData) group <- pData(mouseData)$diet #Running wrench with defaults W <- wrench( cntsMatrix, condition=group ) compositionalFactors <- W$ccf normalizationFactors <- W$nf #Introducing the above normalization factors for the most # commonly used tools is shown below. #If using metagenomeSeq normalizedObject <- mouseData normFactors(normalizedObject) <- normalizationFactors #If using edgeR, we must pass in the compositional factors require(edgeR) edgerobj <- DGEList( counts=cntsMatrix, group = as.matrix(group), norm.factors=compositionalFactors ) #If using DESeq/DESeq2 require(DESeq2) deseq.obj <- DESeqDataSetFromMatrix(countData = cntsMatrix, DataFrame(group), ~ group ) DESeq2::sizeFactors(deseq.obj) <- normalizationFactors
#Obtain counts matrix and some group information require(metagenomeSeq) data(mouseData) cntsMatrix <- MRcounts(mouseData) group <- pData(mouseData)$diet #Running wrench with defaults W <- wrench( cntsMatrix, condition=group ) compositionalFactors <- W$ccf normalizationFactors <- W$nf #Introducing the above normalization factors for the most # commonly used tools is shown below. #If using metagenomeSeq normalizedObject <- mouseData normFactors(normalizedObject) <- normalizationFactors #If using edgeR, we must pass in the compositional factors require(edgeR) edgerobj <- DGEList( counts=cntsMatrix, group = as.matrix(group), norm.factors=compositionalFactors ) #If using DESeq/DESeq2 require(DESeq2) deseq.obj <- DESeqDataSetFromMatrix(countData = cntsMatrix, DataFrame(group), ~ group ) DESeq2::sizeFactors(deseq.obj) <- normalizationFactors