K-nearest neighbors:

We read in input.scone.csv, which is our file modified (and renamed) from the get.marker.names() function. The K-nearest neighbor generation is derived from the Fast Nearest Neighbors (FNN) R package, within our function Fnn(), which takes as input the “input markers” to be used, along with the concatenated data previously generated, and the desired k. We advise the default selection to the total number of cells in the dataset divided by 100, as has been optimized on existing mass cytometry datasets. The output of this function is a matrix of each cell and the identity of its k-nearest neighbors, in terms of its row number in the dataset used here as input.

library(Sconify)
# Markers from the user-generated excel file
marker.file <- system.file('extdata', 'markers.csv', package = "Sconify")
markers <- ParseMarkers(marker.file)

# How to convert your excel sheet into vector of static and functional markers
markers
## $input
##  [1] "CD3(Cd110)Di"           "CD3(Cd111)Di"           "CD3(Cd112)Di"          
##  [4] "CD235-61-7-15(In113)Di" "CD3(Cd114)Di"           "CD45(In115)Di"         
##  [7] "CD19(Nd142)Di"          "CD22(Nd143)Di"          "IgD(Nd145)Di"          
## [10] "CD79b(Nd146)Di"         "CD20(Sm147)Di"          "CD34(Nd148)Di"         
## [13] "CD179a(Sm149)Di"        "CD72(Eu151)Di"          "IgM(Eu153)Di"          
## [16] "Kappa(Sm154)Di"         "CD10(Gd156)Di"          "Lambda(Gd157)Di"       
## [19] "CD24(Dy161)Di"          "TdT(Dy163)Di"           "Rag1(Dy164)Di"         
## [22] "PreBCR(Ho165)Di"        "CD43(Er167)Di"          "CD38(Er168)Di"         
## [25] "CD40(Er170)Di"          "CD33(Yb173)Di"          "HLA-DR(Yb174)Di"       
## 
## $functional
##  [1] "pCrkL(Lu175)Di"  "pCREB(Yb176)Di"  "pBTK(Yb171)Di"   "pS6(Yb172)Di"   
##  [5] "cPARP(La139)Di"  "pPLCg2(Pr141)Di" "pSrc(Nd144)Di"   "Ki67(Sm152)Di"  
##  [9] "pErk12(Gd155)Di" "pSTAT3(Gd158)Di" "pAKT(Tb159)Di"   "pBLNK(Gd160)Di" 
## [13] "pP38(Tm169)Di"   "pSTAT5(Nd150)Di" "pSyk(Dy162)Di"   "tIkBa(Er166)Di"
# Get the particular markers to be used as knn and knn statistics input
input.markers <- markers[[1]]
funct.markers <- markers[[2]]

# Selection of the k. See "Finding Ideal K" vignette
k <- 30

# The built-in scone functions
wand.nn <- Fnn(cell.df = wand.combined, input.markers = input.markers, k = k)
# Cell identity is in rows, k-nearest neighbors are columns
# List of 2 includes the cell identity of each nn, 
#   and the euclidean distance between
#   itself and the cell of interest

