Contents

1 Introduction

The CompensAID package provides quality control for evaluating the potential presence of reference errors in flow cytometry data. It requires pre-processed Flow Cytometry Standard (FCS) files as input. The tool automatically calculates the standard deviation (SD) and median fluorescence intensity (MFI) in the secondary channel for both the negative and positive populations defined in the primary channel. Next, the positive population is divided into equal segments, and for each segment, the Secondary Stain Index (SSI) is computed. This process allows for a more detailed assessment of marker combinations affected by reference errors.

2 Installing CompensAID

Not evaluated, run to install packages

# Install BiocManager if not already installed
if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

# Install devtools if not already installed
if (!requireNamespace("devtools", quietly = TRUE))
    install.packages("devtools")

# Install ggcyto development version from GitHub - otherwise does not work with ggplot2
devtools::install_github("RGLab/ggcyto", ref = "devel")

# Install CompensAID from GitHub
devtools::install_github("Olsman/CompensAID")

# Install ggplot2 from Bioconductor - older version otherwise it ggcyto gives errors
devtools::install_version("ggplot2", version = "3.5.2", repos = "http://cran.us.r-project.org")

# Load libraries
library(ggplot2)
library(ggcyto)
library(CompensAID)

3 Pre-processing

The CompensAID package operates on pre-processed FlowFrame objects to evaluate the presence of reference errors. Below is an example pipeline for preparing your data: 1. Removal of margin events (conventional FCS files) 2. Compensation (conventional FCS files) 3. Logicle transformation 4. Quality control using PeacoQC Note: Doublets and debris are not removed in this pipeline. In this example, scatter channels are excluded entirely from the FCS file.

# Read flowFrame
file <- flowCore::read.FCS(system.file("extdata", "raw.68983.fcs", package = "CompensAID"))

# Remove margins (for conventional only)
file <- PeacoQC::RemoveMargins(ff = file,
                               channels = colnames(flowCore::exprs(file)),
                               output = "frame")

# Compensate data (for conventional only)
file <- flowCore::compensate(file, file@description$SPILL)

# Transform data 
file <- flowCore::transform(file,
                            flowCore::estimateLogicle(file, colnames(flowCore::keyword(file)$SPILL)))

# Run PeacoQC
file <- PeacoQC::PeacoQC(ff = file,
                       channels = colnames(flowCore::exprs(file)),
                       save_fcs = FALSE,
                       plot = FALSE,
                       report = FALSE,
                       output_directory = NULL)
file <- file$FinalFF

# Retain only markers
file <- file[, names(flowCore::markernames(file))]

4 Final data

The data file is processed as shown above

# Pre-processed data 
file <- flowCore::read.FCS(system.file("extdata", "68983.fcs", package = "CompensAID"))

5 CompensAID: run tool

The CompensAID tool is used here with its default parameter settings. A segment.value of four, which defines the number of segments within the positive population, was found to be suitable for both conventional and spectral flow cytometry data. The minimum number of events per segment is set to 50 by default. If you observe segments with very negative Secondary Stain Index (SSI) values and few events, this could indicate a false-positive finding by the tool

# CompensAID
res.compensaid <- CompensAID::CompensAID(ff = file)
## Importing sample: /tmp/RtmpX8h1UB/Rinst3673477be9b1be/CompensAID/extdata/68983.fcs

6 Secondary Stain Index Matrix

Visualize the output of the CompensAID tool.

matrix <- CompensAID::PlotMatrix(output = res.compensaid)
print(matrix)

7 Identify flagged marker combinations

Obtain marker combinations with a Secondary Stain Index < -1.

# Check SSI < -1
index <- which(res.compensaid[["matrix"]] < -1, arr.ind = TRUE)

