File: read_vcfs_as_granges.R

package info (click to toggle)
r-bioc-mutationalpatterns 3.16.0%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 5,360 kB
  • sloc: sh: 8; makefile: 2
file content (323 lines) | stat: -rw-r--r-- 12,478 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
#' Read VCF files into a GRangesList
#'
#' This function reads Variant Call Format (VCF) files into a GRanges object
#' and combines them in a GRangesList.  In addition to loading the files, this
#' function applies the same seqlevel style to the GRanges objects as the
#' reference genome passed in the 'genome' parameter.
#' By default only reads in snv variants.
#'
#' @param vcf_files Character vector of VCF file names
#' @param sample_names Character vector of sample names
#' @param genome BSgenome reference genome object
#' @param group Selector for a seqlevel group.  All seqlevels outside
#'              of this group will be removed.  Possible values:
#'              * 'all' for all chromosomes;
#'              * 'auto' for autosomal chromosomes;
#'              * 'sex' for sex chromosomes;
#'              * 'auto+sex' for autosomal + sex chromosomes (default);
#'              * 'circular' for circular chromosomes;
#'              * 'none' for no filtering, which results in keeping all
#'                seqlevels from the VCF file.
#' @param type The mutation type that will be loaded. All other variants will be filtered out.
#'              Possible values:
#'              * 'snv' (default)
#'              * 'indel'
#'              * 'dbs'
#'              * 'mbs'
#'              * 'all'
#'    When you use 'all', no filtering or merging of variants is performed.
#' @param change_seqnames Boolean. Whether to change the seqlevelsStyle of the
#'   vcf to that of the BSgenome object. (default = TRUE)
#' @param predefined_dbs_mbs Boolean. Whether DBS and MBS variants have been
#'    predefined in your vcf. This function by default assumes that DBS and MBS
#'    variants are present in the vcf as SNVs, which are positioned next to each
#'    other. If your DBS/MBS variants are called separately you should set this
#'    argument to TRUE. (default = FALSE)
#' @param remove_duplicate_variants Boolean. Whether duplicate variants are
#'   removed. This is based on genomic coordinates and does not take the
#'   alternative bases into account. It is generally recommended to keep this
#'   on. Turning this off can result in warnings in plot_rainfall. When a
#'   duplicate SNV is identified as part of a DBS, only a single DBS, instead of
#'   a duplicate DBS will be formed. (default = TRUE)
#'
#' @return A GRangesList containing the GRanges obtained from 'vcf_files'
#'
#' @importFrom magrittr %>%
#'
#' @examples
#' ## The example data set consists of three colon samples, three intestine
#' ## samples and three liver samples.  So, to map each file to its appropriate
#' ## sample name, we create a vector containing the sample names:
#' sample_names <- c(
#'   "colon1", "colon2", "colon3",
#'   "intestine1", "intestine2", "intestine3",
#'   "liver1", "liver2", "liver3"
#' )
#'
#' ## We assemble a list of files we want to load.  These files match the
#' ## sample names defined above.
#' vcf_files <- list.files(system.file("extdata",
#'   package = "MutationalPatterns"
#' ),
#' pattern = "sample.vcf", full.names = TRUE
#' )
#'
#' ## Get a reference genome BSgenome object.
#' ref_genome <- "BSgenome.Hsapiens.UCSC.hg19"
#' library("BSgenome")
#' library(ref_genome, character.only = TRUE)
#'
#' ## This function loads the files as GRanges objects.
#' ## For backwards compatability reasons it only loads SNVs by default
#' vcfs <- read_vcfs_as_granges(vcf_files, sample_names, ref_genome)
#'
#' ## To load all variant types use:
#' vcfs <- read_vcfs_as_granges(vcf_files, sample_names, ref_genome, type = "all")
#'
#' ## Loading only indels can be done like this.
#'
#' ## Select data containing indels.
#' vcf_fnames <- list.files(system.file("extdata", package = "MutationalPatterns"),
#'   pattern = "blood.*vcf", full.names = TRUE
#' )
#' sample_names <- c("AC", "ACC55", "BCH")
#'
#' ## Read data and select only the indels.
#' ## Other mutation types can be read in the same way.
#' read_vcfs_as_granges(vcf_fnames, sample_names, ref_genome, type = "indel")
#' @export
read_vcfs_as_granges <- function(vcf_files,
                                 sample_names,
                                 genome,
                                 group = c("auto+sex", "auto", "sex", "circular", "all", "none"),
                                 type = c("snv", "indel", "dbs", "mbs", "all"),
                                 change_seqnames = TRUE,
                                 predefined_dbs_mbs = FALSE,
                                 remove_duplicate_variants = TRUE) {

  # Match arguments
  type <- match.arg(type)
  group <- match.arg(group)

  # Check sample names
  if (length(vcf_files) != length(sample_names)) {
    stop("Please provide the same number of sample names as VCF files", call. = FALSE)
  }

  # Get the reference genome
  tryCatch(
    error = function(cnd) {
      stop("Please provide the name of a BSgenome object.", call. = FALSE)
    },
    {
      genome <- BSgenome::getBSgenome(genome)
    }
  )

  # Read vcfs
  grl <- purrr::map(vcf_files, .read_single_vcf_as_grange, genome, group, change_seqnames, remove_duplicate_variants) %>%
    GenomicRanges::GRangesList()

  # Filter for mutation type
  if (type != "all") {
    grl <- get_mut_type(grl, type, predefined_dbs_mbs)
  }

  # Set the provided names for the samples.
  names(grl) <- sample_names

  return(grl)
}

