profile_qtl: QTL profiling

Description Usage Arguments Value Author(s) References Examples

View source: R/profile_qtl.R

Description

Generates the score-based genome-wide profile conditional to the selected QTL.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
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, ...)

Arguments

data

an object of class qtlpoly.data.

model

an object of class qtlpoly.model containing the QTL to be profiled.

d.sint

a d value to subtract from logarithm of p-value (LOP-d) for support interval calculation, e.g. d=1.5 (default) represents approximate 95% support interval.

polygenes

if TRUE all QTL but the one being tested are treated as a single polygenic effect, if FALSE (default) all QTL effect variances have to estimated.

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 NULL, no file is produced.

verbose

if TRUE (default), current progress is shown; if FALSE, no output is produced.

x

an object of class qtlpoly.profile to be printed.

pheno.col

a numeric vector with the phenotype column numbers to be plotted; if NULL, all phenotypes from 'data' will be included.

sint

whether "upper" or "lower" support intervals should be printed; if NULL (default), only QTL peak information will be printed.

...

currently ignored

Value

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.

Author(s)

Guilherme da Silva Pereira, gdasilv@ncsu.edu

References

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.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
  
  # 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)
  

qtlpoly documentation built on Jan. 12, 2022, 5:06 p.m.