xCell2Analysis: Perform Cell Type Enrichment Analysis

View source: R/xCell2Analysis.R

xCell2AnalysisR Documentation

Perform Cell Type Enrichment Analysis

Description

This function estimate the relative enrichment of cell types in a bulk gene expression mixture. The analysis leverages gene signatures from a pre-trained xCell2Object to compute enrichment scores for each cell type. It also applies linear transformation and spillover correction to refine the enrichment scores.

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 matrix (genes in rows, samples in columns). The input should use the same gene annotation system as the reference object.

xcell2object

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

minSharedGenes

Minimum fraction of shared genes required between the mixture and the reference object (default: 0.9). If the shared fraction falls below this threshold, the function will stop with an error or warning, as accurate analysis depends on sufficient overlap between the mixture and reference genes.

rawScores

A Boolean indicating whether to return raw enrichment scores (default: FALSE). Raw enrichment scores are computed directly from gene rankings without linear transformation or spillover correction.

spillover

A Boolean to enable spillover correction on the enrichment scores (default: TRUE). Spillover occurs when gene expression patterns overlap between closely related cell types, potentially inflating enrichment scores. Correcting for spillover enhances the specificity of enrichment scores, particularly for related cell types. The strength of this correction can be adjusted using the spilloverAlpha parameter.

spilloverAlpha

A numeric value controlling the strength of spillover correction (default: 0.5). Lower values apply weaker correction, while higher values apply stronger correction. An alpha value of 0.5 is suitable for most cases, but users may tune this parameter based on the similarity of cell types in their reference.

BPPARAM

A BiocParallelParam instance that determines the parallelization strategy (more in "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.

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.

Author(s)

Almog Angel and Dvir Aran

See Also

xCell2Train, for generating 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
data(mix_demo, package = "xCell2")

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

# Example using parallel processing with MulticoreParam
library(BiocParallel)
parallel_param <- MulticoreParam(workers = 2) # Adjust workers as needed
xcell2_res_parallel <- xCell2::xCell2Analysis(
  mix = mix_demo, 
  xcell2object = DICE_demo.xCell2Ref, 
  BPPARAM = parallel_param
)


AlmogAngel/xCell2 documentation built on Dec. 19, 2024, 8:21 a.m.