1 Introduction

scRepertoire is designed to take filter contig outputs from the 10x Genomics Cell Ranger pipeline, processes that data to assign clonotype based on two TCR or Ig chains, and analyze the clonotype dynamics. The latter can be separated into 1) clonotype-only analysis functions, such as unique clonotypes or clonal space quantification and 2) interaction with mRNA expression data using Seurat, SingleCellExperiment or Monocle 3 packages.

1.1 Loading Libraries

suppressMessages(library(scRepertoire))

1.2 Loading and Processing Contig Data

1.2.1 What data to load into scRepertoire?

scRepertoire functions using the filtered_contig_annotations.csv output from the 10x Genomics Cell Ranger. This file is located in the ./outs/ directory of the VDJ alignment folder. To generate a list of contigs to use for scRepertoire, just 1) load the filtered_contig_annotations.csv for each of your samples and then 2) make a list. For example:

S1 <- read.csv(".../Sample1/outs/filtered_contig_annotations.csv")
S2 <- read.csv(".../Sample2/outs/filtered_contig_annotations.csv")
S3 <- read.csv(".../Sample3/outs/filtered_contig_annotations.csv")
S4 <- read.csv(".../Sample4/outs/filtered_contig_annotations.csv")

contig_list <- list(S1, S2, S3, S4)

1.2.2 Multiplexed Experiment?

It is now easy to create the contig list from a multiplexed experiment by first generating a single-cell RNA object (either Seurat or Single Cell Experiment), loading the filtered contig file and then using createHTOContigList(). This function will return a list, separated by the group.by variable(s).

This function depends on the match of barcodes between the single-cell object and contigs. If there is a prefix or different suffix added to the barcode, this will result in no contigs recovered. As of right now, it is recommended you do this step before the integrated, as integration workflows commonly alter the barcodes. There is a multi.run variable that can be used on the integrated object, but it assumes you have modified the barcodes with the Seurat pipeline (automatic addition of _# to end) and your contig list is in the same order.

contigs <- read.csv(".../outs/filtered_contig_annotations.csv")

contig.list <- createHTOContigList(contigs, Seurat.Obj, group.by = "HTO_maxID")

scRepertoire comes with a data set derived from T cells derived from three patients with renal clear cell carcinoma in order to demonstrate the functionality of the R package. More information on the data set can be found at the corresponding manuscript. The samples consist of paired peripheral-blood and tumor-infiltrating runs, effectively creating 6 distinct runs for T cell receptor (TCR) enrichment. We can preview the elements in the list by using the head function and looking at the first contig annotation.

data("contig_list") #the data built into scRepertoire

head(contig_list[[1]])
##            barcode is_cell                   contig_id high_confidence length
## 1 AAACCTGAGAGCTGGT    TRUE AAACCTGAGAGCTGGT-1_contig_1            TRUE    705
## 2 AAACCTGAGAGCTGGT    TRUE AAACCTGAGAGCTGGT-1_contig_2            TRUE    502
## 3 AAACCTGAGCATCATC    TRUE AAACCTGAGCATCATC-1_contig_1            TRUE    693
## 4 AAACCTGAGCATCATC    TRUE AAACCTGAGCATCATC-1_contig_2            TRUE    567
## 5 AAACCTGAGCATCATC    TRUE AAACCTGAGCATCATC-1_contig_5            TRUE    361
## 6 AAACCTGAGTGGTCCC    TRUE AAACCTGAGTGGTCCC-1_contig_1            TRUE    593
##   chain   v_gene d_gene  j_gene c_gene full_length productive             cdr3
## 1   TRB TRBV20-1  TRBD1 TRBJ1-5  TRBC1        TRUE       TRUE CSASMGPVVSNQPQHF
## 2   TRB     None   None TRBJ1-5  TRBC1       FALSE       None             None
## 3   TRB  TRBV5-1  TRBD2 TRBJ2-2  TRBC2        TRUE       TRUE  CASSWSGAGDGELFF
## 4   TRA TRAV12-1   None  TRAJ37   TRAC        TRUE       TRUE CVVNDEGSSNTGKLIF
## 5   TRB     None   None TRBJ1-5  TRBC1       FALSE       None             None
## 6   TRB  TRBV7-9  TRBD1 TRBJ2-5  TRBC2        TRUE       TRUE CASSPSEGGRQETQYF
##                                            cdr3_nt reads umis raw_clonotype_id
## 1 TGCAGTGCTAGCATGGGACCGGTAGTGAGCAATCAGCCCCAGCATTTT 16718    6      clonotype96
## 2                                             None  6706    3      clonotype96
## 3    TGCGCCAGCAGCTGGTCAGGAGCGGGAGACGGGGAGCTGTTTTTT 26719   11      clonotype97
## 4 TGTGTGGTGAACGATGAAGGCTCTAGCAACACAGGCAAACTAATCTTT 18297    6      clonotype97
## 5                                             None   882    4      clonotype97
## 6 TGTGCCAGCAGCCCCTCCGAAGGGGGGAGACAAGAGACCCAGTACTTC 11218    6      clonotype98
##          raw_consensus_id
## 1 clonotype96_consensus_1
## 2                    None
## 3 clonotype97_consensus_2
## 4 clonotype97_consensus_1
## 5                    None
## 6 clonotype98_consensus_1

Some workflows will have the additional labeling of the standard barcode. Before we proceed, we can use the function stripBarcode() in order to avoid any labeling issues down the line. Importantly, stripBarcode() is for removing prefixes on barcodes that have resulted from other pipelines.

No need for stripBarcode function, if the barcodes look like:
* AAACGGGAGATGGCGT-1
* AAACGGGAGATGGCGT

In terms of using stripBarcode(), please think about the following parameters.

column
* The column in which the barcodes are present

connector
* The character that is connecting the barcode with the prefix

num_connects
* The levels of barcode prefix, where X_X_AAACGGGAGATGGCGT-1 == 3, X_AAACGGGAGATGGCGT-1 = 2.

for (i in seq_along(contig_list)) {
    contig_list[[i]] <- stripBarcode(contig_list[[i]], column = 1, connector = "_", num_connects = 3)
}

1.3 Loading from other pipelines

New to this version of the package is support for AIRR, WAT3R, and TRUST4 single-cell outputs. By indicating the directory or directories the outputs from these pipelines are located, loadContigs() will generate a contig list that is compatible with the rest of scRepertoire. In addition, although not necessary, loadContigs() can be used for 10X filtered_contig_annotation.csv files.

contig.output <- c("./Sample1/outs/", "./Sample2/outs/", "./Sample3/outs/")
contig.list <- loadContigs(dir = contig.output, 
                           format = "TRUST4")

2 Combining the Contigs

As the output of CellRanger are quantifications of both the TCRA and TCRB chains, the next step is to create a single list object with the TCR genes (comprised of the VDJC genes) and CDR3 sequences by cell barcode. This is performed using the combineTCR(), where the input is the stripped contig_list. There is also the relabeling of the barcodes by sample and ID information to prevent duplicates.

cells
* T-AB - T cells, alpha-beta TCR
* T-GD - T cells, gamma-delta TCR

removeNA
* TRUE - this is a stringent filter to remove any cell barcode with an NA value in at least one of the chains
* FALSE - the default setting to include and incorporate cells with 1 NA value

removeMulti
* TRUE - this is a stringent filter to remove any cell barcode with more than 2 immune receptor chains
* FALSE - the default setting to include and incorporate cells with > 2 chains

filterMulti
* TRUE - Isolated the top 2 expressed chains in cell barcodes with multiple chains
* FALSE - the default setting to include and incorporate cells with > 2 chains

combined <- combineTCR(contig_list, 
                        samples = c("PY", "PY", "PX", "PX", "PZ","PZ"), 
                        ID = c("P", "T", "P", "T", "P", "T"), 
                        cells ="T-AB")

The output of combineTCR() will be a list of contig data frames that will be reduced to the reads associated with a single cell barcode. It will also combine the multiple reads into clonotype calls by either the nucleotide sequence (CTnt), amino acid sequence (CTaa), the VDJC gene sequence (CTgene) or the combination of the nucleotide and gene sequence (CTstrict). The analogous function for B cells, combineBCR() functions similarly with 2 major caveats: 1) Each barcode can only have a maximum of 2 sequences, if greater exists, the 2 with the highest reads are selected. 2) The strict definition of clonotype (CTstrict) is based on the V gene and >85% normalized Levenshtein distance of the nucleotide sequence. The Levenshtein distance is calculated across all BCR sequences recovered, regardless of the run.

