Expanded STR algorithm (exSTRa)

Rick Tankard

2019-01-18

exSTRa is an R package for the detection of repeat expansions with Illumina next-generation sequencing (NGS) data of cases and controls. exSTRa supports both whole-genome and whole-exome sequencing (WGS and WES). Only paired-end data is supported. This package implements the algorithm as described in:

Rick M. Tankard, Mark F. Bennett, Peter Degorski, Martin B. Delatycki, Paul J. Lockhart, Melanie Bahlo. Detecting Expansions of Tandem Repeats in Cohorts Sequenced with Short-Read Sequencing Data. American Journal of Human Genetics, 103(6):858-873, 2018. https://doi.org/10.1016/j.ajhg.2018.10.015

A table of repeat expansion disorders for the human genome reference hg19 is included. Other references may be provided in the future, or these may be created by the user. The use of GRCh37, that does not include the ‘chr’ prefix in chromosome names, has not yet been tested, but this is intended for the near future.

For general repeat loci, exSTRa can read from the UCSC Simple Repeats Track; either from the UCSC table browser or downloads page (although we recommend not providing the entire file due to current computational limitations). Text files without a header line are assumed to be from the Simple Repeats track, such that you may easily provide a subset of the file. By default, only STRs with 2 to 6 bp length repeat motifs will be analysed, with others filtered.

Data preparation

This section describes the steps required to process your own data. If you wish to try exSTRa with our example data, you may skip to the next section, “Using the R exSTRa package”.

Standard NGS steps

Data should be aligned to the hg19 reference genome. We recommend this is performed in local alignment mode, whicih allows more bases to be soft-clipped, and hence more tolerant of repeat expansions. We have had success finding expansions in data aligned in a non-local mode.

Optionally, duplicate marking may be performed. We have tested this with Picard tools. The effect of omitting this duplicate marking step has not been assessed.

Both BAM and CRAM formats are supported. The resulting BAM/CRAM file should be sorted and indexed, with each representing exactly one individual.

Repeat content collection with Perl package

At present, BAM/CRAM files are summarized by an external Perl package, Bio::STR::exSTRa. Future versions of exSTRa may include this within the R package. For installation instructions, see the package page https://github.com/bahlolab/Bio-STR-exSTRa.

The README.md file of Bio::STR::exSTRa, describes the steps to process the BAM/CRAM files. This will produce a tab-delimited file that is ready for analysis with the R exSTRa package. The R exSTRa package includes an example of this output; assuming exSTRa is installed, this can be found on your system with

system.file("extdata", "HiSeqXTen_WGS_PCR_2.txt", package = "exSTRa")

or online at https://github.com/bahlolab/exSTRa/blob/master/inst/extdata/HiSeqXTen_WGS_PCR_2.txt.

Using the R exSTRa package

The following describes the analysis of the WGS_PCR_2 data from the exSTRa publication described at the top of this vignette. WGS_PCR_2 is a data set of 16 case and 2 control individuals, sequenced to 60x and 30x depth (16 (14 cases) and 2 (case) samples, respectively) on an Illumina HiSeq X Ten system. A PCR library preparation protocol with Illumina TruSeq Nano was used. Data was aligned to the hg19 human genome reference sequence with Bowtie2 in very sensitive local mode, with maximum insert size 1000. The data underwent PCR duplicate marking. Repeat content of reads was obtained with Bio::STR::exSTRa.

Firstly the exSTRa package and dependencies should be installed:

  # install.packages("devtools") # if devtools is not already installed
  # install.packages("data.table") # if data.table is not already installed
    devtools::install_github("bahlolab/exSTRa")

We load required packages. The data.table package should be loaded first, as exSTRa must redefine it’s copy() function as a generic function.

library(data.table) # should be loaded before exSTRa
library(exSTRa)
#> 
#> Attaching package: 'exSTRa'
#> The following object is masked from 'package:data.table':
#> 
#>     copy

Data output from Bio::STR::exSTRa is read, along with the specification of repeat expansion loci. This results in an exstra_score object.

