Title: | Pattern Support for 'qtl2' Package |
---|---|
Description: | Routines in 'qtl2' to study allele patterns in quantitative trait loci (QTL) mapping over a chromosome. Useful in crosses with more than two alleles to identify how sets of alleles, genetically different strands at the same locus, have different response levels. Plots show profiles over a chromosome. Can handle multiple traits together. See <https://github.com/byandell/qtl2pattern>. |
Authors: | Brian S Yandell [aut, cre] |
Maintainer: | Brian S Yandell <[email protected]> |
License: | GPL-3 |
Version: | 1.2.2 |
Built: | 2025-02-09 05:16:12 UTC |
Source: | https://github.com/byandell-sysgen/qtl2pattern |
Create table of alleles for various model fits.
Plot alleles for haplotype, diplotype and top patterns and genome position.
allele1( probD, phe_df = NULL, cov_mx = NULL, map = NULL, K_chr = NULL, patterns = NULL, alt = NULL, blups = FALSE, ... ) ggplot_allele1( object, scan1_object = NULL, map = NULL, pos = NULL, trim = TRUE, legend.position = "none", ... ) ## S3 method for class 'allele1' autoplot(object, ...)
allele1( probD, phe_df = NULL, cov_mx = NULL, map = NULL, K_chr = NULL, patterns = NULL, alt = NULL, blups = FALSE, ... ) ggplot_allele1( object, scan1_object = NULL, map = NULL, pos = NULL, trim = TRUE, legend.position = "none", ... ) ## S3 method for class 'allele1' autoplot(object, ...)
probD |
object of class |
phe_df |
data frame with one phenotype |
cov_mx |
covariate matrix |
map |
Genome map (required if |
K_chr |
kinship matrix |
patterns |
data frame of pattern information |
alt |
Haplotype allele letter(s) for alternative to reference. |
blups |
Create BLUPs if |
... |
Other parameters ignored. |
object |
Object of class |
scan1_object |
Optional object of class |
pos |
Genome position in Mbp (supercedes |
trim |
If |
legend.position |
Legend position (default is |
Table with allele effects across sources.
object of class ggplot
Create a function that will connect to a database of genotype probability information and return a list with 'probs' object and a 'map' object.
create_probs_query_func(dbfile, method_val = "fst", probdir_val = "genoprob")
create_probs_query_func(dbfile, method_val = "fst", probdir_val = "genoprob")
dbfile |
Name of database file |
method_val |
either |
probdir_val |
name of probability directory (default |
Note that this function assumes that probdir_val
has a file with the
physical map with positions in Mbp and other files with genotype probabilities.
See read_probs
for details on how probabilities are read.
See create_variant_query_func
for original idea.
Function with six arguments, 'chr', 'start',
'end', 'allele', 'method' and 'probdir'. It returns a list with 'probs' and 'map' objects
spanning the region specified by the first three arguments.
The 'probs' element should be either a 'calc_genoprob' or 'fst_genoprob' object
(see fst_genoprob
).
dirpath <- "https://raw.githubusercontent.com/rqtl/qtl2data/master/DOex" create_qv <- function(dirpath) { # Download SNP info for DOex from web via RDS. # snpinfo is referenced internally in the created function. 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) function(chr, start, end) { if(chr != "2") return(NULL) if(start < 96.5) start <- 96.5 if(end > 98.5) end <- 98.5 if(start >= end) return(NULL) dplyr::filter(snpinfo, .data$pos >= start, .data$pos <= end) } } query_variants <- create_qv(dirpath) create_qg <- function(dirpath) { # Download Gene info for DOex from web via RDS # gene_tbl is referenced internally in the created function. tmpfile <- tempfile() download.file(file.path(dirpath, "c2_genes.rds"), tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) function(chr, start, end) { if(chr != "2") return(NULL) if(start < 96.5) start <- 96.5 if(end > 98.5) end <- 98.5 if(start >= end) return(NULL) dplyr::filter(gene_tbl, .data$end >= start, .data$start <= end) } } query_genes <- create_qg(dirpath) # Examples for probs require either FST or RDS storage of data.
dirpath <- "https://raw.githubusercontent.com/rqtl/qtl2data/master/DOex" create_qv <- function(dirpath) { # Download SNP info for DOex from web via RDS. # snpinfo is referenced internally in the created function. 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) function(chr, start, end) { if(chr != "2") return(NULL) if(start < 96.5) start <- 96.5 if(end > 98.5) end <- 98.5 if(start >= end) return(NULL) dplyr::filter(snpinfo, .data$pos >= start, .data$pos <= end) } } query_variants <- create_qv(dirpath) create_qg <- function(dirpath) { # Download Gene info for DOex from web via RDS # gene_tbl is referenced internally in the created function. tmpfile <- tempfile() download.file(file.path(dirpath, "c2_genes.rds"), tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) function(chr, start, end) { if(chr != "2") return(NULL) if(start < 96.5) start <- 96.5 if(end > 98.5) end <- 98.5 if(start >= end) return(NULL) dplyr::filter(gene_tbl, .data$end >= start, .data$start <= end) } } query_genes <- create_qg(dirpath) # Examples for probs require either FST or RDS storage of data.
Match up exon start,stop,strand with genes. Use query_genes
to find features; see create_gene_query_func
.
Returns table of gene and its exons.
Uses gene_exon
to plot genes, exons, mRNA with SNPs.
gene_exon( top_snps_tbl, feature_tbl = query_genes(chr_id, range_Mbp[1], range_Mbp[2]) ) ## S3 method for class 'gene_exon' summary(object, gene_name = NULL, top_snps_tbl = NULL, extra = 0.005, ...) ## S3 method for class 'gene_exon' subset(x, gene_val, ...) ggplot_gene_exon( object, top_snps_tbl = NULL, plot_now = TRUE, genes = unique(object$gene), ... ) ## S3 method for class 'gene_exon' autoplot(object, ...)
gene_exon( top_snps_tbl, feature_tbl = query_genes(chr_id, range_Mbp[1], range_Mbp[2]) ) ## S3 method for class 'gene_exon' summary(object, gene_name = NULL, top_snps_tbl = NULL, extra = 0.005, ...) ## S3 method for class 'gene_exon' subset(x, gene_val, ...) ggplot_gene_exon( object, top_snps_tbl = NULL, plot_now = TRUE, genes = unique(object$gene), ... ) ## S3 method for class 'gene_exon' autoplot(object, ...)
top_snps_tbl |
table from |
feature_tbl |
table of features from |
object |
Object of class |
gene_name |
name of gene as character string |
extra |
extra region beyond gene for SNPs (in Mbp) |
... |
arguments passed along to |
x |
Object of class |
gene_val |
Name of gene from object |
plot_now |
plot now if TRUE |
genes |
Names of genes in |
tbl of exon and gene features
tbl of summary
list of ggplots (see gene_exon
)
Brian S Yandell, [email protected]
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) # Collect top SNPs top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) # Download Gene info for DOex from web via RDS tmpfile <- tempfile() download.file(file.path(dirpath, "c2_genes.rds"), tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) # Get Gene exon information. out <- gene_exon(top_snps_tbl, gene_tbl) summary(out, gene = out$gene[1])
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) # Collect top SNPs top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) # Download Gene info for DOex from web via RDS tmpfile <- tempfile() download.file(file.path(dirpath, "c2_genes.rds"), tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) # Get Gene exon information. out <- gene_exon(top_snps_tbl, gene_tbl) summary(out, gene = out$gene[1])
Collapse genoprob according to pattern
genoprob_to_patternprob(probs1, sdp, alleles = FALSE)
genoprob_to_patternprob(probs1, sdp, alleles = FALSE)
probs1 |
object of class |
sdp |
SNP distribution pattern |
alleles |
use allele string if |
object of class calc_genoprob
Brian S Yandell, [email protected]
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) # Convert genotype probabilities to pattern probabilities for pattern 1. pattern_pr <- genoprob_to_patternprob(pr, 7, TRUE) str(pr) str(pattern_pr)
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) # Convert genotype probabilities to pattern probabilities for pattern 1. pattern_pr <- genoprob_to_patternprob(pr, 7, TRUE) str(pr) str(pattern_pr)
Find features that overlap with SNPs
get_feature_snp(snp_tbl, feature_tbl, extend = 0.005)
get_feature_snp(snp_tbl, feature_tbl, extend = 0.005)
snp_tbl |
tbl of SNPs from |
feature_tbl |
tbl of feature information from |
extend |
extend region for SNPs in Mbp (default 0.005) |
tbl of features covering SNPs
Brian S Yandell, [email protected]
Internal routine to find features that overlap with SNPs
get_gene_snp( snp_tbl, feature_tbl, feature_snp = get_feature_snp(snp_tbl, feature_tbl, 0) )
get_gene_snp( snp_tbl, feature_tbl, feature_snp = get_feature_snp(snp_tbl, feature_tbl, 0) )
snp_tbl |
tbl of SNPs from |
feature_tbl |
tbl of feature information from |
feature_snp |
tbl of feature information from |
tbl of genes covering SNPs
Brian S Yandell, [email protected]
Figure out gene locations to make room for gene names. Written original by Dan Gatti 2013-02-13
get.gene.locations( locs, xlim, text_size = 3, str_rect = c("iW", "i"), n_rows = 10, plot_width = 6, ... )
get.gene.locations( locs, xlim, text_size = 3, str_rect = c("iW", "i"), n_rows = 10, plot_width = 6, ... )
locs |
tbl of gene information |
xlim |
X axis limits |
text_size |
size of text (default 3) |
str_rect |
character spacing on left and right of rectangles (default c("iW","i")) |
n_rows |
desired number of rows (default 10) |
plot_width |
width of default plot window (in inches) |
... |
additional parameters (not used) |
list object used by ggplot_feature_tbl
Brian S Yandell, [email protected] Daniel Gatti, [email protected]
https://github.com/dmgatti/DOQTL/blob/master/R/gene.plot.R
Merge all SNPs in small region with LOD peaks across multiple phenotype.
ggplot_merge_feature(object, pheno, plot_by = c("pattern", "consequence"), ...) ## S3 method for class 'merge_feature' autoplot(object, ...) merge_feature( top_snps_tbl, snpinfo, out_lmm_snps, drop = 1.5, dropchar = 0, exons = gene_exon(top_snps_tbl) ) ## S3 method for class 'merge_feature' summary(object, sum_type = c("SNP type", "pattern"), ...)
ggplot_merge_feature(object, pheno, plot_by = c("pattern", "consequence"), ...) ## S3 method for class 'merge_feature' autoplot(object, ...) merge_feature( top_snps_tbl, snpinfo, out_lmm_snps, drop = 1.5, dropchar = 0, exons = gene_exon(top_snps_tbl) ) ## S3 method for class 'merge_feature' summary(object, sum_type = c("SNP type", "pattern"), ...)
object |
of class |
pheno |
name of phenotype to be plotted |
plot_by |
element to plot by (one of |
... |
other arguments not used |
top_snps_tbl |
tbl from |
snpinfo |
SNP information table |
out_lmm_snps |
tbl from |
drop |
include LOD scores within |
dropchar |
number of characters to drop on phenames |
exons |
table from |
sum_type |
one of |
ggplot2 object
tbl with added information on genes and exons
table summary
Brian S Yandell, [email protected]
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) # Collect top SNPs top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) summary(top_snps_tbl) # Download Gene info for DOex from web via RDS tmpfile <- tempfile() download.file(file.path(dirpath, "c2_genes.rds"), tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) out <- merge_feature(top_snps_tbl, snpinfo, scan_snppr, exons = gene_tbl) summary(out, "pattern")
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) # Collect top SNPs top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) summary(top_snps_tbl) # Download Gene info for DOex from web via RDS tmpfile <- tempfile() download.file(file.path(dirpath, "c2_genes.rds"), tmpfile, quiet=TRUE) gene_tbl <- readRDS(tmpfile) unlink(tmpfile) out <- merge_feature(top_snps_tbl, snpinfo, scan_snppr, exons = gene_tbl) summary(out, "pattern")
Plot scan pattern usign ggplot2
Genome scan by pattern set
ggplot_scan1pattern( object, map, plot_type = c("lod", "coef", "coef_and_lod"), patterns = object$patterns$founders, columns = 1:3, min_lod = 3, lodcolumn = seq_along(patterns), facet = "pheno", ... ) ## S3 method for class 'scan1pattern' autoplot(object, ...) scan1pattern( probs1, phe, K = NULL, covar = NULL, map, patterns, condense_patterns = TRUE, blups = FALSE, do_scans = TRUE ) ## S3 method for class 'scan1pattern' summary(object, map, ...)
ggplot_scan1pattern( object, map, plot_type = c("lod", "coef", "coef_and_lod"), patterns = object$patterns$founders, columns = 1:3, min_lod = 3, lodcolumn = seq_along(patterns), facet = "pheno", ... ) ## S3 method for class 'scan1pattern' autoplot(object, ...) scan1pattern( probs1, phe, K = NULL, covar = NULL, map, patterns, condense_patterns = TRUE, blups = FALSE, do_scans = TRUE ) ## S3 method for class 'scan1pattern' summary(object, map, ...)
object |
object of class |
map |
genome map |
plot_type |
type of plot from |
patterns |
data frame of pattern information |
columns |
columns used for coef plot |
min_lod |
minimum LOD peak for contrast to be retained |
lodcolumn |
columns used for scan1 plot (default all |
facet |
Plot facets if multiple phenotypes and patterns provided (default = |
... |
additional parameters passed on to other methods |
probs1 |
object of class |
phe |
data frame with one phenotype |
K |
kinship matrix |
covar |
covariate matrix |
condense_patterns |
remove snp_action from contrasts if TRUE |
blups |
Create BLUPs if |
do_scans |
Do scans if |
object of class ggplot
List containing:
patterns Data frame of summary for top patterns (column founders
has pattern)
dip_set Diplotype sets for contrasts
group Group for each founder pattern
scan Object of class scan1
.
coef Object of class listof_scan1coef
. See package 'qtl2ggplot'.
Brian S Yandell, [email protected]
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) top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) # Summarize to find top patterns patterns <- dplyr::arrange(summary(top_snps_tbl), dplyr::desc(max_lod)) # Scan using patterns. scan_pat <- scan1pattern(pr, DOex$pheno, map = DOex$pmap, patterns = patterns) # Summary of scan1pattern. summary(scan_pat, DOex$pmap)
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) top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) # Summarize to find top patterns patterns <- dplyr::arrange(summary(top_snps_tbl), dplyr::desc(max_lod)) # Scan using patterns. scan_pat <- scan1pattern(pr, DOex$pheno, map = DOex$pmap, patterns = patterns) # Summary of scan1pattern. summary(scan_pat, DOex$pmap)
Extract pattern of diplotypes
Extract pattern of haplotypes
pattern_diplos(sdp, haplos, diplos, cont = NULL) pattern_haplos(sdp, haplos)
pattern_diplos(sdp, haplos, diplos, cont = NULL) pattern_haplos(sdp, haplos)
sdp |
vector of sdp from |
haplos |
vector of haplotype names |
diplos |
vector of diplotype names |
cont |
vector of types of contrasts ( |
matrix of diplotype patterns
matrix of haplotype patterns
Brian S Yandell, [email protected]
Turn genotype probabilities into labels
pattern_label(genos, allele = TRUE) pattern_sdp(label, sdp = NULL, geno_names = sort(unique(label)))
pattern_label(genos, allele = TRUE) pattern_sdp(label, sdp = NULL, geno_names = sort(unique(label)))
genos |
matrix of genotype probabilities at locus |
allele |
Driver has alleles if |
label |
character string from |
sdp |
SNP distribution pattern for plot colors |
geno_names |
unique genotype names (alleles or allele pairs) |
character vector of genotype names.
Read fast database with format fst
. Use first column
of database (must be named 'ind') as rownames if desired. R/qtl2 routines assume data frames have
rownames to use to align individuals.
read_fast(datapath, columns = NULL, rownames = TRUE)
read_fast(datapath, columns = NULL, rownames = TRUE)
datapath |
character string path to database |
columns |
names or indexes for columns to be extracted |
rownames |
use first column of rownames if |
extracted data frame with appropriate rows and columns.
Read object from file stored according to method.
read_probs( chr = NULL, start_val = NULL, end_val = NULL, datapath, allele = TRUE, method, probdir = "genoprob" )
read_probs( chr = NULL, start_val = NULL, end_val = NULL, datapath, allele = TRUE, method, probdir = "genoprob" )
chr |
vector of chromosome identifiers |
start_val , end_val
|
start and end values in Mbp |
datapath |
name of folder with Derived Data |
allele |
read haplotype allele probabilities (if |
method |
method of genoprob storage |
probdir |
genotype probability directory (default |
list with probs
= large object of class calc_genoprob
and map
= physical map for selected chr
Brian S Yandell, [email protected]
Convert strain distribution pattern (sdp) to letter pattern.
sdp_to_pattern(sdp, haplos, symmetric = TRUE) sdp_to_logical(sdp, haplos, symmetric = TRUE)
sdp_to_pattern(sdp, haplos, symmetric = TRUE) sdp_to_logical(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]
dirpath <- "https://raw.githubusercontent.com/rqtl/qtl2data/master/DOex" # 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) # Extract strain distribution pattern. sdp <- snpinfo$sdp # Find out how many alleles. nallele <- ceiling(log2(max(sdp))) out <- sdp_to_pattern(sdp, LETTERS[seq_len(nallele)]) # Show most frequent patterns. head(rev(sort(c(table(out)))))
dirpath <- "https://raw.githubusercontent.com/rqtl/qtl2data/master/DOex" # 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) # Extract strain distribution pattern. sdp <- snpinfo$sdp # Find out how many alleles. nallele <- ceiling(log2(max(sdp))) out <- sdp_to_pattern(sdp, LETTERS[seq_len(nallele)]) # Show most frequent patterns. head(rev(sort(c(table(out)))))
Convert SNP info to map
snpinfo_to_map(snpinfo)
snpinfo_to_map(snpinfo)
snpinfo |
Data frame with SNP information with the following
columns (the last three are generally derived from with
|
map as list of vectors of marker positions.
Collapse genoprob according to pattern
snpprob_collapse( snpprobs, action = c("additive", "add+dom", "non-add", "recessive", "dominant", "basic") )
snpprob_collapse( snpprobs, action = c("additive", "add+dom", "non-add", "recessive", "dominant", "basic") )
snpprobs |
object of class |
action |
SNP gene action type |
object of class calc_genoprob
Brian S Yandell, [email protected]
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) dim(snppr[[1]]) dim(snpprob_collapse(snppr, "additive")[[1]])
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) dim(snppr[[1]]) dim(snpprob_collapse(snppr, "additive")[[1]])
Summary of features with SNP information
## S3 method for class 'feature_snp' summary(object, ...)
## S3 method for class 'feature_snp' summary(object, ...)
object |
tbl of feature information from |
... |
additional parameters ignored |
tbl of feature summaries by type
Brian S Yandell, [email protected]
Show count min and max of features by type
Plot genes as rectangles followed by names. Stagger genes for easy reading. Written original by Dan Gatti 2013-02-13
## S3 method for class 'feature_tbl' summary(object, major = TRUE, ...) ## S3 method for class 'feature_tbl' subset(x, start_val = 0, stop_val = max(x$stop), ...) ggplot_feature_tbl( object, rect_col = "grey70", strand_col = c(`-` = "#1b9e77", `+` = "#d95f02"), type_col = c(gene = "black", pseudogene = "#1b9e77", other = "#d95f02"), text_size = 3, xlim = NULL, snp_pos = top_snps_tbl$pos, snp_lod = top_snps_tbl$lod, top_snps_tbl = NULL, snp_col = "grey70", extend = 0.005, ... ) ## S3 method for class 'feature_tbl' autoplot(object, ...)
## S3 method for class 'feature_tbl' summary(object, major = TRUE, ...) ## S3 method for class 'feature_tbl' subset(x, start_val = 0, stop_val = max(x$stop), ...) ggplot_feature_tbl( object, rect_col = "grey70", strand_col = c(`-` = "#1b9e77", `+` = "#d95f02"), type_col = c(gene = "black", pseudogene = "#1b9e77", other = "#d95f02"), text_size = 3, xlim = NULL, snp_pos = top_snps_tbl$pos, snp_lod = top_snps_tbl$lod, top_snps_tbl = NULL, snp_col = "grey70", extend = 0.005, ... ) ## S3 method for class 'feature_tbl' autoplot(object, ...)
object |
tbl of gene information from |
major |
if |
... |
additional arguments (not used) |
x |
tbl of feature information from |
start_val , stop_val
|
start and stop positions for subset |
rect_col |
fill color of rectangle (default "grey70") |
strand_col |
edge color of rectangle by strand from |
type_col |
color of type from |
text_size |
size of text (default 3) |
xlim |
horizontal axis limits (default is range of features) |
snp_pos |
position of SNPs in bp if used (default NULL) |
snp_lod |
LOD of SNPs (for color plotting) |
top_snps_tbl |
table from |
snp_col |
color of SNP vertical lines (default "grey70") |
extend |
extend region for SNPs in bp (default 0.005) |
tbl of feature summaries by type
tbl of feature summaries by type
data frame of gene information (invisible)
Brian S Yandell, [email protected]
Brian S Yandell, [email protected] Daniel Gatti, [email protected]
https://github.com/dmgatti/DOQTL/blob/master/R/gene.plot.R
Summary of genes overlapping SNPs
## S3 method for class 'gene_snp' summary(object, ...)
## S3 method for class 'gene_snp' summary(object, ...)
object |
tbl of feature information from |
... |
additional parameters ignored |
tbl of feature summaries by type
Brian S Yandell, [email protected]
Separate fine mapping scans by allele pattern.
top_snps_pattern( scan1_output, snpinfo, drop = 1.5, show_all_snps = TRUE, haplos ) ## S3 method for class 'top_snps_pattern' summary(object, sum_type = c("range", "best", "peak"), ...) ## S3 method for class 'top_snps_pattern' subset(x, start_val = 0, end_val = max(x$pos), pheno = NULL, ...)
top_snps_pattern( scan1_output, snpinfo, drop = 1.5, show_all_snps = TRUE, haplos ) ## S3 method for class 'top_snps_pattern' summary(object, sum_type = c("range", "best", "peak"), ...) ## S3 method for class 'top_snps_pattern' subset(x, start_val = 0, end_val = max(x$pos), pheno = NULL, ...)
scan1_output |
output of linear mixed model for |
snpinfo |
Data frame with SNP information with the following
columns (the last three are generally derived from with
|
drop |
include all SNPs within |
show_all_snps |
show all SNPs if |
haplos |
optional argument identify codes for haplotypes |
object |
object of class |
sum_type |
type of summary (one of "range","best") |
... |
additional parameters ignored |
x |
tbl of feature information from |
start_val , end_val
|
start and end positions for subset |
pheno |
phenotype name(s) for subset |
table of top_snps at maximum lod for pattern
table summary
subset of x
Brian S Yandell, [email protected]
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) # Collect top SNPs top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) summary(top_snps_tbl)
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) # Collect top SNPs top_snps_tbl <- top_snps_pattern(scan_snppr, snpinfo) summary(top_snps_tbl)