Title: | Interactive Graphics for QTL Experiments |
---|---|
Description: | Web-based interactive charts (using D3.js) for the analysis of experimental crosses to identify genetic loci (quantitative trait loci, QTL) contributing to variation in quantitative traits. Broman (2015) <doi:10.1534/genetics.114.172742>. |
Authors: | Karl W Broman [aut, cre] , Michael Bostock [ctb, cph] (d3.js library in htmlwidgets/lib, https://d3js.org), jQuery Foundation [cph] (jQuery library in htmlwidgets/lib, https://jquery.com), jQuery contributors [ctb] (jQuery library in htmlwidgets/lib; see https://github.com/jquery/jquery/blob/master/AUTHORS.txt), jQuery UI contributors [ctb] (jQuery UI library in htmlwidgets/lib; see https://jqueryui.com/about/) |
Maintainer: | Karl W Broman <[email protected]> |
License: | GPL-3 | file LICENSE |
Version: | 0.17-2 |
Built: | 2024-10-25 03:00:43 UTC |
Source: | https://github.com/kbroman/qtlcharts |
A QTL is a *quantitative trait locus*: a genetic locus that contributes to variation in a quantitative trait. The goal of [R/qtlcharts](https://kbroman.org/qtlcharts/) is to provide interactive data visualizations for QTL analyses, and to make these visualizations available from R. It is a companion to the [R/qtl](https://rqtl.org) package.
Vignettes online at the [R/qtlcharts website](https://kbroman.org/qtlcharts/):
- [User guide](https://kbroman.org/qtlcharts/assets/vignettes/userGuide.html) - [Developer guide](https://kbroman.org/qtlcharts/assets/vignettes/develGuide.html) - [Use with R Markdown](https://kbroman.org/qtlcharts/assets/vignettes/Rmarkdown.html) \[[Rmd source](https://github.com/kbroman/qtlcharts/blob/gh-pages/assets/vignettes/Rmarkdown.Rmd)\] - [List of chart customization options](https://kbroman.org/qtlcharts/assets/vignettes/chartOpts.html)
Combine multiple runs of estQTLeffects by applying cbind to each component
cbindQTLeffects(..., labels)
cbindQTLeffects(..., labels)
... |
Results of [estQTLeffects()] |
labels |
Vector of labels to use in the combination. |
list of matrices; each component corresponds to a position in the genome and is a matrix with phenotypes x effects
[estQTLeffects()]
library(qtl) data(fake.f2) fake.f2 <- calc.genoprob(fake.f2) sex <- fake.f2$pheno$sex eff.fem <- estQTLeffects(fake.f2[,sex==0], pheno.col=1) eff.mal <- estQTLeffects(fake.f2[,sex==1], pheno.col=1) eff <- cbindQTLeffects(eff.fem, eff.mal, labels=c("female", "male"))
library(qtl) data(fake.f2) fake.f2 <- calc.genoprob(fake.f2) sex <- fake.f2$pheno$sex eff.fem <- estQTLeffects(fake.f2[,sex==0], pheno.col=1) eff.mal <- estQTLeffects(fake.f2[,sex==1], pheno.col=1) eff <- cbindQTLeffects(eff.fem, eff.mal, labels=c("female", "male"))
Calculates the effects of QTL at each position across the genome using Haley-Knott regression, much like [qtl::effectscan()], but considering multiple phenotypes and not plotting the results
estQTLeffects(cross, pheno.col = 1, what = c("means", "effects"))
estQTLeffects(cross, pheno.col = 1, what = c("means", "effects"))
cross |
(Optional) Object of class '"cross"', see [qtl::read.cross()]. |
pheno.col |
Phenotype columns in cross object. |
what |
Indicates whether to calculate phenotype averages for each genotype group or to turn these into additive and dominance effects. |
One should first run [qtl::calc.genoprob()]; if not, it is run with the default arguments.
The estimated effects will be poorly estimated in the case of selective genotyping, as Haley-Knott regression performs poorly in this case.
list of matrices; each component corresponds to a position in the genome and is a matrix with phenotypes x effects
[iplotMScanone()], [qtl::effectscan()] [cbindQTLeffects()]
data(grav) library(qtl) grav <- reduce2grid(calc.genoprob(grav, step=1)) out <- estQTLeffects(grav, phe=seq(1, nphe(grav), by=5))
data(grav) library(qtl) grav <- reduce2grid(calc.genoprob(grav, step=1)) out <- estQTLeffects(grav, phe=seq(1, nphe(grav), by=5))
An anonymized set of gene expression values, for 100 genes all influenced by a common locus, plus a vector of genotypes for the 491 individuals.
data(geneExpr)
data(geneExpr)
A list containing a matrix 'expr' with the gene expression data plus a vector 'genotype' with the genotypes.
Karl W Broman, 2013-05-16
data(geneExpr) # heat map of correlation matrix, linked to scatterplots iplotCorr(geneExpr$expr, geneExpr$genotype, reorder=TRUE)
data(geneExpr) # heat map of correlation matrix, linked to scatterplots iplotCorr(geneExpr$expr, geneExpr$genotype, reorder=TRUE)
Data from a QTL experiment on gravitropism in Arabidopsis, with data on 162 recombinant inbred lines (Ler x Cvi). The outcome is the root tip angle (in degrees) at two-minute increments over eight hours.
data(grav)
data(grav)
An object of class '"cross"'; see [qtl::read.cross()].
Mouse Phenome Database, <https://phenome.jax.org/projects/Moore1b>
Moore CR, Johnson LS, Kwak I-Y, Livny M, Broman KW, Spalding EP (2013) High-throughput computer vision introduces the time axis to a quantitative trait map of a plant growth response. Genetics 195:1077-1086 ([PubMed](https://pubmed.ncbi.nlm.nih.gov/23979570/))
data(grav) times <- attr(grav, "time") phe <- grav$pheno iplotCurves(phe, times, phe[,c(61,121)], phe[,c(121,181)], chartOpts=list(curves_xlab="Time (hours)", curves_ylab="Root tip angle (degrees)", scat1_xlab="Angle at 2 hrs", scat1_ylab="Angle at 4 hrs", scat2_xlab="Angle at 4 hrs", scat2_ylab="Angle at 6 hrs"))
data(grav) times <- attr(grav, "time") phe <- grav$pheno iplotCurves(phe, times, phe[,c(61,121)], phe[,c(121,181)], chartOpts=list(curves_xlab="Time (hours)", curves_ylab="Root tip angle (degrees)", scat1_xlab="Angle at 2 hrs", scat1_ylab="Angle at 4 hrs", scat2_xlab="Angle at 4 hrs", scat2_ylab="Angle at 6 hrs"))
Creates an interactive graph for a large set of box plots (rendered as lines connecting the quantiles), linked to underlying histograms.
iboxplot( dat, qu = c(0.001, 0.01, 0.1, 0.25), orderByMedian = TRUE, breaks = 251, chartOpts = NULL, digits = 5 )
iboxplot( dat, qu = c(0.001, 0.01, 0.1, 0.25), orderByMedian = TRUE, breaks = 251, chartOpts = NULL, digits = 5 )
dat |
Data matrix (individuals x variables) |
qu |
Quantiles to plot (All with 0 < qu < 0.5) |
orderByMedian |
If TRUE, reorder individuals by their median |
breaks |
Number of bins in the histograms, or a vector of locations of the breakpoints between bins (as in [graphics::hist()]) |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotCorr()], [scat2scat()]
n.ind <- 500 n.gene <- 10000 expr <- matrix(rnorm(n.ind * n.gene, (1:n.ind)/n.ind*3), ncol=n.gene) dimnames(expr) <- list(paste0("ind", 1:n.ind), paste0("gene", 1:n.gene)) iboxplot(expr, chartOpts=list(xlab="Mice", ylab="Gene expression"))
n.ind <- 500 n.gene <- 10000 expr <- matrix(rnorm(n.ind * n.gene, (1:n.ind)/n.ind*3), ncol=n.gene) dimnames(expr) <- list(paste0("ind", 1:n.ind), paste0("gene", 1:n.gene)) iboxplot(expr, chartOpts=list(xlab="Mice", ylab="Gene expression"))
Output and render functions for using R/qtlcharts widgets within Shiny applications and interactive Rmd documents.
iboxplot_output(outputId, width = "100%", height = "900") iboxplot_render(expr, env = parent.frame(), quoted = FALSE) idotplot_output(outputId, width = "100%", height = "530") idotplot_render(expr, env = parent.frame(), quoted = FALSE) iheatmap_output(outputId, width = "100%", height = "1000") iheatmap_render(expr, env = parent.frame(), quoted = FALSE) ipleiotropy_output(outputId, width = "100%", height = "580") ipleiotropy_render(expr, env = parent.frame(), quoted = FALSE) iplot_output(outputId, width = "100%", height = "580") iplot_render(expr, env = parent.frame(), quoted = FALSE) iplotCorr_output(outputId, width = "100%", height = "1000") iplotCorr_render(expr, env = parent.frame(), quoted = FALSE) iplotCurves_output(outputId, width = "100%", height = "1000") iplotCurves_render(expr, env = parent.frame(), quoted = FALSE) iplotMScanone_output(outputId, width = "100%", height = "580") iplotMScanone_render(expr, env = parent.frame(), quoted = FALSE) iplotMap_output(outputId, width = "100%", height = "680") iplotMap_render(expr, env = parent.frame(), quoted = FALSE) iplotRF_output(outputId, width = "100%", height = "1000") iplotRF_render(expr, env = parent.frame(), quoted = FALSE) iplotScanone_output(outputId, width = "100%", height = "580") iplotScanone_render(expr, env = parent.frame(), quoted = FALSE) iplotScantwo_output(outputId, width = "100%", height = "1000") iplotScantwo_render(expr, env = parent.frame(), quoted = FALSE) itriplot_output(outputId, width = "100%", height = "530") itriplot_render(expr, env = parent.frame(), quoted = FALSE) scat2scat_output(outputId, width = "100%", height = "530") scat2scat_render(expr, env = parent.frame(), quoted = FALSE)
iboxplot_output(outputId, width = "100%", height = "900") iboxplot_render(expr, env = parent.frame(), quoted = FALSE) idotplot_output(outputId, width = "100%", height = "530") idotplot_render(expr, env = parent.frame(), quoted = FALSE) iheatmap_output(outputId, width = "100%", height = "1000") iheatmap_render(expr, env = parent.frame(), quoted = FALSE) ipleiotropy_output(outputId, width = "100%", height = "580") ipleiotropy_render(expr, env = parent.frame(), quoted = FALSE) iplot_output(outputId, width = "100%", height = "580") iplot_render(expr, env = parent.frame(), quoted = FALSE) iplotCorr_output(outputId, width = "100%", height = "1000") iplotCorr_render(expr, env = parent.frame(), quoted = FALSE) iplotCurves_output(outputId, width = "100%", height = "1000") iplotCurves_render(expr, env = parent.frame(), quoted = FALSE) iplotMScanone_output(outputId, width = "100%", height = "580") iplotMScanone_render(expr, env = parent.frame(), quoted = FALSE) iplotMap_output(outputId, width = "100%", height = "680") iplotMap_render(expr, env = parent.frame(), quoted = FALSE) iplotRF_output(outputId, width = "100%", height = "1000") iplotRF_render(expr, env = parent.frame(), quoted = FALSE) iplotScanone_output(outputId, width = "100%", height = "580") iplotScanone_render(expr, env = parent.frame(), quoted = FALSE) iplotScantwo_output(outputId, width = "100%", height = "1000") iplotScantwo_render(expr, env = parent.frame(), quoted = FALSE) itriplot_output(outputId, width = "100%", height = "530") itriplot_render(expr, env = parent.frame(), quoted = FALSE) scat2scat_output(outputId, width = "100%", height = "530") scat2scat_render(expr, env = parent.frame(), quoted = FALSE)
outputId |
output variable to read from |
width , height
|
Must be a valid CSS unit (like '"100%"', '"400px"', '"auto"') or a number, which will be coerced to a string and have '"px"' appended. |
expr |
An expression that generates a networkD3 graph |
env |
The environment in which to evaluate 'expr'. |
quoted |
Is 'expr' a quoted expression (with 'quote()')? This is useful if you want to save an expression in a variable. |
Creates an interactive graph of phenotypes vs genotypes at a marker.
idotplot(x, y, indID = NULL, group = NULL, chartOpts = NULL, digits = 5)
idotplot(x, y, indID = NULL, group = NULL, chartOpts = NULL, digits = 5)
x |
Vector of groups of individuals (e.g., a genotype) |
y |
Numeric vector (e.g., a phenotype) |
indID |
Optional vector of character strings, shown with tool tips |
group |
Optional vector of categories for coloring points |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplot()], [iplotPXG()]
n <- 100 g <- sample(LETTERS[1:3], n, replace=TRUE) y <- rnorm(n, match(g, LETTERS[1:3])*10, 5) idotplot(g, y)
n <- 100 g <- sample(LETTERS[1:3], n, replace=TRUE) y <- rnorm(n, match(g, LETTERS[1:3])*10, 5) idotplot(g, y)
Creates an interactive heatmap, with each cell linked to plots of horizontal and vertical slices
iheatmap(z, x = NULL, y = NULL, chartOpts = NULL, digits = 5)
iheatmap(z, x = NULL, y = NULL, chartOpts = NULL, digits = 5)
z |
Numeric matrix of dim 'length(x)' x 'length(y)' |
x |
Vector of numeric values for the x-axis |
y |
Vector of numeric values for the y-axis |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
By default, the z-axis limits are from '-max(abs(z))' to 'max(abs(z))', and negative cells are colored blue to white which positive cells are colored white to red.
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotCorr()]
n <- 101 x <- y <- seq(-2, 2, len=n) z <- matrix(ncol=n, nrow=n) for(i in seq(along=x)) for(j in seq(along=y)) z[i,j] <- x[i]*y[j]*exp(-x[i]^2 - y[j]^2) iheatmap(z, x, y)
n <- 101 x <- y <- seq(-2, 2, len=n) z <- matrix(ncol=n, nrow=n) for(i in seq(along=x)) for(j in seq(along=y)) z[i,j] <- x[i]*y[j]*exp(-x[i]^2 - y[j]^2) iheatmap(z, x, y)
Creates an interactive graph of a scatterplot of two phenotypes, plus optionally the LOD curves for the two traits along one chromosome, with a slider for selecting the locations of two QTL which are then indicated on the LOD curves and the corresponding genotypes used to color the points in the scatterplot.
ipleiotropy( cross, scanoneOutput = NULL, pheno.col = 1:2, lodcolumn = 1:2, chr = NULL, interval = NULL, fillgenoArgs = NULL, chartOpts = NULL, digits = 5 )
ipleiotropy( cross, scanoneOutput = NULL, pheno.col = 1:2, lodcolumn = 1:2, chr = NULL, interval = NULL, fillgenoArgs = NULL, chartOpts = NULL, digits = 5 )
cross |
(Optional) Object of class '"cross"', see [qtl::read.cross()]. |
scanoneOutput |
(Optional) object of class '"scanone"', as output from [qtl::scanone()]. |
pheno.col |
Vector indicating two phenotype column in cross object; either numeric or character strings (the latter being the phenotype column names). |
lodcolumn |
Vector of two numeric values indicating LOD score columns to plot. |
chr |
A single chromosome ID, as a character string. |
interval |
A numeric vector of length 2, defining an interval that indicates what portion of the chromosome should be included. |
fillgenoArgs |
List of named arguments to pass to [qtl::fill.geno()], if needed. |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
[qtl::fill.geno()] is used to impute missing genotypes. In this case, arguments to [qtl::fill.geno()] are passed as a list, for example 'fillgenoArgs=list(method="argmax", error.prob=0.002, map.function="c-f")'.
Individual IDs (viewable when hovering over a point in the scatterplot of the two phenotypes) are taken from the input 'cross' object, using the [qtl::getid()] function in R/qtl.
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotScanone()], [iplotMScanone()], [iplotPXG()]
library(qtl) data(fake.bc) fake.bc <- calc.genoprob(fake.bc[5,], step=1) # select chr 5 out <- scanone(fake.bc, method="hk", pheno.col=1:2) ipleiotropy(fake.bc, out) # omit the LOD curves ipleiotropy(fake.bc)
library(qtl) data(fake.bc) fake.bc <- calc.genoprob(fake.bc[5,], step=1) # select chr 5 out <- scanone(fake.bc, method="hk", pheno.col=1:2) ipleiotropy(fake.bc, out) # omit the LOD curves ipleiotropy(fake.bc)
Creates an interactive scatterplot.
iplot(x, y = NULL, group = NULL, indID = NULL, chartOpts = NULL, digits = 5)
iplot(x, y = NULL, group = NULL, indID = NULL, chartOpts = NULL, digits = 5)
x |
Numeric vector of x values |
y |
Numeric vector of y values |
group |
Optional vector of categories for coloring the points |
indID |
Optional vector of character strings, shown with tool tips |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotCorr()], [iplotCurves()], [itriplot()], [idotplot()], [iplotPXG()]
x <- rnorm(100) grp <- sample(1:3, 100, replace=TRUE) y <- x*grp + rnorm(100) iplot(x, y, grp)
x <- rnorm(100) grp <- sample(1:3, 100, replace=TRUE) y <- x*grp + rnorm(100) iplot(x, y, grp)
Creates an interactive graph with an image of a correlation matrix linked to underlying scatterplots.
iplotCorr( mat, group = NULL, rows = NULL, cols = NULL, reorder = FALSE, corr = NULL, scatterplots = TRUE, chartOpts = NULL, digits = 5 )
iplotCorr( mat, group = NULL, rows = NULL, cols = NULL, reorder = FALSE, corr = NULL, scatterplots = TRUE, chartOpts = NULL, digits = 5 )
mat |
Data matrix (individuals x variables) |
group |
Optional vector of groups of individuals (e.g., a genotype) |
rows |
Selected rows of the correlation matrix to include in the image. Ignored if 'corr' is provided. |
cols |
Selected columns of the correlation matrix to include in the image. Ignored if 'corr' is provided. |
reorder |
If TRUE, reorder the variables by clustering. Ignored if 'corr' is provided as a subset of the overall correlation matrix |
corr |
Correlation matrix (optional). |
scatterplots |
If ‘FALSE', don’t have the heat map be linked to scatterplots. |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
'corr' may be provided as a subset of the overall correlation matrix for the columns of 'mat'. In this case, the 'reorder', 'rows' and 'cols' arguments are ignored. The row and column names of 'corr' must match the names of some subset of columns of 'mat'.
Individual IDs are taken from 'rownames(mat)'; they must match 'names(group)'.
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iheatmap()], [scat2scat()], [iplotCurves()]
data(geneExpr) iplotCorr(geneExpr$expr, geneExpr$genotype, reorder=TRUE, chartOpts=list(cortitle="Correlation matrix", scattitle="Scatterplot")) # use Spearman's correlation corr <- cor(geneExpr$expr, method="spearman", use="pairwise.complete.obs") # order by hierarchical clustering o <- hclust(as.dist(1-corr))$order iplotCorr(geneExpr$expr[,o], geneExpr$genotype, corr=corr[o,o], chartOpts=list(cortitle="Spearman correlation", scattitle="Scatterplot"))
data(geneExpr) iplotCorr(geneExpr$expr, geneExpr$genotype, reorder=TRUE, chartOpts=list(cortitle="Correlation matrix", scattitle="Scatterplot")) # use Spearman's correlation corr <- cor(geneExpr$expr, method="spearman", use="pairwise.complete.obs") # order by hierarchical clustering o <- hclust(as.dist(1-corr))$order iplotCorr(geneExpr$expr[,o], geneExpr$genotype, corr=corr[o,o], chartOpts=list(cortitle="Spearman correlation", scattitle="Scatterplot"))
Creates an interactive graph with a panel having a number of curves (say, a phenotype measured over time) linked to one or two (or no) scatter plots (say, of the first vs middle and middle vs last times).
iplotCurves( curveMatrix, times = NULL, scatter1 = NULL, scatter2 = NULL, group = NULL, chartOpts = NULL, digits = 5 )
iplotCurves( curveMatrix, times = NULL, scatter1 = NULL, scatter2 = NULL, group = NULL, chartOpts = NULL, digits = 5 )
curveMatrix |
Matrix (dim n_ind x n_times) with outcomes |
times |
Vector (length n_times) with time points for the columns of curveMatrix |
scatter1 |
Matrix (dim n_ind x 2) with data for the first scatterplot |
scatter2 |
Matrix (dim n_ind x 2) with data for the second scatterplot |
group |
Optional vector of groups of individuals (e.g., a genotype) |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotCorr()], [iplot()], [scat2scat()]
# random growth curves, based on some data times <- 1:16 n <- 100 start <- rnorm(n, 5.2, 0.8) slope1to5 <- rnorm(n, 2.6, 0.5) slope5to16 <- rnorm(n, 0.24 + 0.09*slope1to5, 0.195) y <- matrix(ncol=16, nrow=n) y[,1] <- start for(j in 2:5) y[,j] <- y[,j-1] + slope1to5 for(j in 6:16) y[,j] <- y[,j-1] + slope5to16 y <- y + rnorm(prod(dim(y)), 0, 0.35) iplotCurves(y, times, y[,c(1,5)], y[,c(5,16)], chartOpts=list(curves_xlab="Time", curves_ylab="Size", scat1_xlab="Size at T=1", scat1_ylab="Size at T=5", scat2_xlab="Size at T=5", scat2_ylab="Size at T=16"))
# random growth curves, based on some data times <- 1:16 n <- 100 start <- rnorm(n, 5.2, 0.8) slope1to5 <- rnorm(n, 2.6, 0.5) slope5to16 <- rnorm(n, 0.24 + 0.09*slope1to5, 0.195) y <- matrix(ncol=16, nrow=n) y[,1] <- start for(j in 2:5) y[,j] <- y[,j-1] + slope1to5 for(j in 6:16) y[,j] <- y[,j-1] + slope5to16 y <- y + rnorm(prod(dim(y)), 0, 0.35) iplotCurves(y, times, y[,c(1,5)], y[,c(5,16)], chartOpts=list(curves_xlab="Time", curves_ylab="Size", scat1_xlab="Size at T=1", scat1_ylab="Size at T=5", scat2_xlab="Size at T=5", scat2_ylab="Size at T=16"))
Creates an interactive graph of a genetic marker map.
iplotMap( map, chr = NULL, shift = FALSE, horizontal = FALSE, chartOpts = NULL, digits = 5 )
iplotMap( map, chr = NULL, shift = FALSE, horizontal = FALSE, chartOpts = NULL, digits = 5 )
map |
Object of class '"map"', a list with each component being a vector of marker positions. You can also provide an object of class '"cross"', in which case the map is extracted with [qtl::pull.map()]. |
chr |
(Optional) Vector indicating the chromosomes to plot. |
shift |
If TRUE, shift each chromsome so that the initial marker is at position 0. |
horizontal |
If TRUE, have chromosomes arranged horizontally |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotScanone()], [iplotPXG()]
library(qtl) data(hyper) map <- pull.map(hyper) iplotMap(map, shift=TRUE)
library(qtl) data(hyper) map <- pull.map(hyper) iplotMap(map, shift=TRUE)
Creates an interactive graph of a set of single-QTL genome scans, as calculated by [qtl::scanone()]. If 'cross' or 'effects' are provided, LOD curves will be linked to a panel with estimated QTL effects.
iplotMScanone( scanoneOutput, cross = NULL, lodcolumn = NULL, pheno.col = NULL, times = NULL, effects = NULL, chr = NULL, chartOpts = NULL, digits = 5 )
iplotMScanone( scanoneOutput, cross = NULL, lodcolumn = NULL, pheno.col = NULL, times = NULL, effects = NULL, chr = NULL, chartOpts = NULL, digits = 5 )
scanoneOutput |
Object of class '"scanone"', as output from [qtl::scanone()]. |
cross |
(Optional) Object of class '"cross"', see [qtl::read.cross()]. |
lodcolumn |
Numeric value indicating LOD score column to plot. |
pheno.col |
(Optional) Phenotype column in cross object. |
times |
(Optional) Vector (length equal to the number of LOD score columns) with quantitative values to which the different LOD score columns correspond (times of measurements, or something like age or dose). These need to be ordered and equally-spaced. If omitted, the names of the columns in 'scanoneOutput' are used and treated as qualitative. |
effects |
(Optional) Estimated QTL effects, as obtained with [estQTLeffects()]. |
chr |
(Optional) Optional vector indicating the chromosomes for which LOD scores should be calculated. This should be a vector of character strings referring to chromosomes by name; numeric values are converted to strings. Refer to chromosomes with a preceding - to have all chromosomes but those considered. A logical (TRUE/FALSE) vector may also be used. |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
If 'cross' is provided, Haley-Knott regression is used to estimate QTL effects at each pseudomarker.
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotScanone()]
data(grav) library(qtl) grav <- calc.genoprob(grav, step=1) grav <- reduce2grid(grav) # we're going to subset the phenotypes phecol <- seq(1, nphe(grav), by=5) # the times were saved as an attributed times <- attr(grav, "time")[phecol] # genome scan out <- scanone(grav, phe=phecol, method="hk") # plot with qualitative labels on y-axis iplotMScanone(out) # plot with quantitative y-axis iplotMScanone(out, times=times) # estimate QTL effect for each time point at each genomic position eff <- estQTLeffects(grav, phe=seq(1, nphe(grav), by=5), what="effects") # plot with QTL effects included (and with quantitative y-axis) iplotMScanone(out, effects=eff, times=times, chartOpts=list(eff_ylab="QTL effect", eff_xlab="Time (hrs)"))
data(grav) library(qtl) grav <- calc.genoprob(grav, step=1) grav <- reduce2grid(grav) # we're going to subset the phenotypes phecol <- seq(1, nphe(grav), by=5) # the times were saved as an attributed times <- attr(grav, "time")[phecol] # genome scan out <- scanone(grav, phe=phecol, method="hk") # plot with qualitative labels on y-axis iplotMScanone(out) # plot with quantitative y-axis iplotMScanone(out, times=times) # estimate QTL effect for each time point at each genomic position eff <- estQTLeffects(grav, phe=seq(1, nphe(grav), by=5), what="effects") # plot with QTL effects included (and with quantitative y-axis) iplotMScanone(out, effects=eff, times=times, chartOpts=list(eff_ylab="QTL effect", eff_xlab="Time (hrs)"))
Creates an interactive graph of phenotypes vs genotypes at a marker.
iplotPXG( cross, marker, pheno.col = 1, chartOpts = NULL, fillgenoArgs = NULL, digits = 5 )
iplotPXG( cross, marker, pheno.col = 1, chartOpts = NULL, fillgenoArgs = NULL, digits = 5 )
cross |
Object of class '"cross"', see [qtl::read.cross()]. |
marker |
Character string with marker name. |
pheno.col |
Phenotype column in cross object. |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
fillgenoArgs |
List of named arguments to pass to [qtl::fill.geno()], if needed. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
The function [qtl::fill.geno()] is used to impute missing genotypes, with arguments passed as a list, for example 'fillgenoArgs=list(method="argmax", error.prob=0.002, map.function="c-f")'.
Individual IDs (viewable when hovering over a point) are taken from the input 'cross' object, using the [qtl::getid()] function in R/qtl.
By default, points are colored blue and pink according to whether the marker genotype is observed or inferred, respectively.
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[idotplot()], [iplot()], [iplotScanone()], [iplotMap()]
library(qtl) data(hyper) marker <- sample(markernames(hyper), 1) iplotPXG(hyper, marker) # different colors iplotPXG(hyper, marker, chartOpts=list(pointcolor=c("black", "gray")))
library(qtl) data(hyper) marker <- sample(markernames(hyper), 1) iplotPXG(hyper, marker) # different colors iplotPXG(hyper, marker, chartOpts=list(pointcolor=c("black", "gray")))
Creates an interactive graph of estimated recombination fractions and LOD scores for all pairs of markers.
iplotRF(cross, chr = NULL, chartOpts = NULL, digits = 5)
iplotRF(cross, chr = NULL, chartOpts = NULL, digits = 5)
cross |
Object of class '"cross"', see [qtl::read.cross()]. |
chr |
Optional vector indicating chromosomes to include. This should be a vector of character strings referring to chromosomes by name; numeric values are converted to strings. Refer to chromosomes with a preceding '-' to have all chromosomes but those considered. A logical (TRUE/FALSE) vector may also be used. |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
The usual 'height' and 'width' options in 'chartOpts' are ignored in this plot. Instead, you may provide 'pixelPerCell' (number of pixels per cell in the heat map), 'chrGap' (gap in pixels between chromosomes in the heat map), 'cellHeight' (height in pixels of each cell in the cross-tabulation), 'cellWidth' (width in pixels of each cell in the cross-tabulation), and 'hbot' (height in pixels of the lower panels showing cross-sections of the heat map)
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[qtl::est.rf()], [qtl::plotRF()]
library(qtl) data(fake.f2) fake.f2 <- est.rf(fake.f2) iplotRF(fake.f2)
library(qtl) data(fake.f2) fake.f2 <- est.rf(fake.f2) iplotRF(fake.f2)
Creates an interactive graph of a single-QTL genome scan, as calculated by [qtl::scanone()]. If 'cross' is provided, the LOD curves are linked to a phenotype x genotype plot for a marker: Click on a marker on the LOD curve and see the corresponding phenotype x genotype plot.
iplotScanone( scanoneOutput, cross = NULL, lodcolumn = 1, pheno.col = 1, chr = NULL, pxgtype = c("ci", "raw"), fillgenoArgs = NULL, chartOpts = NULL, digits = 5 )
iplotScanone( scanoneOutput, cross = NULL, lodcolumn = 1, pheno.col = 1, chr = NULL, pxgtype = c("ci", "raw"), fillgenoArgs = NULL, chartOpts = NULL, digits = 5 )
scanoneOutput |
Object of class '"scanone"', as output from [qtl::scanone()]. |
cross |
(Optional) Object of class '"cross"', see [qtl::read.cross()]. |
lodcolumn |
Numeric value indicating LOD score column to plot. |
pheno.col |
(Optional) Phenotype column in cross object. |
chr |
(Optional) Vector indicating the chromosomes for which LOD scores should be calculated. This should be a vector of character strings referring to chromosomes by name; numeric values are converted to strings. Refer to chromosomes with a preceding - to have all chromosomes but those considered. A logical (TRUE/FALSE) vector may also be used. |
pxgtype |
If phenotype x genotype plot is to be shown, should
it be with means |
fillgenoArgs |
List of named arguments to pass to [qtl::fill.geno()], if needed. |
chartOpts |
A list of options for configuring the chart (see the coffeescript code). Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
If 'cross' is provided, [qtl::fill.geno()] is used to impute missing genotypes. In this case, arguments to [qtl::fill.geno()] are passed as a list, for example 'fillgenoArgs=list(method="argmax", error.prob=0.002, map.function="c-f")'.
With 'pxgtype="raw"', individual IDs (viewable when hovering over a point in the phenotype-by-genotype plot) are taken from the input 'cross' object, using the [qtl::getid()] function in R/qtl.
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotMScanone()], [iplotPXG()], [iplotMap()]
library(qtl) data(hyper) hyper <- calc.genoprob(hyper, step=1) out <- scanone(hyper) # iplotScanone with no effects iplotScanone(out, chr=c(1, 4, 6, 7, 15)) # iplotScanone with CIs iplotScanone(out, hyper, chr=c(1, 4, 6, 7, 15)) # iplotScanone with raw phe x gen iplotScanone(out, hyper, chr=c(1, 4, 6, 7, 15), pxgtype='raw')
library(qtl) data(hyper) hyper <- calc.genoprob(hyper, step=1) out <- scanone(hyper) # iplotScanone with no effects iplotScanone(out, chr=c(1, 4, 6, 7, 15)) # iplotScanone with CIs iplotScanone(out, hyper, chr=c(1, 4, 6, 7, 15)) # iplotScanone with raw phe x gen iplotScanone(out, hyper, chr=c(1, 4, 6, 7, 15), pxgtype='raw')
Creates an interactive plot of the results of [qtl::scantwo()], for a two-dimensional, two-QTL genome scan.
iplotScantwo( scantwoOutput, cross = NULL, lodcolumn = 1, pheno.col = 1, chr = NULL, chartOpts = NULL, digits = 5 )
iplotScantwo( scantwoOutput, cross = NULL, lodcolumn = 1, pheno.col = 1, chr = NULL, chartOpts = NULL, digits = 5 )
scantwoOutput |
Output of [qtl::scantwo()] |
cross |
(Optional) Object of class '"cross"', see [qtl::read.cross()]. |
lodcolumn |
Numeric value indicating LOD score column to plot. |
pheno.col |
(Optional) Phenotype column in cross object. |
chr |
(Optional) Optional vector indicating the chromosomes for which LOD scores should be calculated. This should be a vector of character strings referring to chromosomes by name; numeric values are converted to strings. Refer to chromosomes with a preceding - to have all chromosomes but those considered. A logical (TRUE/FALSE) vector may also be used. |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
The estimated QTL effects, and the genotypes in the phenotype x genotype plot, in the right-hand panels, are derived following a single imputation to fill in missing data, and so are a bit crude.
Note that the usual 'height' and 'width' options in 'chartOpts' are ignored here. Instead, you may provide 'pixelPerCell' (number of pixels per cell in the heat map), 'chrGap' (gaps between chr in heat map), 'wright' (width in pixels of right panels), and 'hbot' (height in pixels of each of the lower panels)
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotScanone()]
library(qtl) data(fake.f2) fake.f2 <- calc.genoprob(fake.f2, step=5) out <- scantwo(fake.f2, method="hk", verbose=FALSE) iplotScantwo(out, fake.f2)
library(qtl) data(fake.f2) fake.f2 <- calc.genoprob(fake.f2, step=5) out <- scantwo(fake.f2, method="hk", verbose=FALSE) iplotScantwo(out, fake.f2)
Creates an interactive graph of trinomial probabilities, represented as points within an equilateral triangle.
itriplot(p, indID = NULL, group = NULL, chartOpts = NULL, digits = 5)
itriplot(p, indID = NULL, group = NULL, chartOpts = NULL, digits = 5)
p |
Matrix of trinomial probabilities (n x 3); each row should sum to 1. |
indID |
Optional vector of character strings, shown with tool tips |
group |
Optional vector of categories for coloring the points |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplot()], [iplotPXG()], [idotplot()]
n <- 100 p <- matrix(runif(3*n), ncol=3) p <- p / colSums(p) g <- sample(1:3, n, replace=TRUE) itriplot(p, group=g)
n <- 100 p <- matrix(runif(3*n), ncol=3) p <- p / colSums(p) g <- sample(1:3, n, replace=TRUE) itriplot(p, group=g)
print the installed version of R/qtlcharts
qtlchartsversion()
qtlchartsversion()
Character string with version number
qtlchartsversion()
qtlchartsversion()
A pair of linked scatterplots, where each point the first scatterplot corresponds to a scatter of points in the second scatterplot. The first scatterplot corresponds to a pair of summary measures for a larger dataset.
scat2scat(scat1data, scat2data, group = NULL, chartOpts = NULL, digits = 5)
scat2scat(scat1data, scat2data, group = NULL, chartOpts = NULL, digits = 5)
scat1data |
Matrix with two columns; rownames are used as identifiers. Can have an optional third column with categories for coloring points in the first scatterplot (to be used if 'group' is not provided). |
scat2data |
List of matrices each with at least two columns, to be shown in the second scatterplot. The components of the list correspond to the rows in 'scat1dat'. An option third column can contain categories. Row names identify individual points. |
group |
Categories for coloring points in the first scatterplot; length should be the number of rows in 'scat1data'. |
chartOpts |
A list of options for configuring the chart. Each element must be named using the corresponding option. |
digits |
Round data to this number of significant digits before passing to the chart function. (Use NULL to not round.) |
An object of class 'htmlwidget' that will intelligently print itself into HTML in a variety of contexts including the R console, within R Markdown documents, and within Shiny output bindings.
[iplotCorr()]
# simulate some data p <- 500 n <- 300 SD <- runif(p, 1, 5) r <- runif(p, -1, 1) scat2 <- vector("list", p) for(i in 1:p) scat2[[i]] <- matrix(rnorm(2*n), ncol=2) %*% chol(SD[i]^2*matrix(c(1, r[i], r[i], 1), ncol=2)) scat1 <- cbind(SD=SD, r=r) # plot it scat2scat(scat1, scat2)
# simulate some data p <- 500 n <- 300 SD <- runif(p, 1, 5) r <- runif(p, -1, 1) scat2 <- vector("list", p) for(i in 1:p) scat2[[i]] <- matrix(rnorm(2*n), ncol=2) %*% chol(SD[i]^2*matrix(c(1, r[i], r[i], 1), ncol=2)) scat1 <- cbind(SD=SD, r=r) # plot it scat2scat(scat1, scat2)
Set the default screen size as a global option.
setScreenSize(size = c("normal", "small", "large"), height, width)
setScreenSize(size = c("normal", "small", "large"), height, width)
size |
Character vector representing screen size (normal, small, large). Ignored if height and width are provided. |
height |
(Optional) Height in pixels |
width |
(Optional) Width in pixels |
Used to set a global option, 'qtlchartsScreenSize', that contains the maximum height and maximum width for a chart in the browser.
'"small"', '"normal"', and '"large"' correspond to 600x900, 700x1000, and 1200x1600, for height x width, respectively.
None.
setScreenSize("large")
setScreenSize("large")