R/QFeatures-processing.R

Defines functions sweepSE .applyTransformation

##' @title QFeatures processing
##'
##' @description
##'
##' This manual page describes common quantitative proteomics data
##' processing methods using [QFeatures] objects. In the following
##' functions, if `object` is of class `QFeatures`, and optional assay
##' index or name `i` can be specified to define the assay (by name of
##' index) on which to operate.
##'
##' The following functions are currently available:
##'
##' - `logTransform(object, base = 2, i, pc = 0)` log-transforms (with
##'   an optional pseudocount offset) the assay(s).
##'
##' - `normalize(object, method, i)` normalises the assay(s) according
##'   to `method` (see Details).
##'
##' - `scaleTransform(object, center = TRUE, scale = TRUE, i)` applies
##'   [base::scale()] to `SummarizedExperiment` and `QFeatures`
##'   objects.
##'
##' - `sweep(x, MARGIN, STATS, FUN = "-", check.margin = TRUE, ...)`
##'   sweeps out array summaries from `SummarizedExperiment` and
##'   `QFeatures` objects. See [base::sweep()] for details.
##'
##' See the *Processing* vignette for examples.
##'
##' @details
##'
##' The `method` parameter in `normalize` can be one of `"sum"`,
##' `"max"`, `"center.mean"`, `"center.median"`, `"div.mean"`,
##' `"div.median"`, `"diff.median"`, `"quantiles`", `"quantiles.robust`"
##' or `"vsn"`. The [MsCoreUtils::normalizeMethods()] function returns
##' a vector of available normalisation methods.
##'
##' - For `"sum"` and `"max"`, each feature's intensity is divided by
##'   the maximum or the sum of the feature respectively. These two
##'   methods are applied along the features (rows).
##'
##' - `"center.mean"` and `"center.median"` center the respective
##'   sample (column) intensities by subtracting the respective column
##'   means or medians. `"div.mean"` and `"div.median"` divide by the
##'   column means or medians. These are equivalent to `sweep`ing the
##'   column means (medians) along `MARGIN = 2` with `FUN = "-"` (for
##'   `"center.*"`) or `FUN = "/"` (for `"div.*"`).
##'
##' - `"diff.median"` centers all samples (columns) so that they all
##'   match the grand median by subtracting the respective columns
##'   medians differences to the grand median.
##'
##' - Using `"quantiles"` or `"quantiles.robust"` applies (robust) quantile
##'   normalisation, as implemented in [preprocessCore::normalize.quantiles()]
##'   and [preprocessCore::normalize.quantiles.robust()]. `"vsn"` uses the
##'   [vsn::vsn2()] function.  Note that the latter also glog-transforms the
##'   intensities.  See respective manuals for more details and function
##'   arguments.
##'
##' For further details and examples about normalisation, see
##' [MsCoreUtils::normalize_matrix()].
##'
##' @param object An object of class `QFeatures` or `SummarizedExperiment`.
##'
##' @param x An object of class `QFeatures` or `SummarizedExperiment`
##'     in `sweep`.
##'
##' @param base `numeric(1)` providing the base with respect to which
##'     logarithms are computed. Defaults is 2.
##'
##' @param pc `numeric(1)` with a pseudocount to add to the
##'     quantitative data. Useful when (true) 0 are present in the
##'     data. Default is 0 (no effect).
##'
##' @param center `logical(1)` (default is `TRUE`) value or
##'     numeric-alike vector of length equal to the number of columns
##'     of `object`. See [base::scale()] for details.
##'
##' @param scale `logical(1)` (default is `TRUE`) or a numeric-alike
##'     vector of length equal to the number of columns of
##'     `object`. See [base::scale()] for details.
##'
##' @param method `character(1)` defining the normalisation method to
##'     apply. See Details.
##'
##' @param i A numeric vector or a character vector giving the index or the
##'     name, respectively, of the assay(s) to be processed.
##'
##' @param name A `character(1)` naming the new assay name. Defaults
##'     are `logAssay` for `logTransform`, `scaledAssay` for
##'     `scaleTranform` and `normAssay` for `normalize`.
##'
##' @param MARGIN As in [base::sweep()], a vector of indices giving the
##'     extent(s) of `x` which correspond to `STATS`.
##'
##' @param STATS As in [base::sweep()], the summary statistic which is
##'     to be swept out.
##'
##' @param FUN As in [base::sweep()], the function to be used to carry
##'     out the sweep.
##'
##' @param check.margin As in [base::sweep()], a `logical`.  If `TRUE`
##'     (the default), warn if the length or dimensions of `STATS` do
##'     not match the specified dimensions of `x`.  Set to `FALSE` for
##'     a small speed gain when you know that dimensions match.
##'
##' @param ... Additional parameters passed to inner functions.
##'
##' @return An processed object of the same class as `x` or `object`.
##'
##' @aliases logTransform logTransform,SummarizedExperiment-method logTransform,QFeatures-method
##' @aliases scaleTransform scaleTransform,SummarizedExperiment-method scaleTransform,QFeatures-method
##' @aliases normalize normalize,SummarizedExperiment-method normalize,QFeatures-method
##' @aliases sweep sweep,SummarizedExperiment-method sweep,QFeatures-method
##' @aliases normalizeMethods
##'
##' @name QFeatures-processing
##'
##' @rdname QFeatures-processing
##'
##' @examples
##'
##' MsCoreUtils::normalizeMethods()
NULL

