[Bioc-devel] segfault when using RleList in DataFrames
Leonard Goldstein
goldstein.leonard at gene.com
Mon Dec 7 18:17:05 CET 2015
Thanks Hervé!
Best wishes
Leonard
On Mon, Dec 7, 2015 at 12:08 AM, Hervé Pagès <hpages at fredhutch.org> wrote:
> Hi Leonard,
>
> Thanks for the report. This segfault was due to a long-standing bug
> in the "window" method for Rle objects when the object has length 0:
>
> > window(Rle(), 1, 0) # you might need to try many times before
> # you get the segfault
> *** caught segfault ***
> address 0x12f69d74, cause 'memory not mapped'
>
> Traceback:
> 1: .Call(.NAME, ..., PACKAGE = PACKAGE)
> 2: .Call2("Rle_getStartEndRunAndOffset", x, start, end, PACKAGE =
> "IRanges")
> 3: getStartEndRunAndOffset(x, start(solved_SEW), end(solved_SEW))
> 4: .local(x, ...)
> 5: window(Rle(), 1, 0)
> 6: window(Rle(), 1, 0)
>
> Possible actions:
> 1: abort (with core dump, if enabled)
> 2: normal R exit
> 3: exit R without saving workspace
> 4: exit R saving workspace
> Selection: 2
> Save workspace image? [y/n/c]:
>
> It is now fixed in release (S4Vectors 0.8.4 and IRanges 2.4.5) and
> devel (S4Vectors 0.9.12 and IRanges 2.5.10).
>
> Cheers,
> H.
>
>
>
> On 12/05/2015 03:29 PM, Leonard Goldstein wrote:
>>
>> Hi all,
>>
>> I ran into problems when using an RleList as column in a DataFrame
>> object (see example below).
>>
>> Many thanks in advance for your help.
>>
>> Leonard
>>
>>
>>> sessionInfo()
>>
>> R version 3.2.2 (2015-08-14)
>> Platform: x86_64-pc-linux-gnu (64-bit)
>> Running under: Red Hat Enterprise Linux Server release 6.6 (Santiago)
>>
>> 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] stats4 parallel stats graphics grDevices utils datasets
>> [8] methods base
>>
>> other attached packages:
>> [1] IRanges_2.4.4 S4Vectors_0.8.3 BiocGenerics_0.16.1
>>>
>>>
>>> df <- DataFrame(ID = 1:3)
>>> x <- RleList(IntegerList(vector("list", 3)))
>>> df$rle <- x
>>> df
>>
>> DataFrame with 3 rows and 2 columns
>>
>> *** caught segfault ***
>> address 0x9b00364, cause 'memory not mapped'
>>
>> Traceback:
>> 1: .Call(.NAME, ..., PACKAGE = PACKAGE)
>> 2: .Call2("Rle_getStartEndRunAndOffset", x, start, end, PACKAGE =
>> "S4Vectors")
>> 3: S4Vectors:::getStartEndRunAndOffset(x, start(solved_SEW),
>> end(solved_SEW))
>> 4: .local(x, ...)
>> 5: window(x, start = 1L, width = n)
>> 6: window(x, start = 1L, width = n)
>> 7: .local(x, ...)
>> 8: head(x, 3)
>> 9: head(x, 3)
>>
>> Possible actions:
>> 1: abort (with core dump, if enabled)
>> 2: normal R exit
>> 3: exit R without saving workspace
>> 4: exit R saving workspace
>> Selection:
>>
>> _______________________________________________
>> Bioc-devel at r-project.org mailing list
>> https://stat.ethz.ch/mailman/listinfo/bioc-devel
>>
>
> --
> 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