\name{ScanVcfParam-class} \docType{class} \alias{ScanVcfParam} \alias{ScanVcfParam-class} \alias{ScanVcfParam,missing-method} \alias{ScanVcfParam,ANY-method} \alias{vcfFixed} \alias{vcfInfo} \alias{vcfGeno} \alias{vcfTrimEmpty} \alias{vcfWhich} \title{Parameters for scanning VCF files} \description{ Use \code{ScanVcfParam()} to create a parameter object influencing which records and fields are imported from a VCF file. Record parsing is based on genomic coordinates and requires a Tabix index file. Individual VCF elements can be specified in the \sQuote{fixed}, \sQuote{info} and \sQuote{geno} arguments. } \usage{ ScanVcfParam(fixed=character(), info=character(), geno=character(), trimEmpty=TRUE, which, ...) ## Accessors vcfFixed(object) vcfInfo(object) vcfGeno(object) vcfTrimEmpty(object) vcfWhich(object) } \arguments{ \item{fixed}{A character() vector of fixed fields to be returned. Possible values are ALT, QUAL and FILTER. The CHROM, POS, ID and REF fields are needed to create the \code{GRanges} of variant locations. Because these are essential fields there is no option to request or omit them. If not specified, all fields are returned; if \code{fixed=NA} only REF is returned. } \item{info}{A character() vector of \sQuote{INFO} fields (see \code{\link{scanVcfHeader}}) to be returned. If not specified, all fields are returned; if \code{info=NA} no fields are returned. } \item{geno}{A character() vector of \sQuote{GENO} fields (see \code{\link{scanVcfHeader}}) to be returned. If not specified, all fields are returned; if \code{geno=NA} no fields are returned. } \item{trimEmpty}{A logical(1) indicating whether \sQuote{GENO} fields with no values should be returned. } \item{which}{An object, for which a method is defined (see usage, above), describing the sequences and ranges to be queried. Variants whose \code{POS} lies in the interval(s) \code{[start, end]} are returned. } \item{object}{An instance of class \code{ScanVcfParam}. } \item{\dots}{Arguments passed to methods. } } \section{Objects from the Class}{ Objects can be created by calls of the form \code{ScanVcfParam()}. } \section{Slots}{ \describe{ \item{\code{which}:}{Object of class \code{"RangesList"} indicating which reference sequence and coordinate variants must overlap. } \item{\code{fixed}:}{Object of class \code{"character"} indicating fixed fields to be returned. } \item{\code{info}:}{Object of class \code{"character"} indicating portions of \sQuote{INFO} to be returned. } \item{\code{geno}:}{Object of class \code{"character"} indicating portions of \sQuote{GENO} to be returned. } \item{\code{trimEmpty}:}{Object of class \code{"logical"} indicating whether empty \sQuote{GENO} fields are to be returned. } } } \section{Functions and methods}{ See 'Usage' for details on invocation. Constructor: \describe{ \item{ScanVcfParam:}{Returns a \code{ScanVcfParam} object. The \code{which} argument to the constructor can be one of several types, as documented above.} } Accessors: \describe{ \item{vcfFixed, vcfInfo, vcfGeno, vcfTrimEmpty, vcfWhich:}{ Return the corresponding field from \code{object}. } } Methods: \describe{ \item{show}{Compactly display the object. } } } \author{ Martin Morgan and Valerie Obenchain } \seealso{ \code{\link{readVcf}} } \examples{ ScanVcfParam() ## -------------------------------------- ## 'which' argument ## -------------------------------------- ## To subset on genomic coordinates, create a GRanges, RangedData or ## RangesList with the ranges of interest. This object is supplied ## to ScanVcfParam() as the 'which' argument. which <- RangesList(seq1=IRanges(1000, 2000), seq2=IRanges(c(100, 1000), c(1000, 2000))) ScanVcfParam(which=which) ## -------------------------------------- ## 'fixed', 'info' and 'geno' arguments ## -------------------------------------- ## This param specifies the return of the "GT" 'geno' field and the ## subset of ranges in 'which'. All 'fixed' and 'info' fields will be ## returned. ScanVcfParam(geno="GT", which=which) ## Here two 'fixed' and one 'geno' field are specified ScanVcfParam(fixed=c("ALT", "QUAL"), geno="GT", info=NA) ## Return only the 'fixed' fields ScanVcfParam(geno=NA, info=NA) } \keyword{classes}