Package 'OmicCircos'

Title: High-quality circular visualization of omics data
Description: OmicCircos is an R application and package for generating high-quality circular plots for omics data.
Authors: Ying Hu <[email protected]> Chunhua Yan <[email protected]>
Maintainer: Ying Hu <[email protected]>
License: GPL-2
Version: 1.43.0
Built: 2024-07-03 06:02:21 UTC
Source: https://github.com/bioc/OmicCircos

Help Index


OmicCircos package

Description

OmicCircos is for generating high-quality circular plots for omics data.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>

References

OmicCircos: an R package for simple and circular visualization of omics data. Cancer Inform. 2014 Jan 16;13:13-20. doi: 10.4137/CIN.S13495. eCollection 2014. PMID: 24526832 [PubMed] PMCID: PMC3921174


draw circular

Description

This is the main function of OmicCircos to draw circular plots.

Usage

circos(mapping=mapping, xc=400, yc=400, R=400, W=W,  cir="", type="n", 
    col.v=3, B=FALSE, print.chr.lab=TRUE, col.bar=FALSE, col.bar.po = "topleft",       
    cluster=FALSE, order=NULL, scale=FALSE, cutoff = "n", zoom="", cex=1, lwd=1, 
    col=rainbow(10, alpha=0.8)[7], side="")

Arguments

mapping

data frame or matrix containing mapping information and values. Column 1: the segment or chromosome ID; column 2: the position; column 3: the position label (optional) or the value and additional columns are the values. such as gene expression and copy number. Missing values are allowed and will be ignored.

xc

integer, the circle center x coordinate

yc

integer, the circle center y coordinate

R

integer, the circle radius

W

integer, the circle width

cir

genome reference name (hg19, mm10 ...) or data frame from segAnglePo function or data frame from user's mapping data

.

type

the type is one of

  • "arc": arcs with variable radii

  • "arc2": arcs with the fixed radius

  • "b": bar charts

  • "b2": bar charts (opposite side of cutoff value)

  • "b3": bar charts with the same height

  • "box": box plots

  • "chr": plots of chromosomes or segments

  • "chr2": plots of chromosomes or segments of partial genome

  • "ci95": 95% confidence interval lines

  • "h": histograms

  • "heatmap": heatmaps

  • "heatmap2": heatmaps with genomic coordinates

  • "hightlight.link": link lines for zoom

  • "hist": polygons for multiple samples

  • "hl": highlight

  • "l": lines

  • "label": gene labels or text annotations

  • "label2": gene labels or text annotations with the same circumference coordinate

  • "lh": horizontal lines

  • "link.pg": link polygons based on Bezier curve

  • "link": link lines based on Bezier curve

  • "link2": link lines with smaller intra-chromosome arcs

  • "ls": lines in stair steps

  • "ml": multiple lines (for more than 1 samples)

  • "ml2": multiple horizontal lines

  • "ml3": multiple lines in stair steps

  • "ms": multiple points

  • "quant75": 75% quantile lines

  • "quant90": 90% quantile lines

  • "s": dots

  • "s2": dots with the fixed radius

  • "s.sd": dots proportional to standard deviation

  • "ss": dot sizes proportional to the values

  • "sv": dot sizes proportional to the variances

col.v

column number. The column value will be displayed. If type=heatmap, the number is as the first column.

B

logical: draw background?

print.chr.lab

logical: draw chromosomes or segment labels?

col.bar

logical: draw col.bar? It is for type=heatmap.

col.bar.po

draw col.bar position, e.g. topleft, bottomright.

cluster

logical: cluster and draw Dendrogram at left coner? It is for type=heatmap only.

order

vector: chromosome or segment order

scale

logical: draw scale?

cutoff

numeric: for multiple samples

zoom

vector containing six values: start chromosome, end chromosome, start position, end position, start angle and end angle

lwd

numeric: line width

cex

numeric: fond or point sizes

col

character or vector: color names

side

character (in or out): for type=label(2) only

...

...

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>

References

OmicCircos: an R package for simple and circular visualization of omics data. Cancer Inform. 2014 Jan 16;13:13-20. doi: 10.4137/CIN.S13495. eCollection 2014. PMID: 24526832 [PubMed] PMCID: PMC3921174

Examples

library(OmicCircos);
options(stringsAsFactors = FALSE);

set.seed(1234);

