impute: Quantitative proteomics data imputation

imputeR Documentation

Quantitative proteomics data imputation

Description

The impute method performs data imputation on QFeatures and SummarizedExperiment instance using a variety of methods.

Users should proceed with care when imputing data and take precautions to assure that the imputation produce valid results, in particular with naive imputations such as replacing missing values with 0.

See MsCoreUtils::impute_matrix() for details on the different imputation methods available and strategies.

Usage

impute

## S4 method for signature 'SummarizedExperiment'
impute(object, method, ...)

## S4 method for signature 'QFeatures'
impute(object, method, ..., i, name = "imputedAssay")

Arguments

object

A SummarizedExperiment or QFeatures object with missing values to be imputed.

method

character(1) defining the imputation method. See imputeMethods() for available ones. See MsCoreUtils::impute_matrix() for details.

...

Additional parameters passed to the inner imputation function. See MsCoreUtils::impute_matrix() for details.

i

A logical(1) or a character(1) that defines which element of the QFeatures instance to impute. It cannot be missing and must be of length one.

name

A character(1) naming the new assay name. Default is imputedAssay.

Format

An object of class standardGeneric of length 1.

Examples

MsCoreUtils::imputeMethods()

data(se_na2)
## table of missing values along the rows (proteins)
table(rowData(se_na2)$nNA)
## table of missing values along the columns (samples)
colData(se_na2)$nNA

## non-random missing values
notna <- which(!rowData(se_na2)$randna)
length(notna)
notna

impute(se_na2, method = "min")

if (require("imputeLCMD")) {
  impute(se_na2, method = "QRILC")
  impute(se_na2, method = "MinDet")
}

if (require("norm"))
  impute(se_na2, method = "MLE")

impute(se_na2, method = "mixed",
       randna = rowData(se_na2)$randna,
       mar = "knn", mnar = "QRILC")

## neighbour averaging
x <- se_na2[1:4, 1:6]
assay(x)[1, 1] <- NA ## min value
assay(x)[2, 3] <- NA ## average
assay(x)[3, 1:2] <- NA ## min value and average
## 4th row: no imputation
assay(x)

assay(impute(x, "nbavg"))

rformassspectrometry/QFeatures documentation built on Sept. 18, 2024, 10:39 p.m.