Title: | Data Visualization for QTL Experiments |
---|---|
Description: | Functions to plot QTL (quantitative trait loci) analysis results and related diagnostics. Part of 'qtl2', an upgrade of the 'qtl' package to better handle high-dimensional data and complex cross designs. |
Authors: | Brian S Yandell [aut, cre], Karl W Broman [aut] |
Maintainer: | Brian S Yandell <[email protected]> |
License: | GPL-3 |
Version: | 1.2.4 |
Built: | 2024-11-11 06:02:05 UTC |
Source: | https://github.com/byandell-sysgen/qtl2ggplot |
Set up col, pattern and group for plotting.
color_patterns_get(scan1ggdata, col, palette = NULL)
color_patterns_get(scan1ggdata, col, palette = NULL)
scan1ggdata |
data frame to be used for plotting |
col |
Color for |
palette |
for colors (default |
list of colors
and shapes
.
Set up col, pattern, shape and group for plotting.
color_patterns_pheno( scan1ggdata, lod, pattern, col, shape, patterns, facet = NULL )
color_patterns_pheno( scan1ggdata, lod, pattern, col, shape, patterns, facet = NULL )
scan1ggdata |
data frame to be used for plotting |
lod |
matrix of LOD scores by position and pheno |
pattern |
allele pattern of form |
col |
Color for |
shape |
Shape for |
patterns |
Connect SDP patterns: one of |
facet |
use |
data frame scan1ggdata
with additional objects.
Set up colors for patterns or points
color_patterns_set( scan1output, snpinfo, patterns, col, pattern, show_all_snps, col_hilit, drop_hilit, maxlod )
color_patterns_set( scan1output, snpinfo, patterns, col, pattern, show_all_snps, col_hilit, drop_hilit, maxlod )
scan1output |
output of linear mixed model for |
snpinfo |
Data frame with snp information |
patterns |
Connect SDP patterns: one of |
col |
Color of other points, or colors for patterns |
pattern |
allele pattern as of form |
show_all_snps |
show all SNPs if |
col_hilit |
Color of highlighted points |
drop_hilit |
SNPs with LOD score within this amount of the maximum SNP association will be highlighted. |
maxlod |
Maximum LOD for drop of |
list of col
and pattern
.
Plot estimated QTL effects along a chromosomes.
ggplot_coef( object, map, columns = NULL, col = NULL, scan1_output = NULL, gap = 25, ylim = NULL, bgcolor = "gray90", altbgcolor = "gray85", ylab = "QTL effects", xlim = NULL, ... ) ggplot_coefCC(object, map, colors = qtl2::CCcolors, ...) ## S3 method for class 'scan1coef' autoplot(object, ...)
ggplot_coef( object, map, columns = NULL, col = NULL, scan1_output = NULL, gap = 25, ylim = NULL, bgcolor = "gray90", altbgcolor = "gray85", ylab = "QTL effects", xlim = NULL, ... ) ggplot_coefCC(object, map, colors = qtl2::CCcolors, ...) ## S3 method for class 'scan1coef' autoplot(object, ...)
object |
Estimated QTL effects ("coefficients") as obtained from
|
map |
A list of vectors of marker positions, as produced by
|
columns |
Vector of columns to plot |
col |
Vector of colors, same length as |
scan1_output |
If provided, we make a two-panel plot with coefficients on top and LOD scores below. Should have just one LOD score column; if multiple, only the first is used. |
gap |
Gap between chromosomes. |
ylim |
y-axis limits. If |
bgcolor |
Background color for the plot. |
altbgcolor |
Background color for alternate chromosomes. |
ylab |
y-axis label |
xlim |
x-axis limits. If |
... |
Additional graphics parameters. |
colors |
Colors to use for plotting. |
ggplot_coefCC()
is the same as ggplot_coef()
, but forcing
columns=1:8
and using the Collaborative Cross colors,
CCcolors
.
object of class ggplot
.
# read data iron <- qtl2::read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- qtl2::insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- qtl2::calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno[,1] covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) # calculate coefficients for chromosome 7 coef <- qtl2::scan1coef(probs[,7], pheno, addcovar=covar) # plot QTL effects ggplot2::autoplot(coef, map[7], columns=1:3)
# read data iron <- qtl2::read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- qtl2::insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- qtl2::calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno[,1] covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) # calculate coefficients for chromosome 7 coef <- qtl2::scan1coef(probs[,7], pheno, addcovar=covar) # plot QTL effects ggplot2::autoplot(coef, map[7], columns=1:3)
Plot gene locations for a genomic interval, as rectangles with gene symbol (and arrow indicating strand/direction) below.
ggplot_genes( object, xlim = NULL, minrow = 4, padding = 0.2, colors = c("black", "red3", "green4", "blue3", "orange"), ... ) ## S3 method for class 'genes' autoplot(object, ...)
ggplot_genes( object, xlim = NULL, minrow = 4, padding = 0.2, colors = c("black", "red3", "green4", "blue3", "orange"), ... ) ## S3 method for class 'genes' autoplot(object, ...)
object |
Object of class |
xlim |
x-axis limits (in Mbp) |
minrow |
Minimum number of rows of object |
padding |
Proportion to pad with white space around the object |
colors |
Vectors of colors, used sequentially and then re-used. |
... |
Optional arguments passed to |
None.
filename <- file.path("https://raw.githubusercontent.com/rqtl", "qtl2data/master/DOex", "c2_genes.rds") tmpfile <- tempfile() download.file(filename, tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) ggplot_genes(gene_tbl)
filename <- file.path("https://raw.githubusercontent.com/rqtl", "qtl2data/master/DOex", "c2_genes.rds") tmpfile <- tempfile() download.file(filename, tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) ggplot_genes(gene_tbl)
Plot genes at positions
ggplot_genes_internal( start, end, strand, rect_top, rect_bottom, colors, space, y, dir_symbol, name, xlim, xlab = "Position (Mbp)", ylab = "", bgcolor = "gray92", xat = NULL, legend.position = "none", vlines = NULL, ... )
ggplot_genes_internal( start, end, strand, rect_top, rect_bottom, colors, space, y, dir_symbol, name, xlim, xlab = "Position (Mbp)", ylab = "", bgcolor = "gray92", xat = NULL, legend.position = "none", vlines = NULL, ... )
start , end , strand , rect_top , rect_bottom , colors , space , y , dir_symbol , name , xlim
|
usual parameters |
legend.position , vlines , xlab , ylab , bgcolor , xat
|
hidden parameters |
... |
Additional graphics parameters. |
object of class ggplot
.
Plot object of class listof_scan1coeff
, which is a list of objects of class scan1coef
.
ggplot_listof_scan1coef( object, map, columns = NULL, col = NULL, scan1_output = NULL, facet = "pattern", ... ) ## S3 method for class 'listof_scan1coef' autoplot(object, ...)
ggplot_listof_scan1coef( object, map, columns = NULL, col = NULL, scan1_output = NULL, facet = "pattern", ... ) ## S3 method for class 'listof_scan1coef' autoplot(object, ...)
object |
object of class |
map |
A list of vectors of marker positions, as produced by
|
columns |
Vector of columns to plot |
col |
Vector of colors, same length as |
scan1_output |
If provided, we make a two-panel plot with coefficients on top and LOD scores below. Should have just one LOD score column; if multiple, only the first is used. |
facet |
Plot facets if multiple phenotypes and group provided (default = |
... |
arguments for |
pattern |
Use phenotype names as pattern. |
object of class ggplot
Brian S Yandell, [email protected]
Plot one individual's genome-wide genotypes
ggplot_onegeno( geno, map, ind = 1, chr = NULL, col = NULL, shift = FALSE, chrwidth = 0.5, ... )
ggplot_onegeno( geno, map, ind = 1, chr = NULL, col = NULL, shift = FALSE, chrwidth = 0.5, ... )
geno |
Imputed phase-known genotypes, as a list of matrices
(as produced by |
map |
Marker map (a list of vectors of marker positions). |
ind |
Individual to plot, either a numeric index or an ID (can be a vector). |
chr |
Selected chromosomes to plot; a vector of character strings. |
col |
Vector of colors for the different genotypes. |
shift |
If TRUE, shift the chromosomes so they all start at 0. |
chrwidth |
Total width of rectangles for each chromosome, as a fraction of the distance between them. |
... |
Additional graphics parameters |
object of class ggplot
.
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # inferred genotypes geno <- maxmarg(probs) # plot the inferred genotypes for the first individual ggplot_onegeno(geno, map, shift = TRUE) # plot the inferred genotypes for the first four individuals ggplot_onegeno(geno, map, ind=1:4)
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # inferred genotypes geno <- maxmarg(probs) # plot the inferred genotypes for the first individual ggplot_onegeno(geno, map, shift = TRUE) # plot the inferred genotypes for the first four individuals ggplot_onegeno(geno, map, ind=1:4)
Plot QTL peak locations (possibly with intervals) for multiple traits.
ggplot_peaks( peaks, map, chr = NULL, tick_height = 0.3, gap = 25, bgcolor = "gray90", altbgcolor = "gray85", ... )
ggplot_peaks( peaks, map, chr = NULL, tick_height = 0.3, gap = 25, bgcolor = "gray90", altbgcolor = "gray85", ... )
peaks |
Data frame such as that produced by
|
map |
Marker map, used to get chromosome lengths (and start and end positions). |
chr |
Selected chromosomes to plot; a vector of character strings. |
tick_height |
Height of tick marks at the peaks (a number between 0 and 1). |
gap |
Gap between chromosomes. |
bgcolor |
Background color for the plot. |
altbgcolor |
Background color for alternate chromosomes. |
... |
Additional graphics parameters |
None.
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) Xcovar <- get_x_covar(iron) # perform genome scan out <- scan1(probs, pheno, addcovar=covar, Xcovar=Xcovar) # find peaks above lod=3.5 (and calculate 1.5-LOD support intervals) peaks <- find_peaks(out, map, threshold=3.5, drop=1.5) # color peaks above 6 red; only show chromosomes with peaks plot_peaks(peaks, map) peaks$col <- (peaks$lod > 6) ggplot_peaks(peaks, map[names(map) %in% peaks$chr], col = c("blue","red"), legend.title = "LOD > 6")
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) Xcovar <- get_x_covar(iron) # perform genome scan out <- scan1(probs, pheno, addcovar=covar, Xcovar=Xcovar) # find peaks above lod=3.5 (and calculate 1.5-LOD support intervals) peaks <- find_peaks(out, map, threshold=3.5, drop=1.5) # color peaks above 6 red; only show chromosomes with peaks plot_peaks(peaks, map) peaks$col <- (peaks$lod > 6) ggplot_peaks(peaks, map[names(map) %in% peaks$chr], col = c("blue","red"), legend.title = "LOD > 6")
Plot phenotype vs genotype for a single putative QTL and a single phenotype.
ggplot_pxg( geno, pheno, sort = TRUE, SEmult = NULL, pooledSD = TRUE, jitter = 0.2, bgcolor = "gray90", seg_width = 0.4, seg_lwd = 2, seg_col = "black", hlines = NULL, hlines_col = "white", hlines_lty = 1, hlines_lwd = 1, vlines_col = "gray80", vlines_lty = 1, vlines_lwd = 3, force_labels = TRUE, alternate_labels = FALSE, omit_points = FALSE, ... ) mean_pxg(geno, pheno, dataframe = NULL)
ggplot_pxg( geno, pheno, sort = TRUE, SEmult = NULL, pooledSD = TRUE, jitter = 0.2, bgcolor = "gray90", seg_width = 0.4, seg_lwd = 2, seg_col = "black", hlines = NULL, hlines_col = "white", hlines_lty = 1, hlines_lwd = 1, vlines_col = "gray80", vlines_lty = 1, vlines_lwd = 3, force_labels = TRUE, alternate_labels = FALSE, omit_points = FALSE, ... ) mean_pxg(geno, pheno, dataframe = NULL)
geno |
Vector of genotypes, as produced by
|
pheno |
Vector of phenotypes. |
sort |
If TRUE, sort genotypes from largest to smallest. |
SEmult |
If specified, interval estimates of the within-group
averages will be displayed, as |
pooledSD |
If TRUE and |
jitter |
Amount to jitter the points horizontally, if a vector of length > 0, it is taken to be the actual jitter amounts (with values between -0.5 and 0.5). |
bgcolor |
Background color for the plot. |
seg_width |
Width of segments at the estimated within-group averages |
seg_lwd |
Line width used to plot estimated within-group averages |
seg_col |
Line color used to plot estimated within-group averages |
hlines |
Locations of horizontal grid lines. |
hlines_col |
Color of horizontal grid lines |
hlines_lty |
Line type of horizontal grid lines |
hlines_lwd |
Line width of horizontal grid lines |
vlines_col |
Color of vertical grid lines |
vlines_lty |
Line type of vertical grid lines |
vlines_lwd |
Line width of vertical grid lines |
force_labels |
If TRUE, force all genotype labels to be shown. |
alternate_labels |
If TRUE, place genotype labels in two rows |
omit_points |
If TRUE, omit the points, just plotting the averages (and, potentially, the +/- SE intervals). |
... |
Additional graphics parameters, passed to |
dataframe |
Supplied data frame, or constructed from |
object of class ggplot
.
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # inferred genotype at a 28.6 cM on chr 16 geno <- maxmarg(probs, map, chr=16, pos=28.6, return_char=TRUE) # plot phenotype vs genotype ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver))) # include +/- 2 SE intervals ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver)), SEmult=2) # plot just the means ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver)), omit_points=TRUE) # plot just the means +/- 2 SEs ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver)), omit_points=TRUE, SEmult=2)
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # inferred genotype at a 28.6 cM on chr 16 geno <- maxmarg(probs, map, chr=16, pos=28.6, return_char=TRUE) # plot phenotype vs genotype ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver))) # include +/- 2 SE intervals ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver)), SEmult=2) # plot just the means ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver)), omit_points=TRUE) # plot just the means +/- 2 SEs ggplot_pxg(geno, log10(iron$pheno[,1]), ylab=expression(log[10](Liver)), omit_points=TRUE, SEmult=2)
Plot LOD curves for a genome scan
Plot LOD curves for a genome scan
ggplot_scan1( object, map, lodcolumn = 1, chr = NULL, gap = 25, bgcolor = "gray90", altbgcolor = "gray85", ... ) ## S3 method for class 'scan1' autoplot(object, ...) ggplot_scan1_internal( map, lod, gap = 25, col = NULL, shape = NULL, pattern = NULL, facet = NULL, patterns = c("none", "all", "hilit"), chrName = "Chr", ... )
ggplot_scan1( object, map, lodcolumn = 1, chr = NULL, gap = 25, bgcolor = "gray90", altbgcolor = "gray85", ... ) ## S3 method for class 'scan1' autoplot(object, ...) ggplot_scan1_internal( map, lod, gap = 25, col = NULL, shape = NULL, pattern = NULL, facet = NULL, patterns = c("none", "all", "hilit"), chrName = "Chr", ... )
object |
Output of |
map |
Map of pseudomarker locations. |
lodcolumn |
LOD score column to plot (a numeric index, or a character string for a column name). One or more value(s) allowed. |
chr |
Selected chromosomes to plot; a vector of character strings. |
gap |
Gap between chromosomes. |
bgcolor |
Background color for the plot. |
altbgcolor |
Background color for alternate chromosomes. |
... |
Additional graphics parameters. |
lod |
Matrix of lod (or other) values. |
col |
Colors for points or lines, with labels. |
shape |
Shapes for points. |
pattern |
Use to group values for plotting (default = |
facet |
Plot facets if multiple phenotypes and group provided (default = |
patterns |
Connect SDP patterns: one of |
chrName |
Add prefix chromosome name (default |
None.
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) Xcovar <- get_x_covar(iron) # perform genome scan out <- scan1(probs, pheno, addcovar=covar, Xcovar=Xcovar) # plot the results for selected chromosomes chr <- c(2,7,8,9,15,16) ggplot_scan1(out, map, lodcolumn=1:2, chr=chr, col=c("darkslateblue","violetred"), legend.position=c(0.1,0.9)) # plot just one chromosome ggplot_scan1(out, map, chr=8, lodcolumn=1:2, col=c("darkblue","violetred")) # can also use autoplot from ggplot2 # lodcolumn can also be a column name library(ggplot2) autoplot(out, map, chr=8, lodcolumn=c("liver","spleen"), col=c("darkblue","violetred"))
# load qtl2 package for data and genoprob calculation library(qtl2) # read data iron <- read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) Xcovar <- get_x_covar(iron) # perform genome scan out <- scan1(probs, pheno, addcovar=covar, Xcovar=Xcovar) # plot the results for selected chromosomes chr <- c(2,7,8,9,15,16) ggplot_scan1(out, map, lodcolumn=1:2, chr=chr, col=c("darkslateblue","violetred"), legend.position=c(0.1,0.9)) # plot just one chromosome ggplot_scan1(out, map, chr=8, lodcolumn=1:2, col=c("darkblue","violetred")) # can also use autoplot from ggplot2 # lodcolumn can also be a column name library(ggplot2) autoplot(out, map, chr=8, lodcolumn=c("liver","spleen"), col=c("darkblue","violetred"))
Plot SNP associations, with possible expansion from distinct snps to all snps.
ggplot_snpasso( scan1output, snpinfo, genes = NULL, lodcolumn = 1, show_all_snps = TRUE, drop_hilit = NA, col_hilit = "violetred", col = "darkslateblue", ylim = NULL, gap = 25, minlod = 0, bgcolor = "gray90", altbgcolor = "gray85", ... )
ggplot_snpasso( scan1output, snpinfo, genes = NULL, lodcolumn = 1, show_all_snps = TRUE, drop_hilit = NA, col_hilit = "violetred", col = "darkslateblue", ylim = NULL, gap = 25, minlod = 0, bgcolor = "gray90", altbgcolor = "gray85", ... )
scan1output |
Output of |
snpinfo |
Data frame with SNP information with the following
columns (the last three are generally derived from with
|
genes |
Optional data frame containing gene information for the region, with columns 'start' and 'stop' in Mbp, 'strand' (as '"-"', '"+"', or 'NA'), and 'Name'. If included, a two-panel plot is produced, with SNP associations above and gene locations below. |
lodcolumn |
LOD score column to plot (a numeric index, or a character string for a column name). One or more value(s) allowed. |
show_all_snps |
If TRUE, expand to show all SNPs. |
drop_hilit |
SNPs with LOD score within this amount of the maximum SNP association will be highlighted. |
col_hilit |
Color of highlighted points |
col |
Color of other points |
ylim |
y-axis limits |
gap |
Gap between chromosomes. |
minlod |
Minimum LOD to display. (Mostly for GWAS, in which case using 'minlod=1' will greatly increase the plotting speed, since the vast majority of points would be omittted. |
bgcolor |
Background color for the plot. |
altbgcolor |
Background color for alternate chromosomes. |
... |
Additional graphics parameters. |
object of class ggplot
.
A number of graphics parameters can be passed via '...'. For example, 'bgcolor' to control the background color and 'altbgcolor' to control the background color on alternate chromosomes. 'cex' for character expansion for the points (default 0.5), 'pch' for the plotting character for the points (default 16), and 'ylim' for y-axis limits.
dirpath <- "https://raw.githubusercontent.com/rqtl/qtl2data/master/DOex" # Read DOex example cross from 'qtl2data' DOex <- subset(qtl2::read_cross2(file.path(dirpath, "DOex.zip")), chr = "2") # Download genotype probabilities tmpfile <- tempfile() download.file(file.path(dirpath, "DOex_genoprobs_2.rds"), tmpfile, quiet=TRUE) pr <- readRDS(tmpfile) unlink(tmpfile) # Download SNP info for DOex from web and read as RDS. tmpfile <- tempfile() download.file(file.path(dirpath, "c2_snpinfo.rds"), tmpfile, quiet=TRUE) snpinfo <- readRDS(tmpfile) unlink(tmpfile) snpinfo <- dplyr::rename(snpinfo, pos = pos_Mbp) # Convert to SNP probabilities snpinfo <- qtl2::index_snps(DOex$pmap, snpinfo) snppr <- qtl2::genoprob_to_snpprob(pr, snpinfo) # Scan SNPs. scan_snppr <- qtl2::scan1(snppr, DOex$pheno) # plot results ggplot_snpasso(scan_snppr, snpinfo, show_all_snps=FALSE, patterns="all", drop_hilit=1.5) # can also just type autoplot() if ggplot2 attached library(ggplot2) # plot just subset of distinct SNPs autoplot(scan_snppr, snpinfo, show_all_snps=FALSE, drop_hilit=1.5) # highlight SDP patterns in SNPs; connect with lines. autoplot(scan_snppr, snpinfo, patterns="all",drop_hilit=4) # query function for finding genes in region gene_dbfile <- system.file("extdata", "mouse_genes_small.sqlite", package="qtl2") query_genes <- qtl2::create_gene_query_func(gene_dbfile) genes <- query_genes(2, 97, 98) # plot SNP association results with gene locations autoplot(scan_snppr, snpinfo, patterns="hilit", drop_hilit=1.5, genes=genes)
dirpath <- "https://raw.githubusercontent.com/rqtl/qtl2data/master/DOex" # Read DOex example cross from 'qtl2data' DOex <- subset(qtl2::read_cross2(file.path(dirpath, "DOex.zip")), chr = "2") # Download genotype probabilities tmpfile <- tempfile() download.file(file.path(dirpath, "DOex_genoprobs_2.rds"), tmpfile, quiet=TRUE) pr <- readRDS(tmpfile) unlink(tmpfile) # Download SNP info for DOex from web and read as RDS. tmpfile <- tempfile() download.file(file.path(dirpath, "c2_snpinfo.rds"), tmpfile, quiet=TRUE) snpinfo <- readRDS(tmpfile) unlink(tmpfile) snpinfo <- dplyr::rename(snpinfo, pos = pos_Mbp) # Convert to SNP probabilities snpinfo <- qtl2::index_snps(DOex$pmap, snpinfo) snppr <- qtl2::genoprob_to_snpprob(pr, snpinfo) # Scan SNPs. scan_snppr <- qtl2::scan1(snppr, DOex$pheno) # plot results ggplot_snpasso(scan_snppr, snpinfo, show_all_snps=FALSE, patterns="all", drop_hilit=1.5) # can also just type autoplot() if ggplot2 attached library(ggplot2) # plot just subset of distinct SNPs autoplot(scan_snppr, snpinfo, show_all_snps=FALSE, drop_hilit=1.5) # highlight SDP patterns in SNPs; connect with lines. autoplot(scan_snppr, snpinfo, patterns="all",drop_hilit=4) # query function for finding genes in region gene_dbfile <- system.file("extdata", "mouse_genes_small.sqlite", package="qtl2") query_genes <- qtl2::create_gene_query_func(gene_dbfile) genes <- query_genes(2, 97, 98) # plot SNP association results with gene locations autoplot(scan_snppr, snpinfo, patterns="hilit", drop_hilit=1.5, genes=genes)
Create a list of scan1coef objects using scan1coef
.
Summary of object of class listof_scan1coef
, which is a list of objects of class scan1coef
.
Summary of object of class listof_scan1coef
, which is a list of objects of class scan1coef
.
Subset of object of class listof_scan1coef
, which is a list of objects of class scan1coef
.
listof_scan1coef( probs, phe, K = NULL, covar = NULL, blups = FALSE, center = FALSE, ... ) summary_listof_scan1coef( object, scan1_object, map, coef_names = dimnames(object[[1]])[[2]], center = TRUE, ... ) ## S3 method for class 'listof_scan1coef' summary(object, ...) summary_scan1coef(object, scan1_object, map, ...) ## S3 method for class 'scan1coef' summary(object, ...) subset_listof_scan1coef(x, elements, ...) ## S3 method for class 'listof_scan1coef' subset(x, ...) ## S3 method for class 'listof_scan1coef' x[...]
listof_scan1coef( probs, phe, K = NULL, covar = NULL, blups = FALSE, center = FALSE, ... ) summary_listof_scan1coef( object, scan1_object, map, coef_names = dimnames(object[[1]])[[2]], center = TRUE, ... ) ## S3 method for class 'listof_scan1coef' summary(object, ...) summary_scan1coef(object, scan1_object, map, ...) ## S3 method for class 'scan1coef' summary(object, ...) subset_listof_scan1coef(x, elements, ...) ## S3 method for class 'listof_scan1coef' subset(x, ...) ## S3 method for class 'listof_scan1coef' x[...]
probs |
genotype probabilities object for one chromosome from |
phe |
data frame of phenotypes |
K |
list of length 1 with kinship matrix |
covar |
matrix of covariates |
blups |
Create BLUPs if |
center |
center coefficients if |
... |
ignored |
object |
object of class |
scan1_object |
object from |
map |
A list of vectors of marker positions, as produced by
|
coef_names |
names of effect coefficients (default is all coefficient names) |
x |
object of class |
elements |
indexes or names of list elements in x |
object of class listof_scan1coef
Brian S Yandell, [email protected]
# read data iron <- qtl2::read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- qtl2::insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- qtl2::calc_genoprob(iron, map, error_prob=0.002) # Ensure that covariates have names attribute covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) # Calculate scan1coef on all phenotypes, # returning a list of \code{\link{scan1coef}} objects out <- listof_scan1coef(probs[,7], iron$pheno, addcovar = covar, center = TRUE) # Plot coefficients for all phenotypes ggplot2::autoplot(out, map[7], columns = 1:3) # Summary of coefficients at scan peak scan_pr <- qtl2::scan1(probs[,7], iron$pheno) summary(out, scan_pr, map[7])
# read data iron <- qtl2::read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- qtl2::insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- qtl2::calc_genoprob(iron, map, error_prob=0.002) # Ensure that covariates have names attribute covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) # Calculate scan1coef on all phenotypes, # returning a list of \code{\link{scan1coef}} objects out <- listof_scan1coef(probs[,7], iron$pheno, addcovar = covar, center = TRUE) # Plot coefficients for all phenotypes ggplot2::autoplot(out, map[7], columns = 1:3) # Summary of coefficients at scan peak scan_pr <- qtl2::scan1(probs[,7], iron$pheno) summary(out, scan_pr, map[7])
Convert strain distribution pattern (sdp) to letter pattern. Taken from package 'qtl2pattern' for internal use here.
sdp_to_pattern(sdp, haplos, symmetric = TRUE)
sdp_to_pattern(sdp, haplos, symmetric = TRUE)
sdp |
vector of sdp values |
haplos |
letter codes for haplotypes (required) |
symmetric |
make patterns symmetric if |
vector of letter patterns
Brian S Yandell, [email protected]
Summary of scan1 object
summary_scan1( object, map, snpinfo = NULL, lodcolumn = seq_len(ncol(object)), chr = names(map), sum_type = c("common", "best"), drop = 1.5, show_all_snps = TRUE, ... ) ## S3 method for class 'scan1' summary(object, ...)
summary_scan1( object, map, snpinfo = NULL, lodcolumn = seq_len(ncol(object)), chr = names(map), sum_type = c("common", "best"), drop = 1.5, show_all_snps = TRUE, ... ) ## S3 method for class 'scan1' summary(object, ...)
object |
object from |
map |
A list of vectors of marker positions, as produced by
|
snpinfo |
Data frame with SNP information with the following
columns (the last three are generally derived from with
|
lodcolumn |
one or more lod columns |
chr |
one or more chromosome IDs |
sum_type |
type of summary |
drop |
LOD drop from maximum |
show_all_snps |
show all SNPs if |
... |
other arguments not used |
tbl summary
Brian S Yandell, [email protected]
# read data iron <- qtl2::read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- qtl2::insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- qtl2::calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) Xcovar <- qtl2::get_x_covar(iron) # perform genome scan out <- qtl2::scan1(probs, pheno, addcovar=covar, Xcovar=Xcovar) # summary summary(out, map)
# read data iron <- qtl2::read_cross2(system.file("extdata", "iron.zip", package="qtl2")) # insert pseudomarkers into map map <- qtl2::insert_pseudomarkers(iron$gmap, step=1) # calculate genotype probabilities probs <- qtl2::calc_genoprob(iron, map, error_prob=0.002) # grab phenotypes and covariates; ensure that covariates have names attribute pheno <- iron$pheno covar <- match(iron$covar$sex, c("f", "m")) # make numeric names(covar) <- rownames(iron$covar) Xcovar <- qtl2::get_x_covar(iron) # perform genome scan out <- qtl2::scan1(probs, pheno, addcovar=covar, Xcovar=Xcovar) # summary summary(out, map)