3 Other Processing Functions

3.1 Adding Additional Variables

What if there are more variables to add than just sample and ID? We can add them by using the addVariable() function. All we need is the name of the variable you’d like to add and the specific character or numeric values (variables). As an example, here we add the batches in which the samples were processed and sequenced.

example <- addVariable(combined, 
                      name = "batch", 
                      variables = c("b1", "b1", "b2", "b2", "b2", "b2"))

example[[1]][1:5,ncol(example[[1]])] # This is showing the first 5 values of the new column added
## [1] "b1" "b1" "b1" "b1" "b1"

3.2 Subsetting Contigs

Likewise we can remove specific list elements after combineTCR() using the subsetContig() function. In order to subset, we need to identify the vector we would like to use for subsetting (name) and also the variable values to subset (variables). Below you can see us isolate just the 4 sequencing results from PX and PY.

subset <- subsetContig(combined, 
                       name = "sample", 
                      variables = c("PX", "PY"))

4 Visualizing Contigs

cloneCall
* “gene” - use the VDJC genes comprising the TCR/Ig
* “nt” - use the nucleotide sequence of the CDR3 region
* “aa” - use the amino acid sequence of the CDR3 region
* “strict” - use the VDJC genes comprising the TCR + the nucleotide sequence of the CDR3 region. This is the proper definition of clonotype. For combineBCR() strict refers to the edit distance clusters + Vgene of the Ig.

Important to note, that the clonotype is called using essentially the combination of genes or nt/aa CDR3 sequences for both loci. As of this implementation of scRepertoire, clonotype calling is not incorporating small variations within the CDR3 sequences. As such the gene approach will be the most sensitive, while the use of nt or aa moderately so, and the most specific for clonotypes being strict. Additionally, the clonotype call is trying to incorporate both loci, i.e, both TCRA and TCRB chains and if a single cell barcode has multiple sequences identified (i.e., 2 TCRA chains expressed in one cell). Using the 10x approach, there is a subset of barcodes that only return one of the immune receptor chains, the unreturned chain is assigned an NA value.