# Indices
str(wand.nn[[1]])
##  int [1:1000, 1:30] 224 638 885 336 226 486 274 499 90 189 ...
wand.nn[[1]][1:20, 1:10]
##       [,1] [,2] [,3] [,4] [,5] [,6] [,7] [,8] [,9] [,10]
##  [1,]  224  662  903  777  857   97  995  423  759   450
##  [2,]  638  479  779  558  689  826  613  814  890   310
##  [3,]  885  865  165  271  941  453  456  185  386   422
##  [4,]  336  322  195  813  480  472  426  572  464   525
##  [5,]  226  635  688   23  952  119  847  689  423   335
##  [6,]  486  707   76  108  717  188  949  628  869   531
##  [7,]  274  732   76  256  486  949  381  833   60   785
##  [8,]  499   97  477   71  133  483  498  335    1   777
##  [9,]   90  417  754  919  948  367  282  427  940   152
## [10,]  189  485  339  198  850  246  308  567   50    12
## [11,]  550  182  986  211  720  287  369  742  256   819
## [12,]  822  567  734  155  178  485   50  651  682   189
## [13,]  369  869  628  258  918  751  924  774  486   748
## [14,]  413  869  735  283  369  550  374  545  268    77
## [15,]  968  429   58  634  563   32   29  131  814   584
## [16,]  361  475  158   28  755   62  281  796  707   904
## [17,]  553    6  486  532  622  303  599   35   47   818
## [18,]  976  548  157    1  510  559  260  277  337   952
## [19,]  486  780  774  258  374   13  949  875  751   707
## [20,]  647  550  158  707  281   83  475  981  112    28
# Distance
str(wand.nn[[2]])
##  num [1:1000, 1:30] 2.79 3.15 4.68 3.41 3.42 ...
wand.nn[[2]][1:20, 1:10]
##           [,1]     [,2]     [,3]     [,4]     [,5]     [,6]     [,7]     [,8]
##  [1,] 2.789880 3.055646 3.132630 3.165508 3.252539 3.326228 3.361317 3.367862
##  [2,] 3.148988 3.259271 3.285041 3.339632 3.387981 3.426064 3.466001 3.557338
##  [3,] 4.682335 4.897369 4.966901 4.999660 5.039738 5.046384 5.056376 5.080678
##  [4,] 3.413191 3.932257 3.983984 4.103935 4.231706 4.375022 4.487338 4.490640
##  [5,] 3.422646 3.630493 3.733845 3.978432 3.993102 4.009889 4.012787 4.093523
##  [6,] 2.314570 3.058306 3.192741 3.211358 3.282797 3.304618 3.323551 3.338076
##  [7,] 2.355356 2.908774 2.919602 2.931802 2.977181 3.087312 3.093901 3.129033
##  [8,] 2.808429 3.262643 3.579560 3.805075 3.820019 3.855877 3.929064 3.949969
##  [9,] 4.741335 4.751820 4.903699 4.918202 4.938301 5.049838 5.110778 5.178555
## [10,] 4.376914 4.381185 4.453382 4.556114 4.574178 4.590487 4.769049 4.790982
## [11,] 2.715516 3.031190 3.071348 3.181328 3.194773 3.252136 3.287926 3.291665
## [12,] 3.939728 4.022139 4.126250 4.184798 4.201195 4.271325 4.335344 4.385860
## [13,] 2.866192 2.883953 2.918729 2.920141 3.038580 3.053570 3.063518 3.064875
## [14,] 2.501066 2.749921 2.937049 2.994278 2.995655 3.025830 3.037645 3.049261
## [15,] 3.395493 3.514509 3.544002 3.652081 3.968964 3.986126 3.987274 3.998168
## [16,] 3.107219 3.339896 3.593604 3.609493 3.673527 3.675341 3.694217 3.706624
## [17,] 4.008753 4.190589 4.403804 4.451535 4.502277 4.552564 4.580753 4.603309
## [18,] 3.064093 3.278002 3.458017 3.550795 3.611620 3.653345 3.700618 3.732034
## [19,] 3.350017 3.367601 3.408984 3.460197 3.493728 3.498565 3.556694 3.568768
## [20,] 3.205045 3.301887 3.354461 3.396483 3.488228 3.498656 3.529567 3.567202
##           [,9]    [,10]
##  [1,] 3.372474 3.377001
##  [2,] 3.569509 3.578715
##  [3,] 5.085931 5.087400
##  [4,] 4.493344 4.552603
##  [5,] 4.099559 4.115680
##  [6,] 3.388791 3.422527
##  [7,] 3.211902 3.252272
##  [8,] 4.022563 4.057465
##  [9,] 5.336097 5.384740
## [10,] 4.798983 4.907949
## [11,] 3.318469 3.367012
## [12,] 4.560742 4.592491
## [13,] 3.180808 3.225316
## [14,] 3.158268 3.218674
## [15,] 4.116360 4.191713
## [16,] 3.797831 3.852141
## [17,] 4.639045 4.661527
## [18,] 3.774769 3.785342
## [19,] 3.573260 3.610366
## [20,] 3.583937 3.618509

Finding scone values:

This function iterates through each KNN, and performs a series of calculations. The first is fold change values for each maker per KNN, where the user chooses whether this will be based on medians or means. The second is a statistical test, where the user chooses t test or Mann-Whitney U test. I prefer the latter, because it does not assume any properties of the distributions. Of note, the p values are adjusted for false discovery rate, and therefore are called q values in the output of this function. The user also inputs a threshold parameter (default 0.05), where the fold change values will only be shown if the corresponding statistical test returns a q value below said threshold. Finally, the “multiple.donor.compare” option, if set to TRUE will perform a t test based on the mean per-marker values of each donor. This is to allow the user to make comparisons across replicates or multiple donors if that is relevant to the user’s biological questions. This function returns a matrix of cells by computed values (change and statistical test results, labeled either marker.change or marker.qvalue). This matrix is intermediate, as it gets concatenated with the original input matrix in the post-processing step (see the relevant vignette). We show the code and the output below. See the post-processing vignette, where we show how this gets combined with the input data, and additional analysis is performed.

