seqBED2GDS | R Documentation |
Conversion between PLINK BED format and SeqArray GDS format.
seqBED2GDS(bed.fn, fam.fn, bim.fn, out.gdsfn, compress.geno="LZMA_RA",
compress.annotation="LZMA_RA", chr.conv=TRUE, include.pheno=TRUE,
optimize=TRUE, digest=TRUE, parallel=FALSE, verbose=TRUE)
seqGDS2BED(gdsfile, out.fn, write.rsid=c("auto", "annot_id", "chr_pos_ref_alt"),
multi.row=FALSE, verbose=TRUE)
bed.fn |
the file name of PLINK binary file, genotype information |
fam.fn |
the file name of first six columns of |
bim.fn |
the file name of extended MAP file with 6 columns, variant
information; if missing, determine the file name using |
gdsfile |
character (a GDS file name), or
a |
out.gdsfn |
the file name, output a file of SeqArray format |
out.fn |
the file name of PLINK binary format without extended names |
compress.geno |
the compression method for "genotype"; optional
values are defined in the function |
compress.annotation |
the compression method for the GDS variables,
except "genotype"; optional values are defined in the function
|
chr.conv |
if |
include.pheno |
if |
optimize |
if |
digest |
a logical value (TRUE/FALSE) or a character ("md5", "sha1", "sha256", "sha384" or "sha512"); add hash codes to the GDS file if TRUE or a digest algorithm is specified |
parallel |
|
write.rsid |
|
multi.row |
if |
verbose |
if |
Return the file name of SeqArray file with an absolute path.
Xiuwen Zheng
seqSNP2GDS
, seqVCF2GDS
library(SNPRelate)
# PLINK BED files
bed.fn <- system.file("extdata", "plinkhapmap.bed.gz", package="SNPRelate")
fam.fn <- system.file("extdata", "plinkhapmap.fam.gz", package="SNPRelate")
bim.fn <- system.file("extdata", "plinkhapmap.bim.gz", package="SNPRelate")
# convert bed to gds
seqBED2GDS(bed.fn, fam.fn, bim.fn, "tmp.gds")
seqSummary("tmp.gds")
# convert gds to bed
gdsfn <- seqExampleFileName("gds")
seqGDS2BED(gdsfn, "plink")
# remove the temporary file
unlink(c("tmp.gds", "plink.fam", "plink.bim", "plink.bed"), force=TRUE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.