4.1 Quantify Clonotypes

The first function to explore the clonotypes is quantContig() to return the total or relative numbers of unique clonotypes.

scale
* TRUE - relative percent of unique clonotypes scaled by total size of the size of the clonotype repertoire
* FALSE - Report the total number of unique clonotypes

Important to note, that the clonotype is called using essentially the combination of genes or nt/aa CDR3 sequences for both loci. As of this implementation of scRepertoire, clonotype calling is not incorporating small variations within the CDR3 sequences. As such the gene approach will be the most sensitive, while the use of nt or aa moderately so, and the most specific for clonotypes being strict. Additionally, the clonotype call is trying to incorporate both loci, i.e, both TCRA and TCRB chains and if a single cell barcode has multiple sequences identified (i.e., 2 TCRA chains expressed in one cell). Using the 10x approach, there is a subset of barcodes that only return one of the immune receptor chains, the unreturned chain is assigned an NA value.

4.2 Quantify Clonotypes

The first function to explore the clonotypes is quantContig() to return the total or relative numbers of unique clonotypes.

scale
* TRUE - relative percent of unique clonotypes scaled by total size of the size of the clonotype repertoire
* FALSE - Report the total number of unique clonotypes

chain
+ “both” for combined chain visualization
+ “TRA”, “TRB”, “TRD”, “TRG”, “IGH” or “IGL” to select single chain

quantContig(combined, 
            cloneCall="strict", 
            chain = "both", 
            scale = TRUE)

Within each of the general analysis functions, there is the ability to export the data frame used to create the visualization. To get the exported values, use exportTable = TRUE. It will return the data frame used to make the graph, instead of the visual output.

quantContig_output <- quantContig(combined, 
                                  cloneCall="strict", 
                                  scale = TRUE, 
                                  exportTable = TRUE)
quantContig_output
##   contigs values total   scaled
## 1    2712   PY_P  3208 84.53865
## 2    1585   PY_T  3119 50.81757
## 3     823   PX_P  1068 77.05993
## 4     918   PX_T  1678 54.70799
## 5    1143   PZ_P  1434 79.70711
## 6     768   PZ_T  2768 27.74566

The other option here is to be able to define the visualization by data classes. Here we used the combineTCR() to define the ID variable as part of the naming structure. We can the group.by to specifically use a column in the data set to organize the visualization.

quantContig(combined, 
            cloneCall="gene", 
            group.by = "ID", 
            scale = TRUE)

4.3 Clonotype Abundance

We can also examine the relative distribution of clonotypes by abundance. Here abundanceContig() will produce a line graph with a total number of clonotypes by the number of instances within the sample or run. Like above, we can also group.by this by vectors within the contig object using the group.by variable in the function.

abundanceContig(combined, 
                cloneCall = "gene", 
                scale = FALSE)

abundanceContig(combined, 
                cloneCall = "gene", 
                group.by = "ID", 
                scale = FALSE)

As you can see the peripheral blood sample derived from patient 1 is a relative extreme outlier. Another method to examine the relative abundance is to look at the density by using the scale call in the function.

abundanceContig(combined, 
                group.by = "ID", 
                scale = TRUE)

4.4 Length of Clonotypes

We can look at the length distribution of the CDR3 sequences by calling the lengtheContig() function. Importantly, unlike the other basic visualizations, the cloneCall can only be “nt” or “aa”. Due to the method of calling clonotypes as outlined above, the length should reveal a multimodal curve, this is a product of using the NA for the unreturned chain sequence and multiple chains within a single barcode.

chain
+ “both” for combined chain visualization
+ “TRA”, “TRB”, “TRD”, “TRG”, “IGH” or “IGL” to select single chain

lengthContig(combined, 
             cloneCall="aa", 
             chain = "both") 

lengthContig(combined, 
             cloneCall="aa", 
             chain = "TRA") 

4.5 Compare Clonotypes

We can also look at clonotypes between samples and changes in dynamics by using the compareClonotypes() function.

samples
* Can be used to isolate specific samples based on the name of the list element

graph
* “alluvial” - graph imaged below
* “area” - graph by area of the respective clonotype

number
* The top number of clonotypes to graph, this will be calculated based on the frequency of the individual sample. This can also be left blank.

We can also look at clonotypes between samples and changes in dynamics by using the compareClonotypes() function.

samples
* Can be used to isolate specific samples based on the name of the list element

graph
* “alluvial” - graph imaged below
* “area” - graph by area of the respective clonotype

number
* The top number of clonotypes to graph, this will be calculated based on the frequency of the individual sample. This can also be left blank.

clonotypes
* Can be used to isolate specific clonotype sequences, ensure the call matches the sequences you would like to visualize.

compareClonotypes(combined, 
                  numbers = 10, 
                  samples = c("PX_P", "PX_T"), 
                  cloneCall="aa", 
                  graph = "alluvial")

4.6 Visualize Gene Usage

Last of the basic analysis visualizations is the relative usage of genes of the TCR or BCR, using vizGenes().

gene
+ “V”
+ “D”
+ “J”
+ “C”

chain + “TRB”
+ “TRA”
+ “TRG”
+ “TRD” + “IGH”
+ “IGL”

plot
+ “bar” for a bar chart
+ “heatmap” for a heatmap

