| NullModel-class {podkat} | R Documentation |
NullModelS4 class for storing null models for later usage with the
assocTest method
Objects of this class are created by calling nullModel.
The following slots are defined for NullModel objects:
type:type of model
residuals:residuals of linear model; for type
“bernoulli”, this is simply the trait vector (see
nullModel-methods for details)
model.matrix:model matrix of the (generalized) linear model trained for the covariates (if any)
inv.matrix:pre-computed inverse of some matrix needed for computing the null distribution; only used for types “logistic” and “linear”
P0sqrt:pre-computed square root of matrix P_0
(see Subsections 9.1 and 9.5 of the package vignette);
needed for computing the null distribution in case the small
sample correction is used for a logistic model; computed only
if nullModel is called with adjExact=TRUE.
coefficients:coefficients of (generalized) linear model trained for the covariates (if any)
na.omit:indices of samples omitted from (generalized) linear model because of missing values in target or covariates
n.cases:for binary traits (types “logistic” and “bernoulli”), the number of cases, i.e. the number of 1's in the trait vector
variance:for continuous traits (type “linear”), this is a single numeric value with the variance of residuals of the linear model; for logistic models with binary traits (type “logistic”), this is a vector with variances of the per-sample Bernoulli distributions; for later use of the exact mixture-of-Bernoulli test (type “bernoulli”), this is the variance of the Bernoulli distribution
prob:for logistic models with binary traits (type “logistic”), this is a vector with probabilities of the per-sample Bernoulli distributions; for later use of the exact mixture-of-Bernoulli test (type “bernoulli”), this is the probability of the Bernoulli distribution
type.resampling:which resampling algorithm was used
res.resampling:matrix with residuals sampled under the null hypothesis (if any)
res.resampling.adj:matrix with residuals sampled under the null hypothesis for the purpose of higher moment correction (if any; only used for logistic models with small sample correction)
call:the matched call with which the object was created
This class serves as the general interface for storing the necessary
phenotype information for a later association test. Objects of this
class should only be created by the nullModel function.
Direct modification of object slots is strongly discouraged!
signature(object="NullModel"):
displays basic information about the null model, such as,
the type of the model and the numbers of covariates.
signature(object="NullModel"):
returns the residuals slot.
signature(object="NullModel"):
returns the names of samples in the null model.
signature(object="NullModel"):
returns the coefficients slot.
signature(x="NullModel"):
returns the number of samples that was used to train the null model.
For a NullModel object x and an index vector i
that is a permutation of 1:length(x),
x[i] returns a new NullModel object in which the samples
have been rearranged according to the permutation i. This is
meant for applications in which the order of the samples in a subsequent
association test is different from the order of the samples when the
null model was trained/created.
Ulrich Bodenhofer bodenhofer@bioinf.jku.at
http://www.bioinf.jku.at/software/podkat
## read phenotype data from CSV file (continuous trait + covariates)
phenoFile <- system.file("examples/example1lin.csv", package="podkat")
pheno <-read.table(phenoFile, header=TRUE, sep=",")
## train null model with all covariates in data frame 'pheno'
model <- nullModel(y ~ ., pheno)
model
length(model)
residuals(model)
## read phenotype data from CSV file (binary trait + covariates)
phenoFile <- system.file("examples/example1log.csv", package="podkat")
pheno <-read.table(phenoFile, header=TRUE, sep=",")
## train null model with all covariates in data frame 'pheno'
model <- nullModel(y ~ ., pheno)
model
length(model)
residuals(model)
## "train" simple Bernoulli model on a subset of 100 samples
model <- nullModel(y ~ 0, pheno[1:100, ])
model
length(model)
residuals(model)