Contents

1 Install and Load the Package

To install cytofkit package, start R and run the following codes on R console:

source("https://bioconductor.org/biocLite.R")
biocLite("cytofkit")

Notes: cytofkit GUI is dependent on XQuartz windowing system (X Windows) on Mac (OS X > 10.7). Install XQuartz from http://xquartz.macosforge.org.

Load the Package:

library("cytofkit") 

Read the package description:

?"cytofkit-package"

2 Options for Using cytofkit Package

cytofkit provides three ways to employ the workforce of this package:

2.1 Run with GUI

The easiest way to use cytofkit package is through the GUI. The GUI provides all main options of cytofkit on a visual interface. To launch the GUI, load the package and type the following command:

cytofkit_GUI()  

The interface will appear like below, you can click the information button ! to check the explanation for each entry and customize your own analysis.

cytofkit GUI

cytofkit GUI

Start your analysis as simply as following:

  • Choose the input fcs files from the directory where you store FCS data;
  • Select the markers from the auto-generated marker list;
  • Choose the directory where to save your output;
  • Give a project name as a prefix for the names of result files;
  • Select a data merging method if you have multiple FCS files;
  • Select your clustering method(s), visualization method(s), progression estimation method(s)

Then submit it, that’s all.

Depends on the size of your data, it will take some time to run the analysis. Once done, a window will pop up, showing you the path where the results have been stored, and asking you if open the shiny web APP. If YES, the shiny APP will be deployed locally and opened in your default web browser. Among the saved results, a special R data object with suffix of .RData is for loading the results into the shiny APP. Choose the .RData file on the shiny APP then submit it, your journey of exploring the results starts.

2.2 Run with the Core Function

Cytofkit provides a core function cytofkit() to drive the analysis pipeline of mass cytometry data. Users only need to define several key parameters to start their analysis automatically. One simple example of running cytofkit using the core function is like this:

set.seed(100)
dir <- system.file('extdata',package='cytofkit')
file <- list.files(dir ,pattern='.fcs$', full=TRUE)
parameters <- list.files(dir, pattern='.txt$', full=TRUE)
res <- cytofkit(fcsFiles = file, 
                markers = parameters, 
                projectName = 'cytofkit_test',
                transformMethod = "cytofAsinh", 
                mergeMethod = "ceil",
                fixedNum = 500,                                    ## set at 500 for faster run
                dimReductionMethod = "tsne",
                clusterMethods = c("Rphenograph", "ClusterX"),    ## accept multiple methods
                visualizationMethods = c("tsne", "pca"),          ## accept multiple methods
                progressionMethod = "isomap",
                clusterSampleSize = 500,
                resultDir = getwd(),
                saveResults = TRUE, 
                saveObject = TRUE)

You can customize the parameters for your own need, run ?cytofkit to get information of all the parameters for cytofkit. As running with GUI, once the analysis is done, the results will be saved under resultDir automatically.

2.3 Run with Commands (Step-by-Step)

You can make use of the functions exported from cytofkit to make your analysis more flexible and fit your own need. Here we use a sample data for demo:

2.3.1 Pre-processing

## Loading the FCS data:  
dir <- system.file('extdata',package='cytofkit')
file <- list.files(dir ,pattern='.fcs$', full=TRUE)
paraFile <- list.files(dir, pattern='.txt$', full=TRUE)
parameters <- as.character(read.table(paraFile, header = TRUE)[,1])

## File name
file
## [1] "/tmp/RtmpCY9HqB/Rinst35864554e43/cytofkit/extdata/130515_C2_stim_CD19-.fcs"
## parameters
parameters
##  [1] "(Sm152)Di<Vd2>"    "(Eu153)Di<CD107a>" "(Sm154)Di<CD3>"   
##  [4] "(Gd155)Di<CD152>"  "(Gd156)Di<CD19>"   "(Gd157)Di<TIM3>"  
##  [7] "(Gd158)Di<CD56>"   "(Tb159)Di<IL10>"   "(Gd160)Di<CD28>"  
## [10] "(Dy161)Di<CD38>"   "(Dy162)Di<IL4>"    "(Dy163)Di<CD127>"
## Extract the expression matrix with transformation
data_transformed <- cytof_exprsExtract(fcsFile = file, 
                                       comp = FALSE, 
                                       transformMethod = "cytofAsinh")