Last of the basic analysis visualizations is the relative usage of genes of the TCR or BCR, using vizGenes().

gene
+ “V”
+ “D”
+ “J”
+ “C”

chain + “TRB”
+ “TRA”
+ “TRG”
+ “TRD” + “IGH”
+ “IGL”

plot
+ “bar” for a bar chart
+ “heatmap” for a heatmap

y.axis
Variable to separate the counts along the y-axis. If not specified the default is to separate by the list element names.

scale
+ TRUE to scale the graph by number of genes per sample
+ FALSE to report raw numbers

order
+ “gene” to order by gene name
+ “variance” to order by variance between the separate variable categories

vizGenes(combined, gene = "V", 
         chain = "TRB", 
         plot = "bar", 
         scale = TRUE)

We can also use vizGenes() to look at the usage of genes in a single chain. So for example say we are interesting in the difference in TRB V and J usage between tumor and peripheral blood samples - we can easily take a look at this using the following code:

#Peripheral Blood
vizGenes(combined[c(1,3,5)], 
         gene = "V", 
         chain = "TRB", 
         y.axis = "J", 
         plot = "heatmap", 
         scale = TRUE, 
         order = "variance")

#Tumor Infiltrating
vizGenes(combined[c(2,4,6)], 
         gene = "V", 
         chain = "TRB", 
         y.axis = "J", 
         plot = "heatmap", 
         scale = TRUE, 
         order = "variance")


5 More Advanced Clonal Analysis

After we have completed the basic processing and summary functions in scRepertoire, we can begin to explore the clonotypes of the single-cell data in more detail.

5.1 Clonal Space Homeostasis

By examining the clonal space, we are effectively looking at the relative space occupied by clones at specific proportions. Another way to think about this would be thinking of the total immune receptor sequencing run as a measuring cup. In this cup, we will fill liquids of different viscosity - or different number of clonal proportions. Clonal space homeostasis is asking what percentage of the cup is filled by clones in distinct proportions (or liquids of different viscosity, to extend the analogy). The proportional cutpoints are set under the cloneType variable in the function and can be adjusted, at baseline the bins are as follows:

cloneTypes
* Rare = .0001
* Small = .001
* Medium = .01
* Large = .1
* Hyperexpanded = 1

clonalHomeostasis(combined, cloneCall = "gene")

clonalHomeostasis(combined, cloneCall = "aa")

5.2 Clonal Proportion

Like clonal space homeostasis above, clonal proportion acts to place clones into separate bins. The key difference is instead of looking at the relative proportion of the clone to the total, the clonalProportion() function will rank the clones by total number and place them into bins.

The split represents ranking of clonotypes by copy or frequency of occurrence, meaning 1:10 are the top 10 clonotypes in each sample. The default bins are under the split variable in the function and can be adjusted, but at baseline they are as follows.

split
* 10
* 100
* 1000
* 10000 * 30000 * 100000

clonalProportion(combined, cloneCall = "gene") 

clonalProportion(combined, cloneCall = "nt") 

5.3 Overlap Analysis

If you are interested in measures of similarity between the samples loaded into scRepertoire, using clonalOverlap() can assist in the visualization. Three methods currently can be performed in clonalOverlap() 1) overlap coefficient, 2) Morisita index, or 3) Jaccard index. The former is looking at the overlap of clonotypes scaled to the length of unique clonotypes in the smaller sample. The Morisita index is more complex, it is an ecological measure of the dispersion of individuals within a population, incorporating the size of the population. The Jaccard Similarity Index is very similar to the overlap coefficient - instead of using the length of the smaller sample, the denominator for the Jaccard Index is the union of the two comparisons, leading to a much smaller number.

clonalOverlap(combined, 
              cloneCall = "strict", 
              method = "morisita")

Another addition to scRepertoire is the ability to cluster the samples by the clone size distribution using clonesizeDistribution() adapted from the powerTCR R package. Please read and cite the respective citation if using this function to analyze the similarities of sample clone size distributions. In this function, method refers to the method for which the hierarchical clustering will be based on.

clonesizeDistribution(combined, 
                      cloneCall = "strict", 
                      method="ward.D2")

## NULL

5.4 Diversity Analysis

Diversity can also be measured for samples or by other variables. Diversity is calculated using five metrics: 1) Shannon, 2) inverse Simpson, 3) Chao1, and 4) Abundance-based Coverage Estimator (ACE), and 5) inverse Pielou. With the former two generally used to estimate baseline diversity and Chao/ACE indices used to estimate the richness of the samples. New implementation of this function include downsampling with 100 boot straps (n.boots) using the minimum number of unique clonotypes, as a more robust diversity estimate. To get a table of the bootstrapped values, set return.boots = TRUE.

clonalDiversity(combined, 
                cloneCall = "gene", 
                group.by = "sample", 
                x.axis = "ID", 
                n.boots = 100)

5.5 Scatter Compare

A number of users requested a visualization from the work of Wu, et al 2020, PMID: 32103181 that allows for the direct comparison of clonotypes between 2 samples. scatterClonotype() will organize two samples from the combineTCR/BCR product, count the relative clonotypes, and produce a scatter plot comparing the two samples.

x.axis and y.axis
* Names of the list element to place on the x-axis and y-axis - so for example “PX_P” and “PX_T”

dot.size
* “total” to display the total number of clones between the the x- and y-axis
* Name of the list element to use for size calculation