wand.scone <- SconeValues(nn.matrix = wand.nn, 
                      cell.data = wand.combined, 
                      scone.markers = funct.markers, 
                      unstim = "basal")

wand.scone
## # A tibble: 1,000 x 34
##    `pCrkL(Lu175)Di… `pCREB(Yb176)Di… `pBTK(Yb171)Di.… `pS6(Yb172)Di.I…
##               <dbl>            <dbl>            <dbl>            <dbl>
##  1            0.977            0.984                1            0.610
##  2            0.932            0.984                1            0.887
##  3            0.925            0.984                1            0.670
##  4            0.967            0.984                1            0.624
##  5            0.991            0.984                1            0.608
##  6            0.977            0.984                1            0.910
##  7            0.925            0.984                1            0.961
##  8            0.925            0.984                1            0.851
##  9            0.973            0.984                1            0.515
## 10            0.959            0.997                1            0.630
## # … with 990 more rows, and 30 more variables:
## #   `cPARP(La139)Di.IL7.qvalue` <dbl>, `pPLCg2(Pr141)Di.IL7.qvalue` <dbl>,
## #   `pSrc(Nd144)Di.IL7.qvalue` <dbl>, `Ki67(Sm152)Di.IL7.qvalue` <dbl>,
## #   `pErk12(Gd155)Di.IL7.qvalue` <dbl>, `pSTAT3(Gd158)Di.IL7.qvalue` <dbl>,
## #   `pAKT(Tb159)Di.IL7.qvalue` <dbl>, `pBLNK(Gd160)Di.IL7.qvalue` <dbl>,
## #   `pP38(Tm169)Di.IL7.qvalue` <dbl>, `pSTAT5(Nd150)Di.IL7.qvalue` <dbl>,
## #   `pSyk(Dy162)Di.IL7.qvalue` <dbl>, `tIkBa(Er166)Di.IL7.qvalue` <dbl>,
## #   `pCrkL(Lu175)Di.IL7.change` <dbl>, `pCREB(Yb176)Di.IL7.change` <dbl>,
## #   `pBTK(Yb171)Di.IL7.change` <dbl>, `pS6(Yb172)Di.IL7.change` <dbl>,
## #   `cPARP(La139)Di.IL7.change` <dbl>, `pPLCg2(Pr141)Di.IL7.change` <dbl>,
## #   `pSrc(Nd144)Di.IL7.change` <dbl>, `Ki67(Sm152)Di.IL7.change` <dbl>,
## #   `pErk12(Gd155)Di.IL7.change` <dbl>, `pSTAT3(Gd158)Di.IL7.change` <dbl>,
## #   `pAKT(Tb159)Di.IL7.change` <dbl>, `pBLNK(Gd160)Di.IL7.change` <dbl>,
## #   `pP38(Tm169)Di.IL7.change` <dbl>, `pSTAT5(Nd150)Di.IL7.change` <dbl>,
## #   `pSyk(Dy162)Di.IL7.change` <dbl>, `tIkBa(Er166)Di.IL7.change` <dbl>,
## #   IL7.fraction.cond.2 <dbl>, density <dbl>

For programmers: performing additional per-KNN statistics

If one wants to export KNN data to perform other statistics not available in this package, then I provide a function that produces a list of each cell identity in the original input data matrix, and a matrix of all cells x features of its KNN.

