If FALSE (the default), then edges with the same global
edge id are merged into a single row. Use keep.dup.edges=TRUE
if this merging should not be performed.
Details
TODO
Value
TODO
Author(s)
H. Pages
See Also
This man page is part of the SplicingGraphs package.
Please see ?`SplicingGraphs-package` for an overview of the
package and for an index of its man pages.
Examples
## ---------------------------------------------------------------------
## 1. Make SplicingGraphs object 'sg' from toy gene model (see
## '?SplicingGraphs')
## ---------------------------------------------------------------------
example(SplicingGraphs) # create SplicingGraphs object 'sg'
sg
## 'sg' has 1 element per gene and 'names(sg)' gives the gene ids.
names(sg)
## ---------------------------------------------------------------------
## 2. Basic usage
## ---------------------------------------------------------------------
sgedges(sg["geneD"])
sgnodes(sg["geneD"])
outdeg(sg["geneD"])
indeg(sg["geneD"])
## ---------------------------------------------------------------------
## 3. Sanity checks
## ---------------------------------------------------------------------
check_way1_vs_way2 <- function(res1, res2)
{
edges1 <- res1[res1$ex_or_in != "", ] # remove artificial edges
edges2 <- mcols(unlist(res2, use.names=FALSE))
stopifnot(identical(edges1, edges2))
}
for (i in seq_along(sg)) {
sgi <- sg[i]
## After removal of the artificial edges, the edges returned
## by 'sgedges()' should be the same as those returned
## by 'sgedgesByGene()' on a SplicingGraphs object of length 1.
check_way1_vs_way2(
sgedges(sgi),
sgedgesByGene(sgi))
## After removal of the artificial edges, the edges returned
## by 'sgedges( , keep.dup.edges=TRUE)' should be the same as
## those returned by 'sgedgesByGene( , keep.dup.edges=TRUE)' or by
## 'sgedgesByTranscript()' on a SplicingGraphs object of length 1.
res1 <- DataFrame(sgedges(sgi, keep.dup.edges=TRUE))
check_way1_vs_way2(
res1,
sgedgesByGene(sgi, keep.dup.edges=TRUE))
check_way1_vs_way2(
res1,
sgedgesByTranscript(sgi))
}
Results
R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)
R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.
R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.
Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.
> library(SplicingGraphs)
Loading required package: GenomicFeatures
Loading required package: BiocGenerics
Loading required package: parallel
Attaching package: 'BiocGenerics'
The following objects are masked from 'package:parallel':
clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
clusterExport, clusterMap, parApply, parCapply, parLapply,
parLapplyLB, parRapply, parSapply, parSapplyLB
The following objects are masked from 'package:stats':
IQR, mad, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
rbind, rownames, sapply, setdiff, sort, table, tapply, union,
unique, unsplit
Loading required package: S4Vectors
Loading required package: stats4
Attaching package: 'S4Vectors'
The following objects are masked from 'package:base':
colMeans, colSums, expand.grid, rowMeans, rowSums
Loading required package: IRanges
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: AnnotationDbi
Loading required package: Biobase
Welcome to Bioconductor
Vignettes contain introductory material; view with
'browseVignettes()'. To cite Bioconductor, see
'citation("Biobase")', and for packages 'citation("pkgname")'.
Loading required package: GenomicAlignments
Loading required package: SummarizedExperiment
Loading required package: Biostrings
Loading required package: XVector
Loading required package: Rsamtools
Loading required package: Rgraphviz
Loading required package: graph
Attaching package: 'graph'
The following object is masked from 'package:Biostrings':
complement
Loading required package: grid
Attaching package: 'Rgraphviz'
The following objects are masked from 'package:IRanges':
from, to
The following objects are masked from 'package:S4Vectors':
from, to
Warning messages:
1: replacing previous import 'IRanges::from' by 'Rgraphviz::from' when loading 'SplicingGraphs'
2: replacing previous import 'IRanges::to' by 'Rgraphviz::to' when loading 'SplicingGraphs'
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/SplicingGraphs/sgedges-methods.Rd_%03d_medium.png", width=480, height=480)
> ### Name: sgedges-methods
> ### Title: Extract the edges (and nodes) of a splicing graph
> ### Aliases: sgedges-methods sgedges sgedges,SplicingGraphs-method sgnodes
> ### sgnodes,SplicingGraphs-method sgnodes,IntegerList-method
> ### sgnodes,data.frame-method sgnodes,DataFrame-method outdeg
> ### outdeg,ANY-method outdeg,DataFrame-method indeg indeg,ANY-method
> ### indeg,DataFrame-method
>
> ### ** Examples
>
> ## ---------------------------------------------------------------------
> ## 1. Make SplicingGraphs object 'sg' from toy gene model (see
> ## '?SplicingGraphs')
> ## ---------------------------------------------------------------------
> example(SplicingGraphs) # create SplicingGraphs object 'sg'
SplcnG> ## ---------------------------------------------------------------------
SplcnG> ## 1. Load a toy gene model as a TxDb object
SplcnG> ## ---------------------------------------------------------------------
SplcnG>
SplcnG> library(GenomicFeatures)
SplcnG> suppressWarnings(
SplcnG+ toy_genes_txdb <- makeTxDbFromGFF(toy_genes_gff())
SplcnG+ )
Import genomic features from the file as a GRanges object ... OK
Prepare the 'metadata' data frame ... OK
Make the TxDb object ... OK
SplcnG> ## ---------------------------------------------------------------------
SplcnG> ## 2. Compute all the splicing graphs (1 graph per gene) and return them
SplcnG> ## in a SplicingGraphs object
SplcnG> ## ---------------------------------------------------------------------
SplcnG>
SplcnG> ## Extract the exons grouped by transcript:
SplcnG> ex_by_tx <- exonsBy(toy_genes_txdb, by="tx", use.names=TRUE)
SplcnG> ## Extract the transcripts grouped by gene:
SplcnG> tx_by_gn <- transcriptsBy(toy_genes_txdb, by="gene")
SplcnG> sg <- SplicingGraphs(ex_by_tx, tx_by_gn)
SplcnG> sg
SplicingGraphs object with 5 gene(s) and 13 transcript(s)
SplcnG> ## Alternatively 'sg' can be constructed directly from the TxDb
SplcnG> ## object:
SplcnG> sg2 <- SplicingGraphs(toy_genes_txdb) # same as 'sg'
SplcnG> sg2
SplicingGraphs object with 5 gene(s) and 13 transcript(s)
SplcnG> ## Note that because SplicingGraphs objects have a slot that is an
SplcnG> ## environment (for caching the bubbles), they cannot be compared with
SplcnG> ## 'identical()' (will always return FALSE). 'all.equal()' should be
SplcnG> ## used instead:
SplcnG> stopifnot(isTRUE(all.equal(sg2, sg)))
SplcnG> ## 'sg' has 1 element per gene and 'names(sg)' gives the gene ids:
SplcnG> length(sg)
[1] 5
SplcnG> names(sg)
[1] "geneA" "geneB" "geneC" "geneD" "geneE"
SplcnG> ## ---------------------------------------------------------------------
SplcnG> ## 3. Basic manipulation of a SplicingGraphs object
SplcnG> ## ---------------------------------------------------------------------
SplcnG>
SplcnG> ## Basic accessors:
SplcnG> seqnames(sg)
geneA geneB geneC geneD geneE
chrX chrX chrX chrX chrX
Levels: chrX
SplcnG> strand(sg)
geneA geneB geneC geneD geneE
+ - + + +
Levels: + - *
SplcnG> seqinfo(sg)
Seqinfo object with 1 sequence from an unspecified genome; no seqlengths:
seqnames seqlengths isCircular genome
chrX NA NA <NA>
SplcnG> ## Number of transcripts per gene:
SplcnG> elementNROWS(sg)
geneA geneB geneC geneD geneE
2 2 3 4 2
SplcnG> ## The transcripts of a given gene can be extracted with [[. The result
SplcnG> ## is an *unnamed* GRangesList object containing the exons grouped by
SplcnG> ## transcript:
SplcnG> sg[["geneD"]]
GRangesList object of length 4:
[[1]]
GRanges object with 2 ranges and 5 metadata columns:
seqnames ranges strand | exon_id exon_name exon_rank start_SSid
<Rle> <IRanges> <Rle> | <integer> <character> <integer> <integer>
[1] chrX [601, 630] + | 10 Dx2 1 1
[2] chrX [666, 675] + | 12 Dx4 2 5
end_SSid
<integer>
[1] 3
[2] 6
[[2]]
GRanges object with 2 ranges and 5 metadata columns:
seqnames ranges strand | exon_id exon_name exon_rank start_SSid
[1] chrX [601, 620] + | 9 Dx1 1 1
[2] chrX [651, 700] + | 11 Dx3 2 4
end_SSid
[1] 2
[2] 8
[[3]]
GRanges object with 3 ranges and 5 metadata columns:
seqnames ranges strand | exon_id exon_name exon_rank start_SSid
[1] chrX [601, 620] + | 9 Dx1 1 1
[2] chrX [666, 675] + | 12 Dx4 2 5
[3] chrX [691, 700] + | 13 Dx5 3 7
end_SSid
[1] 2
[2] 6
[3] 8
...
<1 more element>
-------
seqinfo: 1 sequence from an unspecified genome; no seqlengths
SplcnG> ## See '?plotTranscripts' for how to plot those transcripts.
SplcnG>
SplcnG> ## The transcripts of all the genes can be extracted with unlist(). The
SplcnG> ## result is a *named* GRangesList object containing the exons grouped
SplcnG> ## by transcript. The names on the object are the gene ids:
SplcnG> ex_by_tx <- unlist(sg)
SplcnG> ex_by_tx
GRangesList object of length 13:
$geneA
GRanges object with 1 range and 5 metadata columns:
seqnames ranges strand | exon_id exon_name exon_rank start_SSid
<Rle> <IRanges> <Rle> | <integer> <character> <integer> <integer>
[1] chrX [11, 50] + | 2 Ax2 1 1
end_SSid
<integer>
[1] 3
$geneA
GRanges object with 2 ranges and 5 metadata columns:
seqnames ranges strand | exon_id exon_name exon_rank start_SSid
[1] chrX [11, 40] + | 1 Ax1 1 1
[2] chrX [71, 100] + | 3 Ax3 2 4
end_SSid
[1] 2
[2] 5
$geneB
GRanges object with 2 ranges and 5 metadata columns:
seqnames ranges strand | exon_id exon_name exon_rank start_SSid
[1] chrX [251, 300] - | 23 Bx1 1 3
[2] chrX [201, 230] - | 20 Bx2 2 6
end_SSid
[1] 1
[2] 4
...
<10 more elements>
-------
seqinfo: 1 sequence from an unspecified genome; no seqlengths
> sg
SplicingGraphs object with 5 gene(s) and 13 transcript(s)
>
> ## 'sg' has 1 element per gene and 'names(sg)' gives the gene ids.
> names(sg)
[1] "geneA" "geneB" "geneC" "geneD" "geneE"
>
> ## ---------------------------------------------------------------------
> ## 2. Basic usage
> ## ---------------------------------------------------------------------
> sgedges(sg["geneD"])
DataFrame with 13 rows and 5 columns
from to sgedge_id ex_or_in tx_id
<character> <character> <character> <factor> <CharacterList>
1 R 1 geneD:R,1 D4,D1,D2,...
2 1 3 geneD:1,3 ex D4,D3
3 3 5 geneD:3,5 in D4
4 5 6 geneD:5,6 ex D4,D2
5 6 L geneD:6,L D4
... ... ... ... ... ...
9 8 L geneD:8,L D1,D2,D3
10 2 5 geneD:2,5 in D2
11 6 7 geneD:6,7 in D2
12 7 8 geneD:7,8 ex D2
13 3 4 geneD:3,4 in D3
> sgnodes(sg["geneD"])
[1] "R" "1" "2" "3" "4" "5" "6" "7" "8" "L"
> outdeg(sg["geneD"])
R 1 2 3 4 5 6 7 8 L
1 2 2 2 1 1 2 1 1 0
> indeg(sg["geneD"])
R 1 2 3 4 5 6 7 8 L
0 1 1 1 2 2 1 1 2 2
>
> ## ---------------------------------------------------------------------
> ## 3. Sanity checks
> ## ---------------------------------------------------------------------
> check_way1_vs_way2 <- function(res1, res2)
+ {
+ edges1 <- res1[res1$ex_or_in != "", ] # remove artificial edges
+ edges2 <- mcols(unlist(res2, use.names=FALSE))
+ stopifnot(identical(edges1, edges2))
+ }
>
> for (i in seq_along(sg)) {
+ sgi <- sg[i]
+ ## After removal of the artificial edges, the edges returned
+ ## by 'sgedges()' should be the same as those returned
+ ## by 'sgedgesByGene()' on a SplicingGraphs object of length 1.
+ check_way1_vs_way2(
+ sgedges(sgi),
+ sgedgesByGene(sgi))
+ ## After removal of the artificial edges, the edges returned
+ ## by 'sgedges( , keep.dup.edges=TRUE)' should be the same as
+ ## those returned by 'sgedgesByGene( , keep.dup.edges=TRUE)' or by
+ ## 'sgedgesByTranscript()' on a SplicingGraphs object of length 1.
+ res1 <- DataFrame(sgedges(sgi, keep.dup.edges=TRUE))
+ check_way1_vs_way2(
+ res1,
+ sgedgesByGene(sgi, keep.dup.edges=TRUE))
+ check_way1_vs_way2(
+ res1,
+ sgedgesByTranscript(sgi))
+ }
>
>
>
>
>
> dev.off()
null device
1
>