TVTBparam-class | R Documentation |
The TVTBparam
class stores recurrent parameters of the TVTB
package in a convenient format.
TVTBparam(
genos, ranges = GRangesList(),
aaf = "AAF", maf = "MAF", vep = "CSQ", bp = SerialParam(),
svp = ScanVcfParam(which = reduce(unlist(ranges))))
genos |
A |
ranges |
A |
aaf |
INFO key suffix used to store the alternate allele frequency (AAF). |
maf |
INFO key suffix used to store the minor allele frequency (MAF). |
vep |
INFO key suffix used to extract the VEP predictions. See |
bp |
A |
svp |
A |
For each suffix stored in the TVTBparam
object, TVTB
may store data in the VCF
object under the INFO keys defined as follows:
Statistics across all samples in the ExpandedVCF
(e.g. "MAF").
Statistics across samples associated with a given level of a given phenotype (e.g. "gender_male_MAF").
Users are recommended to avoid using those INFO keys for other purposes.
A TVTBparam
object that contains recurrent parameters.
In the following code snippets x
is a TVTBparam
object.
genos(x)
, genos(x) <- value
Gets or sets the Genotypes
object stored in the genos
slot.
ranges(x)
, ranges(x) <- value
List of genomic ranges to group variants during analyses and plots.
ref(x)
, ref(x) <- value
Gets or sets the character
vector
that declares homozygote reference genotypes.
het(x)
, het(x) <- value
Gets or sets the character
vector
that declares heterozygote genotypes.
alt(x)
, alt(x) <- value
Gets or sets the character
vector
that declares homozygote alternate genotypes.
carrier(x)
Gets a character
vectors of concatenated
heterozygote and homozygote alternate genotypes.
See also het
and alt
accessors.
aaf(x)
, aaf(x) <- value
Gets or sets the INFO key suffix used to store the alternate allele frequency (AAF).
maf(x)
, maf(x) <- value
Gets or sets the INFO key suffix used to store the minor allele frequency (MAF).
vep(x)
, maf(x) <- value
Gets or sets the INFO key suffix used to extract the VEP predictions.
bp(x)
, bp(x) <- value
Gets or sets the BiocParallel
parameters.
suffix(x)
Gets a named character
vector that declares individual suffixes
used to store the data for each set of genotypes in the INFO field of the
VCF
object.
Names of this vector are ref
, het
, alt
, aaf
,
and maf
.
svp(x)
, svp(x) <- value
Gets or sets the ScanVcfParam
parameters.
Kevin Rue-Albrecht
Genotypes
,
VCF
,
ExpandedVCF
,
addCountGenos-methods
vepInPhenoLevel-methods
,
variantsInSamples-methods
,
and BiocParallelParam
.
# Constructors ----
grl <- GenomicRanges::GRangesList(GenomicRanges::GRanges(
"15", IRanges::IRanges(48413170, 48434757, names = "SLC24A5")
))
tparam <- TVTBparam(Genotypes("0|0", c("0|1", "1|0"), "1|1"), ranges = grl)
# Accessors ----
## The Genotypes object stored in the genos slot of the TVTBparam object
## return by the genos() accessor.
genos(tparam)
## Genomic ranges stored in the TVTBparam object returned by the ranges()
## accessor.
ranges(tparam)
## Individual genotypes can be extracted with ref(), het(), alt() accessors.
ref(tparam)
het(tparam)
alt(tparam)
## Their individual INFO key suffixes can be extracted with suffix() applied to
## the above accessors.
suffix(tparam)
suffix(tparam)["ref"]
suffix(tparam)["het"]
suffix(tparam)["alt"]
suffix(tparam)["aaf"]
suffix(tparam)["maf"]
## Heterozygote, and alternate heterozygote genotypes are
## returned by the carrier() accessor.
carrier(tparam)
## INFO key suffix of alternate/minor allele frequency returned by the aaf()
## and maf() accessors.
aaf(tparam)
maf(tparam)
## INFO key suffix of the VEP predictions returned by the vep() accessor.
vep(tparam)
## BiocParallel parameters
bp(tparam)
## ScanVcfParam parameters
svp(tparam)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.