str_score <- read_score (
  file = system.file("extdata", "HiSeqXTen_WGS_PCR_2.txt", package = "exSTRa"), 
  database = system.file("extdata", "repeat_expansion_disorders_hg19.txt", package = "exSTRa"),
  groups.regex = c(control = "^WGSrpt_0[24]$", case = "")
)

# an exstra_score object:
str_score
#> exstra_score object with 33072 observations of type named ($data),
#>   for 18 samples. ($samples)
#>   Includes associated STR database of 21 loci. ($db)

When reading your own data, system.file(...) should be replaced with a string of the path to your data. The groups.regex vector does not currently impact the detection of repeat expansions, but does affect plotting. groups.regex should contain regular expressions to match sample names. The order of groups.regex matters, such that the first match will determine the grouping; this is why the “case” regular expression is an empty string that matches anything, hence specifying any non-control sample as a case.

Note the read_score() function automatically filters low repeat content scores according to the length of the repeat motif; shorter repeat motifs are more likely to appear by chance and hence have a higher threshold than longer repeat motifs, that are less likely to appear in a sequence by chance.

The loci are specified in str_score$db as a data.table.

str_score$db$locus
#>  [1] "DM1"     "DM2"     "DRPLA"   "EPM1A"   "FRAXA"   "FRAXE"   "FRDA"   
#>  [8] "FTDALS1" "HD"      "HDL2"    "SBMA"    "SCA1"    "SCA10"   "SCA12"  
#> [15] "SCA17"   "SCA2"    "SCA3"    "SCA36"   "SCA6"    "SCA7"    "SCA8"

ECDF plotting

We can plot an empirical distribution function (ECDF) of repeat scores for a chosen locus. This example is of theHuntington disease (HD) locus. Here, all case samples, including those for other loci, are colored red.

plot(str_score["HD"])

We may also choose a locus directly from the plot() function. Here, we additionally set the color of two known HD cases. All other samples are automatically colored a transparent black.

plot(str_score, "HD", sample_col = c("WGSrpt_10" = "red", "WGSrpt_12" = "blue"))

More information on the plot function on exstra_score objects can be accessed with

?`[.exstra_score`

For convenience and speed in this vignette, we only assess expansions for the Huntington disease, spinocerebellar ataxias 2 and 6, and Friedreich’s ataxia loci.

( str_score_four <- str_score[c("HD", "SCA2", "SCA6", "FRDA")] )
#> exstra_score object with 5468 observations of type named ($data),
#>   for 18 samples. ($samples)
#>   Includes associated STR database of 4 loci. ($db)

For creating multiple ECDF curves, with legends, you may use the plot_multi() function. Here, four ECDFs are created both to the R session (plot_type 1) and to a single PDF file example_images/HiSeqXTen_WGS_PCR_2.pdf.

par(mfrow = c(2, 2))
plot_multi(str_score_four, dir = "example_images", 
  prefix = "HiSeqXTen_WGS_PCR_2", plot_types = c(1, 2), alpha_case = 0.2)
#> Warning in dir.create(dirbase, recursive = TRUE): 'example_images' already
#> exists

Testing for significant expansions

Here, we calculate an aggregated T statistic over quantiles as described in the exSTRa paper. By default, p-values are calculated with a simulation procedure. With default settings, this may take several minutes to complete. This creates an exstra_tsum object.

( tsum <- tsum_test(str_score_four) )
#> Working on locus HD
#> Working on locus SCA2
#> Working on locus SCA6
#> Working on locus FRDA
#> exstra_tsum object with 72 T sum statistics ($stats),
#>   with p-values calculated ($stats),
#>   over 4 loci. ($db)
#> 
#>     T sum statistics summary:
#>     exSTRa T := sum of two sample t-tests
#> 
#> Alternative hypotheses: subject sample has a larger allele than background samples.
#> 
#> alpha  Bonferroni unadjusted
#> 0.0001          0          7 
#> 0.001           0          1 
#> 0.01            7          1 
#> 0.05            1          1 
#> 1              64         62 
#> NA              0          0 
#> 
#> Number of samples: 18 
#> Number of loci:    4 
#> Defined p-values:  72 
#> NA p-values:       0 
#> Function arguments: trim = 0.15, min.quant = 0.5, B = 999