## initial values for simulation data 
seg.num     <- 10;
ind.num     <- 20;
seg.po      <- c(20:50);
link.num    <- 10;
link.pg.num <- 4;
## output simulation data
sim.out <- sim.circos(seg=seg.num, po=seg.po, ind=ind.num, link=link.num, 
  link.pg=link.pg.num);

seg.f     <- sim.out$seg.frame;
seg.v     <- sim.out$seg.mapping;
link.v    <- sim.out$seg.link
link.pg.v <- sim.out$seg.link.pg
seg.num   <- length(unique(seg.f[,1]));

## select segments
seg.name <- paste("chr", 1:seg.num, sep="");
db       <- segAnglePo(seg.f, seg=seg.name);

colors   <- rainbow(seg.num, alpha=0.5);

pdffile  <- "OmicCircos4vignette1.pdf";
pdf(pdffile, 8, 8);
par(mar=c(2, 2, 2, 2));
plot(c(1,800), c(1,800), type="n", axes=FALSE, xlab="", ylab="", main="");

circos(R=400, cir=db, type="chr",  col=colors, print.chr.lab=TRUE, W=4, scale=TRUE);
circos(R=360, cir=db, W=40, mapping=seg.v, col.v=3, type="l",   B=TRUE, col=colors[1], lwd=2, scale=TRUE);
circos(R=320, cir=db, W=40, mapping=seg.v, col.v=3, type="ls",  B=FALSE, col=colors[9], lwd=2, scale=TRUE);
circos(R=280, cir=db, W=40, mapping=seg.v, col.v=3, type="lh",  B=TRUE, col=colors[7], lwd=2, scale=TRUE);
circos(R=240, cir=db, W=40, mapping=seg.v, col.v=19, type="ml",  B=FALSE, col=colors, lwd=2, scale=TRUE);
circos(R=200, cir=db, W=40, mapping=seg.v, col.v=19, type="ml2", B=TRUE, col=colors, lwd=2);
circos(R=160, cir=db, W=40, mapping=seg.v, col.v=19, type="ml3", B=FALSE, cutoff=5, lwd=2);
circos(R=150, cir=db, W=40, mapping=link.v, type="link", lwd=2, col=colors[c(1,7)]);
circos(R=150, cir=db, W=40, mapping=link.pg.v, type="link.pg", lwd=2, col=sample(colors,link.pg.num));

dev.off()

## Not run: 
demo(OmicCircos4vignette1)
demo(OmicCircos4vignette2)
demo(OmicCircos4vignette3)
demo(OmicCircos4vignette4)
demo(OmicCircos4vignette5)
demo(OmicCircos4vignette6)
demo(OmicCircos4vignette7)
demo(OmicCircos4vignette8)
demo(OmicCircos4vignette9)
demo(OmicCircos4vignette10)

## End(Not run)

generate circular skeleton data from user's mapping data

Description

This function creates a data frame and converts the segment pointer positions (linear coordinates) into angle values (the angle based coordinates along circumference). In the data frame, column 1 is unique segment or chromosome names; column 2 is the start angle; column 3 is the end angle; column 4 is the accumulative start position; column 5 is the accumulative end position; column 6 is the start position and column 7 is the end position for each segment or chromosome.

Usage

segAnglePo(seg.dat=seg.dat, seg=seg, angle.start=angle.start, angle.end=angle.end);

Arguments

seg.dat

the segment data should be a matrix or a data frame: column 1 is the segment name or chromosome name; column 2 is the segment start; column 3 is the segment end; column 4 is segment name2 (optional); and column 5 is segment description (optional).

seg

vector: segment names (optional)

angle.start

numeric: plot start angle, angle.start=0 by default (optional)

angle.end

numeric: plot end angle, angle.end=360 by default (optional)

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>

References

OmicCircos: an R package for simple and circular visualization of omics data. Cancer Inform. 2014 Jan 16;13:13-20. doi: 10.4137/CIN.S13495. eCollection 2014. PMID: 24526832 [PubMed] PMCID: PMC3921174

Examples

library(OmicCircos);
options(stringsAsFactors = FALSE);
set.seed(1234);

## initial values for simulation data 
seg.num     <- 10;
ind.num     <- 20;
seg.po      <- c(20:50);
link.num    <- 10;
link.pg.num <- 4;
## output simulation data
sim.out <- sim.circos(seg=seg.num, po=seg.po, ind=ind.num, link=link.num, 
  link.pg=link.pg.num);

