ScpModel-Workflow: Modelling single-cell proteomics data

ScpModel-WorkflowR Documentation

Modelling single-cell proteomics data

Description

Function to estimate a linear model for each feature (peptide or protein) of a single-cell proteomics data set.

Usage

scpModelWorkflow(object, formula, i = 1, name = "model", verbose = TRUE)

scpModelFilterPlot(object, name)

Arguments

object

An object that inherits from the SingleCellExperiment class.

formula

A formula object controlling which variables are to be modelled.

i

A logical, numeric or character indicating which assay of object to use as input for modelling. Only a single assay can be provided. Defaults to the first assays.

name

A character(1) providing the name to use to store or retrieve the modelling results. When retrieving a model and name is missing, the name of the first model found in object is used.

verbose

A logical(1) indicating whether to print progress to the console.

Input data

The main input is object that inherits from the SingleCellExperiment class. The quantitative data will be retrieve using assay(object). If object contains multiple assays, you can specify which assay to take as input thanks to the argument i, the function will then assume assay(object, i) as quantification input .

The objective of modelling single-cell proteomics data is to estimate, for each feature (peptide or protein), the effect of known cell annotations on the measured intensities. These annotations may contain biological information such as the cell line, FACS-derived cell type, treatment, etc. We also highly recommend including technical information, such as the MS acquisition run information or the chemical label (in case of multiplexed experiments). These annotation must be available from colData(object). formula specifies which annotations to use during modelling.

Data modelling workflow

The modelling worflow starts with generating a model matrix for each feature given the colData(object) and formula. The model matrix for peptide i, denoted X_i, is adapted to the pattern of missing values (see section below). Then, the functions fits the model matrix against the quantitative data. In other words, the function determines for each feature i (row in the input data) the contribution of each variable in the model. More formally, the general model definition is:

Y_i = \beta_i X^T_{(i)} + \epsilon_i

where Y is the feature by cell quantification matrix, \beta_i contains the estimated coefficients for feature i with as many coefficients as variables to estimate, X^T_{(i)} is the model matrix generated for feature i, and \epsilon is the feature by cell matrix with residuals.

The coefficients are estimated using penalized least squares regression. Next, the function computes the residual matrix and the effect matrices. An effect matrix contains the data that is captured by a given cell annotation. Formally, for each feature i:

\hat{M^f_i} = \hat{\beta^f_i} X^{fT}_{(i)}

where \hat{M^f} is a cell by feature matrix containing the variables associated to annotation f, \hat{\beta^f_i} are the estimated coefficients associated to annotation f and estimated for feature i, and X^{fT}_{(i)} is the model matrix for peptide i containing only the variables to annotation f.

All the results are stored in an ScpModel object which is stored in the object's metadata. Note that multiple models can be estimated for the same object. In that case, provide the name argument to store the results in a separate ScpModel.

Feature filtering

The proportion of missing values for each features is high in single-cell proteomics data. Many features can typically contain more coefficients to estimate than observed values. These features cannot be estimated and will be ignored during further steps. These features are identified by computing the ratio between the number of observed values and the number of coefficients to estimate. We call it the n/p ratio. Once the model is estimated, use scpModelFilterPlot(object) to explore the distribution of n/p ratios across the features. You can also extract the n/p ratio for each feature using scpModelFilterNPRatio(object). By default, any feature that has an n/p ratio lower than 1 is ignored. However, feature with an n/p ratio close to 1 may lead to unreliable outcome because there are not enough observed data. You could consider the n/p ratio as the average number of replicate per coefficient to estimate. Therefore, you may want to increase the n/p threshold. You can do so using scpModelFilter(object) <- npThreshold.

About missing values

The data modelling workflow is designed to take the presence of missing values into account. We highly recommend to not impute the data before modelling. Instead, the modelling approach will ignore missing values and will generate a model matrix using only the observed values for each feature. However, the model matrices for some features may contain highly correlated variables, leading to near singular designs. We include a small ridge penalty to reduce numerical instability associated to correlated variables.

Author(s)

Christophe Vanderaa, Laurent Gatto

See Also

  • ScpModel for functions to extract information from the ScpModel object

  • ScpModel-VarianceAnalysis, ScpModel-DifferentialAnalysis, ScpModel-ComponentAnalysis to explore the model results

  • scpKeepEffect and scpRemoveBatchEffect to perform batch correction for downstream analyses.

Examples


data("leduc_minimal")
leduc_minimal
## Overview of available cell annotations
colData(leduc_minimal)

####---- Model data ----####

f <- ~ 1 + ## intercept
    Channel + Set + ## batch variables
    MedianIntensity +## normalization
    SampleType ## biological variable
leduc_minimal <- scpModelWorkflow(leduc_minimal, formula = f)

####---- n/p feature filtering ----####

## Get n/p ratios
head(scpModelFilterNPRatio(leduc_minimal))

## Plot n/p ratios
scpModelFilterPlot(leduc_minimal)

## Change n/p ratio threshold
scpModelFilterThreshold(leduc_minimal) <- 2
scpModelFilterPlot(leduc_minimal)

UCLouvain-CBIO/scp documentation built on May 5, 2024, 1:17 a.m.