Plotting an exstra_tsum object highlights significant samples, after Bonferroni correction by default.

par(mfrow = c(2, 2))
plot(tsum)

You may manually set the colors each sample will use with a vector of colors, with names the corresponding sample name.

plot_cols <- c(RColorBrewer::brewer.pal(8, "Set2"), RColorBrewer::brewer.pal(8, "Dark2"), "orange", "blue")
names(plot_cols) <- str_score_four$samples[, sample]
plot_cols
#> WGSrpt_02 WGSrpt_04 WGSrpt_05 WGSrpt_07 WGSrpt_08 WGSrpt_09 WGSrpt_10 
#> "#66C2A5" "#FC8D62" "#8DA0CB" "#E78AC3" "#A6D854" "#FFD92F" "#E5C494" 
#> WGSrpt_11 WGSrpt_12 WGSrpt_13 WGSrpt_14 WGSrpt_15 WGSrpt_16 WGSrpt_17 
#> "#B3B3B3" "#1B9E77" "#D95F02" "#7570B3" "#E7298A" "#66A61E" "#E6AB02" 
#> WGSrpt_18 WGSrpt_19 WGSrpt_20 WGSrpt_21 
#> "#A6761D" "#666666"  "orange"    "blue"
# To demonstrate the colours used:
par(mfrow = c(1, 1))
pie(rep(1, length(plot_cols)), col = plot_cols, labels = names(plot_cols), cex = 0.7)

Note that some data sets may not be able to reach significant levels after correction with the default number of simulations (9999). This can be adjusted with the B parameter of tsum_test(), or a less stringent threshold can be used. Bonferroni correction is too severe here, so we can specify Bonferroni correction only on each locus for the number of samples tested.

par(mfrow = c(2, 2))
plot(tsum, sample_col = plot_cols, correction = "samples")

Or Bonferroni correction may be applied for the number of loci tested.

par(mfrow = c(2, 2))
plot(tsum, sample_col = plot_cols, correction = "loci")

You may obtain a data.table of each sample and locus with the p-value, and if it is significant with the correction method applied. Here, the correction method is Bonferroni per locus.