I also provide a function to find the KNN density estimation independently of the rest of the “scone.values” analysis, to save time if density is all the user wants. With this density estimation, one can perform interesting analysis, ranging from understanding phenotypic density changes along a developmental progression (see post-processing vignette for an example), to trying out density-based binning methods (eg. X-shift). Of note, this density is specifically one divided by the aveage distance to k-nearest neighbors. This specific measure is related to the Shannon Entropy estimate of that point on the manifold (https://hal.archives-ouvertes.fr/hal-01068081/document).

I use this metric to avoid the unusual properties of the volume of a sphere as it increases in dimensions (https://en.wikipedia.org/wiki/Volume_of_an_n-ball). This being said, one can modify this vector to be such a density estimation (example http://www.cs.haifa.ac.il/~rita/ml_course/lectures_old/KNN.pdf), by treating the distance to knn as the radius of a n-dimensional sphere and incoroprating said volume accordingly.

An individual with basic programming skills can iterate through these elements to perform the statistics of one’s choosing. Examples would include per-KNN regression and classification, or feature imputation. The additional functionality is shown below, with the example knn.list in the package being the first ten instances:

# Constructs KNN list, computes KNN density estimation
wand.knn.list <- MakeKnnList(cell.data = wand.combined, nn.matrix = wand.nn)
wand.knn.list[[8]]
## # A tibble: 30 x 51
##    `CD3(Cd110)Di` `CD3(Cd111)Di` `CD3(Cd112)Di` `CD235-61-7-15(… `CD3(Cd114)Di`
##             <dbl>          <dbl>          <dbl>            <dbl>          <dbl>
##  1         0.227         -0.185         1.51              -1.48         -0.473 
##  2         0.304         -0.136        -0.260              0.390        -0.0503
##  3        -0.326         -0.189        -0.584             -1.09          1.24  
##  4         0.127         -0.0391       -0.00847           -1.46         -0.175 
##  5         0.0509        -0.254         0.947             -0.608        -0.125 
##  6        -0.167          0.549         1.49              -0.830         0.266 
##  7        -0.0132        -0.487        -0.389             -0.832        -0.0537
##  8         0.689         -0.238        -0.215             -2.99          0.0304
##  9        -0.541          0.434        -0.508             -0.981        -0.361 
## 10        -0.165         -0.0847       -0.159             -1.76         -0.178 
## # … with 20 more rows, and 46 more variables: `CD45(In115)Di` <dbl>,
## #   `CD19(Nd142)Di` <dbl>, `CD22(Nd143)Di` <dbl>, `IgD(Nd145)Di` <dbl>,
## #   `CD79b(Nd146)Di` <dbl>, `CD20(Sm147)Di` <dbl>, `CD34(Nd148)Di` <dbl>,
## #   `CD179a(Sm149)Di` <dbl>, `CD72(Eu151)Di` <dbl>, `IgM(Eu153)Di` <dbl>,
## #   `Kappa(Sm154)Di` <dbl>, `CD10(Gd156)Di` <dbl>, `Lambda(Gd157)Di` <dbl>,
## #   `CD24(Dy161)Di` <dbl>, `TdT(Dy163)Di` <dbl>, `Rag1(Dy164)Di` <dbl>,
## #   `PreBCR(Ho165)Di` <dbl>, `CD43(Er167)Di` <dbl>, `CD38(Er168)Di` <dbl>,
## #   `CD40(Er170)Di` <dbl>, `CD33(Yb173)Di` <dbl>, `HLA-DR(Yb174)Di` <dbl>,
## #   Time <dbl>, Cell_length <dbl>, `cPARP(La139)Di` <dbl>,
## #   `pPLCg2(Pr141)Di` <dbl>, `pSrc(Nd144)Di` <dbl>, `pSTAT5(Nd150)Di` <dbl>,
## #   `Ki67(Sm152)Di` <dbl>, `pErk12(Gd155)Di` <dbl>, `pSTAT3(Gd158)Di` <dbl>,
## #   `pAKT(Tb159)Di` <dbl>, `pBLNK(Gd160)Di` <dbl>, `pSyk(Dy162)Di` <dbl>,
## #   `tIkBa(Er166)Di` <dbl>, `pP38(Tm169)Di` <dbl>, `pBTK(Yb171)Di` <dbl>,
## #   `pS6(Yb172)Di` <dbl>, `pCrkL(Lu175)Di` <dbl>, `pCREB(Yb176)Di` <dbl>,
## #   `DNA1(Ir191)Di` <dbl>, `DNA2(Ir193)Di` <dbl>, `Viability1(Pt195)Di` <dbl>,
## #   `Viability2(Pt196)Di` <dbl>, wanderlust <dbl>, condition <chr>
# Finds the KNN density estimation for each cell, ordered by column, in the 
# original data matrix
wand.knn.density <- GetKnnDe(nn.matrix = wand.nn)
str(wand.knn.density)
##  num [1:1000] 0.289 0.271 0.189 0.216 0.229 ...