rarefyAssay | R Documentation |
rarefyAssay
randomly subsamples counts within a
SummarizedExperiment
object and returns a new
SummarizedExperiment
containing the original assay and the new
subsampled assay.
rarefyAssay(
x,
assay.type = assay_name,
assay_name = "counts",
sample = min_size,
min_size = min(colSums2(assay(x))),
replace = TRUE,
name = "subsampled",
verbose = TRUE,
...
)
## S4 method for signature 'SummarizedExperiment'
rarefyAssay(
x,
assay.type = assay_name,
assay_name = "counts",
sample = min_size,
min_size = min(colSums2(assay(x, assay.type))),
replace = TRUE,
name = "subsampled",
verbose = TRUE,
...
)
x |
|
assay.type |
|
assay_name |
Deprecated. Use |
sample |
|
min_size |
Deprecated. Use |
replace |
|
name |
|
verbose |
|
... |
additional arguments not used |
Although the subsampling approach is highly debated in microbiome research,
we include the rarefyAssay
function because there may be some
instances where it can be useful.
Note that the output of rarefyAssay
is not the equivalent as the
input and any result have to be verified with the original dataset.
Subsampling/Rarefying may undermine downstream analyses and have unintended consequences. Therefore, make sure this normalization is appropriate for your data.
To maintain the reproducibility, please define the seed using set.seed() before implement this function.
rarefyAssay
return x
with subsampled data.
McMurdie PJ, Holmes S. Waste not, want not: why rarefying microbiome data is inadmissible. PLoS computational biology. 2014 Apr 3;10(4):e1003531.
Gloor GB, Macklaim JM, Pawlowsky-Glahn V & Egozcue JJ (2017) Microbiome Datasets Are Compositional: And This Is Not Optional. Frontiers in Microbiology 8: 2224. doi: 10.3389/fmicb.2017.02224
Weiss S, Xu ZZ, Peddada S, Amir A, Bittinger K, Gonzalez A, Lozupone C, Zaneveld JR, Vázquez-Baeza Y, Birmingham A, Hyde ER. Normalization and microbial differential abundance strategies depend upon data characteristics. Microbiome. 2017 Dec;5(1):1-8.
# When samples in TreeSE are less than specified sample, they will be removed.
# If after subsampling features are not present in any of the samples,
# they will be removed.
data(GlobalPatterns)
tse <- GlobalPatterns
set.seed(123)
tse_subsampled <- rarefyAssay(tse, sample = 60000, name = "subsampled")
tse_subsampled
dim(tse)
dim(assay(tse_subsampled, "subsampled"))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.