#' Read a single VCF file into a GRanges object
#'
#' This function reads a Variant Call Format (VCF) file into a GRanges object
#' In addition to loading the files, this
#' function applies the same seqlevel style to the GRanges objects as the
#' reference genome passed in the 'genome' parameter.
#'
#' @param vcf_file A VCF file name
#' @param genome BSgenome object
#' @param group Selector for a seqlevel group.  All seqlevels outside
#'              of this group will be removed.  Possible values:
#'              * 'all' for all chromosomes;
#'              * 'auto' for autosomal chromosomes;
#'              * 'sex' for sex chromosomes;
#'              * 'auto+sex' for autosomal + sex chromosomes (default);
#'              * 'circular' for circular chromosomes;
#'              * 'none' for no filtering, which results in keeping all
#'                seqlevels from the VCF file.
#' @param change_seqnames Boolean. Whether to change the seqlevelStyle of the
#'   vcf to that of the BSgenome object.
#' @param remove_duplicate_variants Boolean. Whether duplicate variants are
#'   removed. This is based on genomic coordinates and does not take the
#'   alternative bases into account. It is generally recommended to keep this
#'   on. Turning this off can result in warnings in plot_rainfall. When a
#'   duplicate SNV is identified as part of a DBS, only a single DBS, instead of
#'   a duplicate DBS will be formed. (default = TRUE)
#' @return A GRanges object
#' @importFrom magrittr %>%
#' @noRd
#'
.read_single_vcf_as_grange <- function(vcf_file, genome, group, change_seqnames, remove_duplicate_variants) {

  # Use VariantAnnotation's readVcf, but only store the
  # GRanges information in memory.  This speeds up the
  # loading significantly.
  # Muffle the warning about duplicate keys.
  withCallingHandlers(
    {
      gr <- GenomicRanges::granges(VariantAnnotation::readVcf(vcf_file))
    },
    warning = function(w) {
      if (grepl("duplicate keys in header will be forced to unique rownames", conditionMessage(w))) {
        invokeRestart("muffleWarning")
      }
    }
  )


  # Throw a warning when a file is empty.
  # Return a empty GR, to prevent errors with changing the seqlevels.
  if (!length(gr)) {
    warning(paste0(
      "There were 0 variants (before filtering) found in the vcf file: ", vcf_file,
      "\nYou might want to remove this sample from your analysis."
    ), call. = FALSE)
    return(gr)
  }

  # Convert to a single chromosome naming standard.
  if (change_seqnames == TRUE) {
    tryCatch(
      error = function(cnd) {
        message(conditionMessage(cnd))
        stop("The seqlevelStyle of the vcf could not be changed to that of the reference.
                     You can run this function with `change_seqnames = F` and `group = 'none'`, 
                     to prevent this error.
                     However, you then have to make sure that the seqnames (chromosome names) are
                     the same between your vcfs and the reference BSgenome object.
                     (The message of the internal error causing this problem is shown above.)",
          call. = FALSE
        )
      },
      {
        GenomeInfoDb::seqlevelsStyle(gr) <- GenomeInfoDb::seqlevelsStyle(genome)[1]
      }
    )
  }
  
  # Change the genome name of the granges
  genome_name <- GenomeInfoDb::genome(genome)[[1]]
  GenomeInfoDb::genome(gr) <- genome_name

  # Filter for variants with the correct seqlevels
  if (group != "none") {
    tryCatch(
      error = function(cnd) {
        message(conditionMessage(cnd))
        stop("The vcf could not be filtered for the specific seqlevels group.
                     You can run this function with `group = 'none'`, to prevent this error.
                     (The message of the internal error causing this problem is shown above.)",
          call. = FALSE
        )
      },
      {
        gr <- .filter_seqlevels(gr, group, genome)
      }
    )
  }

  # Check for duplicate variants
  if (remove_duplicate_variants == TRUE){
    nr_duplicated <- gr %>%
      duplicated() %>%
      sum()
    if (nr_duplicated) {
      warning(paste0(
        "There were ", nr_duplicated, " duplicated variants in vcf file: ",
        vcf_file,
        " They have been filtered out."
      ), call. = FALSE)
      gr <- BiocGenerics::unique(gr)
    }
  }

  return(gr)
}