## If analysing flow cytometry data, you can set comp to TRUE or 
## provide a transformation matrix to apply compensation

## If you have multiple FCS files, expression can be extracted and combined
combined_data_transformed <- cytof_exprsMerge(fcsFiles = file, comp=FALSE,
                                              transformMethod = "cytofAsinh",
                                              mergeMethod = "all")
## change mergeMethod to apply different combination strategy

## Take a look at the extracted expression matrix
head(data_transformed[ ,1:3])
##                        Cell_length<NA> (Rh103)Di<BC103> (Pd104)Di<BC104>
## 130515_C2_stim_CD19-_1        2.824903     0.0023533088         4.072829
## 130515_C2_stim_CD19-_2        2.725596     0.0007761656         4.125760
## 130515_C2_stim_CD19-_3        2.672016    -0.0002698989         3.705151
## 130515_C2_stim_CD19-_4        2.555494     0.0043997057         3.393448
## 130515_C2_stim_CD19-_5        2.644121    -0.0012975977         3.554248
## 130515_C2_stim_CD19-_6        2.800979    -0.0007743693         3.974557

2.3.2 Cell Subset Detection

## use clustering algorithm to detect cell subsets
## to speed up our test here, we only use 100 cells
data_transformed_1k <- data_transformed[1:100, ]

## run PhenoGraph
cluster_PhenoGraph <- cytof_cluster(xdata = data_transformed_1k, method = "Rphenograph")
##   Running PhenoGraph...  Finding nearest neighbors...DONE ~ 0.025 s
##   Compute jaccard coefficient between nearest-neighbor sets...DONE ~ 0.032 s
##   Build undirected graph from the weighted links...DONE ~ 0.014 s
##   Run louvain clustering on the graph ...DONE ~ 0.035 s
##   Return a community class
##   -Modularity value: 0.4697538 
##   -Number of clusters: 4 DONE!
## run ClusterX
data_transformed_1k_tsne <- cytof_dimReduction(data=data_transformed_1k, method = "tsne")
##   Running t-SNE...with seed 42  DONE
cluster_ClusterX <- cytof_cluster(ydata = data_transformed_1k_tsne,  method="ClusterX")
##   Running ClusterX...    Calculate cutoff distance...0.36  
##     Calculate local Density...DONE!
##     Detect nearest neighbour with higher density...DONE!
##     Peak detection...DONE!
##     Cluster assigning...DONE!
##  DONE!
## run DensVM (takes long time, we skip here)
cluster_DensVM <- cytof_cluster(xdata = data_transformed_1k, 
                                ydata = data_transformed_1k_tsne, method = "DensVM")
## run FlowSOM with cluster number 15
cluster_FlowSOM <- cytof_cluster(xdata = data_transformed_1k, method = "FlowSOM", FlowSOM_k = 12)
##   Running FlowSOM...    Building SOM...
##     Meta clustering to 12 clusters...
##  DONE!
## combine data
data_1k_all <- cbind(data_transformed_1k, data_transformed_1k_tsne, 
                     PhenoGraph = cluster_PhenoGraph, ClusterX=cluster_ClusterX, 
                     FlowSOM=cluster_FlowSOM)
data_1k_all <- as.data.frame(data_1k_all)

2.3.3 Cell Subset Visualization and Interpretation

## PhenoGraph plot on tsne
cytof_clusterPlot(data=data_1k_all, xlab="tsne_1", ylab="tsne_2", 
                  cluster="PhenoGraph", sampleLabel = FALSE)

## PhenoGraph cluster heatmap
PhenoGraph_cluster_median <- aggregate(. ~ PhenoGraph, data = data_1k_all, median)
cytof_heatmap(PhenoGraph_cluster_median[, 2:37], baseName = "PhenoGraph Cluster Median")

## ClusterX plot on tsne
cytof_clusterPlot(data=data_1k_all, xlab="tsne_1", ylab="tsne_2", cluster="ClusterX", sampleLabel = FALSE)

