Title: | Interval Mapping for Quantitative Trait Loci Underlying Neighbor Effects |
---|---|
Description: | To enable quantitative trait loci mapping of neighbor effects, this package extends a single-marker regression to interval mapping. The theoretical background of the method is described in Sato et al. (2021) <doi:10.1093/g3journal/jkab017>. |
Authors: | Yasuhiro Sato [aut, cre] , Kazuya Takeda [aut], Atsushi J. Nagano [aut] |
Maintainer: | Yasuhiro Sato <[email protected]> |
License: | GPL-3 |
Version: | 1.1.2 |
Built: | 2024-12-23 06:38:41 UTC |
Source: | CRAN |
A function to calculate self QTL effects for all individuals, with given deviation coefficients and conditional genotype probabilities.
calc_neiprob( genoprobs, a2, d2, contrasts = NULL, smap, scale, grouping = rep(1, nrow(smap)), d2sq0 = FALSE )
calc_neiprob( genoprobs, a2, d2, contrasts = NULL, smap, scale, grouping = rep(1, nrow(smap)), d2sq0 = FALSE )
genoprobs |
Conditional genotype probabilities as taken from |
a2 |
A numeric scalar indicating additive deviation. |
d2 |
A numeric scalar indicating dominance deviation. |
contrasts |
A vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial positions along an x-axis and y-axis, respectively. |
scale |
A numeric scalar indicating the maximum spatial distance between a focal individual and neighbors to define neighbor effects. |
grouping |
An integer vector assigning each individual to a group. This argument can be used when |
d2sq0 |
An option to make AB/AB interaction effects zero. |
A numeric matrix containing individuals x marker elements for neighbor QTL effects.
Yasuhiro Sato ([email protected])
A function to calculate the proportion or ratio of phenotypic variation explained (PVE or RVE) by neighbor effects for a series of neighbor distance (s_seq
) using mixed models.
calc_pve( genoprobs, pheno, smap, s_seq, addcovar = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), fig = TRUE, contrasts = NULL )
calc_pve( genoprobs, pheno, smap, s_seq, addcovar = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), fig = TRUE, contrasts = NULL )
genoprobs |
Conditional genotype probabilities as taken from |
pheno |
A vector of individual phenotypes. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial positions along an x-axis and y-axis, respectively. |
s_seq |
A numeric vector including a set of the maximum spatial distance between a focal individual and neighbors to define neighbor effects. A scalar is also allowed. |
addcovar |
An optional matrix including additional non-genetic covariates. It contains no. of individuals x no. of covariates. |
grouping |
An optional integer vector assigning each individual to a group. This argument can be used when |
response |
An optional argument to select trait types. The |
fig |
TRUE/FALSE to add a figure of Delta PVE or not. |
contrasts |
An optional vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. If |
This function calls linear or logistic mixed models via the gaston
package (Perdry & Dandine-Roulland 2020).
If "quantitative"
is selected, Var_self
or Var_nei
in the output is given by the proportion of phenotypic variation explained (PVE) by neighbor effects as PVEnei =.
If
"binary"
is selected, Var_self
or Var_nei
is given by the ratio of phenotypic variation explained (RVE) by neighbor effects as RVEnei = and p-values are not available.
This is because a logistic mixed model
logistic.mm.aireml()
called via the gaston
package does not provide and log-likelihood (see Chen et al. 2016 for the theory).
A matrix containing the maximum neighbor distance, phenotypic variation explained by neighbor effects, and p-value by a likelihood ratio test.
scale
Maximum neighbor distance given as an argument
Var_self
Proportion or ratio of phenotypic variation explained (PVE or RVE) by self-genotype effects for linear or logistic mixed models, respectively
Var_nei
Proportion or ratio of phenotypic variation explained (PVE or RVE) by neighbor effects for linear or logistic mixed models, respectively
p-value
p-value by a likelihood ratio test between models with or without neighbor effects. Self effects are tested when the scale is zero
Yasuhiro Sato ([email protected])
Perdry H, Dandine-Roulland C (2020) gaston: Genetic Data Handling (QC, GRM, LD, PCA) & Linear Mixed Models. R package version 1.5.6. https://CRAN.R-project.org/package=gaston
Chen H, Wang C, Conomos M. et al. (2016) Control for population structure and relatedness for binary traits in genetic association studies via logistic mixed models. The American Journal of Human Genetics 98: 653-666.
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) s_seq <- quantile(dist(test_smap),c(0.1*(1:10))) test_pve <- calc_pve(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap, s_seq=s_seq, )
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) s_seq <- quantile(dist(test_smap),c(0.1*(1:10))) test_pve <- calc_pve(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap, s_seq=s_seq, )
A function to decompose qtl
's object of conditional genotype probabilities.
decompose_genoprobs(genoprobs, contrasts = NULL)
decompose_genoprobs(genoprobs, contrasts = NULL)
genoprobs |
Conditional genotype probabilities as taken from |
contrasts |
A vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. |
A list of three numeric matrices for genotype probabilities AA, AB, and BB. Each contains elements of individuals x markers.
AA
Homozygote AA probabilities.
AB
Heterozygote AB probabilities for. NA if inbred lines
BB
Homozygote BB probabilities. NA if backcross lines
Yasuhiro Sato ([email protected])
A function to estimate additive and dominance deviation for self and neighbor QTL effects by a simple regression.
eff_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), fig = TRUE, contrasts = NULL )
eff_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), fig = TRUE, contrasts = NULL )
genoprobs |
Conditional genotype probabilities as taken from |
pheno |
A vector of individual phenotypes. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial position along an x-axis and y-axis, respectively. |
scale |
A numeric scalar indicating the maximum spatial distance between a focal individual and neighbors to define neighbor effects. |
addcovar |
An optional matrix including additional non-genetic covariates. It contains no. of individuals x no. of covariates. |
addQTL |
An optional vector containing marker names that are considered covariates. Namely, this option allows composite interval mapping (Jansen 1993). |
grouping |
An optional integer vector assigning each individual to a group. This argument can be used when |
response |
An optional argument to select trait types. The |
fig |
TRUE/FALSE to plot the effects or not. |
contrasts |
An optional vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. If |
Similar to Haley-Knott regression (Haley & Knott 1992), the additive and dominance deviations are approximated by a regression of trait values on conditional genotype probabilities.
The self QTL effects a1
and d1
are estimated in the same way as the qtl
package performs the Haley-Knott regression.
If contrasts = c(TRUE, TRUE, TRUE)
, neighbor QTL effects a1
and d1
are estimated using a quadratic regression; otherwise, the additive neighbor effects are estimated using a linear regression.
See also Sato, Takeda & Nagano (2021) for the rationale behind the approximation.
A matrix of estimated additive and dominance deviation for self and neighbor effects, with the chromosome numbers and positions. The row names correspond to marker names.
chr
Chromosome number
pos
Marker position
a1
Additive deviation for self effects
d1
Dominance deviation for self effects
a2
Additive deviation for neighbor effects
d2
Dominance deviation for neighbor effects
Yasuhiro Sato ([email protected])
Haley CS, Knott SA (1992) A simple regression method for mapping quantitative trait loci in line crosses using flanking markers. Heredity 69:315-324.
Jansen RC (1993) Interval mapping of multiple quantitative trait loci. Genetics 135:205-211.
Sato Y, Takeda K, Nagano AJ (2021) Neighbor QTL: an interval mapping method for quantitative trait loci underlying plant neighborhood effects. G3; Genes|Genomes|Genetics 11:jkab017.
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_eff <- eff_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap, scale=20, fig=TRUE )
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_eff <- eff_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap, scale=20, fig=TRUE )
A function to reshape qtl
's object of conditional genotype probabilities, and to calculate self QTL effects for all individuals with given deviation coefficients and conditional genotype probabilities.
genoprobs2selfprobs(genoprobs, a1, d1, contrasts = NULL)
genoprobs2selfprobs(genoprobs, a1, d1, contrasts = NULL)
genoprobs |
Conditional genotype probabilities as taken from |
a1 |
A numeric scalar indicating additive deviation. |
d1 |
A numeric scalar indicating dominance deviation. |
contrasts |
A vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. |
A numeric matrix containing individuals x marker elements for self QTL effects.
Yasuhiro Sato ([email protected])
A function to get marker information from a genetic map including observed and pseudo markers
get_markers(genoprobs)
get_markers(genoprobs)
genoprobs |
Conditional genotype probabilities as taken from |
A matrix showing the chromosome numbers (the first column) and positions (the second column) for all markers (row names).
Yasuhiro Sato ([email protected])
A function to test interaction terms between one focal marker and the other markers across a genome.
int_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL, intQTL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), contrasts = NULL )
int_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL, intQTL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), contrasts = NULL )
genoprobs |
Conditional genotype probabilities as taken from |
pheno |
A vector of individual phenotypes. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial positions along an x-axis and y-axis, respectively. |
scale |
A numeric scalar indicating the maximum spatial distance between a focal individual and neighbors to define neighbor effects. |
addcovar |
An optional matrix including additional non-genetic covariates. It contains no. of individuals x no. of covariates. |
addQTL |
A vector containing marker names that are considered covariates. This argument is necessary for |
intQTL |
A name of a focal marker to be tested for its epistasis with the other markers in neighbor effects. The marker name must be included by |
grouping |
An optional integer vector assigning each individual to a group. This argument can be used when |
response |
An optional argument to select trait types. The |
contrasts |
An optional vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. If |
This is an optimal function to test two-way interactions between the main neighbor effect of a focal marker given by intQTL
and the others.
All the main neighbor effects are first estimated using eff_neighbor()
, and then a two-way interaction term between the focal marker effect and its counterpart was considered an additional explanatory variable.
LOD score was compared between models with or without the two-way interaction.
A matrix of LOD scores for neighbor epistasis effects, with the chromosome numbers and positions. The row names correspond to marker names.
chr
Chromosome number
pos
Marker position
LOD_int
LOD score for epistasis in neighbor effects between a focal and the other markers
Yasuhiro Sato ([email protected])
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_int <- int_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap,scale=20, addQTL=c("c1_D1M1","c1_D1M2"),intQTL="c1_D1M1" ) plot_nei(test_int, type="int")
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_int <- int_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap,scale=20, addQTL=c("c1_D1M1","c1_D1M2"),intQTL="c1_D1M1" ) plot_nei(test_int, type="int")
An utility function to extract log-likelihood based on AIC of glm.fit()
logLik_glm.fit(...)
logLik_glm.fit(...)
... |
Arguments to be passed to |
Log-likelihood
A function to calculate a Euclidian distance including at least one neighbor for all individuals.
min_dist(smap, grouping = rep(1, nrow(smap)))
min_dist(smap, grouping = rep(1, nrow(smap)))
smap |
A matrix showing a spatial map. The first and second column include spatial points along a x-axis and y-axis, respectively. |
grouping |
A integer vector assigning each individual to a group. This argument can be useful when a "smap" contains different experimental replicates. Default setting means that all individuals are belong to a single group. |
Return a scalar of the minimum Euclidian distance that allows all individuals to have at least one neighbor.
Yasuhiro Sato ([email protected])
A function to calculate neighbor QTL effects between two individuals, with given deviation coefficients and conditional genotype probabilities.
neiprob(i, j, a2, d2, AA, AB, BB, d2sq0 = FALSE)
neiprob(i, j, a2, d2, AA, AB, BB, d2sq0 = FALSE)
i |
ID of a target individual. |
j |
ID of an interacting neighbor. |
a2 |
A numeric scalar indicating additive deviation. |
d2 |
A numeric scalar indicating dominance deviation. |
AA |
An individual x marker matrix of conditional probabilities for AA genotype. |
AB |
An individual x marker matrix of conditional probabilities for AB genotype. Input NA if heterozygotes are absent. |
BB |
An individual x marker matrix of conditional probabilities for BB genotype. Input NA for backcross lines. |
d2sq0 |
An option to make AB/AB interaction effects zero. |
A numeric vector containing each marker effect for individual i.
Yasuhiro Sato ([email protected])
A function to calculate a genome-wide LOD threshold using permutation tests for self or neighbor effects.
perm_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL = NULL, intQTL = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), type = c("neighbor", "self", "int"), times = 99, p_val = 0.05, n_core = 1L, contrasts = NULL )
perm_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL = NULL, intQTL = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), type = c("neighbor", "self", "int"), times = 99, p_val = 0.05, n_core = 1L, contrasts = NULL )
genoprobs |
Conditional genotype probabilities as taken from |
pheno |
A vector of individual phenotypes. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial positions along an x-axis and y-axis, respectively. |
scale |
A numeric scalar indicating the maximum spatial distance between a focal individual and neighbors to define neighbor effects. |
addcovar |
An optional matrix including additional non-genetic covariates. It contains no. of individuals x no. of covariates. |
addQTL |
An optional vector containing marker names that are considered covariates. Namely, this option allows composite interval mapping (Jansen 1993). |
intQTL |
An option when using |
grouping |
An optional integer vector assigning each individual to a group. This argument can be used when |
response |
An optional argument to select trait types. The |
type |
Select |
times |
No. of permutation iterations. Default at 99 times |
p_val |
A vector indicating upper quantiles for permutation LOD scores |
n_core |
No. of cores for a parallel computation. This does not work for Windows OS. Default is a single-core computation. |
contrasts |
An optional vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. If |
LOD thresholds at given quantiles by p-val
Yasuhiro Sato ([email protected])
plot_nei
scan_neighbor
int_neighbor
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_perm <- perm_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap,scale=20, times=3, p_val=c(1.0,0.5) )
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_perm <- perm_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap,scale=20, times=3, p_val=c(1.0,0.5) )
Plot estimated additive and dominance deviation for self or neighbor effects across a genome
plot_eff(res, type = c("neighbor", "self"))
plot_eff(res, type = c("neighbor", "self"))
res |
Output results of |
type |
An option to select |
Yasuhiro Sato ([email protected])
eff_neighbor
Plot LOD curves for a genome scan of self and neighbor QTL effects.
plot_nei(res, type = c("neighbor", "self", "int"), chr = NULL, th = NULL, ...)
plot_nei(res, type = c("neighbor", "self", "int"), chr = NULL, th = NULL, ...)
res |
Output results of |
type |
Plot |
chr |
An optional vector to select chromosome numbers to be plotted. If |
th |
Add genome-wide threshold by user-defined vectors or Bonferroni correction. Default is no thresholds added. |
... |
Arguments to be passed to |
For the type
argument, "int"
can be selected to draw the results of int_neighbor()
.
In this case, the res
object and type
must match, otherwise it returns an error message.
Yasuhiro Sato ([email protected])
scan_neighbor
int_neighbor
perm_neighbor
Genome scan using a QTL model for self and neighbor effects, with possible allowance for additional covariates and non-normal traits. Theoretical background is described in Sato, Takeda & Nagano (2021).
scan_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), contrasts = NULL )
scan_neighbor( genoprobs, pheno, smap, scale, addcovar = NULL, addQTL = NULL, grouping = rep(1, nrow(smap)), response = c("quantitative", "binary"), contrasts = NULL )
genoprobs |
Conditional genotype probabilities as taken from |
pheno |
A vector of individual phenotypes. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial positions along an x-axis and y-axis, respectively. |
scale |
A numeric scalar indicating the maximum spatial distance between a focal individual and neighbors to define neighbor effects. |
addcovar |
An optional matrix including additional non-genetic covariates. It contains no. of individuals x no. of covariates. |
addQTL |
An optional vector containing marker names that are considered covariates. Namely, this option allows composite interval mapping (Jansen 1993). |
grouping |
An optional integer vector assigning each individual to a group. This argument can be used when |
response |
An optional argument to select trait types. The |
contrasts |
An optional vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. If |
This function calculates LOD score after the additive and dominance deviation are estimated using eff_neighbor()
.
As it adopts a stepwise testing from self to neighbor effects, LOD_self
are the same as standard QTL mapping.
Note that the results return 0 LOD scores for covariate markers when using addQTL
option.
A matrix of LOD scores for self and neighbor effects, with the chromosome numbers and positions. The row names correspond to marker names.
chr
Chromosome number
pos
Marker position
LOD_self
LOD score for self effects
LOD_nei
LOD score for neighbor effects
Yasuhiro Sato ([email protected])
Jansen RC (1993) Interval mapping of multiple quantitative trait loci. Genetics 135:205-211.
Sato Y, Takeda K, Nagano AJ (2021) Neighbor QTL: an interval mapping method for quantitative trait loci underlying plant neighborhood effects. G3; Genes|Genomes|Genetics 11:jkab017.
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_scan <- scan_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap, scale=20 ) plot_nei(test_scan)
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) test_scan <- scan_neighbor(genoprobs=test_genoprobs, pheno=test_cross$pheno$phenotype, smap=test_smap, scale=20 ) plot_nei(test_scan)
A function to calculate self QTL effects for an individual, with given deviation coefficients and conditional genotype probabilities.
selfprob(i, a1, d1, AA, AB, BB)
selfprob(i, a1, d1, AA, AB, BB)
i |
ID of a target individual. |
a1 |
A numeric scalar indicating additive deviation. |
d1 |
A numeric scalar indicating dominance deviation. |
AA |
An individual x marker matrix of conditional probabilities for AA genotype. |
AB |
An individual x marker matrix of conditional probabilities for AB genotype. Input NA if heterozygotes are absent. |
BB |
An individual x marker matrix of conditional probabilities for BB genotype. Input NA for backcross lines. |
A numeric vector containing each marker effect for individual i.
Yasuhiro Sato ([email protected])
A function to simulate neighbor effects with given QTL effects, distance scale, and causal markers.
sim_nei_qtl( genoprobs, a2, d2, smap, scale, grouping = rep(1, nrow(smap)), n_QTL = 1, contrasts = NULL )
sim_nei_qtl( genoprobs, a2, d2, smap, scale, grouping = rep(1, nrow(smap)), n_QTL = 1, contrasts = NULL )
genoprobs |
Conditional genotype probabilities as taken from |
a2 |
A numeric scalar indicating additive deviation. |
d2 |
A numeric scalar indicating dominance deviation. |
smap |
A matrix showing a spatial map for individuals. The first and second column include spatial positions along an x-axis and y-axis, respectively. |
scale |
A numeric scalar indicating the maximum spatial distance between a focal individual and neighbors to define neighbor effects. |
grouping |
An integer vector assigning each individual to a group. This argument can be used when |
n_QTL |
A positive integer indicating the number of causal markers. |
contrasts |
An optional vector composed of three TRUE/FALSE values, which represents the presence/absence of specific genotypes as c(TRUE/FALSE, TRUE/FALSE, TRUE/FALSE) = AA, AB, BB. If |
Major genetic effects, a2
and d2
, are allocated to causal loci randomly selected by n_QTL
, while minor polygenic effects (i.e., 1% of a2
) are allocated to the other loci.
A numeric matrix containing individuals x marker elements for neighbor QTL effects.
true_scale
True distance scale of simulated neighbor effects
true_marker
The name(s) of causal markers
nei_y
Simulated neighbor effects standardized to have zero mean and one variance
Yasuhiro Sato ([email protected])
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) nei_eff <- sim_nei_qtl(genoprobs=test_genoprobs, a2=0.5, d2=0.5, smap=test_smap, scale=20, n_QTL=1) test_scan <- scan_neighbor(genoprobs=test_genoprobs, pheno=nei_eff$nei_y, smap=test_smap, scale=20 ) plot_nei(test_scan)
set.seed(1234) test_map <- qtl::sim.map(len=rep(20,5),n.mar=3,include.x=FALSE) test_cross <- qtl::sim.cross(test_map,n.ind=50) test_smap <- cbind(runif(50,1,100),runif(50,1,100)) test_genoprobs <- qtl::calc.genoprob(test_cross,step=2) nei_eff <- sim_nei_qtl(genoprobs=test_genoprobs, a2=0.5, d2=0.5, smap=test_smap, scale=20, n_QTL=1) test_scan <- scan_neighbor(genoprobs=test_genoprobs, pheno=nei_eff$nei_y, smap=test_smap, scale=20 ) plot_nei(test_scan)