Title: | Genomic Breeding Tools: Genetic Variance Prediction and Cross-Validation |
Version: | 1.3.2 |
Description: | The main attribute of 'PopVar' is the prediction of genetic variance in bi-parental populations, from which the package derives its name. 'PopVar' contains a set of functions that use phenotypic and genotypic data from a set of candidate parents to 1) predict the mean, genetic variance, and superior progeny value of all, or a defined set of pairwise bi-parental crosses, and 2) perform cross-validation to estimate genome-wide prediction accuracy of multiple statistical models. More details are available in Mohammadi, Tiede, and Smith (2015, <doi:10.2135/cropsci2015.01.0030>). A dataset 'think_barley.rda' is included for reference and examples. |
License: | GPL-3 |
Encoding: | UTF-8 |
LazyData: | true |
RoxygenNote: | 7.3.1 |
URL: | https://github.com/UMN-BarleyOatSilphium/PopVar |
Depends: | R (≥ 3.5.0) |
Imports: | BGLR, qtl, rrBLUP, stats, utils, methods, parallel |
Suggests: | knitr, rmarkdown |
VignetteBuilder: | knitr |
BugReports: | https://github.com/UMN-BarleyOatSilphium/PopVar/issues |
NeedsCompilation: | no |
Packaged: | 2024-10-20 14:13:00 UTC; jeffneyhart |
Author: | Tyler Tiede [aut],
Jeffrey Neyhart |
Maintainer: | Jeffrey Neyhart <neyhartje@gmail.com> |
Repository: | CRAN |
Date/Publication: | 2024-10-20 22:40:07 UTC |
PopVar: Genomic Breeding Tools: Genetic Variance Prediction and Cross-Validation
Description
The main attribute of 'PopVar' is the prediction of genetic variance in bi-parental populations, from which the package derives its name. 'PopVar' contains a set of functions that use phenotypic and genotypic data from a set of candidate parents to 1) predict the mean, genetic variance, and superior progeny value of all, or a defined set of pairwise bi-parental crosses, and 2) perform cross-validation to estimate genome-wide prediction accuracy of multiple statistical models. More details are available in Mohammadi, Tiede, and Smith (2015, doi:10.2135/cropsci2015.01.0030). A dataset 'think_barley.rda' is included for reference and examples.
Author(s)
Maintainer: Jeffrey Neyhart neyhartje@gmail.com (ORCID) [contributor]
Authors:
Tyler Tiede tyler.tiede7@gmail.com
Other contributors:
See Also
Useful links:
Report bugs at https://github.com/UMN-BarleyOatSilphium/PopVar/issues
Internal functions
Description
Internal functions generally not to be called by the user.
Usage
par_position(crossing.table, par.entries)
par_name(crossing.mat, par.entries)
tails(GEBVs, tail.p)
maf_filt(G)
XValidate_nonInd(
y.CV = NULL,
G.CV = NULL,
models.CV = NULL,
frac.train.CV = NULL,
nCV.iter.CV = NULL,
burnIn.CV = NULL,
nIter.CV = NULL
)
XValidate_Ind(
y.CV = NULL,
G.CV = NULL,
models.CV = NULL,
nFold.CV = NULL,
nFold.CV.reps = NULL,
burnIn.CV = NULL,
nIter.CV = NULL
)
calc_marker_effects(
M,
y.df,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
nIter,
burnIn
)
Arguments
crossing.table |
The crossing table. |
par.entries |
The parent entries. |
crossing.mat |
The crossing matrix. |
GEBVs |
The genomic estimated breeding values. |
tail.p |
The proportion from the population to select. |
G |
The marker genotypes |
y.CV |
The phenotypes for cross-validation. |
G.CV |
The marker genotypes for cross-validation. |
models.CV |
The models for cross-validation. |
frac.train.CV |
The fraction of data to use as training data in cross-validation. |
nCV.iter.CV |
The number of iterations of cross-validation. |
burnIn.CV |
The burn-in number for cross-validation. |
nIter.CV |
The number of iterations for Bayesian models in cross-validation. |
nFold.CV |
The number of folds in k-fold cross-validation. |
nFold.CV.reps |
The number of replications of k-fold cross-validation. |
M |
The marker matrix. |
y.df |
The phenotype data. |
models |
The models to use. |
nIter |
The number of iterations. |
burnIn |
The burn-in rate. |
Predict genetic variance and genetic correlations in multi-parent populations using a deterministic equation.
Description
Predicts the genotypic mean, genetic variance, and usefulness criterion (superior progeny mean) in a set of multi-parent populations using marker effects and a genetic map. If more than two traits are specified, the function will also return predictions of the genetic correlation in the population and the correlated response to selection.
Usage
mppop.predict(
G.in,
y.in,
map.in,
crossing.table,
parents,
n.parents = 4,
tail.p = 0.1,
self.gen = 10,
DH = FALSE,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
n.core = 1,
...
)
mppop_predict2(
M,
y.in,
marker.effects,
map.in,
crossing.table,
parents,
n.parents = 4,
tail.p = 0.1,
self.gen = 10,
DH = FALSE,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
n.core = 1,
...
)
Arguments
G.in |
See |
y.in |
See |
map.in |
See |
crossing.table |
A |
parents |
See |
n.parents |
Integer number of parents per cross. May be 2 or 4. If |
tail.p |
See |
self.gen |
The number of selfing generations in the potential cross. Can be an integer or |
DH |
Indicator if doubled-haploids are to be induced after the number of selfing generations indicated by
|
models |
See |
n.core |
Number of cores for parallelization. Parallelization is supported only on a Linux or Mac OS operating system; if working on a Windows system, the function is executed on a single core. |
... |
Additional arguments to pass depending on the choice of |
M |
A Matrix of marker genotypes of dimensions |
marker.effects |
A data frame of marker effects. The first column should include the marker name and
subsequent columns should include the marker effects. Supercedes |
Details
Predictions are based on the deterministic equations specified by Allier et al. (2019).
In the case of four-way crosses (i.e. 4 parents), the function assumes that the first two parents are mated,
producing a F_1
offspring; then, the next two parents are mated, producing another F_1
offspring.
The two F_1
offspring are then mated and inbreeding or doubled haploid induction (if specified) proceeds
from there. For example, say cross i uses parents P1, P2, P3, and P4. P1 and P2 are first mated,
producing O1; then, P3 and P4 are mated, producing O2; then, O1 and O2 are mated, producing a segregating family.
The mppop.predict
function takes similarly formatted arguments as the pop.predict
function
in the PopVar
package. For the sake of simplicity, we also include the mppop_predict2
function, which
takes arguments in a format more consistent with other genomewide prediction packages/functions.
If you select a model
other than "rrBLUP", you must specify the following additional arguments:
nIter
: Seepop.predict
.burnIn
: Seepop.predict
.
Value
A data.frame
containing predictions of \mu
, V_G
, and \mu_{sp}
for
each trait for each potential multi-parent cross. When multiple traits are provided, the correlated
responses and correlation between all pairs of traits is also returned.
References
Allier, A., L. Moreau, A. Charcosset, S. Teyssèdre, and C. Lehermeier, 2019 Usefulness Criterion and Post-selection Parental Contributions in Multi-parental Crosses: Application to Polygenic Trait Introgression. G3 (Bethesda) 9: 1469–1479. https://doi.org/https://doi.org/10.1534/g3.119.400129
Examples
# Load data
data("think_barley")
# Vector with 8 parents
parents <- sample(y.in_ex$Entry, 8)
# Create a crossing table with four parents per cross
cross_tab <- as.data.frame(t(combn(x = parents, m = 4)))
names(cross_tab) <- c("parent1", "parent2", "parent3", "parent4")
out <- mppop_predict2(M = G.in_ex_mat, y.in = y.in_ex, map.in = map.in_ex,
crossing.table = cross_tab, models = "rrBLUP")
A genome-wide procedure for predicting genetic variance and correlated response in bi-parental breeding populations
Description
pop.predict
uses phenotypic and genotypic data from a set of individuals known as a training population (TP) and a set of candidate parents, which may or may not be included in the TP, to predict the mean (\mu
), genetic variance (V_G), and superior progeny values (\mu
_sp) of the half-diallel, or a defined set of pairwise bi-parental crosses between parents. When multiple traits are provided pop.predict
will also predict the correlated responses and correlation between all pairwise traits. See Mohammadi, Tiede, and Smith (2015) for further details.
NOTE - \code{pop.predict} writes and reads files to disk so it is highly recommended to set your working directory
Usage
pop.predict(
G.in = NULL,
y.in = NULL,
map.in = NULL,
crossing.table = NULL,
parents = NULL,
tail.p = 0.1,
nInd = 200,
map.plot = FALSE,
min.maf = 0.01,
mkr.cutoff = 0.5,
entry.cutoff = 0.5,
remove.dups = TRUE,
impute = "EM",
nSim = 25,
frac.train = 0.6,
nCV.iter = 100,
nFold = NULL,
nFold.reps = 1,
nIter = 12000,
burnIn = 3000,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
return.raw = FALSE,
saveAt = tempdir()
)
Arguments
G.in |
TIP - Set header= |
y.in |
|
map.in |
|
crossing.table |
Optional |
parents |
Optional |
tail.p |
Optional |
nInd |
Optional |
map.plot |
Optional |
min.maf |
Optional |
mkr.cutoff |
Optional |
entry.cutoff |
Optional |
remove.dups |
Optional |
impute |
Options include |
nSim |
Optional |
frac.train |
Optional |
nCV.iter |
Optional |
nFold |
Optional |
nFold.reps |
Optional |
nIter , burnIn |
Optional |
models |
Optional |
return.raw |
Optional |
saveAt |
When using models other than "rrBLUP" (i.e. Bayesian models), this is a path and prefix for saving temporary files
the are produced by the |
Details
pop.predict
can be used to predict the mean (\mu
), genetic variance (V_G), superior progeny values (\mu
_sp
), as well as the predicted correlated response and correlations between all pairwise traits. The methodology and procedure to do so has been described in Bernardo (2014) and Mohammadi, Tiede, and K.P. Smith (2015). Users familiar with genome-wide prediction, association mapping, and/or linkage mapping will be familiar with the
required inputs of pop.predict
. G.in
includes all of the entries (taxa) in the TP as well as additional entries to be considered as parent candidates. Entries included in G.in
that do have a phenotype for any or all traits in y.in
are considered TP entries for those respective traits. G.in
is filtered according to min.maf
, mkr.cutoff
, entry.cutoff
, and remove.dups
;
remaining missing marker data is imputed using the EM algorithm (Poland et al., 2012) when possible, and the marker mean otherwise, both implemented in A.mat
. For each trait, the TP (i.e. entries with phenotype) is used to:
Perform CV to select a regression model. NOTE - Using the model with the highest CV accuracy is expected to result in the most accurate marker effect estimates (Bernardo, 2014). This expectation, however, is yet to be empirically validated and the user is encouraged to investigate the various models in order to make an educated decision about which one to ultimately use.
Estimate marker effects using the model resulting in the highest CV accuracy
Models include ridge regression BLUP implemented in mixed.solve
(Endelman, 2011) and BayesA, BayesB, BayesC\pi
, Bayesian lasso (BL), and Bayesian ridge regression (BRR) implemented in BGLR
(de los Compos and Rodriguez, 2014).
Information from the map.in
is then used to simulate chromosomal recombination expected in a recombinant inbred line (i.e. F-infinity) (Broman et al., 2003) population (size=nInd
). A function then converts the recombined chromosomal segments of the generic RIL population to the chromosomal segments of the population's respective parents and GEBVs of the simulated progeny are calculated.
The simulation and conversion process is repeated s times, where s = nSim
, to calculate dispersion statistics for \mu
and V_G; the remainder of the values in the predictions
output are means of the s simulations. During each iteration the correlation (r) and correlated response of each pairwise combination of traits is also calculated and their mean across n simulations is returned.
The correlated response of trait.B when predicting trait.A is the mean of trait.B for the (\mu
_sp
) of trait.A, and vice-versa; a correlated response for the bottom tail.p
and upper 1-tail.p
is returned for each trait.
A dataset \code{\link{think_barley.rda}} is provided as an example of the proper formatting of input files and also for users to become familiar with \code{pop.predict}.
Value
A list
containing:
-
predictions
Alist
of dataframes containing predictions of (\mu
), (V_G), and (\mu
_sp). When multiple traits are provided the correlated responses and correlation between all pairwise traits is also included. More specifically, for a given trait pair the correlated response of the secondary trait with both the high and low superior progeny of the primary trait is returned since the favorable values cannot be known byPopVar
. -
preds.per.sim
If return.raw isTRUE
then adataframe
containing the results of each simulation is returned. This is useful for calculating dispersion statistics for traits not provided in the standardpredictions
dataframe. -
CVs
Adataframe
of CV results for each trait/model combination specified. -
models.chosen
Amatrix
listing the statistical model chosen for each trait. -
markers.removed
Avector
of markers removed during filtering for MAF and missing data. -
entries.removed
Avector
of entries removed during filtering for missing data and duplicate entries.
References
Bernardo, R. 2014. Genomewide Selection of Parental Inbreds: Classes of Loci and Virtual Biparental Populations. Crop Sci. 55:2586-2595. Broman, K. W., H. Wu, S. Sen and G.A. Churchill. 2003. R/qtl: QTL mapping in experimental crosses. Bioinformatics 19:889-890. Endelman, J. B. 2011. Ridge regression and other kernels for genomic selection with R package rrBLUP. Plant Genome 4:250-255. doi: 10.3835/plantgenome2011.08.0024 Gustavo de los Campos and Paulino Perez Rodriguez, (2014). BGLR: Bayesian Generalized Linear Regression. R package version 1.0.3. http://CRAN.R-project.org/package=BGLR Mohammadi M., T. Tiede, and K.P. Smith. 2015. PopVar: A genome-wide procedure for predicting genetic variance and correlated response in bi-parental breeding populations. Crop Sci. \emph{Accepted}. Munoz-Amatriain, M., M. J. Moscou, P. R. Bhat, J. T. Svensson, J. Bartos, P. Suchankova, H. Simkova, T. R. Endo, R. D. Fenton, S. Lonardi, A. M. Castillo, S. Chao, L. Cistue, A. Cuesta-Marcos, K. L. Forrest, M. J. Hayden, P. M. Hayes, R. D. Horsley, K. Makoto, D. Moody, K. Sato, M. P. Valles, B. B. H. Wulff, G. J. Muehlbauer, J. Dolezel, and T. J. Close. 2011 An improved consensus linkage map of barley based on flow-sorted chromosomes and single nucleotide polymorphism markers. Plant Gen. 4:238-249. Poland, J., J. Endelman, J. Dawson, J. Rutkoski, S. Wu, Y. Manes, S. Dreisigacker, J. Crossa, H. Sanches-Villeda, M. Sorrells, and J.-L. Jannink. 2012. Genomic Selection in Wheat Breeding using Genotyping-by-Sequencing. Plant Genome 5:103-113.
Examples
## Not run:
# Load data
data("think_barley")
## The following examples only use the model 'rrBLUP' for the sake of testing. Functions
## BGLR package write temporary files to the disk.
##
## Ex. 1 - Predict a defined set of crosses
## This example uses CV method 1 (see Details of x.val() function)
ex1.out <- pop.predict(G.in = G.in_ex, y.in = y.in_ex,
map.in = map.in_ex, crossing.table = cross.tab_ex,
nSim=5, nCV.iter=2, models = "rrBLUP")
ex1.out$predictions ## Predicted parameters
ex1.out$CVs ## CV results
## Ex. 2 - Use only rrBLUP and Bayesian lasso (BL) models
ex3.out <- pop.predict(G.in = G.in_ex, y.in = y.in_ex,
map.in = map.in_ex, crossing.table = cross.tab_ex,
models = c("rrBLUP"), nSim=5, nCV.iter=10)
## Ex. 3 - Same as Ex. 3, but return all raw SNP and prediction data for each simulated population
ex4.out <- pop.predict(G.in = G.in_ex, y.in = y.in_ex,
map.in = map.in_ex, crossing.table = cross.tab_ex,
models = c("rrBLUP"), nSim=5, nCV.iter=2, return.raw = TRUE)
## End(Not run)
Predict genetic variance and genetic correlations in bi-parental populations using a deterministic model
Description
Generates predictions of the genetic variance and genetic correlation in bi-parental populations using a set of deterministic equations instead of simulations.
Usage
pop.predict2(
G.in,
y.in,
map.in,
crossing.table,
parents,
tail.p = 0.1,
self.gen = Inf,
DH = FALSE,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
...
)
pop_predict2(
M,
y.in,
marker.effects,
map.in,
crossing.table,
parents,
tail.p = 0.1,
self.gen = Inf,
DH = FALSE,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
...
)
Arguments
G.in |
See |
y.in |
See |
map.in |
See |
crossing.table |
See |
parents |
See |
tail.p |
See |
self.gen |
The number of selfing generations in the potential cross. Can be an integer or |
DH |
Indicator if doubled-haploids are to be induced after the number of selfing generations indicated by
|
models |
See |
... |
Additional arguments to pass depending on the choice of |
M |
A Matrix of marker genotypes of dimensions |
marker.effects |
A data frame of marker effects. The first column should include the marker name and
subsequent columns should include the marker effects. Supercedes |
Details
Predictions are based on the deterministic equations specified by Zhong and Jannink (2007), Allier et al. (2019), and Neyhart et al. (2019).
If you select a model
other than "rrBLUP", you must specify the following additional arguments:
nIter
: Seepop.predict
.burnIn
: Seepop.predict
.
Value
A data.frame
containing predictions of \mu
, V_G
, and \mu_{sp}
for
each trait for each potential bi-parental cross. When multiple traits are provided, the correlated
responses and correlation between all pairs of traits is also returned.
Functions
-
pop_predict2()
:
References
Zhong, S., and J.-L. Jannink, 2007 Using quantitative trait loci results to discriminate among crosses on the basis of their progeny mean and variance. Genetics 177: 567–576. https://doi.org/10.1534/ genetics.107.075358
Allier, A., L. Moreau, A. Charcosset, S. Teyssèdre, and C. Lehermeier, 2019 Usefulness Criterion and Post-selection Parental Contributions in Multi-parental Crosses: Application to Polygenic Trait Introgression. G3 9: 1469–1479. doi: 10.1534/g3.119.400129
Neyhart, J.L., A.J. Lorenz, and K.P. Smith, 2019 Multi-trait Improvement by Predicting Genetic Correlations in Breeding Crosses. G3 9: 3153-3165. doi: 10.1534/g3.119.400406
Examples
# Load data
data("think_barley")
# Use example data to make predictions
out <- pop.predict2(G.in = G.in_ex_imputed, y.in = y.in_ex, map.in = map.in_ex,
crossing.table = cross.tab_ex)
# Provide a vector of parents to predict all possible crosses (some parents
# have missing phenotypic data)
out <- pop.predict2(G.in = G.in_ex_imputed, y.in = y.in_ex, map.in = map.in_ex,
parents = y.in_ex$Entry[1:5])
# Make predictions for 5 crosses with various levels of inbreeding
out_list <- lapply(X = 1:10, FUN = function(self.gen) {
out <- pop.predict2(G.in = G.in_ex_imputed, y.in = y.in_ex, map.in = map.in_ex,
crossing.table = cross.tab_ex[1:5,], self.gen = self.gen)
out$self.gen <- self.gen
out })
# Plot predictions of grain yield genetic variance over levels of inbreeding
dat <- do.call("rbind", lapply(out_list, subset, trait == "Yield"))
plot(pred_varG ~ self.gen, data = dat, type = "b",
subset = parent1 == parent1[1] & parent2 == parent2[1])
# Load data
data("think_barley")
# Use example data to make predictions
out <- pop_predict2(M = G.in_ex_mat, y.in = y.in_ex, map.in = map.in_ex,
crossing.table = cross.tab_ex)
# Provide a vector of parents to predict all possible crosses (some parents
# have missing phenotypic data)
out <- pop_predict2(M = G.in_ex_mat, y.in = y.in_ex, map.in = map.in_ex,
parents = y.in_ex$Entry[1:10])
An example barley dataset
Description
A sample dataset, previously described in Sallam et al. (2014) is provided as an example of the proper formatting of input files and also for users to become familiar with PopVar
; the think_barley
dataset is useful in demonstrating both pop.predict
and x.val
. Note that a number of entries are missing data for one or both traits,
which is representative of a real breeding scenario where phenotypic data may not be available for all parent candidates.
Format
The names of the example files are:
- G.in_ex
A set of 245 barley lines genotyped with 742 SNP markers
- G.in_ex_mat
A n x p matrix of n = 245 barley lines genotyped with p = 742 SNP markers
- G.in_ex_imputed
A n x p matrix of n = 245 barley lines and p = 742 imputed SNP marker genotypes
- y.in_ex
Phenotypes of four traits for a portion of the 245 barley lines, Fusarium head blight (FHB), deoxynivalenol (DON) in ppm, grain yield in bushels/acre, and plant height in cm.
- map.in_ex
Genetic map (i.e. chromosome assignment and genetic distance (cM) between markers) of the 742 SNP markers based on Munoz-Amatriain et al., 2011
- cross.tab_ex
A table of user-defined crosses
References
Sallam, A.H., J.B. Endelman, J-L. Jannink, and K.P. Smith. 2015. Assessing Genomic Selection Prediction Accuracy in a Dynamic Barley Breeding Population. Plant Gen. 8(1)
Estimate genome-wide prediction accuracy using cross-validation
Description
x.val
performs cross-validation (CV) to estimate the accuracy of genome-wide prediction (otherwise known as genomic selection) for a specific training population (TP), i.e. a set of individuals for which phenotypic and genotypic data is available. Cross-validation can be conducted via one of two methods within x.val
, see Details
for more information.
NOTE - \code{x.val}, specifically \code{\link[BGLR]{BGLR}} writes and reads files to disk so it is highly recommended to set your working directory
Usage
x.val(
G.in = NULL,
y.in = NULL,
min.maf = 0.01,
mkr.cutoff = 0.5,
entry.cutoff = 0.5,
remove.dups = TRUE,
impute = "EM",
frac.train = 0.6,
nCV.iter = 100,
nFold = NULL,
nFold.reps = 1,
return.estimates = FALSE,
CV.burnIn = 750,
CV.nIter = 1500,
models = c("rrBLUP", "BayesA", "BayesB", "BayesC", "BL", "BRR"),
saveAt = tempdir()
)
Arguments
G.in |
TIP - Set header= |
y.in |
|
min.maf |
Optional |
mkr.cutoff |
Optional |
entry.cutoff |
Optional |
remove.dups |
Optional |
impute |
Options include |
frac.train |
Optional |
nCV.iter |
Optional |
nFold |
Optional |
nFold.reps |
Optional |
return.estimates |
Optional |
CV.burnIn |
Optional |
CV.nIter |
Optional |
models |
Optional |
saveAt |
When using models other than "rrBLUP" (i.e. Bayesian models), this is a path and prefix for saving temporary files
the are produced by the |
Details
Two CV methods are available within PopVar
:
-
CV method 1
: During each iteration a training (i.e. model training) set will be randomly sampled from the TP of sizeN*(frac.train)
, where N is the size of the TP, and the remainder of the TP is assigned to the validation set. The accuracies of individual models are expressed as average Pearson's correlation coefficient (r) between the genome estimated breeding value (GEBV) and observed phenotypic values in the validation set across allnCV.iter
iterations. Due to its amendibility to various TP sizes, CV method 1 is the default CV method inpop.predict
. -
CV method 2
:nFold
independent validation sets are sampled from the TP and predicted by the remainder. For example, ifnFold = 10
the TP will be split into 10 equal sets, each containing1/10
-th of the TP, which will be predicted by the remaining9/10
-ths of the TP. The accuracies of individual models are expressed as the average (r) between the GEBV and observed phenotypic values in the validation set across allnFold
folds. The process can be repeatednFold.reps
times withnFold
new independent sets being sampled each replication, in which case the reported prediction accuracies are averages across all folds and replications.
Value
A list containing:
-
CVs
Adataframe
of CV results for each trait/model combination specified If
return.estimates
isTRUE
the additional items will be returned:-
models.used
Alist
of the models chosen to estimate marker effects for each trait -
mkr.effects
Avector
of marker effect estimates for each trait generated by the respective prediction model used -
betas
Alist
of beta values for each trait generated by the respective prediction model used
-
Examples
## The following examples only use the model 'rrBLUP' for the sake of testing. Functions
## BGLR package write temporary files to the disk.
## CV using method 1 with 25 iterations
CV.mthd1 <- x.val(G.in = G.in_ex, y.in = y.in_ex, nCV.iter = 25, models = "rrBLUP")
CV.mthd1$CVs
## CV using method 2 with 5 folds and 3 replications
x.val(G.in = G.in_ex, y.in = y.in_ex, nFold = 5, nFold.reps = 3, models = "rrBLUP")