VcfFilterRules-class | R Documentation |
The VcfFilterRules
class can stores multiple types of filters
applicable to various slots of VCF
objects.
All arguments must be VcfFixedRules
, VcfInfoRules
,
VcfVepRules
, VcfFilterRules
of FilterRules
objects.
In the following code snippets x
is a VcfFilterRules
object.
active(x)
, active(x)<-
Get or set the active state of each filter rule in x
.
Inherited from FilterRules
vep(x)
, vep(x)<-
Gets or sets the INFO key where the Ensembl VEP predictions to use for filtering are stored.
type(x)
Gets the type of each filter stored in a VcfFilterRules
object.
Read-only.
VcfFilterRules(...)
constructs an VcfFilterRules
object from
VcfFixedRules
, VcfInfoRules
, VcfVepRules
,
and VcfFilterRules
objects in ...
.
In the code snippets below, x
is a VcfFilterRules
object.
x[i, drop = TRUE]
: Subsets the filter rules using the
same interface as for Vector
.
If all filter rules are of the same type and drop=TRUE
(default),
the resulting object is re-typed to the most specialised class,
if possible. In other words, if all remaining filter rules are of
type "vep"
, the object will be type as VcfVepRules
.
x[[i]]
: Extracts an expression or function via the same interface
as for List
.
x[i] <- value
: Replaces a filter rule by one of any valid class
(VcfFixedRules
, VcfInfoRules
, VcfVepRules
, or
VcfFilterRules
).
The active state(s), name(s), and type(s) (if applicable)
are transferred from value
.
x[[i]] <- value
: The same interface as for
List
. The default active state for new
rules is TRUE
.
In the following code snippets x
is an object of class
VcfFilterRules
,
while values
and ...
are objects from any of the
classes VcfFixedRules
, VcfInfoRules
, VcfVepRules
,
or VcfFilterRules
:
append(x, values, after = length(x))
:
Appends the values onto x
at the index given by after
.
c(x, ...,)
:
Concatenates the filters objects in ...
onto the end of x
.
As described in the S4Vectors
documentation:
eval(expr, envir, enclos)
Evaluates each active rule in a VcfFilterRules
instance
(passed as the expr
argument)
in their respective context of a VCF
object
(passed as the envir
argument).
evalSeparately(expr, envir, enclos)
:
subsetByFilter(x, filter)
summary(object)
See eval,FilterRules,ANY-method
for details.
Kevin Rue-Albrecht
FilterRules
,
VcfFixedRules
,
VcfInfoRules
,
VcfVepRules
,
and VCF
.
# Constructors ----
fixedR <- VcfFixedRules(list(
pass = expression(FILTER == "PASS"),
qual = expression(QUAL > 20)
))
fixedR
infoR <- VcfInfoRules(list(
common = expression(MAF > 0.1), # minor allele frequency
present = expression(ALT + HET > 0) # count of non-REF homozygotes
))
# ...is synonym to...
infoR <- VcfInfoRules(list(
common = expression(MAF > 0.1), # minor allele frequency
present = expression(ALT > 0 | HET > 0)
))
infoR
vepR <- VcfVepRules(list(
missense = expression(Consequence %in% c("missense_variant")),
CADD = expression(CADD_PHRED > 15)
))
vepR
vcfRules <- VcfFilterRules(fixedR, infoR, vepR)
vcfRules
# Accessors ----
## Type of each filter stored in the VcfFilterRules object
type(vcfRules)
# Example data ----
# VCF file
vcfFile <- system.file("extdata", "moderate.vcf", package = "TVTB")
# TVTB parameters
tparam <- TVTBparam(Genotypes("0|0", c("0|1", "1|0"), "1|1"))
# Pre-process variants
vcf <- VariantAnnotation::readVcf(vcfFile, param = tparam)
vcf <- VariantAnnotation::expand(vcf, row.names = TRUE)
vcf <- addOverallFrequencies(vcf, tparam)
# Applying filters to VCF objects ----
## Evaluate filters
eval(vcfRules, vcf)
## Evaluate filters separately
as.data.frame(evalSeparately(vcfRules, vcf))
# Interestingly, the only common missense variant has a lower CADD score
## Deactivate the CADD score filter
active(vcfRules)["CADD"] <- FALSE
## Subset VCF by filters (except CADD, deactivated above)
subsetByFilter(vcf, vcfRules)
# Subsetting and Replacement ----
v123 <- vcfRules[1:3]
# Extract the expression
v5expr <- vcfRules[[5]]
# Subset the object
v5obj <- vcfRules[5]
# Replace the expression (active reset to TRUE, original name retained)
v123[[2]] <- v5expr
# Replace the rule (active state and name transferred from v5obj)
v123[2] <- v5obj
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.