graph
* “proportion” for the relative proportion the clonotype represents across all clonotypes
* “count” for the total count of clonotypes by sample

scatterClonotype(combined, cloneCall ="gene", 
                 x.axis = "PY_P", 
                 y.axis = "PY_T",
                 dot.size = "total",
                 graph = "proportion")


6 Interacting with mRNA Expression

As mentioned previously, this data set is derived from work performed in the laboratory of Weizhou Zhang. For the purposes of the vignette, we have randomly sampled 100 cells from the the fully integrated Seurat object to minimize the size of the package. A full version of the Seurat object is available via GitHub under the Getting Data header. We will use both Seurat and SingleCellExperiment (SCE) with scater to perform the further visualizations in tandem.

library(Seurat)
library(scater)
screp_example <- get(data("screp_example"))
sce <- suppressMessages(UpdateSeuratObject(screp_example))
sce <- as.SingleCellExperiment(screp_example)

#Seurat Format
DimPlot(screp_example)

##Single Cell Experiment Format
plotUMAP(sce, colour_by = "seurat_clusters")

Here you can see we have 12 total clusters (C1-12), which we have labeled as such for simplicity. We can also get a little more granular information on the number of cells by using the table() function.

table(screp_example$seurat_clusters)
## 
##  0  1  2  3  4  5  6  7  8  9 10 11 
## 15 20  8 11 11 12 11  4  4  2  2  0

Next we can take the clonotypic information and attach it to our Seurat object using the combineExpression() function. Importantly, the major requirement for the attachment is matching contig cell barcodes and barcodes in the row names of the meta data of the Seurat or SCE object. If these do not match, the attachment will fail. Based on ease, we suggest you make the changes to the Seurat object row names.

We can call (cloneCall) the 4 variations of clonotypes: 1) VDJC genes, 2) CDR3 amino acid sequence, 3) CDR3 nucleotide sequence, or 4) VDJC genes and CDR3 nucleotide sequence. The attaching function will also calculate the frequency of the clonotype based on the group.by variable. If blank, group.by will calculate frequencies of clonotypes by individual run, but because we have 6 samples of paired peripheral and tumor T cells, we are actually going to use the group.by variable to call “sample” in order to calculate frequencies across both the peripheral blood and tumor T cells of the same patient.

In order to categorize the frequency, we have the variable proportion which if TRUE allows for the relative proportion or when FALSE will use absolute frequency to define clonotype groups cloneTypes acts as a bin to place labels. As a default, cloneTypes is set to equal cloneTypes=c(Rare = 1e-4, Small = 0.001, Medium = 0.01, Large = 0.1, Hyperexpanded = 1). However, below you can see an example of using total frequency as expansion assignments.

screp_example <- combineExpression(combined, 
                                   screp_example, 
                                   cloneCall="gene", 
                                   group.by = "sample", 
                                   proportion = FALSE, 
                                   cloneTypes=c(Single=1, Small=5, Medium=20, Large=100, Hyperexpanded=500))

sce <- combineExpression(combined, 
                          sce, 
                          cloneCall="gene", 
                          group.by = "sample", 
                          proportion = TRUE)

6.1 Combining both TCR and BCR

If you have performed TCR/BCR enrichment or want to add info for both gamma-delta and alpha-beta T cells, now we can just make a single list and use combineExpression(). Major note if there are duplicate barcodes (if a cell has both Ig and TCR), the immune receptor information will not be added. As an anecdote, the testing data we used to improve this function had 5-6% of barcode overlap.

#This is an example of the process, will not be evaluated during knit
TCR <- combineTCR(...)
BCR <- combineBCR(...)
list.receptors <- c(TCR, BCR)


seurat <- combineExpression(list.receptors, 
                            seurat, 
                            cloneCall="gene", 
                            proportion = TRUE)

We first want to look at the distribution of peripheral versus tumor T cells. We can use the same color scheme as the rest of the scRepertoire package by calling the object colorblind_vector using the following hex codes.

colorblind_vector <- colorRampPalette(rev(c("#0D0887FF", "#47039FFF", 
              "#7301A8FF", "#9C179EFF", "#BD3786FF", "#D8576BFF",
              "#ED7953FF","#FA9E3BFF", "#FDC926FF", "#F0F921FF")))

DimPlot(screp_example, group.by = "Type") + NoLegend() +
    scale_color_manual(values=colorblind_vector(2))

Now we can look at the distribution of the clonotype bins by first ordering the clonoType as a factor, this prevents the coloring from being in alphabetical order. Next we use the DimPlot() function call in Seurat with our scale_color_manual additional layer.

slot(screp_example, "meta.data")$cloneType <- factor(slot(screp_example, "meta.data")$cloneType, 
                levels = c("Hyperexpanded (100 < X <= 500)", 
                           "Large (20 < X <= 100)", 
                           "Medium (5 < X <= 20)", 
                           "Small (1 < X <= 5)", 
                           "Single (0 < X <= 1)", NA))

DimPlot(screp_example, group.by = "cloneType") 

6.2 clonalOverlay

