CNVfilteR 1.0.4
Many tools for germline copy number variant (CNV) detection from NGS data have been developed. Usually, these tools were designed for different input data like WGS, WES or panel data, and their performance may depend on the CNV size. Available benchmarks show that all these tools obtain false positives, sometimes reaching a very high number of them.
With the aim of reducing the number of false positives, CNVfilteR identifies those germline CNVs that can be discarded. This task is performed by using the germline single nucleotide variant (SNV) calls that are usually obtained in common NGS pipelines. As VCF field interpretation is key when working with these files, CNVfilteR specifically supports VCFs produced by VarScan2, Strelka/Strelka2, freeBayes, HaplotypeCaller (GATK), and UnifiedGenotyper (GATK). Additionally, results can be plotted using the functions provided by the R/Bioconductor packages karyoploteR and CopyNumberPlots.
CNVfilteR is a Bioconductor package and to install it we have to use BiocManager.
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("CNVfilteR")
We can also install the package from github to get the latest devel version, but beware that it might be incompatible with the release version of Bioconductor.
BiocManager::install("jpuntomarcos/CNVfilteR")
Below we show a full example that covers the usual steps: CNVs data loading, variants loading, identifying false postives and plotting the results.
First, we can load some CNV tool results:
library(CNVfilteR)
cnvs.file <- system.file("extdata", "DECoN.CNVcalls.csv", package = "CNVfilteR", mustWork = TRUE)
cnvs.gr <- loadCNVcalls(cnvs.file = cnvs.file, chr.column = "Chromosome", start.column = "Start", end.column = "End", cnv.column = "CNV.type", sample.column = "Sample", genome = "hg19")
Then, we load the SNVs stored in a couple of VCF files.
vcf.files <- c(system.file("extdata", "variants.sample1.vcf.gz", package = "CNVfilteR", mustWork = TRUE),
system.file("extdata", "variants.sample2.vcf.gz", package = "CNVfilteR", mustWork = TRUE))
vcfs <- loadVCFs(vcf.files, cnvs.gr = cnvs.gr)
## Scanning file /tmp/RtmpqbMEKO/Rinst32983005b074/CNVfilteR/extdata/variants.sample1.vcf.gz...
## VarScan2 was found as source in the VCF metadata, RD will be used as ref allele depth field, AD will be used as alt allele depth field.
## Scanning file /tmp/RtmpuwuY98/variants.sample2.vcf.gz...
## VarScan2 was found as source in the VCF metadata, RD will be used as ref allele depth field, AD will be used as alt allele depth field.
We observe that the function recognized VarScan2 as the source, so fields
were selected and allele frequency consequently. Now we can call
filterCNVs()
to identify those CNVs that can be discarded.
results <- filterCNVs(cnvs.gr, vcfs)
names(results)
## [1] "cnvs" "variantsForEachCNV" "filterParameters"
And we can check those CNVs that can be filtered out:
results$cnvs[results$cnvs$filter == TRUE]
## GRanges object with 5 ranges and 9 metadata columns:
## seqnames ranges strand | cnv sample filter
## <Rle> <IRanges> <Rle> | <character> <character> <character>
## 3 chr2 48025751-48028294 * | duplication sample1 TRUE
## 16 chr17 41243453-41247939 * | duplication sample1 TRUE
## 17 chr17 41251793-41256973 * | duplication sample1 TRUE
## 19 chr13 32900637-32929425 * | deletion sample2 TRUE
## 20 chr17 59870959-59938900 * | deletion sample2 TRUE
## n.total.variants n.hm.variants n.ht.discard.CNV n.ht.confirm.CNV
## <character> <character> <character> <character>
## 3 5 0 3 2
## 16 2 0 2 0
## 17 1 0 1 0
## 19 10 4 6
## 20 1 0 1
## score cnv.id
## <character> <character>
## 3 2.539241834223 3
## 16 1.99927691434002 16
## 17 0.994339688718521 17
## 19 19
## 20 20
## -------
## seqinfo: 24 sequences from hg19 genome
As an example, we can observe that the CNV with cnv.id
=3 has 4 variants
matching it: 2 in favor of discarding it, two against discarding it. If we want
to know more about the variants matching a certain CNV we can do:
results$variantsForEachCNV[["3"]]
## seqnames start end width strand ref alt ref.support alt.support
## 1 chr2 48026019 48026019 1 * G C 516 521
## 2 chr2 48027019 48027019 1 * G C 1528 964
## 3 chr2 48027182 48027182 1 * G A 1506 971
## 4 chr2 48027434 48027434 1 * A T 1593 1462
## 5 chr2 48027763 48027763 1 * G A 900 863
## alt.freq total.depth indel type score
## 1 50.2411 1037 FALSE ht 0.9999596
## 2 38.6838 2492 FALSE ht -0.3169985
## 3 39.2006 2477 FALSE ht -0.1417051
## 4 47.8560 3055 FALSE ht 0.9981889
## 5 48.9507 1763 FALSE ht 0.9997969
Two variants are close to the default expected heterozygous frequency,
0.5, so they obtain a positive score. The other two variants are not
so clearly close to the default expected duplication value, 0.33, so
they obtain a low negative score. All these default values and others can
be modified in the filterCNVs()
function.
Finally, we may be interested in plotting the results. For example, we can
plot easily the duplication CNV with cnv.id
=3 and all the variants matching
it.
plotVariantsForCNV(results, "3")
We can do the same to plot the deletion CNV with cnv.id
=19, where all
variants discard the CNV except one homozygous variant that does not give us any
information for supporting or discarding the CNV:
plotVariantsForCNV(results, "19")
On the opposite side, we can observe those CNVs that cannot be discarded:
results$cnvs[results$cnvs$filter != TRUE]
## GRanges object with 15 ranges and 9 metadata columns:
## seqnames ranges strand | cnv sample filter
## <Rle> <IRanges> <Rle> | <character> <character> <character>
## 1 chr2 47641409-47641557 * | duplication sample1
## 2 chr2 47698105-47698201 * | duplication sample1
## 4 chr3 10091059-10091189 * | duplication sample1
## 5 chr3 14219967-14220068 * | duplication sample1
## 6 chr3 37042447-37042544 * | duplication sample1
## .. ... ... ... . ... ... ...
## 12 chr13 32906410-32907524 * | duplication sample1
## 13 chr13 32945094-32945237 * | duplication sample1
## 14 chr13 32953455-32969070 * | duplication sample1
## 15 chr17 41209070-41215390 * | duplication sample1
## 18 chr17 41267744-41276113 * | duplication sample1
## n.total.variants n.hm.variants n.ht.discard.CNV n.ht.confirm.CNV
## <character> <character> <character> <character>
## 1 0 0
## 2 0 0
## 4 0 0
## 5 0 0
## 6 1 0 1 0
## .. ... ... ... ...
## 12 0 0
## 13 0 0
## 14 1 0 0 1
## 15 0 0
## 18 0 0
## score cnv.id
## <character> <character>
## 1 1
## 2 2
## 4 4
## 5 5
## 6 0.372384685846936 6
## .. ... ...
## 12 12
## 13 13
## 14 -0.996334798387489 14
## 15 15
## 18 18
## -------
## seqinfo: 24 sequences from hg19 genome
For example, the CNV with cnv.id
=14 has one variant matching it. If we get
the variant info, we see that the variant has an allele frequency very close to
the default expected duplication value, 0.66.
results$variantsForEachCNV[["14"]]
## seqnames start end width strand ref alt ref.support alt.support
## 1 chr13 32968607 32968607 1 * A G 14 25
## alt.freq total.depth indel type score
## 1 64.1026 39 FALSE ht -0.9963348
So, no evidence was found for discarding the CNV. We can also plot the CNV and the variant:
plotVariantsForCNV(results, "3")
CNVfilteR functions expect germline CNVs calls to be a
GRanges
object with a few specificic metadata columns:
You can create this object yourself, but
CNVfilter provides the proper function to do this,
loadCNVcalls()
. This function can interpret any tsv o csv file by indicating
which columns store the information. For example, in the following code, the
chr.column
column is stored in the “Chromosome” column of the cnvs.file
.
cnvs.file <- system.file("extdata", "DECoN.CNVcalls.csv", package = "CNVfilteR", mustWork = TRUE)
cnvs.gr <- loadCNVcalls(cnvs.file = cnvs.file, chr.column = "Chromosome", start.column = "Start", end.column = "End", cnv.column = "CNV.type", sample.column = "Sample", genome = "hg19")
loadCNVcalls()
can interpret different types of CNVs. Among other options,
separator can be selected using the sep
parameter (defaults to \t),
and first lines can be skipped using the skip
parameter (defaults to 0). Also,
the value used in cnv.column
to store the CNV type can be modified
using the deletion
and duplication
parameters (defaults to “deletion” and
“duplication”, respectively). If, for example, our cnv.column
uses “CN1” and
“CN3” for deletion and duplication (respectively), we should indicate
the function to use these values:
cnvs.gr.2 <- loadCNVcalls(cnvs.file = cnvs.file.2, deletion = "CN1", duplication = "CN3", chr.column = "Chromosome", start.column = "Start", end.column = "End", cnv.column = "CNV.type", sample.column = "Sample")
Some CNV tools generate results where the CNV location is stored in a single
column with the format chr:start-end (i.e. 1:538001-540000). In this
case, we can call loadCNVcalls()
using the coord.column
instead of the
chr.column
, start.column
and end.column
columns.
Common NGS pipelines produce germline variant calling (SNVs or INDELs)
in a VCF file. However, VCF interpretation is
challenging due to the flexibility provided by the VCF format definition.
It is not straightforward to interpret correctly the fields in the VCF file
and compute the allele frequency. loadVCFs()
interprets automatically
VCFs produced by VarScan2, Strelka/Strelka2, freeBayes, HaplotypeCaller (GATK),
and UnifiedGenotyper (GATK).
In the following example the function recognizes VarScan2 as the source.
vcf.files <- c(system.file("extdata", "variants.sample1.vcf.gz", package = "CNVfilteR", mustWork = TRUE),
system.file("extdata", "variants.sample2.vcf.gz", package = "CNVfilteR", mustWork = TRUE))
vcfs <- loadVCFs(vcf.files, cnvs.gr = cnvs.gr)
## Scanning file /tmp/RtmpqbMEKO/Rinst32983005b074/CNVfilteR/extdata/variants.sample1.vcf.gz...
## VarScan2 was found as source in the VCF metadata, RD will be used as ref allele depth field, AD will be used as alt allele depth field.
## Scanning file /tmp/RtmpuwuY98/variants.sample2.vcf.gz...
## VarScan2 was found as source in the VCF metadata, RD will be used as ref allele depth field, AD will be used as alt allele depth field.
We can also load the VCF file spicifying how to interpret it, which can be useful if the VCF was generated by a caller not supported by CNVfilteR. For example we can specify the ref/alt fields:
vcfs <- loadVCFs(vcf.files, cnvs.gr = cnvs.gr, vcf.source = "MyCaller", ref.support.field = "RD", alt.support.field = "AD")
## Scanning file /tmp/RtmpqbMEKO/Rinst32983005b074/CNVfilteR/extdata/variants.sample1.vcf.gz...
## VCF source MyCaller is not supported, but ref.support.field/alt.support.field were provided.
## Scanning file /tmp/RtmpuwuY98/variants.sample2.vcf.gz...
## VCF source MyCaller is not supported, but ref.support.field/alt.support.field were provided.
Alternatively, we can set the list.support.field
parameter so that field
will be loaded assuming that it is a list in this order: reference allele,
alternative allele. As an example:
vcf.file3 <- c(system.file("extdata", "variants.sample3.vcf", package = "CNVfilteR", mustWork = TRUE))
vcfs3 <- loadVCFs(vcf.file3, cnvs.gr = cnvs.gr, vcf.source = "MyCaller", list.support.field = "AD")
## Scanning file /tmp/RtmpuwuY98/variants.sample3.vcf.gz...
## VCF source MyCaller is not supported, but list.support.field was provided.
As VCF will be used for determining if a CNV can be filtered out, it is
recommended to use a VCF file free of artifacts to improve the accuracy
when identifying those CNVs. loadVCFs
will load only
variants with a minimum total depth
greater than min.total.depth
(defaults to 30). Additionally, we can exclude
those regions that have already known artifacts with the parameter
regions.to.exclude
. In this example, we exclude PMS2, PRSS1, and
FANCD2 genes because they are pseudogenes with alignments artifacts:
regions.to.exclude <- GRanges(seqnames = c("chr3","chr7", "chr7"), ranges = IRanges(c(10068098, 6012870, 142457319), c(10143614, 6048756, 142460923)))
vcfs4 <- loadVCFs(vcf.files, cnvs.gr = cnvs.gr, regions.to.exclude = regions.to.exclude)
## Scanning file /tmp/RtmpqbMEKO/Rinst32983005b074/CNVfilteR/extdata/variants.sample1.vcf.gz...
## VarScan2 was found as source in the VCF metadata, RD will be used as ref allele depth field, AD will be used as alt allele depth field.
## Scanning file /tmp/RtmpuwuY98/variants.sample2.vcf.gz...
## VarScan2 was found as source in the VCF metadata, RD will be used as ref allele depth field, AD will be used as alt allele depth field.
Also, the parameter exclude.indels
indicates whether to exclude INDELs when
loading the variants. TRUE is the default and recommended value given
that INDELs allele frequency varies differently than SNVs. Including
INDELs may allow the algorithm to identify more CNVs to discard with a greater
risk of identifying them wrongly.
The function loadVCFs()
also adapts to different needs. If sample.names
parameter is
not provided, the sample names included in the VCF itself will be used.
Both single-sample and multi-sample VCFs are accepted, but when
multi-sample VCFs are used, sample.names
parameter must be NULL.
If VCF is not compressed with bgzip, the function compresses it and generates the .gz file. If .tbi file does not exist for a given VCF file, the function also generates it. All files are generated in a temporary folder.
See loadVCFs()
documentation to see other parameters info.
The task of identifying false positives is performed by the filterCNVs()
function. It checks all the variants (SNVs and optionally INDELs) matching
each CNV present in cnvs.gr
to identify those CNVs that can be filtered. It
returns an S3 object with 3 elements: cnvs
, variantsForEachCNV
and
filterParameters
:
results <- filterCNVs(cnvs.gr, vcfs)
tail(results$cnvs)
## GRanges object with 6 ranges and 9 metadata columns:
## seqnames ranges strand | cnv sample filter
## <Rle> <IRanges> <Rle> | <character> <character> <character>
## 15 chr17 41209070-41215390 * | duplication sample1
## 16 chr17 41243453-41247939 * | duplication sample1 TRUE
## 17 chr17 41251793-41256973 * | duplication sample1 TRUE
## 18 chr17 41267744-41276113 * | duplication sample1
## 19 chr13 32900637-32929425 * | deletion sample2 TRUE
## 20 chr17 59870959-59938900 * | deletion sample2 TRUE
## n.total.variants n.hm.variants n.ht.discard.CNV n.ht.confirm.CNV
## <character> <character> <character> <character>
## 15 0 0
## 16 2 0 2 0
## 17 1 0 1 0
## 18 0 0
## 19 10 4 6
## 20 1 0 1
## score cnv.id
## <character> <character>
## 15 15
## 16 1.99927691434002 16
## 17 0.994339688718521 17
## 18 18
## 19 19
## 20 20
## -------
## seqinfo: 24 sequences from hg19 genome
Observe that those CNVs that can be filtered have the value TRUE in the
column filter
. CNVfilteR employs two
different strategies for identifying those CNVs:
ht.deletions.threshold
% of heterozygous variants matching
the CNV. Default ht.deletions.threshold
value is 15, so 15% is required.score
is >=
dup.threshold.score
after computing all
heterozygous variants matching that CNV. Default dup.threshold.score
value
is 0.5. How the score is computed for each variant is explained in
the next section.The scoring model for determining whether a certain duplication CNV can be discarded is based on the allele frequency for each heterozygous variant matching that CNV:
expected.ht.mean
). So, a variant with an allele frequency close to 50%
gives us evidence of the non-existence of a duplication CNV, so the CNV could
be discarded.expected.dup.ht.mean1
) when the variant is not in the same allele
than the duplication CNV, and 66.6% (expected.dup.ht.mean2
) when the variant
is in the same allele than
the duplication CNV call. So, a variant with an allele frequency close to
33.3% or 66.6% gives us evidence of the existence of duplication CNV.Although we can expect that most of the variants are close to the expected means
(33.3%, 50%, and 66.6%), many of them can be far from any expected mean. The
scoring model implemented in the filterCNVs()
function measures
the evidence - for discarding a certain CNV - using
a scoring model. The scoring model is based on the fuzzy logic, where elements
can have any value between 1 (True) and 0 (False). Following this idea,
each variant will be scored with a value between 0 and 1 depending on
how close is the allele frequency to the nearest expected mean.
The total score
is computed among all the variants matching the CNV. If the
score
is greater than the dup.threshold.score
, the CNV can be discarded.
A common way of applying the fuzzy logic is using the sigmoid function. CNVfilteR uses the sigmoid function implemented in the pracma package, which is defined as \[ \begin{aligned} y = 1/(1 + e^{-c1(x−c2)}) \end{aligned} \]
The scoring model is built on 6 sigmoids defined on 6 different intervals. The
c1 parameter is 2 by default (sigmoid.c1
), and the c2 parameter
is defined for the 6 sigmoids (sigmoid.c2.vector
).
expected.dup.ht.mean1
], c2=28expected.dup.ht.mean1
, sigmoid.int1
], c2=38.3sigmoid.int1
, expected.ht.mean
], c2=44.7expected.ht.mean
, sigmoid.int2
], c2=55.3sigmoid.int2
, expected.dup.ht.mean2
], c`=61.3expected.dup.ht.mean2
, 80], c2=71.3Where sigmoid.int1
is the mean between expected.dup.ht.mean1
and
expected.ht.mean
, and sigmoid.int2
is the mean between
expected.dup.ht.mean2
and expected.ht.mean
.
The scoring model can be plotted using the plotScoringModel()
function.
p <- results$filterParameters
plotScoringModel(expected.ht.mean = p$expected.ht.mean,
expected.dup.ht.mean1 = p$expected.dup.ht.mean1,
expected.dup.ht.mean2 = p$expected.dup.ht.mean2,
sigmoid.c1 = p$sigmoid.c1,
sigmoid.c2.vector = p$sigmoid.c2.vector)
And the scoring model can be modified when calling the filterCNVs()
function.
Let’s see how the model changes when we modify the parameter sigmoid.c1
(1 instead of 2):
plotScoringModel(expected.ht.mean = p$expected.ht.mean,
expected.dup.ht.mean1 = p$expected.dup.ht.mean1,
expected.dup.ht.mean2 = p$expected.dup.ht.mean2,
sigmoid.c1 = 1,
sigmoid.c2.vector = p$sigmoid.c2.vector)
We can also modify the sigmoid.c2.vector
parameter for each sigmoid function. For example, to
make the central sigmoids narrower:
plotScoringModel(expected.ht.mean = p$expected.ht.mean,
expected.dup.ht.mean1 = p$expected.dup.ht.mean1,
expected.dup.ht.mean2 = p$expected.dup.ht.mean2,
sigmoid.c1 = p$sigmoid.c1,
sigmoid.c2.vector = c(28, 38.3, 46.7, 53.3, 61.3, 71.3))
We can plot easily a certain CNV and its matching variants. For example, the
duplication CNV with cnv.id
=17 can be plotted as follows:
plotVariantsForCNV(results, "17")
Some parameters can be customized, like points.cex
and points.pch
. It is
also possible to plot all CNVs in a global schema where all the chromosomes are
plotted:
cnvs.file <- system.file("extdata", "DECoN.CNVcalls.2.csv",
package = "CNVfilteR", mustWork = TRUE)
cnvs.gr.2 <- loadCNVcalls(cnvs.file = cnvs.file, chr.column = "Chromosome",
start.column = "Start", end.column = "End",
cnv.column = "CNV.type", sample.column = "Sample",
genome = "hg19")
plotAllCNVs(cnvs.gr.2)
Note that if a CNV is too small, it will not be visible when calling
plotAllCNVs()
.
sessionInfo()
## R version 3.6.2 (2019-12-12)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.10-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.10-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] stats4 parallel stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] BSgenome.Hsapiens.UCSC.hg19.masked_1.3.99
## [2] BSgenome.Hsapiens.UCSC.hg19_1.4.0
## [3] BSgenome_1.54.0
## [4] rtracklayer_1.46.0
## [5] Biostrings_2.54.0
## [6] XVector_0.26.0
## [7] GenomicRanges_1.38.0
## [8] GenomeInfoDb_1.22.0
## [9] IRanges_2.20.2
## [10] S4Vectors_0.24.3
## [11] BiocGenerics_0.32.0
## [12] CNVfilteR_1.0.4
## [13] knitr_1.27.2
## [14] BiocStyle_2.14.4
##
## loaded via a namespace (and not attached):
## [1] colorspace_1.4-1 biovizBase_1.34.1
## [3] htmlTable_1.13.3 base64enc_0.1-3
## [5] dichromat_2.0-0 rstudioapi_0.10
## [7] bit64_0.9-7 AnnotationDbi_1.48.0
## [9] splines_3.6.2 Formula_1.2-3
## [11] Rsamtools_2.2.1 cluster_2.1.0
## [13] dbplyr_1.4.2 png_0.1-7
## [15] BiocManager_1.30.10 compiler_3.6.2
## [17] httr_1.4.1 backports_1.1.5
## [19] assertthat_0.2.1 Matrix_1.2-18
## [21] lazyeval_0.2.2 exomeCopy_1.32.0
## [23] acepack_1.4.1 htmltools_0.4.0
## [25] prettyunits_1.1.1 tools_3.6.2
## [27] gtable_0.3.0 glue_1.3.1
## [29] GenomeInfoDbData_1.2.2 dplyr_0.8.3
## [31] rappdirs_0.3.1 Rcpp_1.0.3
## [33] Biobase_2.46.0 vctrs_0.2.2
## [35] xfun_0.12 stringr_1.4.0
## [37] lifecycle_0.1.0 ensembldb_2.10.2
## [39] XML_3.99-0.3 zlibbioc_1.32.0
## [41] scales_1.1.0 VariantAnnotation_1.32.0
## [43] karyoploteR_1.12.4 hms_0.5.3
## [45] ProtGenerics_1.18.0 SummarizedExperiment_1.16.1
## [47] rhdf5_2.30.1 AnnotationFilter_1.10.0
## [49] RColorBrewer_1.1-2 yaml_2.2.0
## [51] curl_4.3 memoise_1.1.0
## [53] gridExtra_2.3 ggplot2_3.2.1
## [55] biomaRt_2.42.0 rpart_4.1-15
## [57] latticeExtra_0.6-29 stringi_1.4.5
## [59] RSQLite_2.2.0 CopyNumberPlots_1.2.3
## [61] cn.mops_1.32.0 checkmate_1.9.4
## [63] GenomicFeatures_1.38.1 BiocParallel_1.20.1
## [65] rlang_0.4.4 pkgconfig_2.0.3
## [67] matrixStats_0.55.0 bitops_1.0-6
## [69] pracma_2.2.9 evaluate_0.14
## [71] lattice_0.20-38 purrr_0.3.3
## [73] Rhdf5lib_1.8.0 GenomicAlignments_1.22.1
## [75] htmlwidgets_1.5.1 bit_1.1-15.1
## [77] tidyselect_1.0.0 magrittr_1.5
## [79] bookdown_0.17 R6_2.4.1
## [81] magick_2.3 Hmisc_4.3-0
## [83] DelayedArray_0.12.2 DBI_1.1.0
## [85] pillar_1.4.3 foreign_0.8-75
## [87] survival_3.1-8 RCurl_1.98-1.1
## [89] nnet_7.3-12 tibble_2.1.3
## [91] crayon_1.3.4 BiocFileCache_1.10.2
## [93] rmarkdown_2.1 bamsignals_1.18.0
## [95] jpeg_0.1-8.1 progress_1.2.2
## [97] grid_3.6.2 data.table_1.12.8
## [99] blob_1.2.1 digest_0.6.23
## [101] regioneR_1.18.1 bezier_1.1.2
## [103] openssl_1.4.1 munsell_0.5.0
## [105] askpass_1.1