lmebreed {lme4breeding}R Documentation

Fit mixed-effects models incorporating relationships

Description

Fit linear or generalized linear mixed models incorporating the effects of relationships.

Usage

lmebreed(formula, data, family = NULL, REML = TRUE,
           relmat = list(), addmat = list(), control = list(),
           start = NULL, verbose = TRUE, subset, weights,
           na.action, offset, contrasts = NULL, model = TRUE,
           x = TRUE, dateWarning=TRUE, returnParams=FALSE,
           rotation=FALSE, coefOutRotation=8, ...)

Arguments

relmat

an optional named list of relationship matrices (not the inverse). Internally the Cholesky decomposition of those matrices will be computed. The names of the elements must correspond to the names of grouping factors for random-effects terms in the formula argument.

addmat

an optional named list of customized incidence matrices. The names of the elements must correspond to the names of grouping factors for random-effects terms in the formula argument. Depending on the use-case the element in the list may be a single matrix or a list of matrices. Please see examples and vignettes to learn how to use it.

dateWarning

an logical value indicating if you want to be warned when a new version of lme4breeding is available on CRAN.

returnParams

an logical value indicating if you want to only get the incidence matrices of the model.

rotation

an logical value indicating if you want to compute the eigen decomposition of the relationship matrix to rotate y and X and accelerate the computation. See details.

coefOutRotation

a numeric value denoting the inter-quantile outlier coefficient to be used in the rotation of the response when using the eigen decomposition to avoid overshooting.

formula

as in lmer

data

as in lmer

family

as in glmer

REML

as in lmer

control

as in lmer

start

as in lmer

verbose

as in lmer

subset

as in lmer

weights

as in lmer

na.action

as in lmer

offset

as in lmer

contrasts

as in lmer

model

as in lmer

x

as in lmer

...

as in lmer

Details

All arguments to this function are the same as those to the function lmer except relmat and addmat which must be named lists. Each name must correspond to the name of a grouping factor in a random-effects term in the formula. The observed levels of that factor must be contained in the rownames and columnames of the relmat. Each relmat is the relationship matrix restricted to the observed levels and applied to the model matrix for that term. The incidence matrices in the addmat argument must match the dimensions of the final fit (pay attention to missing data in responses).

It is important to remember that when you use the relmat argument you are providing the square root of the relationship matrix and to recover the correct BLUPs for those effects you need to use the ranef function which internally multiple those BLUPs the square root of the relationship matrix one more time to recover the correct BLUPs.

The argument rotation applies the eigen decomposition proposed by Lee and Van der Werf in 2016 and makes the genetic evaluation totally sparse leading to incredible gains in speed compared to the classical approach. Internally, the eigen decomposition UDU' is carried in the relationship matrix. The U matrix is then taken to the n x n level (n being the number of records), and post-multiplied by a matrix of records presence (n x n) using the element-wise multiplication of two matrices (Schur product). By default is not activated since this may not provide the exact same variance components than other software due to numerical reasons. If you would like to obtain the exact same variance components than other software you will have to keep rotation=FALSE. This will slow down considerably the speed. Normally when the rotation is activated and variance components differ slightly with other software they will still provide highly similar estimates at the speed of hundreds or thousands of times faster. Please consider this.

When using the optimizer argument inside the lmerControl keep in mind that the number of iterations is called differently depending on the optimizer. For Nelder_Mead, bobyqa and nlminbwrap is called "maxfun" whereas for nloptwrap is called "maxeval". This should be passed inside a list in the optCtrl argument. For example:

lmebreed(... , control = lmerControl( optimizer="Nelder_Mead", optCtrl=list(maxfun=100) ), ... )

Example Datasets

The package has been equiped with several datasets to learn how to use the lme4breeding package:

* DT_halfdiallel, DT_fulldiallel and DT_mohring datasets have examples to fit half and full diallel designs.

* DT_h2 to calculate heritability

* DT_cornhybrids and DT_technow datasets to perform genomic prediction in hybrid single crosses

* DT_wheat dataset to do genomic prediction in single crosses in species displaying only additive effects.

* DT_cpdata dataset to fit genomic prediction models within a biparental population coming from 2 highly heterozygous parents including additive, dominance and epistatic effects.

* DT_polyploid to fit genomic prediction and GWAS analysis in polyploids.

* DT_gryphon data contains an example of an animal model including pedigree information.

* DT_btdata dataset contains an animal (birds) model.

* DT_legendre simulated dataset for random regression model.

* DT_sleepstudy dataset to know how to translate lme4 models to sommer models.

* DT_ige dataset to show how to fit indirect genetic effect models.

Value

a lmebreed object.

Author(s)

Giovanny Covarrubias-Pazaran

References

Giovanny Covarrubias-Pazaran (2024). lme4breeding: enabling genetic evaluation in the age of genomic data. To be submitted to Bioinformatics.

Douglas Bates, Martin Maechler, Ben Bolker, Steve Walker (2015). Fitting Linear Mixed-Effects Models Using lme4. Journal of Statistical Software, 67(1), 1-48.

Lee & Van der Werf (2016). MTG2: an efficient algorithm for multivariate linear mixed model analysis based on genomic information. Bioinformatics, 32(9), 1420-1422.

Examples


data(DT_example)
DT <- DT_example
A <- A_example

ansMain <- lmebreed(Yield ~ Env + (1|Name),
                        relmat = list(Name = A ),
                        data=DT)
vc <- VarCorr(ansMain); print(vc,comp=c("Variance"))

BLUP <- ranef(ansMain, condVar=TRUE)$Name
SEs <- attr(BLUP, which="postVar")[,,]


[Package lme4breeding version 1.0.30 Index]