Using the dimensional reduction graphs as a reference, we can also generate an overlay of the position of clonally expanded cells using clonalOverlay(). Select the reduction for the visualization, default is “PCA” and the freq.cutpoint or lowest clonal frequency or proportion to generate the contour plot. We can modify the contours by selecting the number of bins or number of contours drawn. clonalOverlay() can be used to look across all cells or faceted by a meta data variable using facet. As we facet, the overall dimensional reduction will be maintained, while the contour plots will adjust based on the facet variable. Coloring of the dot plot is based on the active identity of the single-cell object. This visualization was authored by Dr. Francesco Mazziotta from Johns Hopkins and inspired by Drs. Carmona and Andreatta and their work with ProjectTIL.

clonalOverlay(screp_example, 
              reduction = "umap", 
              freq.cutpoint = 1, 
              bins = 10, 
              facet = "Patient") + 
              guides(color = "none")

6.3 clonalNetwork

Similar to clonalOverlay(), we can look at the network interaction of clonotypes shared between clusters along the single-cell dimensional reduction using clonalNetwork(). This function shows the relative proportion of clones that come from the starting node, with the ending node indicated by the arrow.

Filtering for clones can be accomplished using 3 methods:
filter.clones
* Select a number to isolate the clones comprising the top X number of cells (filter.clones = 2000)
* Select “min” to make sure all groups are scaled to the size of the minimum group

filter.identity
* For the identity chosen to visualize, show the to and from network connections for a single group

filter.proportion * Remove clones that comprise less than a certain proportion of clones in groups.

#ggraph needs to be loaded due to issues with ggplot
library(ggraph)

#No Identity filter
clonalNetwork(screp_example, 
              reduction = "umap", 
              identity = "seurat_clusters",
              filter.clones = NULL,
              filter.identity = NULL,
              cloneCall = "aa")

#Examining Cluster 2 only
clonalNetwork(screp_example, 
              reduction = "umap", 
              identity = "ident",
              filter.identity = "C2",
              cloneCall = "aa")

We can also look at the clonotypes by calling specific sequences in the highlightClonotypes() below. In order to highlight the clonotypes, we first need to use the cloneCall the type of sequence we will be using and then the specific sequences themselves using sequence. Below you can see the steps to highlight the two most prominent sequences “CAVNGGSQGNLIF_CSAEREDTDTQYF” with a frequency = 482 and “NA_CATSATLRVVAEKLFF” with a frequency = 287.

screp_example <- highlightClonotypes(screp_example, 
                    cloneCall= "aa", 
                    sequence = c("CAVNGGSQGNLIF_CSAEREDTDTQYF", 
                                 "NA_CATSATLRVVAEKLFF"))

Seurat::DimPlot(screp_example, group.by = "highlight")

6.4 occupiedscRepertoire

We can also look at the count of cells by cluster assigned into specific frequency ranges, by using the occupiedscRepertoire() function and selecting the x.axis to display cluster or other variables in the meta data of the single cell object.

  • proportion can be used to looked at a relative level groupings
  • label will still return absolute number of clonotypes
occupiedscRepertoire(screp_example, 
                     x.axis = "seurat_clusters")

occupiedscRepertoire(screp_example, 
                     x.axis = "ident", 
                     proportion = TRUE, 
                     label = FALSE)

6.5 alluvialClonotypes

After the metadata has been modified, we can look at clonotypes across multiple categories using the alluvialClonotypes() function. To understand the basic concepts of this graphing method, I’d highly recommend reading this post, essentially we are able to use the plots to examine the interchange of categorical variables. Because this function will produce a graph with each clonotype arranged by called stratification, this will take some time depending on the size of your total cells. To expedite, we will actually subset the Seurat object before using alluvialClonotypes().

alluvialClonotypes(screp_example, 
                   cloneCall = "aa", 
                   y.axes = c("Patient", "ident", "Type"), 
                   color = "CAVNGGSQGNLIF_CSAEREDTDTQYF") + 
    scale_fill_manual(values = c("grey", colorblind_vector(5)[5]))

alluvialClonotypes(screp_example, 
                   cloneCall = "gene", 
                   y.axes = c("Patient", "ident", "Type"), 
                   color = "seurat_clusters") 

6.6 getCirclize

Like alluvial graphs, we can also visualize the interconnection of clusters using the chord diagrams from the circlize R package. The first step is getting the data frame output to feed into the chordDiagram() function in circlize, which can be done using getCirclize(). This will calculate the relative number of clonotypes shared based on the group.by variable using the product of combineExpression(). Important to note, getCirclize() will create a matrix the size of the group.by variable and then simplify into set of instructions to be read by the circlize R package. The output of which is the total number of unique and shared clones by the group.by variable - this will be greater than the number of unique clonotypes due to the shared clone information.

library(circlize)
library(scales)

circles <- getCirclize(screp_example, 
                       group.by = "seurat_clusters")

#Just assigning the normal colors to each cluster
grid.cols <- hue_pal()(length(unique(screp_example$seurat_clusters)))
names(grid.cols) <- unique(screp_example$seurat_clusters)

#Graphing the chord diagram
chordDiagram(circles, self.link = 1, grid.col = grid.cols)

This can also be used if we wanted to explore just the tumor-specific T cells by just subsetting the single-cell object.

subset <- subset(screp_example, Type == "T")


circles <- getCirclize(subset, group.by = "ident")

grid.cols <- scales::hue_pal()(length(unique(subset@active.ident)))
names(grid.cols) <- levels(subset@active.ident)

