Quick standalone gene/region visualization from BigWig files
Source:R/plot_quick_view.R
plot_quick_view.RdA zero-configuration entry point for visualizing any genomic locus with BigWig signal tracks. Requires only a gene symbol (or genomic coordinates) and one or more BigWig file paths. No database setup, no track connection CSV, no chromatin states — just signal overlaid on the gene model.
Usage
plot_quick_view(
gene = NULL,
region = NULL,
bw_paths,
labels = NULL,
colors = NULL,
genome = "hg38",
txdb = NULL,
organism = NULL,
padding = 5000L,
mirror = FALSE,
signal_style = c("line", "polygon"),
signal_layout = "auto",
cex_cell_label = 1.4,
cex_axis = 1.2,
cex_coord = 1.3,
cex_annotation = 1.1,
cex_gene = 1.2,
cex_title = 1.5,
cex_signal_label = 1.2,
quantile_cap = 0.98,
scale_factor = 1.1,
export = NULL,
width = 10,
height = 8
)Arguments
- gene
Character or NULL. HGNC gene symbol (e.g.,
"INS","GCG"). Exactly one ofgeneorregionmust be provided.- region
List or NULL. Genomic coordinates as
list(chr = "chr11", start = 2159779, end = 2161209). Exactly one ofgeneorregionmust be provided.- bw_paths
Named character vector of BigWig file paths. Names are used as sample labels. If names contain both
"atac"and"rna"(case-insensitive), mirror pairing is attempted.- labels
Character vector or NULL. Override sample labels. If NULL, uses
names(bw_paths).- colors
Character vector or NULL. Override sample colors. If NULL, uses a default colorblind-friendly palette. Must have the same length as
bw_paths.- genome
Character. String name of the genome assembly associated with
bw_paths(e.g."hg38","mm10","rn6","dm6"). The value must match the assembly recorded in the suppliedtxdb; a validator compares the string againstGenomeInfoDb::genome(txdb). Defaults to"hg38"for convenience with built-in fallbacks.- txdb
Character or NULL. TxDb specification in
"package::object"form (e.g."TxDb.Rnorvegicus.UCSC.rn6.refGene:: TxDb.Rnorvegicus.UCSC.rn6.refGene"). If NULL, auto-resolved for the built-in convenience genomes"hg38"and"mm10"; for any other genome,txdbmust be supplied explicitly.- organism
Character or NULL. org.db specification (e.g.
"org.Rn.eg.db"). If NULL, auto-resolved for the built-in convenience genomes"hg38"and"mm10"; for any other genome,organismmust be supplied explicitly.- padding
Integer. Base pairs of padding around gene/region (default 5000).
- mirror
Logical. Enable mirrored ATAC/RNA layout (default FALSE).
- signal_style
Character. Signal rendering style:
"line"(default) draws vertical bars at each position (IGV/UCSC browser style),"polygon"draws filled area charts.- signal_layout
Character. Signal layout mode:
"auto","stacked", or"mirrored".- cex_cell_label
Numeric. Font size for cell type labels (default 1.4).
- cex_axis
Numeric. Font size for axis labels (default 1.2).
- cex_coord
Numeric. Font size for coordinate bar (default 1.3).
- cex_annotation
Numeric. Font size for annotation labels (default 1.1).
- cex_gene
Numeric. Font size for gene model labels (default 1.2).
- cex_title
Numeric. Font size for plot title (default 1.5).
- cex_signal_label
Numeric. Font size for signal indicators (default 1.2).
- quantile_cap
numeric. Percentile for capping extreme signal peaks (default 0.98). Peaks above this percentile are clipped to prevent axis compression.
- scale_factor
numeric. Y-axis headroom multiplier (default 1.1). Values above 1.0 add whitespace above the tallest signal peak.
- export
Character or NULL. File path for export (pdf, eps, png, tiff).
- width
Numeric. Export width in inches (default 10).
- height
Numeric. Export height in inches (default 8).
Value
Invisible NULL. A multi-panel base-R plot is drawn
on the current graphics device (or exported to file if
export is set).
Details
This is the simplest way to use epiRomics for exploratory
visualization. Specify a gene by symbol (e.g.,
"INS") or a region by coordinates. The function
auto-resolves the TxDb annotation database for the
specified genome assembly.
Examples
## plot_quick_view requires a TxDb (Suggests) and a real BigWig file.
## This example confirms that input validation fires correctly.
if (requireNamespace("TxDb.Hsapiens.UCSC.hg38.knownGene", quietly = TRUE)) {
bw_path <- make_example_bigwig()
tryCatch(
plot_quick_view(
region = list(chr = "chr1", start = 1000L, end = 51000L),
bw_paths = c(Sample = bw_path),
genome = "hg38"
),
error = function(e) message(e$message)
)
file.remove(bw_path)
}
#> 2169 genes were dropped because they have exons located on both strands of
#> the same reference sequence or on more than one reference sequence, so cannot
#> be represented by a single genomic range.
#> Use 'single.strand.genes.only=FALSE' to get all the genes in a GRangesList
#> object, or use suppressMessages() to suppress this message.
#> [1] TRUE