[Bioc-devel] Issues with GenomicRanges updates

Hervé Pagès hpages at fredhutch.org
Mon Feb 12 20:45:14 CET 2018


Hi Nathan,

You're right, there was a remaining issue. It should be fixed in
IRanges 2.13.26. Make sure to also get S4Vectors 0.17.32 (IRanges
2.13.26 requires it).

The 2 updated packages should become available via biocLite() in
about 24 hours but you can get them now directly from
git.bioconductor.org or from GitHub:

   https://github.com/Bioconductor/S4Vectors
   https://github.com/Bioconductor/IRanges

Thanks for reporting this and sorry for the inconvenience.

H.

On 02/12/2018 04:28 AM, Nathan Sheffield wrote:
> Hi Herve,
> 
> The updates have indeed solved those issues for that sample -- However, 
> when you try to apply a more complicated function, I am getting the same 
> error. Here's a reproducible example of the error again, this time using 
> the latest GenomicRanges and IRanges packages:
> 
> library(GenomicRanges)
> data("sample_input", package="LOLA")
> as.list(userSets)
> lapply(userSets, length)
> 
> # That works now
> 
> 
> data("sample_universe", package="LOLA")
> 
> f = function(x, userUniverse) {
>          fo = findOverlaps(x, userUniverse)
>          y = userUniverse[unique(subjectHits(fo))]
>          return(y)
>   }
> 
> lapply(userSets, f, userUniverse)
> 
> Error in lapply_CompressedList(X, FUN, ...) :
>    invalid output element of class "GRanges"
> 
> 
> This works under earlier versions. Can you see if this fails in your 
> setup? Or am I doing something strange in the function? Here's my 
> sessionInfo showing I updated to the latest packages:
> 
> -Nathan
> 
>  > sessionInfo()
> R Under development (unstable) (2018-02-04 r74204)
> Platform: x86_64-pc-linux-gnu (64-bit)
> Running under: Debian GNU/Linux 9 (stretch)
> 
> Matrix products: default
> BLAS: /usr/lib/openblas-base/libblas.so.3
> LAPACK: /usr/lib/libopenblasp-r0.2.19.so
> 
> locale:
>   [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C
>   [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8
>   [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=C
>   [7] LC_PAPER=en_US.UTF-8       LC_NAME=C
>   [9] LC_ADDRESS=C               LC_TELEPHONE=C
> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
> 
> attached base packages:
> [1] parallel  stats4    stats     graphics  grDevices utils datasets
> [8] methods   base
> 
> other attached packages:
> [1] GenomicRanges_1.31.19 GenomeInfoDb_1.15.5   IRanges_2.13.25
> [4] S4Vectors_0.17.31     BiocGenerics_0.25.3
> 
> loaded via a namespace (and not attached):
> [1] zlibbioc_1.25.0        compiler_3.5.0 tools_3.5.0
> [4] XVector_0.19.8         GenomeInfoDbData_1.1.0 RCurl_1.95-4.10
> [7] bitops_1.0-6
> 
> On 02/10/2018 10:40 PM, Nathan Sheffield wrote:
>> Hi Herve,
>>
>> Never mind, I see now I am still a day old, looks like I was looking 
>> at your sessionInfo paste and thought it was mine, whoops. I'll give 
>> it another try tomorrow with the new versions.
>>
>> other attached packages:
>> [1] GenomicRanges_1.31.18 GenomeInfoDb_1.15.5 IRanges_2.13.24
>> [4] S4Vectors_0.17.31     BiocGenerics_0.25.3
>>
>> -Nathan
>>
>> On 02/10/2018 07:33 PM, Nathan Sheffield wrote:
>>> According to my `sessionInfo` (see below), those are the versions I 
>>> had been using:
>>>
>>> other attached packages:
>>> [1] GenomicRanges_1.31.19 GenomeInfoDb_1.15.5   IRanges_2.13.25
>>> [4] S4Vectors_0.17.31     BiocGenerics_0.25.3
>>>
>>> And I had pulled them from biocLite...what's going on?
>>>
>>> -Nathan
>>>
>>>
>>> On 02/10/2018 06:19 PM, Hervé Pagès wrote:
>>>> Hi Nathan,
>>>>
>>>> I can't reproduce this with the latest versions of S4Vectors (0.17.31),
>>>> IRanges (2.13.25), and GenomicRanges (1.31.19). Note that these 
>>>> versions
>>>> will only become available via biocLite() tomorrow but you can get them
>>>> directly from git.bioconductor.org.
>>>>
>>>> With these versions, as.list, lapply, and mclapply work for on 
>>>> userSets:
>>>>
>>>>   > library(GenomicRanges)
>>>>   > data("sample_input", package="LOLA")
>>>>   > as.list(userSets)
>>>>   $setA
>>>>   GRanges object with 3142 ranges and 0 metadata columns:
>>>>            seqnames               ranges strand
>>>>               <Rle>            <IRanges> <Rle>
>>>>        [1]     chr1   [ 437151,  438164]      *
>>>>        [2]     chr1   [ 875730,  878363]      *
>>>>        [3]     chr1   [ 933387,  937410]      *
>>>>        [4]     chr1   [ 967966,  970238]      *
>>>>        [5]     chr1   [1016863, 1017439]      *
>>>>        ...      ...                  ...    ...
>>>>     [3138]     chrY [ 9364545,  9364859]      *
>>>>     [3139]     chrY [ 9385471,  9385777]      *
>>>>     [3140]     chrY [14532115, 14533600]      *
>>>>     [3141]     chrY [23696580, 23696878]      *
>>>>     [3142]     chrY [26959489, 26959716]      *
>>>>     -------
>>>>     seqinfo: 69 sequences from an unspecified genome; no seqlengths
>>>>
>>>>   $setB
>>>>   GRanges object with 5831 ranges and 0 metadata columns:
>>>>            seqnames               ranges strand
>>>>               <Rle>            <IRanges> <Rle>
>>>>        [1]     chr1     [ 28735,  29810]      *
>>>>        [2]     chr1     [544738, 546649]      *
>>>>        [3]     chr1     [713984, 714547]      *
>>>>        [4]     chr1     [762416, 763445]      *
>>>>        [5]     chr1     [805198, 805628]      *
>>>>        ...      ...                  ...    ...
>>>>     [5827]     chrY [20508190, 20508452]      *
>>>>     [5828]     chrY [21154603, 21155040]      *
>>>>     [5829]     chrY [21238448, 21240005]      *
>>>>     [5830]     chrY [26979889, 26980116]      *
>>>>     [5831]     chrY [28773315, 28773544]      *
>>>>     -------
>>>>     seqinfo: 69 sequences from an unspecified genome; no seqlengths
>>>>
>>>>   > lapply(userSets, length)
>>>>   $setA
>>>>   [1] 3142
>>>>
>>>>   $setB
>>>>   [1] 5831
>>>>
>>>>   > mclapply(userSets, length)
>>>>   $setA
>>>>   [1] 3142
>>>>
>>>>   $setB
>>>>   [1] 5831
>>>>
>>>> Note that you should not need to call as.list() on a GRangesList object
>>>> before passing it to lapply() or mclapply().
>>>>
>>>> Let me know if the problem persist after you update.
>>>>
>>>> Best,
>>>> H.
>>>>
>>>> > sessionInfo()
>>>> R Under development (unstable) (2017-12-11 r73889)
>>>> Platform: x86_64-pc-linux-gnu (64-bit)
>>>> Running under: Ubuntu 16.04.3 LTS
>>>>
>>>> Matrix products: default
>>>> BLAS: /home/hpages/R/R-3.5.r73889/lib/libRblas.so
>>>> LAPACK: /home/hpages/R/R-3.5.r73889/lib/libRlapack.so
>>>>
>>>> locale:
>>>>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C
>>>>  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8
>>>>  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8
>>>>  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C
>>>>  [9] LC_ADDRESS=C               LC_TELEPHONE=C
>>>> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
>>>>
>>>> attached base packages:
>>>> [1] parallel  stats4    stats     graphics  grDevices utils datasets
>>>> [8] methods   base
>>>>
>>>> other attached packages:
>>>> [1] GenomicRanges_1.31.19 GenomeInfoDb_1.15.5 IRanges_2.13.25
>>>> [4] S4Vectors_0.17.31     BiocGenerics_0.25.3
>>>>
>>>> loaded via a namespace (and not attached):
>>>> [1] zlibbioc_1.25.0        compiler_3.5.0         tools_3.5.0
>>>> [4] XVector_0.19.8         GenomeInfoDbData_1.1.0 RCurl_1.95-4.10
>>>> [7] bitops_1.0-6
>>>>
>>>>
>>>> On 02/10/2018 02:48 PM, Nathan Sheffield wrote:
>>>>> I'm having some issues getting my package LOLA to pass R CMD check 
>>>>> using the updated dev versions of GenomicRanges et al.
>>>>>
>>>>> It seems like any time I try to apply something across a 
>>>>> "CompressedGRangesList" object, it's giving errors when I use 
>>>>> mclapply from parallel. Here's a reproducible example:
>>>>>
>>>>> data("sample_input", package="LOLA")
>>>>> library(parallel)
>>>>> mclapply(userSets, length)
>>>>>
>>>>> (loads packages...)
>>>>>
>>>>> Error in lapply_CompressedList(X, FUN, ...) :
>>>>>    invalid output element of class "GRanges"
>>>>>
>>>>>
>>>>> It works with regular lapply:
>>>>>
>>>>>  > lapply(userSets, length)
>>>>> $setA
>>>>> [1] 3142
>>>>>
>>>>> $setB
>>>>> [1] 5831
>>>>>
>>>>>
>>>>> This is running on the bioconductor docker devel_core2 container, 
>>>>> and I've then gone and updated to the latest dev versions of these 
>>>>> packages with `biocLite()`.
>>>>>
>>>>> I earlier ran into issues using `as.list()` on these same 
>>>>> CompressedGRangesList objects. It used to be that I had to call 
>>>>> as.list when they were just GRangesList objects, but now that's 
>>>>> failing, and so I've had to go take all calls to as.list out of my 
>>>>> code. This has solved that issue (I guess the updates made the 
>>>>> as.list calls unnecessary), but you can see it's still causing errors:
>>>>>
>>>>> as.list(userSets)
>>>>> Error in lapply_CompressedList(X, FUN, ...) :
>>>>>    invalid output element of class "GRanges"
>>>>>
>>>>> -Nathan
>>>>>
>>>>> ```
>>>>>  > sessionInfo()
>>>>> R Under development (unstable) (2018-02-04 r74204)
>>>>> Platform: x86_64-pc-linux-gnu (64-bit)
>>>>> Running under: Debian GNU/Linux 9 (stretch)
>>>>>
>>>>> Matrix products: default
>>>>> BLAS: /usr/lib/openblas-base/libblas.so.3
>>>>> LAPACK: /usr/lib/libopenblasp-r0.2.19.so
>>>>>
>>>>> locale:
>>>>>   [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C
>>>>>   [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8
>>>>>   [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=C
>>>>>   [7] LC_PAPER=en_US.UTF-8       LC_NAME=C
>>>>>   [9] LC_ADDRESS=C               LC_TELEPHONE=C
>>>>> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
>>>>>
>>>>> attached base packages:
>>>>> [1] stats4    parallel  stats     graphics  grDevices utils datasets
>>>>> [8] methods   base
>>>>>
>>>>> other attached packages:
>>>>> [1] GenomicRanges_1.31.18 GenomeInfoDb_1.15.5 IRanges_2.13.24
>>>>> [4] S4Vectors_0.17.31     BiocGenerics_0.25.3
>>>>>
>>>>> loaded via a namespace (and not attached):
>>>>> [1] zlibbioc_1.25.0        compiler_3.5.0 XVector_0.19.8
>>>>> [4] tools_3.5.0            GenomeInfoDbData_1.1.0 RCurl_1.95-4.10
>>>>> [7] bitops_1.0-6
>>>>> ```
>>>>>
>>>>> _______________________________________________
>>>>> Bioc-devel at r-project.org mailing list
>>>>> https://urldefense.proofpoint.com/v2/url?u=https-3A__stat.ethz.ch_mailman_listinfo_bioc-2Ddevel&d=DwIDaQ&c=eRAMFD45gAfqt84VtBcfhQ&r=BK7q3XeAvimeWdGbWY_wJYbW0WYiZvSXAJJKaaPhzWA&m=FrZj-JCIQePYSae83Zo7pEJhLJDjFiRPOmff3fGshM4&s=dO7xndrTgL8lx8igvje5h0wKEFwaH0hwLq04G5UNcRs&e= 
>>>>>
>>>>
>>>
>>> _______________________________________________
>>> Bioc-devel at r-project.org mailing list
>>> https://urldefense.proofpoint.com/v2/url?u=https-3A__stat.ethz.ch_mailman_listinfo_bioc-2Ddevel&d=DwIDaQ&c=eRAMFD45gAfqt84VtBcfhQ&r=BK7q3XeAvimeWdGbWY_wJYbW0WYiZvSXAJJKaaPhzWA&m=EbJhDYrzCeCr-dtYMuctYT_W_fHl7mHAwOejsN1xqh4&s=NOmxvotdUgUeYSKGl9Pc11OVQW9TQttdRk81Fg5JNwE&e= 
>>>
>>
>> _______________________________________________
>> Bioc-devel at r-project.org mailing list
>> https://urldefense.proofpoint.com/v2/url?u=https-3A__stat.ethz.ch_mailman_listinfo_bioc-2Ddevel&d=DwIDaQ&c=eRAMFD45gAfqt84VtBcfhQ&r=BK7q3XeAvimeWdGbWY_wJYbW0WYiZvSXAJJKaaPhzWA&m=EbJhDYrzCeCr-dtYMuctYT_W_fHl7mHAwOejsN1xqh4&s=NOmxvotdUgUeYSKGl9Pc11OVQW9TQttdRk81Fg5JNwE&e= 
>>
> 
> _______________________________________________
> Bioc-devel at r-project.org mailing list
> https://urldefense.proofpoint.com/v2/url?u=https-3A__stat.ethz.ch_mailman_listinfo_bioc-2Ddevel&d=DwIDaQ&c=eRAMFD45gAfqt84VtBcfhQ&r=BK7q3XeAvimeWdGbWY_wJYbW0WYiZvSXAJJKaaPhzWA&m=EbJhDYrzCeCr-dtYMuctYT_W_fHl7mHAwOejsN1xqh4&s=NOmxvotdUgUeYSKGl9Pc11OVQW9TQttdRk81Fg5JNwE&e= 
> 

-- 
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 fredhutch.org
Phone:  (206) 667-5791
Fax:    (206) 667-1319



More information about the Bioc-devel mailing list