Snapshot-class: Class '"Snapshot"'

Description Usage Arguments Methods Accessors Fields Class-Based Methods Author(s) See Also Examples

Description

A Snapshot-class to visualize genomic data from BAM files with zoom and pan functionality.

Usage

1

Arguments

files

A character() or BamFileList specifying the file(s) to be visualized.

range

A GRanges object specifying the range to be visualized.

...

Additional, optional, arguments to be passed to the Snapshot initialize function. Arguments include:

functions:

A SnapshotFunctionList of functions, in addition to built-in ‘fine_coverage’, ‘coarse_coverage’, ‘multifine_coverage’, to be used for visualization.

currentFunction:

character(1) naming the function, from functions to be used for data input and visualization. The default chooses a function based on the scale at which the data is being visualized.

annTrack:

Annotation track. If built-in visualization functions are to be used, annTrack should be a GRanges instance and the first column of its elementMeatdata would be used to annotate the range.

fac:

Character(1) indicating which factor used for grouping the sample files. The factor should be included in the elementMetadata of files, otherwise ignored. Used only to visualize multiple files.

.auto_display:

logical(1) indicating whether the visualization is to be updated when show is invoked.

.debug

logical(1) indicating whether debug messages are to be printed.

Methods

zoom

signature(x = "Snapshot"): Zoom (in or out) the current plot.

pan

signature(x = "Snapshot"): Pan (right or left) the current plot.

togglefun

signature(x = "Snapshot"): Toggle the current functions which imported records are to be immediately evaluated. Note that the active range will be changed to the current active window.

togglep

signature(x = "Snapshot"): Toggle the panning effects.

togglez

signature(x = "Snapshot"): Toggle the zooming effects.

Accessors

show

signature(object = "Snapshot"): Display a Snapshot object.

files

signature(x = "Snapshot"): Get the files field (object of class BamFileList) of a Snapshot object.

functions

signature(x = "Snapshot"): Get the functions field (object of SnapshotFunctionList) of a Snapshot object.

view

signature(x = "Snapshot"): Get the view field (object of SpTrellis) of a Snapshot object.

vrange

signature(x = "Snapshot"): Get the .range field (object of GRanges) of a Snapshot object.

getTrellis

signature(x = "Snapshot"): Get the trellis object, a field of the SpTrellis object.

Fields

.debug:

Object of class function to display messages while in debug mode

.auto_display:

Object of class logical to automatically display the coverage plot.

.range:

Object of class GRanges indicating which ranges of records to be imported from BAM fields.

.zin:

Object of class logical indicating whether the current zooming effect is zoom in.

.pright:

Object of class logical indicating whether the current panning effect is right.

.data:

Object of class data.frame containing coverage a position is represented for each strand and BAM file.

.data_dirty:

Object of class logical indicating whether to re-evaluate the imported records.

.initial_functions:

Object of class SnapshotFunctionList available by the Snapshot object.

.current_function:

Object of class character of the function the imported recorded are currently evaluated and visualized.

annTrack:

Default to NULL if not intended to visualize the annotation track. If default visualization function(s) is intended to be used to plot the annotation, annTrack has to be a GRanges instance.

functions:

Object of class SnapshotFunctionList of customized functions to evaluate and visualize the imported records.

files:

Object of class BamFileList to be imported.

view:

Object of class SpTrellis that is essentially a reference class wrapper of Trellis objects.

Class-Based Methods

display():

Display the current Snapshot object.

pan():

Pan (right or left) the current plot.

zoom():

Zoom (in or out) the current plot.

toggle(zoom, pan, currentFunction):

Toggle zooming, panning effects or the currentFuction in which the imported records are to be evaluated and visualized.

Author(s)

Martin Morgan and Chao-Jen Wong cwon2@fhcrc.org

See Also

SpTrellis

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
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
## example 1: Importing specific ranges of records

file <- system.file("extdata", "SRR002051.chrI-V.bam",
                    package="yeastNagalakshmi")
which <-  GRanges("chrI", IRanges(1, 2e5))
s <- Snapshot(file, range=which)

## methods
zoom(s) # zoom in
## zoom in to a specific region
zoom(s, range=GRanges("chrI", IRanges(7e4, 7e4+8000)))
pan(s)  # pan right
togglez(s) # change effect of zooming
zoom(s) # zoom out
togglep(s) # change effect of panning
pan(s)

## accessors
functions(s)
vrange(s)
show(s)
ignore.strand(s)
view(s) ## extract the spTrellis object
getTrellis(s) ## extract the trellis object

## example 2: ignore strand
s <- Snapshot(file, range=which, ignore.strand=TRUE)

##
## example 3: visualizing annotation track
##

library(GenomicFeatures)

getAnnGR <- function(txdb, which) {
    ex <- exonsBy(txdb, by="gene")
    seqlevels(ex, pruning.mode="coarse") <- seqlevels(which)
    r <- range(ex)
    gr <- unlist(r)
    values(gr)[["gene_id"]] <- rep.int(names(r), times=lengths(r))
    gr
}

txdbFile <- system.file("extdata", "sacCer2_sgdGene.sqlite",
                    package="yeastNagalakshmi")
# txdb <- makeTxDbFromUCSC(genome="sacCer2", tablename="sgdGene")
txdb <- loadDb(txdbFile)
which <-  GRanges("chrI", IRanges(1, 2e5))
gr <- getAnnGR(txdb, which)
## note that the first column of the elementMetadata annotates of the
## range of the elements.
gr

s <- Snapshot(file, range=which, annTrack=gr)
annTrack(s)
## zoom in to an interesting region
zoom(s, range=GRanges("chrI", IRanges(7e4, 7e4+8000)))

togglez(s) ## zoom out
zoom(s)

pan(s)

## example 4, 5, 6: multiple BAM files with 'multicoarse_covarage'
## and 'multifine_coverage' view.

## Resolution does not automatically switch for views of multiple
## files. It is important to note if width(which) < 10,000, use
## multifine_coverage.  Otherwise use multicoarse_coverage
file <- system.file("extdata", "SRR002051.chrI-V.bam",
                    package="yeastNagalakshmi")
which <-  GRanges("chrI", IRanges(1, 2e5))
s <- Snapshot(c(file, file), range=which,
              currentFunction="multicoarse_coverage")

## grouping files and view by 'multicoarse_coverage'
bfiles <- BamFileList(c(a=file, b=file))
values(bfiles) <- DataFrame(sampleGroup=factor(c("normal", "tumor")))
values(bfiles)
s <- Snapshot(bfiles, range=which,
              currentFunction="multicoarse_coverage", fac="sampleGroup")

## grouping files and view by 'multifine_coverage'
which <- GRanges("chrI", IRanges(7e4, 7e4+8000))
s <- Snapshot(bfiles, range=which,
              currentFunction="multifine_coverage", fac="sampleGroup")

ShortRead documentation built on Nov. 8, 2020, 8:02 p.m.