se_normalize | R Documentation |
Normalize SummarizedExperiment data
se_normalize(
se,
method = c("quantile", "jammanorm", "limma_batch_adjust", "TMM", "TMMwsp", "RLE"),
assay_names = NULL,
output_method_prefix = NULL,
output_assay_names = NULL,
genes = NULL,
samples = NULL,
params = list(quantile = list(ties = TRUE), jammanorm = list(controlGenes = NULL,
minimum_mean = 0, controlSamples = NULL, centerGroups = NULL, useMedian = FALSE,
noise_floor = NULL, noise_floor_value = NULL), limma_batch_adjust = list(batch =
NULL, group = NULL), TMM = list(refColumn = NULL, logratioTrim = 0.3, sumTrim = 0.05,
doWeighting = TRUE, Acutoff = NULL), TMMwsp = list(refColumn = NULL, logratioTrim =
0.3, sumTrim = 0.05, doWeighting = TRUE, Acutoff = NULL), RLE = list(refColumn =
NULL, logratioTrim = 0.3,
sumTrim = 0.05, doWeighting = TRUE, Acutoff = NULL)),
normgroup = NULL,
floor = 0,
enforce_norm_floor = TRUE,
output_sep = "_",
override = TRUE,
populate_mcols = TRUE,
verbose = FALSE,
...
)
se |
|
method |
|
assay_names |
|
output_method_prefix |
Consider these arguments: assay_name="counts", method="limma_batch_adjust", output_method_prefix="lba" The assay_name created during normalization will be |
output_assay_names |
Therefore the order of |
genes |
|
samples |
|
params |
|
normgroup |
|
output_sep |
|
override |
|
populate_mcols |
|
verbose |
|
... |
additional arguments are passed to |
This function applies one or more data normalization methods
to an input SummarizedExperiment
object. The normalization is
applied to one or more matrix data stored in assays(se)
,
each one is run independently.
Note that supplying genes
and samples
will apply normalization
to only those genes
and samples
, and this data will be
stored in the full SummarizedExperiment
object se
with
NA
values used to fill any values not present in genes
or samples
.
For example if assay_names
contains two assay names,
and method
contains two methods, the output will include
four normalizations, where each assay name is normalized two ways.
The output assay names will be something like "assay1_method1"
,
"assay1_method2"
, "assay2_method1"
, "assay2_method2"
.
It is not always necessary to normalize data by multiple different
methods, however when two methods are similar and need to be
compared, the SummarizedExperiment
object is a convenient
place to store different normalization results for downstream
comparison. Further, the method se_contrast_stats()
is able
to apply equivalent statistical contrasts to each normalization,
and returns an array of statistical hits which is convenient
for direct comparison of results.
This method calls matrix_normalize()
to perform each normalization
step, see that function description for details on each method.
SummarizedExperiment
object where the normalized output
is added to assays(se)
using the naming format method_assayname
.
Other jamses SE utilities:
make_se_test()
,
se_collapse_by_column()
,
se_collapse_by_row()
,
se_detected_rows()
,
se_rbind()
,
se_to_rowcoldata()
if (jamba::check_pkg_installed("farrisdata")) {
# se_normalize
# suppressPackageStartupMessages(library(SummarizedExperiment))
GeneSE <- farrisdata::farrisGeneSE;
samples <- colnames(GeneSE);
genes <- rownames(GeneSE);
GeneSE <- se_normalize(GeneSE,
genes=genes,
samples=samples,
assay_names=c("raw_counts", "counts"),
method="jammanorm",
params=list(jammanorm=list(minimum_mean=5)))
SummarizedExperiment::mcols(SummarizedExperiment::assays(GeneSE))
names(SummarizedExperiment::assays(GeneSE))
# review normalization factor values
round(digits=3, attr(
SummarizedExperiment::assays(GeneSE)$jammanorm_raw_counts, "nf"))
# the data in "counts" was already normalized
# so the normalization factors are very near 0 as expected
round(digits=3,
attr(SummarizedExperiment::assays(GeneSE)$jammanorm_counts, "nf"))
# note that housekeeper genes are supplied in params
# also this demonstrates output_method_prefix
set.seed(123);
hkgenes <- sample(rownames(GeneSE), 1000)
GeneSE <- se_normalize(GeneSE,
genes=genes,
samples=samples,
assay_names=c("raw_counts"),
method="jammanorm",
output_method_prefix="hkjammanorm",
params=list(jammanorm=list(minimum_mean=5,
controlGenes=hkgenes)))
SummarizedExperiment::mcols(SummarizedExperiment::assays(GeneSE))
# example showing quantile normalization
GeneSE <- se_normalize(GeneSE,
assay_names=c("raw_counts"),
method="quantile")
SummarizedExperiment::mcols(SummarizedExperiment::assays(GeneSE))
# example showing quantile normalization with custom output_assay_names
GeneSE <- se_normalize(GeneSE,
assay_names=c("raw_counts"),
method="quantile",
output_assay_names="newquantile_raw_counts")
SummarizedExperiment::mcols(SummarizedExperiment::assays(GeneSE))
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.