| Covariance | R Documentation |
R6 Class representing a covariance function and data
R6 Class representing a covariance function and data
For the generalised linear mixed model
Y \sim F(μ,σ)
μ = h^-1(Xβ + Zγ)
γ \sim MVN(0,D)
where h is the link function, this class defines Z and D. The covariance is defined by a covariance function, data, and parameters. A new instance can be generated with $new(). The class will generate the relevant matrices Z and D automatically.
**Intitialisation**
A covariance function is specified as an additive formula made up of
components with structure (1|f(j)). The left side of the vertical bar
specifies the covariates in the model that have a random effects structure.
The right side of the vertical bar specify the covariance function 'f' for
that term using variable named in the data 'j'. If there are multiple
covariates on the left side, it is assumed their random effects are
correlated, e.g. (1+x|f(j)). Additive functions are assumed to be
independent, for example, (1|f(j))+(x|f(j)) would create random effects
with zero correlation for the intercept and the parameter on covariate x.
Covariance functions on the right side of the vertical bar are multiplied
together, i.e. (1|f(j)*g(t)).
There are several common functions included for a named variable in data x:
* gr(x): Indicator function (1 parameter)
* fexp(x): Exponential function (2 parameters)
* ar1(x): AR1 function (1 parameter)
* sqexp(x): Squared exponential (1 parameter)
* matern(x): Matern function (2 parameters)
* bessel(x): Modified Bessel function of the 2nd kind (1 parameter)
Parameters are provided to the covariance function as a vector. The parameters in the vector for each function should be provided in the order the covariance functions are written are written. For example, * Formula: '~(1|gr(j))+(1|gr(j*t))'; parameters: 'c(0.25,0.1)' * Formula: '~(1|gr(j)*fexp(t))'; parameters: 'c(0.25,1,0.5)' Note that it is also possible to specify a group membership with two variable alternatively as '(1|gr(j)*gr(t))', for example, but this will require two parameters to be specified, so it is recommended against.
dataData frame with data required to build covariance
formulaCovariance function formula. See 'help(Covariance$new())' for details.
parametersList of lists holding the model parameters. See 'help(Covariance$new())' for details.
ZDesign matrix
DCovariance matrix of the random effects
n()Return the size of the design
Covariance$n()
Scalar
new()Create a new Covariance object
Covariance$new(formula = NULL, data = NULL, parameters = NULL, verbose = TRUE)
formulaFormula describing the covariance function. See Details
dataData frame with data required for constructing the covariance.
parametersList of lists with parameter values for the functions in the model formula. See Details.
verboseLogical whether to provide detailed output.
A Covariance object
df <- nelder(~(cl(5)*t(5)) > ind(5))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = c(0.25,0.7),
data= df)
check()Check if anything has changed and update matrices if so.
Covariance$check(verbose = TRUE)
verboseLogical whether to report if any changes detected.
NULL
df <- nelder(~(cl(5)*t(5)) > ind(5))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
cov$parameters <- list(list(0.01,0.1))
cov$check(verbose=FALSE)
print()Show details of Covariance object
Covariance$print()
...ignored
df <- nelder(~(cl(5)*t(5)) > ind(5))
Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
subset()Keep specified indices and removes the rest
Covariance$subset(index)
indexvector of indices to keep
df <- nelder(~(cl(10)*t(5)) > ind(10))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
cov$subset(1:100)
sampleD()Generate a new D matrix
D is the covariance matrix of the random effects terms in the generalised linear mixed model. This function will return a matrix D for a given set of parameters.
Covariance$sampleD(parameters)
parameterslist of lists, see initialize()
matrix
df <- nelder(~(cl(10)*t(5)) > ind(10))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
cov$sampleD(list(list(0.01,0.1)))
clone()The objects of this class are cloneable with this method.
Covariance$clone(deep = FALSE)
deepWhether to make a deep clone.
## ------------------------------------------------
## Method `Covariance$new`
## ------------------------------------------------
df <- nelder(~(cl(5)*t(5)) > ind(5))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = c(0.25,0.7),
data= df)
## ------------------------------------------------
## Method `Covariance$check`
## ------------------------------------------------
df <- nelder(~(cl(5)*t(5)) > ind(5))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
cov$parameters <- list(list(0.01,0.1))
cov$check(verbose=FALSE)
## ------------------------------------------------
## Method `Covariance$print`
## ------------------------------------------------
df <- nelder(~(cl(5)*t(5)) > ind(5))
Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
## ------------------------------------------------
## Method `Covariance$subset`
## ------------------------------------------------
df <- nelder(~(cl(10)*t(5)) > ind(10))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
cov$subset(1:100)
## ------------------------------------------------
## Method `Covariance$sampleD`
## ------------------------------------------------
df <- nelder(~(cl(10)*t(5)) > ind(10))
cov <- Covariance$new(formula = ~(1|gr(j)*ar1(t)),
parameters = list(list(0.05,0.8)),
data= df)
cov$sampleD(list(list(0.01,0.1)))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.