## ClusterX cluster heatmap
ClusterX_cluster_median <- aggregate(. ~ ClusterX, data = data_1k_all, median)
cytof_heatmap(ClusterX_cluster_median[, 2:37], baseName = "ClusterX Cluster Median")

## FlowSOM plot on tsne
cytof_clusterPlot(data=data_1k_all, xlab="tsne_1", ylab="tsne_2", 
                  cluster="FlowSOM", sampleLabel = FALSE)

## FlowSOM cluster heatmap
FlowSOM_cluster_median <- aggregate(. ~ FlowSOM, data = data_1k_all, median)
cytof_heatmap(FlowSOM_cluster_median[, 2:37], baseName = "FlowSOM Cluster Median")

2.3.4 Inference of Subset Progression

## Inference of PhenoGraph cluster relatedness
PhenoGraph_progression <- cytof_progression(data = data_transformed_1k, 
                                            cluster = cluster_PhenoGraph, 
                                            method="isomap", clusterSampleSize = 50, 
                                            sampleSeed = 5)
##   Running ISOMAP...  DONE
p_d <- data.frame(PhenoGraph_progression$sampleData, 
                  PhenoGraph_progression$progressionData, 
                  cluster = PhenoGraph_progression$sampleCluster, 
                  check.names = FALSE)

## cluster relatedness plot
cytof_clusterPlot(data=p_d, xlab="isomap_1", ylab="isomap_2", 
                  cluster="cluster", sampleLabel = FALSE)

## marker expression profile
markers <- c("(Sm150)Di<GranzymeB>", "(Yb173)Di<Perforin>")

cytof_colorPlot(data=p_d, xlab="isomap_1", ylab="isomap_2", zlab = markers[1], limits = range(p_d[,1:52]))

cytof_colorPlot(data=p_d, xlab="isomap_1", ylab="isomap_2", zlab = markers[2], limits = range(p_d[,1:52]))

cytof_progressionPlot(data=p_d, markers=markers, orderCol="isomap_1", clusterCol = "cluster")

## Inference of ClusterX cluster relatedness
ClusterX_progression <- cytof_progression(data = data_transformed_1k, 
                                          cluster = cluster_ClusterX, 
                                          method="isomap", 
                                          clusterSampleSize = 30, 
                                          sampleSeed = 3)
##   Running ISOMAP...  DONE
c_d <- data.frame(ClusterX_progression$sampleData, 
                  ClusterX_progression$progressionData,
                  cluster=ClusterX_progression$sampleCluster, 
                  check.names = FALSE)

## cluster relatedness plot
cytof_clusterPlot(data=c_d, xlab="isomap_1", ylab="isomap_2", 
                  cluster="cluster", sampleLabel = FALSE)

## marker expression profile
markers <- c("(Sm150)Di<GranzymeB>", "(Yb173)Di<Perforin>")

cytof_colorPlot(data=c_d, xlab="isomap_1", ylab="isomap_2", zlab = markers[1], limits = range(c_d[,1:52]))

cytof_colorPlot(data=c_d, xlab="isomap_1", ylab="isomap_2", zlab = markers[2], limits = range(c_d[,1:52]))

cytof_progressionPlot(data=c_d, markers, orderCol="isomap_1", clusterCol = "cluster")

2.3.5 Post-processing

## save analysis results to FCS file
cytof_addToFCS(data_1k_all, rawFCSdir=dir, analyzedFCSdir="analysed_FCS", 
               transformed_cols = c("tsne_1", "tsne_2"), 
               cluster_cols = c("PhenoGraph", "ClusterX", "FlowSOM"))

In addition, expression values for cells in a specified cluster can be extracted using cytof_clusterMtrx()

## See documentation, this function uses the output of main cytofkit cuntion as its input
?cytof_clusterMtrx

2.3.6 Get package update news

cytofkitNews()

3 Session Information