seg.f     <- sim.out$seg.frame;
seg.v     <- sim.out$seg.mapping;
link.v    <- sim.out$seg.link
link.pg.v <- sim.out$seg.link.pg
seg.num   <- length(unique(seg.f[,1]));

## select segments
seg.name <- paste("chr", 1:seg.num, sep="");
db       <- segAnglePo(seg.f, seg=seg.name);

colors   <- rainbow(seg.num, alpha=0.5);
pdffile  <- "OmicCircos4vignette2.pdf";
pdf(pdffile, 8, 8);
par(mar=c(2, 2, 2, 2));
plot(c(1,800), c(1,800), type="n", axes=FALSE, xlab="", ylab="", main="");

circos(R=400, type="chr", cir=db, col=colors, print.chr.lab=TRUE, W=4, scale=TRUE);
circos(R=360, cir=db, W=40, mapping=seg.v, col.v=8, type="box",   B=TRUE, col=colors[1], lwd=0.1, scale=TRUE);
circos(R=320, cir=db, W=40, mapping=seg.v, col.v=8, type="hist",  B=TRUE, col=colors[3], lwd=0.1, scale=TRUE);
circos(R=280, cir=db, W=40, mapping=seg.v, col.v=8, type="ms",  B=TRUE, col=colors[7], lwd=0.1, scale=TRUE);
circos(R=240, cir=db, W=40, mapping=seg.v, col.v=3, type="h",  B=FALSE,  col=colors[2], lwd=0.1);
circos(R=200, cir=db, W=40, mapping=seg.v, col.v=3, type="s", B=TRUE, col=colors, lwd=0.1);
circos(R=160, cir=db, W=40, mapping=seg.v, col.v=3, type="b", B=FALSE, col=colors, lwd=0.1);
circos(R=150, cir=db, W=40, mapping=link.v, type="link", lwd=2, col=colors[c(1,7)]);
circos(R=150, cir=db, W=40, mapping=link.pg.v, type="link.pg", lwd=2, col=sample(colors,link.pg.num));

dev.off()

~~ Methods for Function segAnglePo ~~

Description

~~ Methods for function segAnglePo ~~

Methods

signature(seg.dat = "data.frame")
signature(seg.dat = "GRanges")

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>

References

OmicCircos: an R package for simple and circular visualization of omics data. Cancer Inform. 2014 Jan 16;13:13-20. doi: 10.4137/CIN.S13495. eCollection 2014. PMID: 24526832 [PubMed] PMCID: PMC3921174


circular data simulation

Description

This function generates data for user to test the circos functions

Usage

sim.circos(seg=10, po=c(20,50), ind=10, link=10, link.pg=10);

Arguments

seg

integer, the segment number. The default is 10.

po

vector, the segment positions. The default is c(20:50)

ind

integer, the number of samples. The default is 10.

link

integer, the number of links. The default is 10.

link.pg

integer, the number of link ploygons. The default is 10.

Value

sim.circos returns a list containing at least the following components:

seg.frame

data.frame, segment data

seg.mapping

data.frame, mapping data

seg.link

data.fame, link data

seg.link.pg

data.frame, link polygon data

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>

References

OmicCircos: an R package for simple and circular visualization of omics data. Cancer Inform. 2014 Jan 16;13:13-20. doi: 10.4137/CIN.S13495. eCollection 2014. PMID: 24526832 [PubMed] PMCID: PMC3921174

Examples

library(OmicCircos);
options(stringsAsFactors = FALSE);
set.seed(1234);

## initial values for simulation data 
seg.num     <- 10;
ind.num     <- 20;
seg.po      <- c(20:50);
link.num    <- 10;
link.pg.num <- 4;
## output simulation data
sim.out <- sim.circos(seg=seg.num, po=seg.po, ind=ind.num, link=link.num, 
  link.pg=link.pg.num);

seg.f     <- sim.out$seg.frame;
seg.v     <- sim.out$seg.mapping;
link.v    <- sim.out$seg.link
link.pg.v <- sim.out$seg.link.pg
seg.num   <- length(unique(seg.f[,1]));

## select segments
seg.name <- paste("chr", 1:seg.num, sep="");
db       <- segAnglePo(seg.f, seg=seg.name);