chordDiagram(circles, 
             self.link = 1, 
             grid.col = grid.cols, 
             directional = 1, 
             direction.type =  "arrows",
             link.arr.type = "big.arrow")

6.7 Clustering Clonotypes

The nucleotide or amino acid sequences of the chains can be used to cluster clonotypes, by examining the edit distance of the sequences. This approach is underlying the combineBCR() function, but now can be applied to the T cell receptors. It will add to the end of each list element a cluster by generating a network connected by the similarity in sequence. This network is directed by the threshold variable, where 0.85 is the normalized mean edit distance.

Clonotype Cluster by distance:

  • TCRA:LD.1 [chain:] + :LD + [number]
  • TCRA.123 [chain:] + [number]

LD denotes if the cluster was called using the normalized Levenshtein distance

If performed over the number of samples, such as the list elements, group.by can used to calculate only the clusters on the setting of patient sample (group.by = “Patient”) or tissue type (group.by = “Type”). This will add the selected group to the beginning of the cluster designation.

sub_combined <- clusterTCR(combined[[2]], 
                           chain = "TRA", 
                           sequence = "aa", 
                           threshold = 0.85, 
                           group.by = NULL)

6.8 StartracDiversity

From the excellent work by Lei Zhang, et al in Lineage tracking reveals dynamic relationships of T cells in colorectal cancer, the authors introduce new methods for looking at clonotypes by cellular origins and cluster identification. Their startrac software has been incorporated into scRepertoire. If you are using the output of this specific function, please cite their excellent work.

In order to use the StartracDiversity() function, you will need to include the product of the combinedSeurat() function. The second requirement is a column header in the meta data of the Seurat object that has tissue of origin. In the example data,type corresponds to the column “Type”, which includes the “P” and “T” classifier. The indices can be subseted for a specific patient or examined overall using the by variable. Importantly, the function uses only the strict definition of clonotype of the VDJC genes and the CDR3 nucleotide sequence.

The indices output includes:
+ expa - Clonal Expansion
+ migr - Cross-tissue Migration
+ tran - State Transition

StartracDiversity(screp_example, 
                  type = "Type", 
                  sample = "Patient", 
                  by = "overall")
## [2022-11-01 18:57:01] initialize Startrac ...
## [2022-11-01 18:57:01] calculate startrac index ...
## [2022-11-01 18:57:01] calculate pairwise index ...
## [2022-11-01 18:57:01] calculate indices of each patient ...
## [2022-11-01 18:57:01] collect result
## [2022-11-01 18:57:01] return

6.9 Clonotype Bias

A new metric proposed by Massimo et al, clonal bias, like Startrac is a clonal metric that seeks to quantify how individual clones are skewed towards a specific cellular compartment or cluster. With the minimal example present in the in-built vignette, this will not be as enlightening as when we use a full single-cell object.

split.by
* Variable used to calculate the baseline frequencies

group.by
* The compartment/variable used for the purpose of the comparison

clonotypeBias(screp_example, 
              cloneCall = "aa", 
              split.by = "Type", 
              group.by = "seurat_clusters",
              n.boots = 20, 
              min.expand =1)

7 Working with clonotypes after clustering

For users that would like greater ability to use the meta data in the Seurat objects to perform the analysis that scRepertoire provides, there is also the option of using the expression2List() function that will take the meta data and output the data as a list by cluster.

combined2 <- expression2List(screp_example, split.by = "seurat_clusters")
combined3 <- expression2List(sce, split.by = "seurat_clusters")

7.1 1. Clonal Diversity

clonalDiversity(combined2, cloneCall = "nt")

clonalDiversity(combined3, cloneCall = "nt")

7.2 2. Clonal Homeostasis

clonalHomeostasis(combined2, cloneCall = "nt")

clonalHomeostasis(combined3, cloneCall = "nt")

7.3 3. Clonal Proportion

clonalProportion(combined2, cloneCall = "nt")

clonalProportion(combined3, cloneCall = "nt")

7.4 4. Clonal Overlap

clonalOverlap(combined2, cloneCall="aa", method="overlap")

clonalOverlap(combined3, cloneCall="aa", method="overlap")


8 Conclusion

This has been a general overview of the capabilities for scRepertoire from the initial processing and visualization to attach to the mRNA expression values in a Seurat object. If you have any questions, comments or suggestions, feel free to visit the github repository or email me.

###Session Info

