xCell2Analysis: Perform Cell Type Enrichment Analysis

View source: R/xCell2Analysis.R

xCell2AnalysisR Documentation

Perform Cell Type Enrichment Analysis

Description

Estimates the relative enrichment of cell types in a bulk gene expression mixture. This function uses gene signatures from a pre-trained xCell2Object to compute enrichment scores, with options for linear transformation and spillover correction to improve specificity.

Usage

xCell2Analysis(
  mix,
  xcell2object,
  minSharedGenes = 0.9,
  rawScores = FALSE,
  spillover = TRUE,
  spilloverAlpha = 0.5,
  BPPARAM = BiocParallel::SerialParam()
)

Arguments

mix

A bulk mixture of gene expression data (genes in rows, samples in columns). The input must use the same gene annotation system as the reference object.

xcell2object

A pre-trained reference object of class xCell2Object, created using xCell2Train. Pre-trained references for common cases are provided within the package.

minSharedGenes

Minimum fraction of shared genes required between the mixture and the reference object (default: 0.9). If the shared fraction is below this threshold, the function stops with an error or warning, as sufficient overlap is necessary for accurate analysis.

rawScores

Logical; if TRUE, returns raw enrichment scores without transformation or spillover correction (default: FALSE).

spillover

Logical; enables spillover correction on enrichment scores (default: TRUE). Spillover occurs when closely related cell types share gene expression patterns, inflating enrichment scores. Correction enhances specificity, particularly for related cell types.

spilloverAlpha

Numeric value controlling spillover correction strength (default: 0.5). Lower values apply weaker correction, while higher values apply stronger correction.

BPPARAM

A BiocParallelParam instance to define parallelization strategy (see "Details"). Default is BiocParallel::SerialParam().

Details

The xCell2Analysis function performs cell type enrichment analysis by leveraging gene signatures from a pre-trained xCell2Object. It computes enrichment scores for each cell type in the provided bulk gene expression mixture (mix), applies linear transformations, and corrects for spillover. Spillover correction addresses the overlap of gene expression patterns between closely related cell types, improving the specificity of the enrichment scores.

## Parallelization with BPPARAM To achieve faster processing by running computations in parallel, xCell2Analysis supports parallelization through the BPPARAM parameter. Users can define a parallelization strategy using BiocParallelParam from the BiocParallel package. For example, MulticoreParam is suitable for multi-core processing on Linux and macOS, while SnowParam or SerialParam are better suited for Windows systems. Refer to the BiocParallel documentation for further guidance on parallelization strategies.

## Relationship with Other Function(s) The pre-trained xCell2Object used in xCell2Analysis is created via the xCell2Train function.

The xCell2Analysis function computes enrichment scores for each cell type using gene signatures from a pre-trained xCell2Object. Linear transformations and spillover corrections refine the results, improving specificity when cell types have overlapping gene expression patterns.

Parallelization with BPPARAM: Computations can be parallelized using the BPPARAM parameter. Supported strategies include:

  • MulticoreParam for multi-core processing (Linux/macOS).

  • SnowParam or SerialParam for Windows systems.

See the BiocParallel documentation.

Relationship with Other Functions: The input reference object (xCell2Object) is created via xCell2Train.

Value

A data frame containing the cell type enrichment for each sample in the input matrix, as estimated by xCell2. Each row corresponds to a cell type, and each column corresponds to a sample.

A data frame containing enrichment scores for each cell type and sample. Rows correspond to cell types and columns to samples.

Author(s)

Almog Angel and Dvir Aran

See Also

xCell2Train, for creating the reference object used in this analysis.

Examples

# For detailed example read xCell2 vignette.

library(xCell2)

# Load "ready to use" xCell2 reference object or generate a new one using `xCell2Train`
data(DICE_demo.xCell2Ref, package = "xCell2")

# Load demo bulk RNA-Seq gene expression mixture

# For detailed examples, see the xCell2 vignette.

library(xCell2)

# Load pre-trained reference object
data(DICE_demo.xCell2Ref, package = "xCell2")

# Load demo bulk gene expression mixture
data(mix_demo, package = "xCell2")

# Perform cell type enrichment analysis
xcell2_res <- xCell2::xCell2Analysis(
  mix = mix_demo, 
  xcell2object = DICE_demo.xCell2Ref
)

# Parallel processing example with BiocParallel
library(BiocParallel)
parallel_param <- MulticoreParam(workers = 2)
xcell2_res_parallel <- xCell2::xCell2Analysis(
  mix = mix_demo, 
  xcell2object = DICE_demo.xCell2Ref, 
  BPPARAM = parallel_param
)


AlmogAngel/xCell2 documentation built on Jan. 24, 2025, 10:39 a.m.