sessionInfo()
## R version 3.4.2 (2017-09-28)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.3 LTS
## 
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.6-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.6-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        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       
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] cytofkit_1.10.0 plyr_1.8.4      ggplot2_2.2.1   BiocStyle_2.6.0
## 
## loaded via a namespace (and not attached):
##   [1] Rtsne_0.13                  VGAM_1.0-4                 
##   [3] minqa_1.2.4                 colorspace_1.3-2           
##   [5] RcppEigen_0.3.3.3.0         class_7.3-14               
##   [7] rprojroot_1.2               htmlTable_1.9              
##   [9] corpcor_1.6.9               base64enc_0.1-3            
##  [11] proxy_0.4-19                MatrixModels_0.4-1         
##  [13] ggrepel_0.7.0               mvtnorm_1.0-6              
##  [15] codetools_0.2-15            splines_3.4.2              
##  [17] doParallel_1.0.11           robustbase_0.92-7          
##  [19] knitr_1.17                  Formula_1.2-2              
##  [21] nloptr_1.0.4                pbkrtest_0.4-7             
##  [23] cluster_2.0.6               graph_1.56.0               
##  [25] shiny_1.0.5                 rrcov_1.4-3                
##  [27] compiler_3.4.2              backports_1.1.1            
##  [29] Matrix_1.2-11               lazyeval_0.2.1             
##  [31] acepack_1.4.1               htmltools_0.3.6            
##  [33] quantreg_5.34               tools_3.4.2                
##  [35] igraph_1.1.2                gtable_0.2.0               
##  [37] RANN_2.5.1                  reshape2_1.4.2             
##  [39] Rcpp_0.12.13                Biobase_2.38.0             
##  [41] RJSONIO_1.3-0               gdata_2.18.0               
##  [43] nlme_3.1-131                iterators_1.0.8            
##  [45] lmtest_0.9-35               laeken_0.4.6               
##  [47] stringr_1.2.0               lme4_1.1-14                
##  [49] mime_0.5                    miniUI_0.1.1               
##  [51] gtools_3.5.0                XML_3.98-1.9               
##  [53] DEoptimR_1.0-8              MASS_7.3-47                
##  [55] zoo_1.8-0                   scales_0.5.0               
##  [57] VIM_4.7.0                   colourpicker_1.0           
##  [59] parallel_3.4.2              SparseM_1.77               
##  [61] RColorBrewer_1.1-2          yaml_2.1.14                
##  [63] curl_3.0                    gridExtra_2.3              
##  [65] rpart_4.1-11                latticeExtra_0.6-28        
##  [67] stringi_1.1.5               pcaPP_1.9-72               
##  [69] foreach_1.4.3               permute_0.9-4              
##  [71] e1071_1.6-8                 flowCore_1.44.0            
##  [73] checkmate_1.8.5             destiny_2.6.0              
##  [75] TTR_0.23-2                  caTools_1.17.1             
##  [77] BiocGenerics_0.24.0         boot_1.3-20                
##  [79] rlang_0.1.2                 pkgconfig_2.0.1            
##  [81] matrixStats_0.52.2          bitops_1.0-6               
##  [83] evaluate_0.10.1             lattice_0.20-35            
##  [85] labeling_0.3                htmlwidgets_0.9            
##  [87] pdist_1.2                   magrittr_1.5               
##  [89] bookdown_0.5                R6_2.2.2                   
##  [91] gplots_3.0.1                Hmisc_4.0-3                
##  [93] foreign_0.8-69              mgcv_1.8-22                
##  [95] xts_0.10-0                  survival_2.41-3            
##  [97] sp_1.2-5                    nnet_7.3-12                
##  [99] FlowSOM_1.10.0              tibble_1.3.4               
## [101] tsne_0.1-3                  car_2.1-5                  
## [103] shinyFiles_0.6.2            KernSmooth_2.23-15         
## [105] rmarkdown_1.6               grid_3.4.2                 
## [107] data.table_1.10.4-3         vegan_2.4-4                
## [109] ConsensusClusterPlus_1.42.0 vcd_1.4-3                  
## [111] digest_0.6.12               xtable_1.8-2               
## [113] httpuv_1.3.5                stats4_3.4.2               
## [115] munsell_0.4.3               smoother_1.1               
## [117] tcltk_3.4.2