phenoRegressor.rrBLUP {GROAN} | R Documentation |
SNP-BLUP or G-BLUP using rrBLUP package
Description
This is a wrapper around rrBLUP
function mixed.solve
.
It can either work with genotypes (in form of a SNP matrix) or with kinships (in form of a covariance
matrix). In the first case the function will implement a SNP-BLUP, in the second a G-BLUP. An error is
returned if both SNPs and covariance matrix are passed.
In rrBLUP terms, genotypes are modeled as random effects (matrix Z), covariances as matrix K, and
extra covariates, if present, as fixed effects (matrix X).
Please note that this function won't work if rrBLUP package is not installed.
Usage
phenoRegressor.rrBLUP(
phenotypes,
genotypes = NULL,
covariances = NULL,
extraCovariates = NULL,
...
)
Arguments
phenotypes |
phenotypes, a numeric array (n x 1), missing values are predicted |
genotypes |
SNP genotypes, one row per phenotype (n), one column per marker (m), values in 0/1/2 for
diploids or 0/1/2/...ploidy for polyploids. Can be NULL if |
covariances |
square matrix (n x n) of covariances. |
extraCovariates |
optional extra covariates set, one row per phenotype (n), one column per covariate (w). If NULL no extra covariates are considered. |
... |
extra parameters are passed to rrBLUP::mixed.solve |
Value
The function returns a list with the following fields:
-
predictions
: an array of (k) predicted phenotypes -
hyperparams
: named vector with the following keys: Vu, Ve, beta, LL -
extradata
: list with information on trained model, coming frommixed.solve
See Also
Other phenoRegressors:
phenoRegressor.BGLR()
,
phenoRegressor.RFR()
,
phenoRegressor.SVR()
,
phenoRegressor.dummy()
,
phenoregressor.BGLR.multikinships()
Examples
## Not run:
#using the GROAN.KI dataset, we regress on the dataset and predict the first ten phenotypes
phenos = GROAN.KI$yield
phenos[1:10] = NA
#calling the regressor with ridge regression BLUP on SNPs and kinship
results.SNP.BLUP = phenoRegressor.rrBLUP(
phenotypes = phenos,
genotypes = GROAN.KI$SNPs,
SE = TRUE, return.Hinv = TRUE #rrBLUP-specific parameters
)
results.G.BLUP = phenoRegressor.rrBLUP(
phenotypes = phenos,
covariances = GROAN.KI$kinship,
SE = TRUE, return.Hinv = TRUE #rrBLUP-specific parameters
)
#examining the predictions
plot(GROAN.KI$yield, results.SNP.BLUP$predictions,
main = '[SNP-BLUP] Train set (black) and test set (red) regressions',
xlab = 'Original phenotypes', ylab = 'Predicted phenotypes')
abline(a=0, b=1)
points(GROAN.KI$yield[1:10], results.SNP.BLUP$predictions[1:10], pch=16, col='red')
plot(GROAN.KI$yield, results.G.BLUP$predictions,
main = '[G-BLUP] Train set (black) and test set (red) regressions',
xlab = 'Original phenotypes', ylab = 'Predicted phenotypes')
abline(a=0, b=1)
points(GROAN.KI$yield[1:10], results.G.BLUP$predictions[1:10], pch=16, col='red')
#printing correlations
correlations = data.frame(
model = 'SNP-BLUP',
test_set_correlations = cor(GROAN.KI$yield[1:10], results.SNP.BLUP$predictions[1:10]),
train_set_correlations = cor(GROAN.KI$yield[-(1:10)], results.SNP.BLUP$predictions[-(1:10)])
)
correlations = rbind(correlations, data.frame(
model = 'G-BLUP',
test_set_correlations = cor(GROAN.KI$yield[1:10], results.G.BLUP$predictions[1:10]),
train_set_correlations = cor(GROAN.KI$yield[-(1:10)], results.G.BLUP$predictions[-(1:10)])
))
print(correlations)
## End(Not run)