VCFArray-classes: VCFArray constructor and coercion methods.

extract_arrayR Documentation

VCFArray constructor and coercion methods.

Description

extract_array: the function to extract data from a VCF file, by taking VCFArraySeed as input. This function is required by the DelayedArray for the seed contract.

VCFArray: The function to convert data entries inside VCF file into the VCFArray instance.

Usage

## S4 method for signature 'VCFArraySeed'
extract_array(x, index)

VCFArray(file, vindex = character(), name = NA, pfix = NULL)

Arguments

x

the VCFArraySeed object

index

in extract_array(), an unnamed list of subscripts as positive integer vectors, one vector per dimension in x. Empty and missing subscripts (represented by integer(0) and NULL list elements, respectively) are allowed. The subscripts can contain duplicated indices. They cannot contain NAs or non-positive values.

file

takes values for charater string (specifying the VCF file path), VcfFile object, and RangedVcfStack object.

vindex

in VCFArray(), the character string specifying the index file path. This argument is required if an remote VCF file is used for the file argument.

name

the data entry from VCF file to be read into VCFArraySeed / VCFArray. For VCFArray. This argument should always be specified.

pfix

the category that the name belongs to. Available values are fixed, info, and info. Can also Check vcfFields(file) for matching name and pfix.

Value

VCFArray class object.

Examples

fl <- system.file("extdata", "chr22.vcf.gz",
                  package="VariantAnnotation")
va <- VCFArray(fl, name = "GT")
va
vcf <- VariantAnnotation::VcfFile(fl)
va1 <- VCFArray(vcf, name = "GT")
va1
all.equal(va, va1)
## Not run: 
## RangedVcfStack class
extdata <- system.file(package = "GenomicFiles", "extdata")
files <- dir(extdata, pattern="^CEUtrio.*bgz$", full=TRUE)[1:2]
names(files) <- sub(".*_([0-9XY]+).*", "\\1", basename(files))
seqinfo <- as(readRDS(file.path(extdata, "seqinfo.rds")), "Seqinfo")
stack <- GenomicFiles::VcfStack(files, seqinfo)
gr <- as(GenomicFiles::seqinfo(stack)[rownames(stack)], "GRanges")
## RangedVcfStack
rgstack <- GenomicFiles::RangedVcfStack(stack, rowRanges = gr)
rgstack
va2 <- VCFArray(rgstack, name = "SB")
va2

## End(Not run)
## coercion
as(va[1:10, ], "array")

Bioconductor/VCFArray documentation built on Oct. 26, 2023, 6:38 a.m.