(ps <- p_values(tsum, correction = "samples"))
#>     locus    sample        tsum      p.value   p.value.sd   B signif
#>  1:  FRDA WGSrpt_02  0.22438025 4.285714e-01 3.723684e-03 999  FALSE
#>  2:  FRDA WGSrpt_04  2.37758148 2.152032e-02 1.091872e-03 999  FALSE
#>  3:  FRDA WGSrpt_05 -1.07940062 8.261692e-01 2.851854e-03 999  FALSE
#>  4:  FRDA WGSrpt_07  1.57425576 8.863927e-02 2.138604e-03 999  FALSE
#>  5:  FRDA WGSrpt_08 -0.79878551 7.549908e-01 3.236465e-03 999  FALSE
#>  6:  FRDA WGSrpt_09 12.15687245 5.560807e-05 5.610840e-05 999   TRUE
#>  7:  FRDA WGSrpt_10 -0.04995430 5.193238e-01 3.759497e-03 999  FALSE
#>  8:  FRDA WGSrpt_11 13.60845349 5.560807e-05 5.610840e-05 999   TRUE
#>  9:  FRDA WGSrpt_12 -0.94879960 7.943613e-01 3.041445e-03 999  FALSE
#> 10:  FRDA WGSrpt_13  0.57407040 3.117945e-01 3.485538e-03 999  FALSE
#> 11:  FRDA WGSrpt_14  1.63576022 7.985319e-02 2.039608e-03 999  FALSE
#> 12:  FRDA WGSrpt_15 -1.03104591 8.154924e-01 2.919056e-03 999  FALSE
#> 13:  FRDA WGSrpt_16 -1.07499268 8.250570e-01 2.859034e-03 999  FALSE
#> 14:  FRDA WGSrpt_17 -0.82994546 7.628872e-01 3.200503e-03 999  FALSE
#> 15:  FRDA WGSrpt_18  0.25955474 4.149475e-01 3.707437e-03 999  FALSE
#> 16:  FRDA WGSrpt_19 -1.00755997 8.095980e-01 2.954567e-03 999  FALSE
#> 17:  FRDA WGSrpt_20 -1.24996212 8.594784e-01 2.615383e-03 999  FALSE
#> 18:  FRDA WGSrpt_21  0.42092705 3.591726e-01 3.609940e-03 999  FALSE
#> 19:    HD WGSrpt_02  0.78651664 2.339432e-01 3.185374e-03 999  FALSE
#> 20:    HD WGSrpt_04 -1.17279494 8.590891e-01 2.618408e-03 999  FALSE
#> 21:    HD WGSrpt_05  0.72372001 2.522382e-01 3.267851e-03 999  FALSE
#> 22:    HD WGSrpt_07 -1.69374058 9.389423e-01 1.802393e-03 999  FALSE
#> 23:    HD WGSrpt_08  0.39030004 3.608964e-01 3.613723e-03 999  FALSE
#> 24:    HD WGSrpt_09  0.42849174 3.488295e-01 3.586176e-03 999  FALSE
#> 25:    HD WGSrpt_10  3.91108897 5.004727e-04 1.682878e-04 999   TRUE
#> 26:    HD WGSrpt_11 -0.60989099 7.144525e-01 3.398816e-03 999  FALSE
#> 27:    HD WGSrpt_12  4.76740655 5.560807e-05 5.610840e-05 999   TRUE
#> 28:    HD WGSrpt_13 -1.23399806 8.708224e-01 2.524133e-03 999  FALSE
#> 29:    HD WGSrpt_14 -0.91600406 8.017016e-01 3.000452e-03 999  FALSE
#> 30:    HD WGSrpt_15  0.95760367 1.901796e-01 2.952910e-03 999  FALSE
#> 31:    HD WGSrpt_16  0.65945956 2.730913e-01 3.352505e-03 999  FALSE
#> 32:    HD WGSrpt_17 -2.02905300 9.656342e-01 1.371784e-03 999  FALSE
#> 33:    HD WGSrpt_18 -0.53305756 6.908747e-01 3.477486e-03 999  FALSE
#> 34:    HD WGSrpt_19  1.38595973 1.017072e-01 2.274349e-03 999  FALSE
#> 35:    HD WGSrpt_20  0.40540384 3.560029e-01 3.602853e-03 999  FALSE
#> 36:    HD WGSrpt_21 -2.55875459 9.893789e-01 7.733569e-04 999  FALSE
#> 37:  SCA2 WGSrpt_02 -1.05875451 8.670967e-01 2.554776e-03 999  FALSE
#> 38:  SCA2 WGSrpt_04 -0.76845003 7.884669e-01 3.073249e-03 999  FALSE
#> 39:  SCA2 WGSrpt_05 -0.30080795 6.205861e-01 3.651311e-03 999  FALSE
#> 40:  SCA2 WGSrpt_07 -2.84199617 9.981093e-01 3.317130e-04 999  FALSE
#> 41:  SCA2 WGSrpt_08  0.46825467 3.037313e-01 3.460267e-03 999  FALSE
#> 42:  SCA2 WGSrpt_09 -2.12973781 9.874326e-01 8.400624e-04 999  FALSE
#> 43:  SCA2 WGSrpt_10  0.74965893 2.093644e-01 3.061356e-03 999  FALSE
#> 44:  SCA2 WGSrpt_11  1.02316649 1.377412e-01 2.593127e-03 999  FALSE
#> 45:  SCA2 WGSrpt_12  0.85434106 1.796697e-01 2.888720e-03 999  FALSE
#> 46:  SCA2 WGSrpt_13 -0.64930020 7.505978e-01 3.255831e-03 999  FALSE
#> 47:  SCA2 WGSrpt_14  0.73265499 2.144803e-01 3.088493e-03 999  FALSE
#> 48:  SCA2 WGSrpt_15 -0.21260246 5.835511e-01 3.709442e-03 999  FALSE
#> 49:  SCA2 WGSrpt_16  1.18776049 1.014291e-01 2.271590e-03 999  FALSE
#> 50:  SCA2 WGSrpt_17 -0.89585180 8.260579e-01 2.852574e-03 999  FALSE
#> 51:  SCA2 WGSrpt_18  8.17346266 5.560807e-05 5.610840e-05 999   TRUE
#> 52:  SCA2 WGSrpt_19  2.39053989 5.950064e-03 5.786767e-04 999  FALSE
#> 53:  SCA2 WGSrpt_20  5.76565693 5.560807e-05 5.610840e-05 999   TRUE
#> 54:  SCA2 WGSrpt_21 -1.52309932 9.457821e-01 1.704716e-03 999  FALSE
#> 55:  SCA6 WGSrpt_02 -3.79937114 9.997776e-01 1.256740e-04 999  FALSE
#> 56:  SCA6 WGSrpt_04 -0.07896286 5.284435e-01 3.756219e-03 999  FALSE
#> 57:  SCA6 WGSrpt_05  8.87692558 5.560807e-05 5.610840e-05 999   TRUE
#> 58:  SCA6 WGSrpt_07  9.00122726 5.560807e-05 5.610840e-05 999   TRUE
#> 59:  SCA6 WGSrpt_08 -0.68122300 7.433687e-01 3.286726e-03 999  FALSE
#> 60:  SCA6 WGSrpt_09 -0.77901515 7.747873e-01 3.143410e-03 999  FALSE
#> 61:  SCA6 WGSrpt_10  0.89980173 1.904020e-01 2.954231e-03 999  FALSE
#> 62:  SCA6 WGSrpt_11 -0.77233141 7.730078e-01 3.152175e-03 999  FALSE
#> 63:  SCA6 WGSrpt_12  0.54601268 2.981149e-01 3.441923e-03 999  FALSE
#> 64:  SCA6 WGSrpt_13 -0.16113762 5.601957e-01 3.734963e-03 999  FALSE
#> 65:  SCA6 WGSrpt_14 -0.78958334 7.772341e-01 3.131225e-03 999  FALSE
#> 66:  SCA6 WGSrpt_15  1.52409423 7.373631e-02 1.966437e-03 999  FALSE
#> 67:  SCA6 WGSrpt_16  0.55114384 2.965579e-01 3.436728e-03 999  FALSE
#> 68:  SCA6 WGSrpt_17  0.88474768 1.936829e-01 2.973532e-03 999  FALSE
#> 69:  SCA6 WGSrpt_18 -0.77629825 7.742868e-01 3.145883e-03 999  FALSE
#> 70:  SCA6 WGSrpt_19 -2.32054394 9.848746e-01 9.200550e-04 999  FALSE
#> 71:  SCA6 WGSrpt_20  1.57876816 6.689651e-02 1.879916e-03 999  FALSE
#> 72:  SCA6 WGSrpt_21 -4.16169382 9.999444e-01 7.970361e-05 999  FALSE
#>     locus    sample        tsum      p.value   p.value.sd   B signif