sessionInfo()
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] scales_1.2.1                circlize_0.4.15            
##  [3] ggraph_2.1.0                scater_1.26.0              
##  [5] scuttle_1.8.0               SingleCellExperiment_1.20.0
##  [7] SummarizedExperiment_1.28.0 Biobase_2.58.0             
##  [9] GenomicRanges_1.50.0        GenomeInfoDb_1.34.0        
## [11] IRanges_2.32.0              S4Vectors_0.36.0           
## [13] BiocGenerics_0.44.0         MatrixGenerics_1.10.0      
## [15] matrixStats_0.62.0          sp_1.5-0                   
## [17] SeuratObject_4.1.2          Seurat_4.2.0               
## [19] scRepertoire_1.8.0          ggplot2_3.3.6              
## [21] BiocStyle_2.26.0           
## 
## loaded via a namespace (and not attached):
##   [1] utf8_1.2.2                reticulate_1.26          
##   [3] tidyselect_1.2.0          htmlwidgets_1.5.4        
##   [5] BiocParallel_1.32.0       grid_4.2.1               
##   [7] Rtsne_0.16                ScaledMatrix_1.6.0       
##   [9] munsell_0.5.0             codetools_0.2-18         
##  [11] ica_1.0-3                 future_1.28.0            
##  [13] miniUI_0.1.1.1            withr_2.5.0              
##  [15] spatstat.random_3.0-0     colorspace_2.0-3         
##  [17] progressr_0.11.0          highr_0.9                
##  [19] knitr_1.40                ggalluvial_0.12.3        
##  [21] ROCR_1.0-11               tensor_1.5               
##  [23] listenv_0.8.0             labeling_0.4.2           
##  [25] GenomeInfoDbData_1.2.9    polyclip_1.10-4          
##  [27] farver_2.1.1              parallelly_1.32.1        
##  [29] vctrs_0.5.0               generics_0.1.3           
##  [31] xfun_0.34                 R6_2.5.1                 
##  [33] doParallel_1.0.17         ggbeeswarm_0.6.0         
##  [35] graphlayouts_0.8.3        rsvd_1.0.5               
##  [37] VGAM_1.1-7                isoband_0.2.6            
##  [39] bitops_1.0-7              spatstat.utils_3.0-1     
##  [41] cachem_1.0.6              DelayedArray_0.24.0      
##  [43] assertthat_0.2.1          promises_1.2.0.1         
##  [45] beeswarm_0.4.0            rgeos_0.5-9              
##  [47] gtable_0.3.1              beachmat_2.14.0          
##  [49] globals_0.16.1            goftest_1.2-3            
##  [51] tidygraph_1.2.2           MatrixModels_0.5-1       
##  [53] rlang_1.0.6               GlobalOptions_0.1.2      
##  [55] splines_4.2.1             lazyeval_0.2.2           
##  [57] spatstat.geom_3.0-3       BiocManager_1.30.19      
##  [59] yaml_2.3.6                reshape2_1.4.4           
##  [61] abind_1.4-5               httpuv_1.6.6             
##  [63] tools_4.2.1               bookdown_0.29            
##  [65] cubature_2.0.4.5          ellipsis_0.3.2           
##  [67] spatstat.core_2.4-4       jquerylib_0.1.4          
##  [69] RColorBrewer_1.1-3        ggridges_0.5.4           
##  [71] Rcpp_1.0.9                plyr_1.8.7               
##  [73] sparseMatrixStats_1.10.0  zlibbioc_1.44.0          
##  [75] purrr_0.3.5               RCurl_1.98-1.9           
##  [77] rpart_4.1.19              deldir_1.0-6             
##  [79] pbapply_1.5-0             viridis_0.6.2            
##  [81] cowplot_1.1.1             zoo_1.8-11               
##  [83] ggrepel_0.9.1             cluster_2.1.4            
##  [85] magrittr_2.0.3            data.table_1.14.4        
##  [87] magick_2.7.3              scattermore_0.8          
##  [89] SparseM_1.81              lmtest_0.9-40            
##  [91] RANN_2.6.1                truncdist_1.0-2          
##  [93] fitdistrplus_1.1-8        gsl_2.1-7.1              
##  [95] patchwork_1.1.2           mime_0.12                
##  [97] evaluate_0.17             xtable_1.8-4             
##  [99] shape_1.4.6               gridExtra_2.3            
## [101] compiler_4.2.1            tibble_3.1.8             
## [103] crayon_1.5.2              KernSmooth_2.23-20       
## [105] htmltools_0.5.3           mgcv_1.8-41              
## [107] later_1.3.0               tidyr_1.2.1              
## [109] powerTCR_1.18.0           DBI_1.1.3                
## [111] tweenr_2.0.2              MASS_7.3-58.1            
## [113] Matrix_1.5-1              permute_0.9-7            
## [115] cli_3.4.1                 parallel_4.2.1           
## [117] evd_2.3-6.1               igraph_1.3.5             
## [119] pkgconfig_2.0.3           plotly_4.10.0            
## [121] spatstat.sparse_3.0-0     foreach_1.5.2            
## [123] vipor_0.4.5               bslib_0.4.0              
## [125] stringdist_0.9.9          XVector_0.38.0           
## [127] stringr_1.4.1             digest_0.6.30            
## [129] sctransform_0.3.5         RcppAnnoy_0.0.20         
## [131] vegan_2.6-4               spatstat.data_3.0-0      
## [133] rmarkdown_2.17            leiden_0.4.3             
## [135] uwot_0.1.14               DelayedMatrixStats_1.20.0
## [137] evmix_2.12                quantreg_5.94            
## [139] shiny_1.7.3               lifecycle_1.0.3          
## [141] nlme_3.1-160              jsonlite_1.8.3           
## [143] BiocNeighbors_1.16.0      viridisLite_0.4.1        
## [145] fansi_1.0.3               pillar_1.8.1             
## [147] lattice_0.20-45           fastmap_1.1.0            
## [149] httr_1.4.4                survival_3.4-0           
## [151] glue_1.6.2                png_0.1-7                
## [153] iterators_1.0.14          ggforce_0.4.1            
## [155] stringi_1.7.8             sass_0.4.2               
## [157] BiocSingular_1.14.0       dplyr_1.0.10             
## [159] irlba_2.3.5.1             future.apply_1.9.1