plotIBDproportions: Plot The Proportion of Pairs IBD

Description Usage Arguments See Also Examples

Description

plotIBDproportions() plots the proportion of pairs IBD for each SNP across the genome.

Usage

1
2
3
4
5
6
7
plotIBDproportions(ibd.proportions, interval = NULL,
  annotation.genes = NULL, annotation.genes.color = NULL,
  highlight.genes = NULL, highlight.genes.labels = TRUE,
  highlight.genes.color = NULL, highlight.genes.alpha = 0.1,
  line.color = NULL, add.rug = TRUE, plot.title = NULL,
  add.legend = TRUE, facet.label = TRUE, facet.scales = "fixed",
  subpop.facet = FALSE)

Arguments

ibd.proportions

A data frame containing the proportion of pairs IBD at each SNP. See the returned Value in getIBDproportion for more details. If multiple subpopulations are specified (column name "subpop") then the proportions for each subpopulation will be plotted, either in a single facet or over multiple facets. See http://docs.ggplot2.org/current/facet_grid.html on faceting. If multiple populations (column name "pop") are specified then the proportions for each population will be plotted on a separate facet, with all subpopulations in a single facet. If there are many populations or subpopulations (>8) it may be better to subset the populations to those of interest before plotting. Genomic locations of annotation genes can be included in the figure and specific genes or regions can be highlighted.

interval

A vector of length 3 containing the genomic locations of a specific region to plot. This vector should contain the chromosome ID, the start of the interval in base-pairs and the end of the interval in base-pairs; in this order respectively. The default is interval=NULL which will plot the proportions over all chromosomes in ibd.proportions.

annotation.genes

A data frame containing information on annotation genes to be included in the figure. This data frame must have at least 5 columns of information:

  1. Chromosome (type "numeric" or "integer")

  2. Gene name (type "character")

  3. Start location of the gene in base-pairs (type "numeric" or "integer")

  4. End location of the gene in base-pairs (type "numeric" or "integer")

  5. Gene strand (+ or -) (type "character")

annotation.genes must contain the following headers chr, name, start, end and strand. This data frame does not have to be in a specific order, however it must contain all of the above information with respective labels. The default is annotation.genes=NULL.

annotation.genes.color

A vector of characters or numeric values containing the two colors representing gene stand (positive (+) or negative (-))

highlight.genes

A data frame containing information of genes or regions to highlight. The data frame must have at least 4 columns of information:

  1. Chromosome (type "numeric" or "integer")

  2. Gene name (type "character")

  3. Start location of the gene in base-pairs (type "numeric" or "integer")

  4. End location of the gene in base-pairs (type "numeric" or "integer")

highlight.genes should contain the following headers chr, name, start and end. This data frame does not have to be in a specific order, however it must contain all of the above information with respective labels. The default is highlight.genes=NULL.

highlight.genes.labels

Logical. Whether to include gene names as labels in the figure. The default is highlight.genes.labels=FALSE.

highlight.genes.color

Character string or numeric value. A single color that will be used to highlight a region/gene. The default is highlight.genes.color=NULL.

highlight.genes.alpha

Numeric. A single value between 0 and 1 indicating the gene color transparency. The default is highlight.genes.alpha=0.1.

line.color

A vector of characters or numeric values denoting the color of lines to be plotted. If there are multiple populations/subpopulations then the number of colors specified should equal the number of unique populations/subpopulations combinations. The default is line.color=NULL which will use isoRelate default colors.

add.rug

Logical. Whether to include SNP positions as a rug in the figure. The default is add.rug=FALSE

plot.title

A character string of a title to be added to the figure The default is plot.title=NULL which does not add a title to the plot.

add.legend

Logical. Whether a legend containing subpopulation information should be plotted. The default is add.legend=FALSE.

facet.label

Logical. Whether to include facet labels if multiple populations/subpopulations (column names "pop" and "subpop") are specified.

facet.scales

A character string of either "fixed", "free", "free_x" or "free_y" specifying the facet axis-scales. The default is facet.scales="fixed"

subpop.facet

Logical. Whether to plot subpopulations in separate facets. The default is subpop.facet=FALSE. If subpop.facet=TRUE and there are multiple populations then subpopulations will **not** be drawn in separate facets.

See Also

getIBDproportion

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
# generate a binary IBD matrix
my_matrix <- getIBDmatrix(ped.genotypes = png_genotypes,
                          ibd.segments = png_ibd)

# calculate the proportion of pairs IBD at each SNP
my_proportion <- getIBDproportion(ped.genotypes = png_genotypes,
                                  ibd.matrix = my_matrix,
                                  groups = NULL)

# plot the proportion of pairs IBD
plotIBDproportions(ibd.proportions = my_proportion,
                   interval = NULL,
                   annotation.genes = NULL,
                   annotation.genes.color = NULL,
                   highlight.genes = NULL,
                   highlight.genes.labels = TRUE,
                   highlight.genes.color = NULL,
                   highlight.genes.alpha = 0.1,
                   add.rug = FALSE,
                   plot.title = "Proportion of pairs IBD in PNG",
                   add.legend = FALSE,
                   line.color = NULL,
                   facet.label = TRUE,
                   facet.scales = "fixed",
                   subpop.facet = FALSE)

# creating a stratification dataset
my_groups <- png_genotypes[[1]][,1:3]
my_groups[1:10,"pid"] <- "a"
my_groups[11:25,"pid"] <- "b"
my_groups[26:38,"pid"] <- "c"

my_proportion <- getIBDproportion(ped.genotypes = png_genotypes,
                                  ibd.matrix = my_matrix,
                                  groups = my_groups)

# plot the proportion of pairs IBD
plotIBDproportions(ibd.proportions = my_proportion,
                   interval = NULL,
                   annotation.genes = NULL,
                   annotation.genes.color = NULL,
                   highlight.genes = NULL,
                   highlight.genes.labels = FALSE,
                   highlight.genes.color = NULL,
                   highlight.genes.alpha = 0.1,
                   line.color = NULL,
                   add.rug = FALSE,
                   plot.title = "Proportion of pairs IBD in PNG - with stratification",
                   add.legend = FALSE,
                   facet.label = TRUE,
                   facet.scales = "fixed",
                   subpop.facet = TRUE)

bahlolab/isoRelate documentation built on May 11, 2019, 5:25 p.m.