R/methods.R
, R/methods_SE.R
quantile_normalise_abundance-methods.Rd
quantile_normalise_abundance() takes as input A `tbl` (with at least three columns for sample, feature and transcript abundance) or `SummarizedExperiment` (more convenient if abstracted to tibble with library(tidySummarizedExperiment)) and Scales transcript abundance compansating for sequencing depth (e.g., with TMM algorithm, Robinson and Oshlack doi.org/10.1186/gb-2010-11-3-r25).
quantile_normalise_abundance(
.data,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "limma_normalize_quantiles",
target_distribution = NULL,
action = "add"
)
# S4 method for class 'spec_tbl_df'
quantile_normalise_abundance(
.data,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "limma_normalize_quantiles",
target_distribution = NULL,
action = "add"
)
# S4 method for class 'tbl_df'
quantile_normalise_abundance(
.data,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "limma_normalize_quantiles",
target_distribution = NULL,
action = "add"
)
# S4 method for class 'tidybulk'
quantile_normalise_abundance(
.data,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "limma_normalize_quantiles",
target_distribution = NULL,
action = "add"
)
# S4 method for class 'SummarizedExperiment'
quantile_normalise_abundance(
.data,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "limma_normalize_quantiles",
target_distribution = NULL,
action = NULL
)
# S4 method for class 'RangedSummarizedExperiment'
quantile_normalise_abundance(
.data,
.sample = NULL,
.transcript = NULL,
.abundance = NULL,
method = "limma_normalize_quantiles",
target_distribution = NULL,
action = NULL
)
A `tbl` (with at least three columns for sample, feature and transcript abundance) or `SummarizedExperiment` (more convenient if abstracted to tibble with library(tidySummarizedExperiment))
The name of the sample column
The name of the transcript/gene column
The name of the transcript/gene abundance column
A character string. Either "limma_normalize_quantiles" for limma::normalizeQuantiles or "preprocesscore_normalize_quantiles_use_target" for preprocessCore::normalize.quantiles.use.target for large-scale datasets.
A numeric vector. If NULL the target distribution will be calculated by preprocessCore. This argument only affects the "preprocesscore_normalize_quantiles_use_target" method.
A character string between "add" (default) and "only". "add" joins the new information to the input tbl (default), "only" return a non-redundant tbl with the just new information.
A tbl object with additional columns with scaled data as `<NAME OF COUNT COLUMN>_scaled`
A tbl object with additional columns with scaled data as `<NAME OF COUNT COLUMN>_scaled`
A tbl object with additional columns with scaled data as `<NAME OF COUNT COLUMN>_scaled`
A tbl object with additional columns with scaled data as `<NAME OF COUNT COLUMN>_scaled`
A `SummarizedExperiment` object
A `SummarizedExperiment` object
`r lifecycle::badge("maturing")`
Tranform the feature abundance across samples so to have the same quantile distribution (using preprocessCore).
Underlying method
If `limma_normalize_quantiles` is chosen
.data |>limma::normalizeQuantiles()
If `preprocesscore_normalize_quantiles_use_target` is chosen
.data |> preprocessCore::normalize.quantiles.use.target( target = preprocessCore::normalize.quantiles.determine.target(.data) )
tidybulk::se_mini |>
quantile_normalise_abundance()
#> class: SummarizedExperiment
#> dim: 527 5
#> metadata(0):
#> assays(2): count count_scaled
#> rownames(527): ABCB4 ABCB9 ... ZNF324 ZNF442
#> rowData names(1): entrez
#> colnames(5): SRR1740034 SRR1740035 SRR1740043 SRR1740058 SRR1740067
#> colData names(5): Cell.type time condition days dead