To obtain only the significant samples, you can either use data.table subsetting:

ps[identity(signif)]
#>    locus    sample      tsum      p.value   p.value.sd   B signif
#> 1:  FRDA WGSrpt_09 12.156872 5.560807e-05 0.0000561084 999   TRUE
#> 2:  FRDA WGSrpt_11 13.608453 5.560807e-05 0.0000561084 999   TRUE
#> 3:    HD WGSrpt_10  3.911089 5.004727e-04 0.0001682878 999   TRUE
#> 4:    HD WGSrpt_12  4.767407 5.560807e-05 0.0000561084 999   TRUE
#> 5:  SCA2 WGSrpt_18  8.173463 5.560807e-05 0.0000561084 999   TRUE
#> 6:  SCA2 WGSrpt_20  5.765657 5.560807e-05 0.0000561084 999   TRUE
#> 7:  SCA6 WGSrpt_05  8.876926 5.560807e-05 0.0000561084 999   TRUE
#> 8:  SCA6 WGSrpt_07  9.001227 5.560807e-05 0.0000561084 999   TRUE

or when retrieving the data.table from p_values():

p_values(tsum, only.signif = TRUE, correction = "samples")
#>    locus    sample      tsum      p.value   p.value.sd   B signif
#> 1:  FRDA WGSrpt_09 12.156872 5.560807e-05 0.0000561084 999   TRUE
#> 2:  FRDA WGSrpt_11 13.608453 5.560807e-05 0.0000561084 999   TRUE
#> 3:    HD WGSrpt_10  3.911089 5.004727e-04 0.0001682878 999   TRUE
#> 4:    HD WGSrpt_12  4.767407 5.560807e-05 0.0000561084 999   TRUE
#> 5:  SCA2 WGSrpt_18  8.173463 5.560807e-05 0.0000561084 999   TRUE
#> 6:  SCA2 WGSrpt_20  5.765657 5.560807e-05 0.0000561084 999   TRUE
#> 7:  SCA6 WGSrpt_05  8.876926 5.560807e-05 0.0000561084 999   TRUE
#> 8:  SCA6 WGSrpt_07  9.001227 5.560807e-05 0.0000561084 999   TRUE

