[BioC] GRangesList replacement method bug
Hervé Pagès
hpages at fhcrc.org
Thu Dec 15 19:46:28 CET 2011
Hi Richard,
On 11-12-14 02:01 PM, Richard Bourgon wrote:
> BioC team:
>
> I encountered a bug with replacement methods for GRangesList objects (see below for code).
>
> Basically, you can set to non-existant elementMetaData columns for a GRanges object, but not if it is a part of a GRangesList. Setting to already existing columns is fine, however.
>
> Richard Bourgon, Ph.D.
> Genentech, Inc.
> Bioinformatics& Computational Biology
> 1 DNA Way, MS 93
> South San Francisco, CA 94080-4990
> (650) 467-2064
>
> ***
>
> gr<- GRanges( "a", IRanges(1,10), x = TRUE )
> grl<- GRangesList( y = gr, z = gr )
> values( grl[["y"]] )$x<- FALSE # Fine for existing column
> values( grl[["y"]] )$w<- FALSE # Fails for new column
>
> Error in .Method(..., deparse.level = deparse.level) :
> number of columns for arg 2 do not match those of first arg
>
> values( gr )$w<- FALSE # No problem with new column here
The error message doesn't help but this error is actually intended. This
is because all the individual GRanges elements in a GRangesList object
must have exactly the same columns (in the same order). For example,
you get the same error if you use the GRangesList() constructor to try
to put together 2 GRanges objects that don't have the same cols:
> GRangesList(GRanges("a", IRanges(1,10), x=TRUE),
GRanges("a", IRanges(1,10)))
Error in .Method(..., deparse.level = deparse.level) :
number of columns for arg 2 do not match those of first arg
or if they have the same columns but not in the same order:
> GRangesList(GRanges("a", IRanges(1,10), x=TRUE, w=2),
GRanges("a", IRanges(1,10), w=3, x=FALSE))
Error in .Method(..., deparse.level = deparse.level) :
column names for arg 2 do not match those of first arg
When you do 'values(grl[[1]])$w <- FALSE', you are trying to add a
column to the 1st GRanges in the GRangesList, but then its columns
would not be the same as those of the other elements anymore.
Note that if 'grl' has only 1 element, it works:
> grl <- GRangesList(GRanges("a", IRanges(1,10)))
> values(grl[[1]])$w <- FALSE
> grl
GRangesList of length 1:
[[1]]
GRanges with 1 range and 1 elementMetadata value:
seqnames ranges strand | w
<Rle> <IRanges> <Rle> | <logical>
[1] a [1, 10] * | FALSE
---
seqlengths:
a
NA
Also note that you can always add, remove, shuffle, rename the columns
of all the elements in a GRangesList object by unlisting/relisting the
object. For example to add a "w" column to all the GRanges in 'grl':
unlisted <- unlist(grl, use.names=FALSE)
values(unlisted)$w <- FALSE
relist(unlisted, grl)
More generally speaking, and maybe the name of the container doesn't
help, but conceptually a GRangesList object is not just a list of
arbitrary/unrelated GRanges objects. All the GRanges objects in the
GRangesList share the same columns and sequence levels.
Hope this help,
H.
> sessionInfo()
>
> R version 2.14.0 (2011-10-31)
> Platform: x86_64-unknown-linux-gnu (64-bit)
>
> locale:
> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C LC_TIME=en_US.UTF-8
> [4] LC_COLLATE=en_US.UTF-8 LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
> [7] LC_PAPER=C LC_NAME=C LC_ADDRESS=C
> [10] LC_TELEPHONE=C LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
>
> attached base packages:
> [1] datasets utils grDevices graphics stats methods base
>
> other attached packages:
> [1] rtracklayer_1.14.4 RCurl_1.8-0 bitops_1.0-4.1 GenomicRanges_1.6.4
> [5] IRanges_1.12.5 gdata_2.8.2 BiocInstaller_1.2.1
>
> loaded via a namespace (and not attached):
> [1] Biobase_2.14.0 Biostrings_2.22.0 BSgenome_1.22.0 gtools_2.6.2 tools_2.14.0
> [6] XML_3.6-2 zlibbioc_1.0.0
>
> _______________________________________________
> Bioconductor mailing list
> Bioconductor at r-project.org
> https://stat.ethz.ch/mailman/listinfo/bioconductor
> Search the archives: http://news.gmane.org/gmane.science.biology.informatics.conductor
--
Hervé Pagès
Program in Computational Biology
Division of Public Health Sciences
Fred Hutchinson Cancer Research Center
1100 Fairview Ave. N, M1-B514
P.O. Box 19024
Seattle, WA 98109-1024
E-mail: hpages at fhcrc.org
Phone: (206) 667-5791
Fax: (206) 667-1319
More information about the Bioconductor
mailing list