This vignette outlines a workflow of detecting nuclear-mitochondrial DNA fusions from Variant Call Format (VCF) using the svaNUMT
package.
This package uses a breakend-centric event notation adopted from the StructuralVariantAnnotation
package. In short, breakends are stored in a GRanges object with strand used to indicate breakpoint orientation. where breakpoints are represented using a partner
field containing the name of the breakend at the other side of the breakend. This notation was chosen as it simplifies the annotations of RTs which are detected at breakend-level.
VCF data is parsed into a VCF
object using the readVCF
function from the Bioconductor package VariantAnnotation
. Simple filters could be applied to a VCF
object to remove unwanted calls. The VCF
object is then converted to a GRanges
object with breakend-centric notations using StructuralVariantAnnotation
. More information about VCF
objects and breakend-centric GRanges object can be found by consulting the vignettes in the corresponding packages with browseVignettes("VariantAnnotation")
and browseVignettes("StructuralVariantAnnotation")
.
library(StructuralVariantAnnotation)
library(VariantAnnotation)
library(svaNUMT)
vcf <- readVcf(system.file("extdata", "chr1_numt_pe_HS25.sv.vcf", package = "svaNUMT"))
gr <- breakpointRanges(vcf)
Note that StructuralVariantAnnotation
requires the GRanges
object to be composed entirely of valid breakpoints. Please consult the vignette of the StructuralVariantAnnotation
package for ensuring breakpoint consistency.
Function svaNUMT
searches for NUMT events by identifying breakends supporting the fusion of nuclear chromosome and mitochondrial genome. svaNUMT
returns identified breakends supporting candidate NUMTs in 2 lists of list of GRanges, grouped by chromosome and insertion sites.
NUMT <- numtDetect(gr, max_ins_dist = 20)
The breakends supporting the insertion sites and the MT sequence are arranged by the order of events. Below is an example of a detected NUMT event, where MT sequence MT:15737-15836
followed by polyadenylation is inserted between chr1:1688363-1688364
.
GRangesList(NU=NUMT$NU$`1`[[1]], MT=NUMT$MT$`1`[[1]])
#> GRangesList object of length 2:
#> $NU
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss1fb_4o 1 1688363 + | NA C
#> gridss1bf_1o 1 1688364 - | NA C
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss1fb_4o C[MT:15737[ 3928.49 PASS gridss1fb_4o
#> gridss1bf_1o ]MT:15836]AAAAAAAAAA.. 3581.13 PASS gridss1bf_1o
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss1fb_4o gridss1fb_4h BND NA 0
#> gridss1bf_1o gridss1bf_1h BND NA AAAAAAAAAAAAA 13
#> HOMLEN
#> <numeric>
#> gridss1fb_4o 0
#> gridss1bf_1o 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
#>
#> $MT
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss1fb_4h MT 15737 - | NA G
#> gridss1bf_1h MT 15836 + | NA A
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss1fb_4h ]1:1688363]G 3928.49 PASS gridss1fb_4h
#> gridss1bf_1h AAAAAAAAAAAAAA[1:168.. 3581.13 PASS gridss1bf_1h
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss1fb_4h gridss1fb_4o BND NA 0
#> gridss1bf_1h gridss1bf_1o BND NA AAAAAAAAAAAAA 13
#> HOMLEN
#> <numeric>
#> gridss1fb_4h 0
#> gridss1bf_1h 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
Below is an example to subset the detected NUMTs by a genomic region given seqnames
, start
, and end
. For region chr1:1000000-3000000
, there are 3 NUMTs detected.
seqnames = 1
start = 1000000
end = 3000000
i <- sapply(NUMT$NU[[seqnames]], function(x)
sum(countOverlaps(x, GRanges(seqnames = seqnames, IRanges(start, end))))>0)
list(NU=NUMT$NU[[seqnames]][i], MT=NUMT$MT[[seqnames]][i])
#> $NU
#> $NU[[1]]
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss1fb_4o 1 1688363 + | NA C
#> gridss1bf_1o 1 1688364 - | NA C
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss1fb_4o C[MT:15737[ 3928.49 PASS gridss1fb_4o
#> gridss1bf_1o ]MT:15836]AAAAAAAAAA.. 3581.13 PASS gridss1bf_1o
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss1fb_4o gridss1fb_4h BND NA 0
#> gridss1bf_1o gridss1bf_1h BND NA AAAAAAAAAAAAA 13
#> HOMLEN
#> <numeric>
#> gridss1fb_4o 0
#> gridss1bf_1o 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
#>
#> $NU[[2]]
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss1fb_5o 1 1791082-1791083 + | NA G
#> gridss1bf_2o 1 1791084 - | NA A
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss1fb_5o G[MT:2592[ 1929.85 PASS gridss1fb_5o
#> gridss1bf_2o ]MT:3592]AAAAAAAAAAAA 2894.91 PASS gridss1bf_2o
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss1fb_5o gridss1fb_5h BND NA 0
#> gridss1bf_2o gridss1bf_2h BND NA AAAAAAAAAAA 11
#> HOMLEN
#> <numeric>
#> gridss1fb_5o 1
#> gridss1bf_2o 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
#>
#> $NU[[3]]
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss2fb_3o 1 2869079 + | NA G
#> gridss2bf_2o 1 2869080 - | NA A
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss2fb_3o G[MT:2786[ 2472.12 PASS gridss2fb_3o
#> gridss2bf_2o ]MT:2985]AAAAAAAAAAA.. 2456.81 PASS gridss2bf_2o
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss2fb_3o gridss2fb_3h BND NA 0
#> gridss2bf_2o gridss2bf_2h BND NA AAAAAAAAAAAAAAA 15
#> HOMLEN
#> <numeric>
#> gridss2fb_3o 0
#> gridss2bf_2o 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
#>
#>
#> $MT
#> $MT[[1]]
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss1fb_4h MT 15737 - | NA G
#> gridss1bf_1h MT 15836 + | NA A
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss1fb_4h ]1:1688363]G 3928.49 PASS gridss1fb_4h
#> gridss1bf_1h AAAAAAAAAAAAAA[1:168.. 3581.13 PASS gridss1bf_1h
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss1fb_4h gridss1fb_4o BND NA 0
#> gridss1bf_1h gridss1bf_1o BND NA AAAAAAAAAAAAA 13
#> HOMLEN
#> <numeric>
#> gridss1fb_4h 0
#> gridss1bf_1h 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
#>
#> $MT[[2]]
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss1fb_5h MT 2592-2593 - | NA G
#> gridss1bf_2h MT 3592 + | NA G
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss1fb_5h ]1:1791082]G 1929.85 PASS gridss1fb_5h
#> gridss1bf_2h GAAAAAAAAAAA[1:17910.. 2894.91 PASS gridss1bf_2h
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss1fb_5h gridss1fb_5o BND NA 0
#> gridss1bf_2h gridss1bf_2o BND NA AAAAAAAAAAA 11
#> HOMLEN
#> <numeric>
#> gridss1fb_5h 1
#> gridss1bf_2h 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
#>
#> $MT[[3]]
#> GRanges object with 2 ranges and 12 metadata columns:
#> seqnames ranges strand | paramRangeID REF
#> <Rle> <IRanges> <Rle> | <factor> <character>
#> gridss2fb_3h MT 2786 - | NA T
#> gridss2bf_2h MT 2985 + | NA C
#> ALT QUAL FILTER sourceId
#> <character> <numeric> <character> <character>
#> gridss2fb_3h ]1:2869079]T 2472.12 PASS gridss2fb_3h
#> gridss2bf_2h CAAAAAAAAAAAAAAA[1:2.. 2456.81 PASS gridss2bf_2h
#> partner svtype svLen insSeq insLen
#> <character> <character> <numeric> <character> <numeric>
#> gridss2fb_3h gridss2fb_3o BND NA 0
#> gridss2bf_2h gridss2bf_2o BND NA AAAAAAAAAAAAAAA 15
#> HOMLEN
#> <numeric>
#> gridss2fb_3h 0
#> gridss2bf_2h 0
#> -------
#> seqinfo: 86 sequences from an unspecified genome
One way of visualising paired breakpoints is by circos plots. Here we use the package circlize
to demonstrate breakpoint visualisation. The bedpe2circos
function takes BEDPE-formatted dataframes (see breakpointgr2bedpe()
) and plotting parameters for the circos.initializeWithIdeogram()
and circos.genomicLink()
functions from circlize
.
To generate a simple circos plot of one candidate NUMT event:
library(circlize)
numt_chr_prefix <- c(NUMT$NU$`1`[[2]], NUMT$MT$`1`[[2]])
GenomeInfoDb::seqlevelsStyle(numt_chr_prefix) <- "UCSC"
pairs <- breakpointgr2pairs(numt_chr_prefix)
pairs
To see supporting breakpoints clearly, we generate the circos plot according to the loci of event.
circos.initializeWithIdeogram(
data.frame(V1=c("chr1", "chrM"),
V2=c(1791073,1),
V3=c(1791093,16571),
V4=c("p15.4",NA),
V5=c("gpos50",NA)), sector.width = c(0.2, 0.8))
#circos.initializeWithIdeogram()
circos.genomicLink(as.data.frame(S4Vectors::first(pairs)),
as.data.frame(S4Vectors::second(pairs)))
circos.clear()
sessionInfo()
#> R version 4.1.1 (2021-08-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.14-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] circlize_0.4.13 svaNUMT_1.0.0
#> [3] StructuralVariantAnnotation_1.10.0 VariantAnnotation_1.40.0
#> [5] Rsamtools_2.10.0 Biostrings_2.62.0
#> [7] XVector_0.34.0 SummarizedExperiment_1.24.0
#> [9] Biobase_2.54.0 MatrixGenerics_1.6.0
#> [11] matrixStats_0.61.0 rtracklayer_1.54.0
#> [13] GenomicRanges_1.46.0 GenomeInfoDb_1.30.0
#> [15] IRanges_2.28.0 S4Vectors_0.32.0
#> [17] BiocGenerics_0.40.0
#>
#> loaded via a namespace (and not attached):
#> [1] httr_1.4.2 sass_0.4.0 bit64_4.0.5
#> [4] jsonlite_1.7.2 bslib_0.3.1 assertthat_0.2.1
#> [7] highr_0.9 BiocFileCache_2.2.0 blob_1.2.2
#> [10] BSgenome_1.62.0 GenomeInfoDbData_1.2.7 yaml_2.2.1
#> [13] progress_1.2.2 pillar_1.6.4 RSQLite_2.2.8
#> [16] lattice_0.20-45 glue_1.4.2 digest_0.6.28
#> [19] colorspace_2.0-2 htmltools_0.5.2 Matrix_1.3-4
#> [22] XML_3.99-0.8 pkgconfig_2.0.3 biomaRt_2.50.0
#> [25] zlibbioc_1.40.0 purrr_0.3.4 BiocParallel_1.28.0
#> [28] tibble_3.1.5 KEGGREST_1.34.0 generics_0.1.1
#> [31] ellipsis_0.3.2 cachem_1.0.6 GenomicFeatures_1.46.0
#> [34] magrittr_2.0.1 crayon_1.4.1 memoise_2.0.0
#> [37] evaluate_0.14 fansi_0.5.0 xml2_1.3.2
#> [40] tools_4.1.1 prettyunits_1.1.1 GlobalOptions_0.1.2
#> [43] hms_1.1.1 BiocIO_1.4.0 lifecycle_1.0.1
#> [46] stringr_1.4.0 DelayedArray_0.20.0 AnnotationDbi_1.56.0
#> [49] compiler_4.1.1 jquerylib_0.1.4 rlang_0.4.12
#> [52] grid_4.1.1 RCurl_1.98-1.5 rappdirs_0.3.3
#> [55] rjson_0.2.20 bitops_1.0-7 rmarkdown_2.11
#> [58] restfulr_0.0.13 curl_4.3.2 DBI_1.1.1
#> [61] R6_2.5.1 GenomicAlignments_1.30.0 knitr_1.36
#> [64] dplyr_1.0.7 fastmap_1.1.0 bit_4.0.4
#> [67] utf8_1.2.2 filelock_1.0.2 shape_1.4.6
#> [70] stringi_1.7.5 parallel_4.1.1 Rcpp_1.0.7
#> [73] vctrs_0.3.8 png_0.1-7 dbplyr_2.1.1
#> [76] tidyselect_1.1.1 xfun_0.27