snpgdsVCF2GDS | R Documentation |
Reformat Variant Call Format (VCF) file(s)
snpgdsVCF2GDS(vcf.fn, out.fn, method=c("biallelic.only", "copy.num.of.ref"),
snpfirstdim=FALSE, compress.annotation="LZMA_RA", compress.geno="",
ref.allele=NULL, ignore.chr.prefix="chr", verbose=TRUE)
vcf.fn |
the file name of VCF format, |
out.fn |
the file name of output GDS |
method |
either "biallelic.only" by default or "copy.num.of.ref", see details |
snpfirstdim |
if TRUE, genotypes are stored in the individual-major mode, (i.e, list all SNPs for the first individual, and then list all SNPs for the second individual, etc) |
compress.annotation |
the compression method for the GDS variables,
except "genotype"; optional values are defined in the function
|
compress.geno |
the compression method for "genotype"; optional
values are defined in the function |
ref.allele |
|
ignore.chr.prefix |
a vector of character, indicating the prefix of chromosome which should be ignored, like "chr"; it is not case-sensitive |
verbose |
if |
GDS – Genomic Data Structures used for storing genetic array-oriented data, and the file format used in the gdsfmt package.
VCF – The Variant Call Format (VCF), which is a generic format for storing DNA polymorphism data such as SNPs, insertions, deletions and structural variants, together with rich annotations.
If there are more than one file names in vcf.fn
,
snpgdsVCF2GDS
will merge all dataset together if they all contain
the same samples. It is useful to combine genetic/genomic data together if
VCF data are divided by chromosomes.
method = "biallelic.only"
: to exact bi-allelic and polymorhpic
SNP data (excluding monomorphic variants);
method = "copy.num.of.ref"
: to extract and store dosage (0, 1, 2)
of the reference allele for all variant sites, including bi-allelic SNPs,
multi-allelic SNPs, indels and structural variants.
Haploid and triploid calls are allowed in the transfer, the variable
snp.id
stores the original the row index of variants, and the
variable snp.rs.id
stores the rs id.
When snp.chromosome
in the GDS file is character, SNPRelate treats
a chromosome as autosome only if it can be converted to a numeric value (
like "1", "22"). It uses "X" and "Y" for non-autosomes instead of numeric
codes. However, some software format chromosomes in VCF files with a prefix
"chr". Users should remove that prefix when importing VCF files by setting
ignore.chr.prefix = "chr"
.
The extended GDS format is implemented in the SeqArray package to support the storage of single nucleotide variation (SNV), insertion/deletion polymorphism (indel) and structural variation calls. It is strongly suggested to use SeqArray for large-scale whole-exome and whole-genome sequencing variant data instead of SNPRelate.
Return the file name of GDS format with an absolute path.
Xiuwen Zheng
The variant call format and VCFtools. Danecek P, Auton A, Abecasis G, Albers CA, Banks E, DePristo MA, Handsaker RE, Lunter G, Marth GT, Sherry ST, McVean G, Durbin R; 1000 Genomes Project Analysis Group. Bioinformatics. 2011 Aug 1;27(15):2156-8. Epub 2011 Jun 7.
http://corearray.sourceforge.net/
snpgdsBED2GDS
# the VCF file
vcf.fn <- system.file("extdata", "sequence.vcf", package="SNPRelate")
cat(readLines(vcf.fn), sep="\n")
snpgdsVCF2GDS(vcf.fn, "test1.gds", method="biallelic.only")
snpgdsSummary("test1.gds")
snpgdsVCF2GDS(vcf.fn, "test2.gds", method="biallelic.only", snpfirstdim=TRUE)
snpgdsSummary("test2.gds")
snpgdsVCF2GDS(vcf.fn, "test3.gds", method="copy.num.of.ref", snpfirstdim=TRUE)
snpgdsSummary("test3.gds")
snpgdsVCF2GDS(vcf.fn, "test4.gds", method="copy.num.of.ref")
snpgdsSummary("test4.gds")
snpgdsVCF2GDS(vcf.fn, "test5.gds", method="copy.num.of.ref",
ref.allele=c("A", "T", "T", "T", "A"))
snpgdsSummary("test5.gds")
# open "test1.gds"
(genofile <- snpgdsOpen("test1.gds"))
read.gdsn(index.gdsn(genofile, "sample.id"))
read.gdsn(index.gdsn(genofile, "snp.rs.id"))
read.gdsn(index.gdsn(genofile, "genotype"))
# close the file
snpgdsClose(genofile)
# open "test2.gds"
(genofile <- snpgdsOpen("test2.gds"))
read.gdsn(index.gdsn(genofile, "sample.id"))
read.gdsn(index.gdsn(genofile, "snp.rs.id"))
read.gdsn(index.gdsn(genofile, "genotype"))
# close the file
snpgdsClose(genofile)
# open "test3.gds"
(genofile <- snpgdsOpen("test3.gds"))
read.gdsn(index.gdsn(genofile, "sample.id"))
read.gdsn(index.gdsn(genofile, "snp.rs.id"))
read.gdsn(index.gdsn(genofile, "genotype"))
# close the file
snpgdsClose(genofile)
# open "test4.gds"
(genofile <- snpgdsOpen("test4.gds"))
read.gdsn(index.gdsn(genofile, "sample.id"))
read.gdsn(index.gdsn(genofile, "snp.rs.id"))
read.gdsn(index.gdsn(genofile, "snp.allele"))
read.gdsn(index.gdsn(genofile, "genotype"))
# close the file
snpgdsClose(genofile)
# open "test5.gds"
(genofile <- snpgdsOpen("test5.gds"))
read.gdsn(index.gdsn(genofile, "sample.id"))
read.gdsn(index.gdsn(genofile, "snp.rs.id"))
read.gdsn(index.gdsn(genofile, "snp.allele"))
read.gdsn(index.gdsn(genofile, "genotype"))
# close the file
snpgdsClose(genofile)
# delete the temporary files
unlink(paste("test", 1:5, ".gds", sep=""), force=TRUE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.