# Obtain marker names
library(dplyr)
## 
## Attaching package: 'dplyr'
## The following objects are masked from 'package:stats':
## 
##     filter, lag
## The following objects are masked from 'package:base':
## 
##     intersect, setdiff, setequal, union
combination <- data.frame(primary.channel = rownames(res.compensaid[["matrix"]])[index[, "col"]],
                          secondary.channel = rownames(res.compensaid[["matrix"]])[index[, "row"]]) |>
  dplyr::mutate(primary.marker = flowCore::markernames(file)[primary.channel],
                secondary.marker = flowCore::markernames(file)[secondary.channel])

8 Dot plots

Visualize the flagged marker combinations.

# Visualize dot plot with SSI scores
dotPlot <- CompensAID::PlotDotSSI(output = res.compensaid,
                                  og = file,
                                  primary = combination$primary.marker[1],
                                  secondary = combination$secondary.marker[1],
                                  showScores = TRUE)
## Warning: `aes_()` was deprecated in ggplot2 3.0.0.
## ℹ Please use tidy evaluation idioms with `aes()`
## ℹ The deprecated feature was likely used in the ggcyto package.
##   Please report the issue to the authors.
## This warning is displayed once per session.
## Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
## generated.
print(dotPlot)

9 Session info

sessionInfo()
## R Under development (unstable) (2026-01-15 r89304)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] dplyr_1.1.4      BiocStyle_2.39.0
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1     farver_2.1.2         S7_0.2.1            
##  [4] bitops_1.0-9         fastmap_1.2.0        XML_3.99-0.20       
##  [7] digest_0.6.39        lifecycle_1.0.5      CompensAID_0.99.6   
## [10] magrittr_2.0.4       compiler_4.6.0       rlang_1.1.7         
## [13] sass_0.4.10          tools_4.6.0          yaml_2.3.12         
## [16] data.table_1.18.0    knitr_1.51           labeling_0.4.3      
## [19] interp_1.1-6         plyr_1.8.9           RColorBrewer_1.1-3  
## [22] abind_1.4-8          KernSmooth_2.23-26   withr_3.0.2         
## [25] purrr_1.2.1          RProtoBufLib_2.23.0  BiocGenerics_0.57.0 
## [28] grid_4.6.0           polyclip_1.10-7      stats4_4.6.0        
## [31] latticeExtra_0.6-31  caTools_1.18.3       ggplot2_4.0.1       
## [34] scales_1.4.0         gtools_3.9.5         MASS_7.3-65         
## [37] tinytex_0.58         dichromat_2.0-0.1    cli_3.6.5           
## [40] rmarkdown_2.30       ncdfFlow_2.57.0      generics_0.1.4      
## [43] otel_0.2.0           rstudioapi_0.18.0    reshape2_1.4.5      
## [46] IDPmisc_1.1.21       cachem_1.1.0         flowCore_2.23.1     
## [49] stringr_1.6.0        BiocManager_1.30.27  matrixStats_1.5.0   
## [52] vctrs_0.7.0          jsonlite_2.0.0       carData_3.0-5       
## [55] cytolib_2.23.0       bookdown_0.46        car_3.1-3           
## [58] S4Vectors_0.49.0     Formula_1.2-5        Rgraphviz_2.55.0    
## [61] ParallelLogger_3.5.1 magick_2.9.0         jpeg_0.1-11         
## [64] flowDensity_1.45.0   jquerylib_0.1.4      hexbin_1.28.5       
## [67] tidyr_1.3.2          glue_1.8.0           stringi_1.8.7       
## [70] gtable_0.3.6         deldir_2.0-4         ggcyto_1.39.1       
## [73] tibble_3.3.1         pillar_1.11.1        htmltools_0.5.9     
## [76] gplots_3.3.0         graph_1.89.1         R6_2.6.1            
## [79] evaluate_1.0.5       flowWorkspace_4.23.1 lattice_0.22-7      
## [82] Biobase_2.71.0       png_0.1-8            backports_1.5.0     
## [85] flowViz_1.75.0       bslib_0.9.0          Rcpp_1.1.1          
## [88] gridExtra_2.3        checkmate_2.3.3      xfun_0.56           
## [91] pkgconfig_2.0.3