VCFHeader-class {VariantAnnotation}R Documentation

VCFHeader instances

Description

The VCFHeader class holds Variant Call Format (VCF) file header information and is produced from a call to scanVcfHeader.

Details

The VCFHeader class is holds header information from a from VCF file.

Slots :

reference

character() vector

sample

character() vector

header

DataFrameList-class

Constructor

VCFHeader(reference = character(), samples = character(), header = DataFrameList(), ...)

Accessors

In the following code snippets x is a VCFHeader object.

samples(x): Returns a character() vector of names of samples.

header(x): Returns all information in the header slot which includes meta, info and geno if present.

meta(x), meta(x)<- value: The getter returns a DataFrameList with possible names of META, SAMPLE and PEDIGREE. The META DataFrame is a catch-all for non-standard header fields. This includes any information represented as simple key-value pairs. The replacement value can be a DataFrame or DataFrameList.

fixed(x), fixed(x)<- value: Returns a DataFrameList of information pertaining to any of ‘REF’, ‘ALT’, ‘FILTER’ and ‘QUAL’. Replacement value must be a DataFrameList with one or more of the following names, ‘QUAL’, ‘FILTER’, ‘REF’ and ‘ALT’.

info(x), info(x)<- value: Gets or sets a DataFrame of ‘INFO’ information. Replacement value must be a DataFrame with 3 columns named ‘Number’, ‘Type’ and ‘Description’.

geno(x), geno(x)<- value: Returns a DataFrame of ‘FORMAT’ information. Replacement value must be a DataFrame with 3 columns named ‘Number’, ‘Type’ and ‘Description’.

reference(x): Returns a character() vector with names of reference sequences. Not relevant for scanVcfHeader.

Arguments

reference

A character() vector of sequences.

sample

A character() vector of sample names.

header

A DataFrameList of parsed header lines (preceeded by “##”) present in the VCF file.

...

Additional arguments passed to methods.

Author(s)

Valerie Obenchain

See Also

scanVcfHeader, DataFrameList

Examples

  fl <- system.file("extdata", "structural.vcf", package="VariantAnnotation")
  hdr <- scanVcfHeader(fl)

  fixed(hdr)
  info(hdr)
  geno(hdr)

[Package VariantAnnotation version 1.24.5 Index]