normalizeSCP: Normalize single-cell proteomics (SCP) data

normalizeSCPR Documentation

Normalize single-cell proteomics (SCP) data

Description

This function normalises an assay in a QFeatures according to the supplied method (see Details). The normalized data is added as a new assay

Usage

normalizeSCP(object, i, name = "normAssay", method, ...)

Arguments

object

An object of class QFeatures.

i

A numeric vector or a character vector giving the index or the name, respectively, of the assay(s) to be processed.

name

A character(1) naming the new assay name. Defaults is are normAssay.

method

character(1) defining the normalisation method to apply. See Details.'

...

Additional parameters passed to MsCoreUtils::normalizeMethods().

Details

The method parameter in normalize can be one of "sum", "max", "center.mean", "center.median", "div.mean", "div.median", "diff.meda", ⁠"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 sweeping 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().

Value

A QFeatures object with an additional assay containing the normalized data.

See Also

QFeatures::normalize for more details about normalize

Examples


data("scp1")
scp1
normalizeSCP(scp1, i = "proteins", name = "normproteins", 
             method = "center.mean")


UCLouvain-CBIO/scp documentation built on Oct. 12, 2024, 2:37 a.m.