writeVcf {VariantAnnotation} | R Documentation |
Write Variant Call Format (VCF) files to disk
## S4 method for signature 'VCF,character' writeVcf(obj, filename, index = FALSE, ...) ## S4 method for signature 'VCF,connection' writeVcf(obj, filename, index = FALSE, ...)
obj |
Object containing data to be written out. At present only accepts VCF. |
filename |
The character() name of the VCF file, or a connection
(e.g., |
index |
Whether to bgzip the output file and generate a tabix index. |
... |
Additional arguments, passed to methods.
|
A VCF file can be written out from data in a VCF
file. More
general methods to write out from other objects may be added in the future.
writeVcf
conforms to the VCF standards on the 1000 Genomes Project
web site (see references). When 'fileformat' is not present in the header
data, format is written out according to the standard at the time of the
Bioconductor release, e.g., for Bioconductor 3.1 the format is VCFv4.2.
Large VCF files (i.e., > 1e5 records) are written out in chunks; VCF
files with < 1e5 records are not chunked. The optimal number of records
per chunk depends on both the number of records and complexity of the
data. Currently writeVcf
determines records per chunk based on the
total number of records only. To override this behavior or experiment with
other values use nchunk
as an integer or NA. An integer value
represents the number of records per chunk reguardless of the size of the VCF;
NA disables all chunking.
writeVcf(vcf, tempfile()) ## default chunking
writeVcf(vcf, tempfile(), nchunk = 1e6) ## chunk by 1e6
writeVcf(vcf, tempfile(), nchunk = NA) ## no chunking
VCF file
Valerie Obenchain and Michael Lawrence
http://vcftools.sourceforge.net/specs.html outlines the VCF specification.
http://samtools.sourceforge.net/mpileup.shtml contains
information on the portion of the specification implemented by
bcftools
.
http://samtools.sourceforge.net/ provides information on
samtools
.
fl <- system.file("extdata", "ex2.vcf", package="VariantAnnotation") out1.vcf <- tempfile() out2.vcf <- tempfile() in1 <- readVcf(fl, "hg19") writeVcf(in1, out1.vcf) in2 <- readVcf(out1.vcf, "hg19") writeVcf(in2, out2.vcf) in3 <- readVcf(out2.vcf, "hg19") stopifnot(all(in2 == in3)) ## write incrementally out3.vcf <- tempfile() con <- file(out3.vcf, open="a") writeVcf(in1[1:2,], con) writeVcf(in1[-(1:2),], con) close(con) readVcf(out3.vcf, "hg19")