Other useful functions

Subsetting

We may subset the exstra_score object by locus or samples. This is performed by x[loci] or x[loci, samples].

For a single locus:

exstra_wgs_pcr_2["HD"]
#> exstra_score object with 1621 observations of type named ($data),
#>   for 18 samples. ($samples)
#>   Includes associated STR database of 1 locus. ($db)

A single sample:

exstra_wgs_pcr_2[, "WGSrpt_10"]
#> exstra_score object with 1971 observations of type named ($data),
#>   for 1 samples. ($samples)
#>   Includes associated STR database of 21 loci. ($db)

The square brackets also allow filtering loci or samples with data.table syntax, on the x$db and x$samples. For example, we can extract all triplet repeats:

exstra_wgs_pcr_2[unit_length == 3]
#> exstra_score object with 24889 observations of type named ($data),
#>   for 18 samples. ($samples)
#>   Includes associated STR database of 16 loci. ($db)
exstra_wgs_pcr_2[unit_length == 3]$db$locus
#>  [1] "DM1"   "DRPLA" "FRAXA" "FRAXE" "FRDA"  "HD"    "HDL2"  "SBMA" 
#>  [9] "SCA1"  "SCA12" "SCA17" "SCA2"  "SCA3"  "SCA6"  "SCA7"  "SCA8"

or extract all coding repeats:

exstra_wgs_pcr_2[gene_region == "coding"]
#> exstra_score object with 15275 observations of type named ($data),
#>   for 18 samples. ($samples)
#>   Includes associated STR database of 9 loci. ($db)
exstra_wgs_pcr_2[gene_region == "coding"]$db$locus
#> [1] "DRPLA" "HD"    "SBMA"  "SCA1"  "SCA17" "SCA2"  "SCA3"  "SCA6"  "SCA7"

Similarly, this may be performed on sample information:

exstra_wgs_pcr_2[, group == "case"]
#> exstra_score object with 28928 observations of type named ($data),
#>   for 16 samples. ($samples)
#>   Includes associated STR database of 21 loci. ($db)

Combining data

For combining exstra_score objects. This may be useful if the BAM/CRAM files are analysed separately, but then are required to be analysed together in R. In the following dummy example, we split up the str_score variable into two samples, then recombine them.

data_08 <- str_score[, "WGSrpt_08"]
data_13 <- str_score[, "WGSrpt_13"]
( combined_data <- rbind_score_list(list(data_08, data_13)) )
#> exstra_score object with 4082 observations of type named ($data),
#>   for 2 samples. ($samples)
#>   Includes associated STR database of 21 loci. ($db)