File: methods-filterVcf.R

package info (click to toggle)
r-bioc-variantannotation 1.10.5-1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 2,172 kB
  • ctags: 109
  • sloc: ansic: 1,088; sh: 4; makefile: 2
file content (150 lines) | stat: -rw-r--r-- 4,854 bytes parent folder | download
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
setMethod("filterVcf", "character",
    function(file, genome, destination, ..., verbose=TRUE,
             index=FALSE, prefilters=FilterRules(),
             filters=FilterRules(), param=ScanVcfParam())
{
   if (file.exists(destination))
        stop(sprintf("file '%s' exists and will not be over-written", 
             destination))

    tbx <- open(TabixFile(file, yieldSize=100000))
    on.exit(close(tbx))
 

    filterVcf(tbx, genome=genome, destination=destination, ...,
              verbose=verbose, index=index, prefilters=prefilters,
              filters=filters, param=param)
})

.unlistScan <- function(...)
    unlist(scanTabix(...), use.names=FALSE)

.prefilter <-
    function(tbxFile, verbose, prefilters, param, ...)
{
    if (verbose)
        message("starting prefilter")
    if (length(vcfWhich(param))) {
        warning("vcfWhich(param) ignored when using prefilter")
        vcfWhich(param) <- GRanges()
    }
    if (!isOpen(tbxFile)) {
        open(tbxFile)
        on.exit(close(tbxFile), add=TRUE)
    }

    prefilteredFilename <- tempfile()
    prefiltered <- file(prefilteredFilename, "w")
    needsClosing <- TRUE
    on.exit(if (needsClosing) close(prefiltered), add=TRUE)
 
    ## copy header
    writeLines(headerTabix(tbxFile)$header, prefiltered)

    ## prefilter
    param <- vcfWhich(param)
    nTotal <- 0L
    while (length(tbxChunk <- .unlistScan(tbxFile, ..., param=param))) {
        if (verbose)
            message("prefiltering ", nTotal <- nTotal + length(tbxChunk),
                    " records")
        tbxChunk <- subsetByFilter(tbxChunk, prefilters)
        writeLines(tbxChunk, prefiltered)
    }
    close(prefiltered)
    needsClosing <- FALSE

    prefilteredFilename
}

.filter <-
    function(tbxFile, genome, destination, verbose, filters, param, ...)
{
    if (verbose)
        message("starting filter")
    if (length(vcfWhich(param))) {
        yieldSize(tbxFile) <- NA_integer_ 
        vcfChunk <- readVcf(tbxFile, genome, ..., param=param)
        if (verbose)
            message("filtering ", nrow(vcfChunk), " records")
        vcfChunk <- subsetByFilter(vcfChunk, filters)
        writeVcf(vcfChunk, filtered)
    } else {
        if (!isOpen(tbxFile)) {
            open(tbxFile)
            on.exit(close(tbxFile), add=TRUE)
        }

        filtered <- file(destination, open="w")
        needsClosing <- TRUE
        on.exit(if (needsClosing) close(filtered), add=TRUE)

        nTotal <- 0L
        while (nrow(vcfChunk <- readVcf(tbxFile, genome, ..., param=param))) {
            if (verbose)
                message("filtering ", nTotal <- nTotal + nrow(vcfChunk),
                        " records")
            vcfChunk <- subsetByFilter(vcfChunk, filters)
            writeVcf(vcfChunk, filtered)
        }
        close(filtered)
        needsClosing <- FALSE
    }
    if (verbose)
        message("completed filtering")
    destination
}

setMethod("filterVcf", "TabixFile",
    function(file, genome, destination, ..., verbose = TRUE,
             index = FALSE,
             prefilters = FilterRules(), filters = FilterRules(),
             param = ScanVcfParam())
{
    if (!isSingleString(destination))
        stop("'destination' must be character(1)")
    if (!isTRUEorFALSE(verbose))
        stop("'verbose' must be TRUE or FALSE")
    if (!isTRUEorFALSE(index))
        stop("'index' must be TRUE or FALSE")

    if (!length(prefilters) && !length(filters))
        stop("no 'prefilters' or 'filters' specified")

    if (length(prefilters)) {
        yieldSize <- yieldSize(file)
        file <- .prefilter(file, verbose, prefilters, param, ...)
        if(verbose){
            message(sprintf("prefiltered to %s", file))
            }
        if (length(filters)) {
            ## TabixFile needs to be bgzipped and indexed
            ## FIXME: all records are read at next stage, so no need to index?
            if (verbose)
                message("prefilter compressing and indexing ", sQuote(file))
            gzFilename <- bgzip(file, overwrite = TRUE)
            indexTabix(gzFilename, format = "vcf4")
            file <- TabixFile(gzFilename, yieldSize=yieldSize)
        } else {
            ## file.rename does not work across file systems, so be expensive
            file.copy(file, destination)
            unlink(file)
        }
      }

    if (length(filters)) {
        file <- .filter(file, genome, destination, verbose, filters,
                        param, ...)
        } # if filters
 
    if (index) {
        if (verbose)
            message("compressing and indexing ", sQuote(file))
        gzFilename <- sprintf("%s.gz", destination)
        gzFilename <- bgzip(file, gzFilename, overwrite = TRUE)
        destination <- indexTabix(gzFilename, format = "vcf")
        unlink(file)
    }

    invisible(destination)
})