colors   <- rainbow(seg.num, alpha=0.5);
pdffile  <- "OmicCircos4vignette3.pdf";
pdf(pdffile, 8, 8);
par(mar=c(2, 2, 2, 2));
plot(c(1,800), c(1,800), type="n", axes=FALSE, xlab="", ylab="", main="");

circos(R=400, type="chr", cir=db, col=colors, print.chr.lab=TRUE, W=4, scale=TRUE);
circos(R=360, cir=db, W=40, mapping=seg.v, col.v=8, type="quant90", B=FALSE, col=colors, lwd=2, scale=TRUE);
circos(R=320, cir=db, W=40, mapping=seg.v, col.v=3, type="sv", B=TRUE, col=colors[7],  scale=TRUE);
circos(R=280, cir=db, W=40, mapping=seg.v, col.v=3, type="ss", B=FALSE, col=colors[3],  scale=TRUE);
circos(R=240, cir=db, W=40, mapping=seg.v, col.v=8, type="heatmap", lwd=3);
circos(R=200, cir=db, W=40, mapping=seg.v, col.v=3, type="s.sd", B=FALSE, col=colors[4]);
circos(R=160, cir=db, W=40, mapping=seg.v, col.v=3, type="ci95", B=TRUE, col=colors[4], lwd=2);
circos(R=150, cir=db, W=40, mapping=link.v, type="link", lwd=2, col=colors[c(1,7)]);
circos(R=150, cir=db, W=40, mapping=link.pg.v, type="link.pg", lwd=2, col=sample(colors,link.pg.num));

the.col1=rainbow(10, alpha=0.5)[3];
highlight <- c(160, 410, 6, 2, 6, 10, the.col1, the.col1);
circos(R=110, cir=db, W=40, mapping=highlight, type="hl", lwd=1);

the.col1=rainbow(10, alpha=0.1)[3];
the.col2=rainbow(10, alpha=0.5)[1];
highlight <- c(160, 410, 3, 12, 3, 20, the.col1, the.col2);
circos(R=110, cir=db, W=40, mapping=highlight, type="hl", lwd=2);

dev.off()

TCGA BRCA expression and cnv association

Description

The p-values of associations between the TCGA Breast Cancer copy number and gene expression data.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


copy number data of TCGA Breast Cancer

Description

Examples of TCGA Breast Cancer DNA copy number data from 500 genes and 60 samples.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


TCGA Breast Cancer gene fusion data.

Description

Examples of TCGA Breast Cancer gene fusion data from 18 fusion proteins

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


TCGA BRCA expression data

Description

Examples of TCGA Breast Cancer expression data from 500 genes and 60 samples

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


TCGA BRCA Sample names and subtypes

Description

Names and subtypes of 60 samples for TCGA Breast Cancer expression and copy number data

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


BRCA PAM50 gene list (hg18)

Description

Breast cancer PAM 50 gene list (hg18).

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


chromosome banding colors

Description

Chromosome banding colors from UCSC Genome Browser

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


human hg18 circumference coordinates

Description

Human hg18 circumference coordinates (angles) are calculated from hg18 chromosome size data using the seqAnglPo function. The chromosome size data stored in UCSC.hg18.chr.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


human hg18 segment data.

Description

Human hg18 chromome size and binding data obtained from UCSC Genome Browser cytogenetics table.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


human hg19 circumference coordinates

Description

Human hg19 circumference coordinates (angles) are calculated from hg19 chromosome size data using the seqAnglPo function. The chromosome size data stored in UCSC.hg19.chr.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


human hg19 segment data

Description

Human hg19 chromome size and binding data obtained from UCSC Genome Browser cytogenetics table.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


mouse mm10 circumference coordinates

Description

Mouse mm10 circumference coordinates (angles) are calculated from mm10 chromosome size data using the seqAnglPo function. The mouse chromosome size data stored in UCSC.mm10.chr.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


mouse mm10 segment data.

Description

Mouse mm10 chromsome size and binding data obtained from UCSC Genome Browser cytogenetics table.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


mouse mm9 circumference coordinates

Description

Mouse mm9 circumference coordinates (angles) are calculated from mm9 chromosome size data using the seqAnglPo function. The mouse chromosome size data stored in UCSC.mm9.chr.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>


mouse mm9 segment data.

Description

Mouse mm9 chromsome size and binding data obtained from UCSC Genome Browser cytogenetics table.

Author(s)

Ying Hu <[email protected]> Chunhua Yan <[email protected]>