#' Filter a GRanges object based on seqlevels
#'
#' This function filters a GRanges object based on a group of seqnames.
#'
#' @param gr GRanges object
#' @param group Selector for a seqlevel group.  All seqlevels outside
#'              of this group will be removed.  Possible values:
#'              * 'auto' for autosomal chromosomes;
#'              * 'sex' for sex chromosomes;
#'              * 'auto+sex' for autosomal + sex chromosomes (default);
#'              * 'circular' for circular chromosomes;
#'              * 'none' for no filtering, which results in keeping all
#'                seqlevels from the VCF file.
#' @param genome BSgenome object
#' @return A GRanges object
#' @noRd
#'
.filter_seqlevels <- function(gr, group, genome) {
  groups <- c()
  # These variables are needed to extract the possible seqlevels
  ref_style <- GenomeInfoDb::seqlevelsStyle(genome)
  ref_organism <- GenomeInfoDb::organism(genome)

  if (group == "auto+sex") {
    groups <- c(
      GenomeInfoDb::extractSeqlevelsByGroup(
        species = ref_organism,
        style = ref_style,
        group = "auto"
      ),
      GenomeInfoDb::extractSeqlevelsByGroup(
        species = ref_organism,
        style = ref_style,
        group = "sex"
      )
    )

    # In some cases, the seqlevelsStyle returns multiple styles.
    # In this case, we need to do a little more work to extract
    # a vector of seqlevels from it.
    groups_names <- names(groups)
    if (!is.null(groups_names)) {
      # The seqlevels in the groups are now duplicated.
      # The following code deduplicates the list items, so that
      # creating a data frame will work as expected.
      unique_names <- unique(groups_names)
      groups <- lapply(unique_names, function(x) groups[groups_names == x])
      groups <- lapply(groups, unlist, recursive = FALSE)

      # In case there are multiple styles applied, we only use the first.
      groups <- unique(as.vector(groups[[1]]))
    }
  }
  else {
    groups <- GenomeInfoDb::extractSeqlevelsByGroup(
      species = ref_organism,
      style = ref_style,
      group = group
    )
    groups <- unique(as.vector(t(groups)))
  }

  # The provided VCF files may not contain all chromosomes that are
  # available in the reference genome.  Therefore, we only take the
  # chromosomes that are actually available in the VCF file,
  # belonging to the filter group.
  groups <- BiocGenerics::intersect(groups, GenomeInfoDb::seqlevels(gr))

  # We use 'pruning.mode = "tidy"' to minimize the deleterious effect
  # on variants, yet, remove all variants that aren't in the filter
  # group.  By default, keepSeqlevels would produce an error.
  gr <- GenomeInfoDb::keepSeqlevels(gr, groups, pruning.mode = "tidy")

  return(gr)
}