Title: | Random-Effect Multiple QTL Mapping in Autopolyploids |
---|---|
Description: | Performs random-effect multiple interval mapping (REMIM) in full-sib families of autopolyploid species based on restricted maximum likelihood (REML) estimation and score statistics, as described in Pereira et al. (2020) <doi:10.1534/genetics.120.303080>. |
Authors: | Guilherme da Silva Pereira [aut] , Marcelo Mollinari [ctb] , Gabriel de Siqueira Gesteira [ctb, cre] , Zhao-Bang Zeng [ctb] , Long Qu [ctb] (R code for variance component tests using score statistics in R/varComp.R), Giovanny Covarrubias-Pazaran [ctb] (C code for fitting mixed models with REML estimation in src/MNR.cpp) |
Maintainer: | Gabriel de Siqueira Gesteira <[email protected]> |
License: | GPL-3 |
Version: | 0.2.4 |
Built: | 2024-11-14 23:29:31 UTC |
Source: | https://github.com/gabrielgesteira/qtlpoly |
A dataset of the B2721 population which derived from a cross between two tetraploid potato varieties: Atlantic × B1829-5.
B2721
B2721
An object of class mappoly.data
from the package mappoly.
Marcelo Mollinari, [email protected]
Mollinari M, Garcia AAF (2019) Linkage analysis and haplotype phasing in experimental autopolyploid populations with high ploidy level using hidden Markov models, G3: Genes|Genomes|Genetics 9 (10): 3297-3314. doi:10.1534/g3.119.400378
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Pereira GS, Mollinari M, Schumann MJ, Clough ME, Zeng ZB, Yencho C (2021) The recombination landscape and multiple QTL mapping in a Solanum tuberosum cv. ‘Atlantic’-derived F_1 population. Heredity. doi:10.1038/s41437-021-00416-x.
library(mappoly) print(B2721)
library(mappoly) print(B2721)
Computes breeding values for each genotyped individual based on multiple QTL models
breeding_values(data, fitted) ## S3 method for class 'qtlpoly.bvalues' plot(x, pheno.col = NULL, ...)
breeding_values(data, fitted) ## S3 method for class 'qtlpoly.bvalues' plot(x, pheno.col = NULL, ...)
data |
an object of class |
fitted |
an object of class |
x |
an object of class |
pheno.col |
a numeric vector with the phenotype column numbers to be plotted; if |
... |
currently ignored |
An object of class qtlpoly.bvalues
which is a list of results
for each trait containing the following components:
pheno.col |
a phenotype column number. |
y.hat |
a column matrix of breeding value for each individual. |
A ggplot2 histogram with the distribution of breeding values.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) #5,7 data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data = data, model = remim.mod, probs = "joint", polygenes = "none") # Predict genotypic values y.hat = breeding_values(data = data, fitted = fitted.mod) plot(y.hat)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) #5,7 data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data = data, model = remim.mod, probs = "joint", polygenes = "none") # Predict genotypic values y.hat = breeding_values(data = data, fitted = fitted.mod) plot(y.hat)
Performs interval mapping using the single-QTL, fixed-effect model proposed by Hackett et al. (2001).
feim( data = data, pheno.col = NULL, w.size = 15, sig.lod = 7, d.sint = 1.5, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.feim' print(x, pheno.col = NULL, sint = NULL, ...)
feim( data = data, pheno.col = NULL, w.size = 15, sig.lod = 7, d.sint = 1.5, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.feim' print(x, pheno.col = NULL, sint = NULL, ...)
data |
an object of class |
pheno.col |
a numeric vector with the phenotype columns to be analyzed; if |
w.size |
a number representing the window size (in centiMorgans) to be avoided on either side of QTL already in the model when looking for a new QTL, e.g. 15 (default). |
sig.lod |
the vector of desired significance LOD thresholds (usually permutation-based) for declaring a QTL for each trait, e.g. 5 (default); if a single value is provided, the same LOD threshold will be applied to all traits. |
d.sint |
a |
plot |
a suffix for the file's name containing plots of every algorithm step, e.g. "remim" (default); if |
verbose |
if |
x |
an object of class |
sint |
whether |
... |
currently ignored |
An object of class qtlpoly.feim
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
LRT |
a vector containing LRT values. |
LOD |
a vector containing LOD scores. |
AdjR2 |
a vector containing adjusted |
qtls |
a data frame with information from the mapped QTL. |
lower |
a data frame with information from the lower support interval of mapped QTL. |
upper |
a data frame with information from the upper support interval of mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Hackett CA, Bradshaw JE, McNicol JW (2001) Interval mapping of quantitative trait loci in autotetraploid species, Genetics 159: 1819-1832.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 5) # Perform remim feim.mod = feim(data = data, sig.lod = 7)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 5) # Perform remim feim.mod = feim(data = data, sig.lod = 7)
Fits alternative multiple QTL models by performing variance component estimation using REML.
fit_model( data, model, probs = "joint", polygenes = "none", keep = TRUE, verbose = TRUE, pheno.col = NULL ) ## S3 method for class 'qtlpoly.fitted' summary(object, pheno.col = NULL, ...)
fit_model( data, model, probs = "joint", polygenes = "none", keep = TRUE, verbose = TRUE, pheno.col = NULL ) ## S3 method for class 'qtlpoly.fitted' summary(object, pheno.col = NULL, ...)
data |
an object of class |
model |
an object of class |
probs |
a character string indicating if either |
polygenes |
a character string indicating if either |
keep |
if |
verbose |
if |
pheno.col |
a numeric vector with the phenotype column numbers to be summarized; if |
object |
an object of class |
... |
currently ignored |
An object of class qtlpoly.fitted
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
fitted |
a sommer object of class |
qtls |
a data frame with information from the mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Covarrubias-Pazaran G (2016) Genome-assisted prediction of quantitative traits using the R package sommer. PLoS ONE 11 (6): 1–15. doi:10.1371/journal.pone.0156744.
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data=data, model=remim.mod, probs="joint", polygenes="none")
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data=data, model=remim.mod, probs="joint", polygenes="none")
Fits alternative multiple QTL models by performing variance component estimation using REML.
fit_model2( data, model, probs = "joint", polygenes = "none", keep = TRUE, verbose = TRUE, pheno.col = NULL )
fit_model2( data, model, probs = "joint", polygenes = "none", keep = TRUE, verbose = TRUE, pheno.col = NULL )
data |
an object of class |
model |
an object of class |
probs |
a character string indicating if either |
polygenes |
a character string indicating if either |
keep |
if |
verbose |
if |
pheno.col |
a numeric vector with the phenotype column numbers to be summarized; if |
An object of class qtlpoly.fitted
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
fitted |
a sommer object of class |
qtls |
a data frame with information from the mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Covarrubias-Pazaran G (2016) Genome-assisted prediction of quantitative traits using the R package sommer. PLoS ONE 11 (6): 1–15. doi:10.1371/journal.pone.0156744.
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data=data, model=remim.mod, probs="joint", polygenes="none")
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data=data, model=remim.mod, probs="joint", polygenes="none")
A dataset of a hypothetical autohexaploid full-sib population containing three homology groups
hexafake
hexafake
An object of class mappoly.data
which contains a
list with the following components:
ploidy level = 6
number individuals = 300
total number of markers = 1500
the names of the individuals
the names of the markers
a vector containing the dosage in
parent P for all n.mrk
markers
a vector containing the dosage in
parent Q for all n.mrk
markers
a vector indicating the chromosome each marker belongs. Zero indicates that the marker was not assigned to any chromosome
Physical position of the markers into the sequence
a matrix containing the dosage for each markers (rows)
for each individual (columns). Missing data are represented by
ploidy_level + 1 = 7
There are no phenotypes in this simulation
There are no phenotypes in this simulation
vector containing p-values for all markers associated to the chi-square test for the expected segregation patterns under Mendelian segregation
Marcelo Mollinari, [email protected]
Mollinari M, Garcia AAF (2019) Linkage analysis and haplotype phasing in experimental autopolyploid populations with high ploidy level using hidden Markov models, G3: Genes|Genomes|Genetics 9 (10): 3297-3314. doi:10.1534/g3.119.400378
library(mappoly) plot(hexafake)
library(mappoly) plot(hexafake)
A real autotetraploid potato map containing 12 homology groups from a tetraploid potato full-sib family (Atlantic x B1829-5).
maps4x
maps4x
An object of class "mappoly.map"
from the package mappoly, which is a list of 12 linkage groups (LGs)
Marcelo Mollinari, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2019) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Mollinari M, Garcia AAF (2019) Linkage analysis and haplotype phasing in experimental autopolyploid populations with high ploidy level using hidden Markov models, G3: Genes|Genomes|Genetics 9 (10): 3297-3314. doi:10.1534/g3.119.400378
Pereira GS, Mollinari M, Schumann MJ, Clough ME, Zeng ZB, Yencho C (2021) The recombination landscape and multiple QTL mapping in a Solanum tuberosum cv. ‘Atlantic’-derived F_1 population. Heredity. doi:10.1038/s41437-021-00416-x.
library(mappoly) plot_map_list(maps4x)
library(mappoly) plot_map_list(maps4x)
A simulated map containing three homology groups of a hypotetical cross between two autohexaploid individuals.
maps6x
maps6x
An object of class "mappoly.map"
from the package mappoly, which is a list of three linkage groups (LGs):
538 markers distributed along 112.2 cM
329 markers distributed along 54.6 cM
443 markers distributed along 98.2 cM
Marcelo Mollinari, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2019) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Mollinari M, Garcia AAF (2019) Linkage analysis and haplotype phasing in experimental autopolyploid populations with high ploidy level using hidden Markov models, G3: Genes|Genomes|Genetics 9 (10): 3297-3314. doi:10.1534/g3.119.400378
library(mappoly) plot_map_list(maps6x)
library(mappoly) plot_map_list(maps6x)
Adds or removes QTL manually from a given model.
modify_qtl( model, pheno.col = NULL, add.qtl = NULL, drop.qtl = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.modify' print(x, pheno.col = NULL, ...)
modify_qtl( model, pheno.col = NULL, add.qtl = NULL, drop.qtl = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.modify' print(x, pheno.col = NULL, ...)
model |
an object of class |
pheno.col |
a phenotype column number whose model will be modified or printed. |
add.qtl |
a marker position number to be added. |
drop.qtl |
a marker position number to be removed. |
verbose |
if |
x |
an object of class |
... |
currently ignored |
An object of class qtlpoly.modify
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Modify model modified.mod = modify_qtl(model = remim.mod, pheno.col = 1, drop.qtl = 18)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Modify model modified.mod = modify_qtl(model = remim.mod, pheno.col = 1, drop.qtl = 18)
Creates a null model (with no QTL) for each trait.
null_model( data, offset.data = NULL, pheno.col = NULL, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.null' print(x, pheno.col = NULL, ...) ## S3 method for class 'qtlpoly.null' print(x, pheno.col = NULL, ...)
null_model( data, offset.data = NULL, pheno.col = NULL, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.null' print(x, pheno.col = NULL, ...) ## S3 method for class 'qtlpoly.null' print(x, pheno.col = NULL, ...)
data |
an object of class |
offset.data |
a data frame with the same dimensions of |
pheno.col |
a numeric vector with the phenotype columns to be analyzed; if |
n.clusters |
number of parallel processes to spawn. |
plot |
a suffix for the file's name containing simple plots of every QTL search round, e.g. "null" (default); if |
verbose |
if |
x |
an object of class |
... |
currently ignored |
An object of class qtlpoly.null
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL ( |
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null models null.mod = null_model(data = data, pheno.col = 1, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null models null.mod = null_model(data = data, pheno.col = 1, n.clusters = 1)
Creates a null model (with no QTL) for each trait.
null_model2( data, offset.data = NULL, pheno.col = NULL, n.clusters = NULL, plot = NULL, verbose = TRUE )
null_model2( data, offset.data = NULL, pheno.col = NULL, n.clusters = NULL, plot = NULL, verbose = TRUE )
data |
an object of class |
offset.data |
a data frame with the same dimensions of |
pheno.col |
a numeric vector with the phenotype columns to be analyzed; if |
n.clusters |
number of parallel processes to spawn. |
plot |
a suffix for the file's name containing simple plots of every QTL search round, e.g. "null" (default); if |
verbose |
if |
An object of class qtlpoly.null
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL ( |
Guilherme da Silva Pereira, [email protected], Gabriel de Siqueira Gesteira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null models null.mod = null_model(data = data, pheno.col = 1, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null models null.mod = null_model(data = data, pheno.col = 1, n.clusters = 1)
Tests each QTL at a time and updates its position (if it changes) or drops the QTL (if non-significant).
optimize_qtl( data, offset.data = NULL, model, sig.bwd = 0.05, score.null = NULL, polygenes = FALSE, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.optimize' print(x, pheno.col = NULL, ...)
optimize_qtl( data, offset.data = NULL, model, sig.bwd = 0.05, score.null = NULL, polygenes = FALSE, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.optimize' print(x, pheno.col = NULL, ...)
data |
an object of class |
offset.data |
a data frame with the same dimensions of |
model |
an object of class |
sig.bwd |
the desired score-based p-value threshold for backward elimination, e.g. 0.0001 (default). |
score.null |
an object of class |
polygenes |
if |
n.clusters |
number of parallel processes to spawn. |
plot |
a suffix for the file's name containing plots of every QTL optimization round, e.g. "optimize" (default); if |
verbose |
if |
x |
an object of class |
pheno.col |
a numeric vector with the phenotype columns to be printed; if |
... |
currently ignored |
An object of class qtlpoly.optimize
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
Zou F, Fine JP, Hu J, Lin DY (2004) An efficient resampling method for assessing genome-wide statistical significance in mapping quantitative trait loci. Genetics 168 (4): 2307-16. doi:10.1534/genetics.104.031427
read_data
, null_model
, search_qtl
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null model null.mod = null_model(data = data, pheno.col = 1,n.clusters = 1) # Perform forward search search.mod = search_qtl(data = data, model = null.mod, w.size = 15, sig.fwd = 0.01, n.clusters = 1) # Optimize model optimize.mod = optimize_qtl(data = data, model = search.mod, sig.bwd = 0.0001, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null model null.mod = null_model(data = data, pheno.col = 1,n.clusters = 1) # Perform forward search search.mod = search_qtl(data = data, model = null.mod, w.size = 15, sig.fwd = 0.01, n.clusters = 1) # Optimize model optimize.mod = optimize_qtl(data = data, model = search.mod, sig.bwd = 0.0001, n.clusters = 1)
Stores maximum LOD scores for a number of permutations of given phenotypes.
permutations( data, offset.data = NULL, pheno.col = NULL, n.sim = 1000, probs = c(0.9, 0.95), n.clusters = NULL, seed = 123, verbose = TRUE ) ## S3 method for class 'qtlpoly.perm' print(x, pheno.col = NULL, probs = c(0.9, 0.95), ...) ## S3 method for class 'qtlpoly.perm' plot(x, pheno.col = NULL, probs = c(0.9, 0.95), ...)
permutations( data, offset.data = NULL, pheno.col = NULL, n.sim = 1000, probs = c(0.9, 0.95), n.clusters = NULL, seed = 123, verbose = TRUE ) ## S3 method for class 'qtlpoly.perm' print(x, pheno.col = NULL, probs = c(0.9, 0.95), ...) ## S3 method for class 'qtlpoly.perm' plot(x, pheno.col = NULL, probs = c(0.9, 0.95), ...)
data |
an object of class |
offset.data |
a subset of the data object to be used in permutation calculations. |
pheno.col |
a numeric vector with the phenotype columns to be analyzed; if |
n.sim |
a number of simulations, e.g. 1000 (default). |
probs |
a vector of probability values in [0, 1] representing the quantiles, e.g. c(0.90, 0.95) for the 90% and 95% quantiles. |
n.clusters |
a number of parallel processes to spawn. |
seed |
an integer for the |
verbose |
if |
x |
an object of class |
... |
currently ignored |
An object of class qtlpoly.perm
which contains a list of results
for each trait with the maximum LOD score per permutation.
LOD score thresholds for given quantiles for each trait.
A ggplot2 histogram with the distribution of ordered maximum LOD scores and thresholds for given quantiles for each trait.
Guilherme da Silva Pereira, [email protected]
Churchill GA, Doerge RW (1994) Empirical threshold values for quantitative trait mapping, Genetics 138: 963-971.
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Perform permutations perm = permutations(data = data, pheno.col = 1, n.sim = 10, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Perform permutations perm = permutations(data = data, pheno.col = 1, n.sim = 10, n.clusters = 1)
A subset of phenotypes from a tetraploid potato full-sib family (Atlantic x B1829-5).
pheno4x
pheno4x
A data frame of phenotypes with 156 named individuals in rows and three named phenotypes in columns, which are:
Foliage maturity evaluated in 2007.
Foliage maturity evaluated in 2008.
Foliage maturity evaluated in 2014.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Pereira GS, Mollinari M, Schumann MJ, Clough ME, Zeng ZB, Yencho C (2021) The recombination landscape and multiple QTL mapping in a Solanum tuberosum cv. ‘Atlantic’-derived F_1 population. Heredity. doi:10.1038/s41437-021-00416-x.
head(pheno4x)
head(pheno4x)
A simulated data set of phenotypes for a hipotetical autohexaploid species map.
pheno6x
pheno6x
A data frame of phenotypes with 300 named individuals in rows and three named phenotypes in columns, which are:
3 QTLs, with heritabilities of 0.20 (LG 1 at 32.03 cM), 0.15 (LG 1 at 95.02 cM) and 0.30 (LG 2 at 40.01 cM).
1 QTL, with heritability of 0.15 (LG 3 at 34.51 cM).
no QTLs.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2019) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
head(pheno6x)
head(pheno6x)
Plots profiled logarithm of score-based P-values (LOP) from individual or combined traits.
plot_profile( data = data, model = model, pheno.col = NULL, sup.int = FALSE, main = NULL, legend = "bottom", ylim = NULL, grid = FALSE )
plot_profile( data = data, model = model, pheno.col = NULL, sup.int = FALSE, main = NULL, legend = "bottom", ylim = NULL, grid = FALSE )
data |
an object of class |
model |
an object of class |
pheno.col |
a numeric vector with the phenotype column numbers to be plotted; if |
sup.int |
if |
main |
a character string with the main title; if |
legend |
legend position (either "bottom", "top", "left" or "right"); if |
ylim |
a numeric value pair supplying the limits of y-axis, e.g. c(0,10); if |
grid |
if |
A ggplot2 with the LOP profiles for each trait.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Plot profile plot_profile(data = data, model = remim.mod, grid = FALSE)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Plot profile plot_profile(data = data, model = remim.mod, grid = FALSE)
Creates a plot where dot sizes and colors represent the QTLs heritabilities and their p-values, respectively.
plot_qtl( data = data, model = model, fitted = fitted, pheno.col = NULL, main = NULL, drop.pheno = TRUE, drop.lgs = TRUE )
plot_qtl( data = data, model = model, fitted = fitted, pheno.col = NULL, main = NULL, drop.pheno = TRUE, drop.lgs = TRUE )
data |
an object of class |
model |
an object of class |
fitted |
an object of class |
pheno.col |
the desired phenotype column numbers to be plotted. The order here specifies the order of plotting (from top to bottom.) |
main |
plot title; if |
drop.pheno |
if |
drop.lgs |
if |
A ggplot2 with dots representing the QTLs.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data, remim.mod, probs="joint", polygenes="none") # Plot QTL plot_qtl(data, remim.mod, fitted.mod)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data, remim.mod, probs="joint", polygenes="none") # Plot QTL plot_qtl(data, remim.mod, fitted.mod)
Creates a plot where colored bars represent the support intervals for QTL peaks (black dots).
plot_sint(data, model, pheno.col = NULL, main = NULL, drop = FALSE)
plot_sint(data, model, pheno.col = NULL, main = NULL, drop = FALSE)
data |
an object of class |
model |
an object of class |
pheno.col |
a numeric vector with the phenotype column numbers to be plotted; if |
main |
a character string with the main title; if |
drop |
if |
A ggplot2 with QTL bars for each linkage group.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Plot support intervals plot_sint(data = data, model = remim.mod)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Plot support intervals plot_sint(data = data, model = remim.mod)
Generates the score-based genome-wide profile conditional to the selected QTL.
profile_qtl( data, model, d.sint = 1.5, polygenes = FALSE, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.profile' print(x, pheno.col = NULL, sint = NULL, ...)
profile_qtl( data, model, d.sint = 1.5, polygenes = FALSE, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.profile' print(x, pheno.col = NULL, sint = NULL, ...)
data |
an object of class |
model |
an object of class |
d.sint |
a |
polygenes |
if |
n.clusters |
number of parallel processes to spawn. |
plot |
a suffix for the file's name containing plots of every QTL profiling round, e.g. "profile" (default); if |
verbose |
if |
x |
an object of class |
pheno.col |
a numeric vector with the phenotype column numbers to be plotted; if |
sint |
whether |
... |
currently ignored |
An object of class qtlpoly.profile
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL. |
lower |
a data frame with information from the lower support interval of mapped QTL. |
upper |
a data frame with information from the upper support interval of mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null model null.mod = null_model(data, pheno.col = 1, n.clusters = 1) # Perform forward search search.mod = search_qtl(data = data, model = null.mod, w.size = 15, sig.fwd = 0.01, n.clusters = 1) # Optimize model optimize.mod = optimize_qtl(data = data, model = search.mod, sig.bwd = 0.0001, n.clusters = 1) # Profile model profile.mod = profile_qtl(data = data, model = optimize.mod, d.sint = 1.5, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null model null.mod = null_model(data, pheno.col = 1, n.clusters = 1) # Perform forward search search.mod = search_qtl(data = data, model = null.mod, w.size = 15, sig.fwd = 0.01, n.clusters = 1) # Optimize model optimize.mod = optimize_qtl(data = data, model = search.mod, sig.bwd = 0.0001, n.clusters = 1) # Profile model profile.mod = profile_qtl(data = data, model = optimize.mod, d.sint = 1.5, n.clusters = 1)
Computes allele specific and allele combination (within-parent) heritable effects from multiple QTL models.
qtl_effects(ploidy = 6, fitted, pheno.col = NULL, verbose = TRUE) ## S3 method for class 'qtlpoly.effects' plot(x, pheno.col = NULL, p1 = "P1", p2 = "P2", ...)
qtl_effects(ploidy = 6, fitted, pheno.col = NULL, verbose = TRUE) ## S3 method for class 'qtlpoly.effects' plot(x, pheno.col = NULL, p1 = "P1", p2 = "P2", ...)
ploidy |
a numeric value of ploidy level of the cross (currently, only 2, 4 or 6). |
fitted |
a fitted multiple QTL model of class |
pheno.col |
a numeric vector with the phenotype column numbers to be plotted; if |
verbose |
if |
x |
an object of class |
p1 |
a character string with the first parent name, e.g. |
p2 |
a character string with the second parent name, e.g. |
... |
currently ignored |
An object of class qtlpoly.effects
which is a list of results
for each containing the following components:
pheno.col |
a phenotype column number. |
y.hat |
a vector with the predicted values. |
A ggplot2 barplot with parental allele and allele combination effects.
Guilherme da Silva Pereira, [email protected], with modifications by Gabriel Gesteira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Kempthorne O (1955) The correlation between relatives in a simple autotetraploid population, Genetics 40: 168-174.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data, model=remim.mod, probs="joint", polygenes="none") # Estimate effects est.effects = qtl_effects(ploidy = 4, fitted = fitted.mod, pheno.col = 1) # Plot results plot(est.effects)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1) # Fit model fitted.mod = fit_model(data, model=remim.mod, probs="joint", polygenes="none") # Estimate effects est.effects = qtl_effects(ploidy = 4, fitted = fitted.mod, pheno.col = 1) # Plot results plot(est.effects)
Reads files in specific formats and creates a qtlpoly.data
object to be used in subsequent analyses.
read_data( ploidy = 6, geno.prob, geno.dose = NULL, double.reduction = FALSE, pheno, weights = NULL, step = 1, verbose = TRUE ) ## S3 method for class 'qtlpoly.data' print(x, detailed = FALSE, ...)
read_data( ploidy = 6, geno.prob, geno.dose = NULL, double.reduction = FALSE, pheno, weights = NULL, step = 1, verbose = TRUE ) ## S3 method for class 'qtlpoly.data' print(x, detailed = FALSE, ...)
ploidy |
a numeric value of ploidy level of the cross. |
geno.prob |
an object of class |
geno.dose |
an object of class |
double.reduction |
if |
pheno |
a data frame of phenotypes (columns) with individual names (rows) identical to individual names in |
weights |
a data frame of phenotype weights (columns) with individual names (rows) identical to individual names in |
step |
a numeric value of step size (in centiMorgans) where tests will be performed, e.g. 1 (default); if |
verbose |
if |
x |
an object of class |
detailed |
if |
... |
currently ignored |
An object of class qtlpoly.data
which is a list containing the following components:
ploidy |
a scalar with ploidy level. |
nlgs |
a scalar with the number of linkage groups. |
nind |
a scalar with the number of individuals. |
nmrk |
a scalar with the number of marker positions. |
nphe |
a scalar with the number of phenotypes. |
lgs.size |
a vector with linkage group sizes. |
cum.size |
a vector with cumulative linkage group sizes. |
lgs.nmrk |
a vector with number of marker positions per linkage group. |
cum.nmrk |
a vector with cumulative number of marker positions per linkage group. |
lgs |
a list with selected marker positions per linkage group. |
lgs.all |
a list with all marker positions per linkage group. |
step |
a scalar with the step size. |
pheno |
a data frame with phenotypes. |
G |
a list of relationship matrices for each marker position. |
Z |
a list of conditional probability matrices for each marker position for genotypes. |
X |
a list of conditional probability matrices for each marker position for alleles. |
Pi |
a matrix of identical-by-descent shared alleles among genotypes. |
Guilherme da Silva Pereira, [email protected], with minor updates by Gabriel de Siqueira Gesteira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1)
Reads files in specific formats and creates a qtlpoly.data
object to be used in subsequent analyses.
read_data2( ploidy = 6, geno.prob, geno.dose = NULL, type = c("genome", "mds", "custom"), double.reduction = FALSE, pheno, weights = NULL, step = 1, verbose = TRUE )
read_data2( ploidy = 6, geno.prob, geno.dose = NULL, type = c("genome", "mds", "custom"), double.reduction = FALSE, pheno, weights = NULL, step = 1, verbose = TRUE )
ploidy |
a numeric value of ploidy level of the cross. |
geno.prob |
an object of class |
geno.dose |
an object of class |
type |
either "genome", "mds", or "custom" from the |
double.reduction |
if |
pheno |
a data frame of phenotypes (columns) with individual names (rows) identical to individual names in |
weights |
a data frame of phenotype weights (columns) with individual names (rows) identical to individual names in |
step |
a numeric value of step size (in centiMorgans) where tests will be performed, e.g. 1 (default); if |
verbose |
if |
An object of class qtlpoly.data
which is a list containing the following components:
ploidy |
a scalar with ploidy level. |
nlgs |
a scalar with the number of linkage groups. |
nind |
a scalar with the number of individuals. |
nmrk |
a scalar with the number of marker positions. |
nphe |
a scalar with the number of phenotypes. |
lgs.size |
a vector with linkage group sizes. |
cum.size |
a vector with cumulative linkage group sizes. |
lgs.nmrk |
a vector with number of marker positions per linkage group. |
cum.nmrk |
a vector with cumulative number of marker positions per linkage group. |
lgs |
a list with selected marker positions per linkage group. |
lgs.all |
a list with all marker positions per linkage group. |
step |
a scalar with the step size. |
pheno |
a data frame with phenotypes. |
G |
a list of relationship matrices for each marker position. |
Z |
a list of conditional probability matrices for each marker position for genotypes. |
X |
a list of conditional probability matrices for each marker position for alleles. |
Pi |
a matrix of identical-by-descent shared alleles among genotypes. |
Guilherme da Silva Pereira, [email protected], Gabriel de Siqueira Gesteira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1)
Automatic function that performs REMIM algorithm using score statistics.
remim( data, pheno.col = NULL, w.size = 15, sig.fwd = 0.01, sig.bwd = 1e-04, score.null = NULL, d.sint = 1.5, polygenes = FALSE, n.clusters = NULL, n.rounds = Inf, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.remim' print(x, pheno.col = NULL, sint = NULL, ...)
remim( data, pheno.col = NULL, w.size = 15, sig.fwd = 0.01, sig.bwd = 1e-04, score.null = NULL, d.sint = 1.5, polygenes = FALSE, n.clusters = NULL, n.rounds = Inf, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.remim' print(x, pheno.col = NULL, sint = NULL, ...)
data |
an object of class |
pheno.col |
a numeric vector with the phenotype columns to be analyzed or printed; if |
w.size |
the window size (in centiMorgans) to avoid on either side of QTL already in the model when looking for a new QTL, e.g. 15 (default). |
sig.fwd |
the desired score-based significance level for forward search, e.g. 0.01 (default). |
sig.bwd |
the desired score-based significance level for backward elimination, e.g. 0.001 (default). |
score.null |
an object of class |
d.sint |
a |
polygenes |
if |
n.clusters |
number of parallel processes to spawn. |
n.rounds |
number of search rounds; if |
plot |
a suffix for the file's name containing plots of every algorithm step, e.g. "remim"; if |
verbose |
if |
x |
an object of class |
sint |
whether |
... |
currently ignored |
An object of class qtlpoly.remim
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL. |
lower |
a data frame with information from the lower support interval of mapped QTL. |
upper |
a data frame with information from the upper support interval of mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Kao CH, Zeng ZB, Teasdale RD (1999) Multiple interval mapping for quantitative trait loci. Genetics 152 (3): 1203–16.
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
Zou F, Fine JP, Hu J, Lin DY (2004) An efficient resampling method for assessing genome-wide statistical significance in mapping quantitative trait loci. Genetics 168 (4): 2307-16. doi:10.1534/genetics.104.031427
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1)
Automatic function that performs REMIM algorithm using score statistics.
remim2( data, pheno.col = NULL, w.size = 15, sig.fwd = 0.01, sig.bwd = 1e-04, score.null = NULL, d.sint = 1.5, polygenes = FALSE, n.clusters = NULL, n.rounds = Inf, plot = NULL, verbose = TRUE )
remim2( data, pheno.col = NULL, w.size = 15, sig.fwd = 0.01, sig.bwd = 1e-04, score.null = NULL, d.sint = 1.5, polygenes = FALSE, n.clusters = NULL, n.rounds = Inf, plot = NULL, verbose = TRUE )
data |
an object of class |
pheno.col |
a numeric vector with the phenotype columns to be analyzed or printed; if |
w.size |
the window size (in centiMorgans) to avoid on either side of QTL already in the model when looking for a new QTL, e.g. 15 (default). |
sig.fwd |
the desired score-based significance level for forward search, e.g. 0.01 (default). |
sig.bwd |
the desired score-based significance level for backward elimination, e.g. 0.001 (default). |
score.null |
an object of class |
d.sint |
a |
polygenes |
if |
n.clusters |
number of parallel processes to spawn. |
n.rounds |
number of search rounds; if |
plot |
a suffix for the file's name containing plots of every algorithm step, e.g. "remim"; if |
verbose |
if |
sint |
whether |
An object of class qtlpoly.remim
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL. |
lower |
a data frame with information from the lower support interval of mapped QTL. |
upper |
a data frame with information from the upper support interval of mapped QTL. |
Guilherme da Silva Pereira, [email protected], Getúlio Caixeta Ferreira, [email protected]
Kao CH, Zeng ZB, Teasdale RD (1999) Multiple interval mapping for quantitative trait loci. Genetics 152 (3): 1203–16.
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
Zou F, Fine JP, Hu J, Lin DY (2004) An efficient resampling method for assessing genome-wide statistical significance in mapping quantitative trait loci. Genetics 168 (4): 2307-16. doi:10.1534/genetics.104.031427
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim2(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Search for QTL remim.mod = remim2(data = data, pheno.col = 1, w.size = 15, sig.fwd = 0.0011493379, sig.bwd = 0.0002284465, d.sint = 1.5, n.clusters = 1)
Searches for QTL and adds them one at a time to a multiple random-effect QTL model based on score statistics.
search_qtl( data, offset.data = NULL, model, w.size = 15, sig.fwd = 0.2, score.null = NULL, polygenes = FALSE, n.rounds = Inf, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.search' print(x, pheno.col = NULL, ...)
search_qtl( data, offset.data = NULL, model, w.size = 15, sig.fwd = 0.2, score.null = NULL, polygenes = FALSE, n.rounds = Inf, n.clusters = NULL, plot = NULL, verbose = TRUE ) ## S3 method for class 'qtlpoly.search' print(x, pheno.col = NULL, ...)
data |
an object of class |
offset.data |
a data frame with the same dimensions of |
model |
an object of class |
w.size |
the window size (in cM) to avoid on either side of QTL already in the model when looking for a new QTL. |
sig.fwd |
the desired score-based p-value threshold for forward search, e.g. 0.01 (default). |
score.null |
an object of class |
polygenes |
if |
n.rounds |
number of search rounds; if |
n.clusters |
number of parallel processes to spawn. |
plot |
a suffix for the file's name containing plots of every QTL search round, e.g. "search" (default); if |
verbose |
if |
x |
an object of class |
pheno.col |
a numeric vector with the phenotype column numbers to be printed; if |
... |
currently ignored |
An object of class qtlpoly.search
which contains a list of results
for each trait with the following components:
pheno.col |
a phenotype column number. |
stat |
a vector containing values from score statistics. |
pval |
a vector containing p-values from score statistics. |
qtls |
a data frame with information from the mapped QTL. |
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
Qu L, Guennel T, Marshall SL (2013) Linear score tests for variance components in linear mixed models and applications to genetic association studies. Biometrics 69 (4): 883–92.
Zou F, Fine JP, Hu J, Lin DY (2004) An efficient resampling method for assessing genome-wide statistical significance in mapping quantitative trait loci. Genetics 168 (4): 2307-16. doi:10.1534/genetics.104.031427
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null model null.mod = null_model(data, pheno.col = 1, n.clusters = 1) # Perform forward search search.mod = search_qtl(data, model = null.mod, w.size = 15, sig.fwd = 0.01, n.clusters = 1)
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Build null model null.mod = null_model(data, pheno.col = 1, n.clusters = 1) # Perform forward search search.mod = search_qtl(data, model = null.mod, w.size = 15, sig.fwd = 0.01, n.clusters = 1)
Simulate new phenotypes with a given number of QTL and creates new object with the same structure of class qtlpoly.data
from an existing genetic map.
simulate_qtl( data, mu = 0, h2.qtl = c(0.3, 0.2, 0.1), var.error = 1, linked = FALSE, n.sim = 1000, missing = TRUE, w.size = 20, seed = 123, verbose = TRUE ) ## S3 method for class 'qtlpoly.simul' print(x, detailed = FALSE, ...)
simulate_qtl( data, mu = 0, h2.qtl = c(0.3, 0.2, 0.1), var.error = 1, linked = FALSE, n.sim = 1000, missing = TRUE, w.size = 20, seed = 123, verbose = TRUE ) ## S3 method for class 'qtlpoly.simul' print(x, detailed = FALSE, ...)
data |
an object of class |
mu |
simulated phenotype mean, e.g. 0 (default). |
h2.qtl |
vector with QTL heritabilities, e.g. |
var.error |
simulated error variance, e.g. 1 (default). |
linked |
if |
n.sim |
number of simulations, e.g. 1000 (default). |
missing |
if |
w.size |
the window size (in centiMorgans) between two (linked) QTL, e.g. 20 (default). |
seed |
integer for the |
verbose |
if |
x |
an object of class |
detailed |
if |
... |
currently ignored |
An object of class qtlpoly.sim
which contains a list of results
with the same structure of class qtlpoly.data
.
Guilherme da Silva Pereira, [email protected]
Pereira GS, Gemenet DC, Mollinari M, Olukolu BA, Wood JC, Mosquera V, Gruneberg WJ, Khan A, Buell CR, Yencho GC, Zeng ZB (2020) Multiple QTL mapping in autopolyploids: a random-effect model approach with application in a hexaploid sweetpotato full-sib population, Genetics 215 (3): 579-595. doi:10.1534/genetics.120.303080.
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Simulate new phenotypes sim.dat = simulate_qtl(data = data, n.sim = 1) sim.dat
# Estimate conditional probabilities using mappoly package library(mappoly) library(qtlpoly) genoprob4x = lapply(maps4x[c(5)], calc_genoprob) data = read_data(ploidy = 4, geno.prob = genoprob4x, pheno = pheno4x, step = 1) # Simulate new phenotypes sim.dat = simulate_qtl(data = data, n.sim = 1) sim.dat