## -------------------------------------------------------
##   Internal
## -------------------------------------------------------

## Internal function that applies a given function .FUN on one or more
## assays given by i.
##
## @param object A QFeatures object
##
## @param i The index of one or multiple assays
##
## @param name The name of the new assays to add. Must have the same
##     length as i
##
## @param .FUN A function or the name (character(1)) of a function
##     that takes an object that inherits from the
##     SummarizedExperiment class and returns an object of the same
##     class.
##
## @param ... Further argument passed to .FUN
.applyTransformation <- function(object, i, name, .FUN, ...) {
    ## Check arguments
    i <- .normIndex(object, i)
    stopifnot(length(i) == length(name))
    if (is.character(.FUN)) .FUN <- get(.FUN)
    ## Create and add new assays with the transformed data
    for (j in seq_along(i)) {
        from <- i[[j]]
        to <- name[[j]]
        object <- addAssay(object, .FUN(object[[from]], ...), to)
        object <- addAssayLinkOneToOne(object, from = from, to = to)
    }
    object
}

## -------------------------------------------------------
##   Transformations
## -------------------------------------------------------

##' @exportMethod logTransform
##' @rdname QFeatures-processing
setMethod("logTransform",
          "SummarizedExperiment",
          function(object, base = 2, pc = 0) {
              assay(object) <- log(assay(object) + pc, base)
              object
          })

##' @rdname QFeatures-processing
setMethod("logTransform",
          "QFeatures",
          function(object, i, name = "logAssay", base = 2, pc = 0) {
              .applyTransformation(object, i, name, logTransform,
                                   base = base, pc = pc)
          })

##' @exportMethod scaleTransform
##' @rdname QFeatures-processing
setMethod("scaleTransform", "SummarizedExperiment",
          function(object, center = TRUE, scale = TRUE) {
              e <- scale(assay(object), center = center, scale = scale)
              attr(e, "scaled:center") <- NULL
              attr(e, "scaled:scale") <- NULL
              assay(object) <- e
              object
          })

##' @rdname QFeatures-processing
setMethod("scaleTransform", "QFeatures",
          function(object, i, name = "scaledAssay", center = TRUE, scale = TRUE) {
              .applyTransformation(object, i, name, scaleTransform,
                                   center = center, scale = scale)
          })

## -------------------------------------------------------
##   Normalisation (normalize)
## -------------------------------------------------------

##' @importFrom BiocGenerics normalize
##' @exportMethod normalize
##' @rdname QFeatures-processing
setMethod("normalize", "SummarizedExperiment",
          function(object,
                   method,
                   ...) {
              e <- MsCoreUtils::normalize_matrix(assay(object), method, ...)
              rownames(e) <- rownames(assay(object))
              colnames(e) <- colnames(assay(object))
              assay(object) <- e
              object
          })

## normalise <- normalize

##' @rdname QFeatures-processing
setMethod("normalize", "QFeatures",
          function(object, i, name = "normAssay", method, ...) {
              .applyTransformation(object, i, name, normalize,
                                   method = method, ...)
          })


## -------------------------------------------------------
##   Sweep
## -------------------------------------------------------

sweepSE <- function(x, MARGIN, STATS, FUN = "-", check.margin = TRUE, ...) {
    e <- base::sweep(assay(x), MARGIN, STATS, FUN, check.margin, ...)
    rownames(e) <- rownames(assay(x))
    colnames(e) <- colnames(assay(x))
    assay(x) <- e
    x
}

##' @exportMethod sweep
##' @rdname QFeatures-processing
setMethod("sweep", "SummarizedExperiment",
          function(x, MARGIN, STATS, FUN = "-", check.margin = TRUE, ...)
              sweepSE(x, MARGIN, STATS, FUN, check.margin, ...))


##' @rdname QFeatures-processing
setMethod("sweep", "QFeatures",
          function(x, MARGIN, STATS, FUN = "-", check.margin = TRUE, ..., i, name = "sweptAssay") {
              .applyTransformation(x, i, name, sweep,
                                   MARGIN = MARGIN, STATS = STATS,
                                   FUN = FUN,
                                   check.margin = check.margin, ...)
          })
lgatto/Features documentation built on Sept. 22, 2024, 7:13 p.m.