[Bioc-devel] Bug in DEXSeq, me thinks

Steve Lianoglou mailinglist.honeypot at gmail.com
Tue Mar 13 10:33:41 CET 2012


And while we're at it, I think there's another place you need to swap
out an lapply for an sapply (or set simplify=FALSE) in the
`divideWork` function.

Specifically, it seems that when `testablegenes` is divisible by mc.cores:

(1) This sapply call returns a matrix back to `rowlist`, when you
really want a list, since you will `do.call(c, rowlist)` later:

rowlist <- sapply(1:mc.cores, function(i){which(fData(ecs)$geneID %in%
subgenes[[i]])})

(2) Ditto for your `stopifnot` check. The `do.call(c, sapply(...))`
explodes because the sapply returns a matrix:

stopifnot(all(as.character(fData(ecs)$geneID)[rows] == do.call(c,
sapply(allecs, function(x){as.character(fData(x)$geneID)}))))

... or ... at least these things are happening for me right now.

I have a weird feeling that I might have hosed something, since I
reckon people would have tripped over this by now since it's easy to
fall into this scenario where:

length(testablegenes) %% mc.cores == 0

when mc.cores == 2

So ... hopefully I'm not just chasing my tail then asking you to chase yours.

Thanks,
-steve




On Tue, Mar 13, 2012 at 5:10 AM, Steve Lianoglou
<mailinglist.honeypot at gmail.com> wrote:
> And then add names back to muhats, of course:
>
> On Mon, Mar 12, 2012 at 8:32 PM, Steve Lianoglou
> <mailinglist.honeypot at gmail.com> wrote:
>> Hi,
>>
>> It just so happened that DEXSeq, for reasons I couldn't explain, kept
>> throwing errors during a call to `estimateDispersions`. I couldn't
>> explain it because I really thought I didn't change much, but I've got
>> some magic combination of samples I'm currently trying to
>> estimateDispersion on that none of them throw any errors while
>> estimate `muhats` (line 245 of core_functions.R):
>>
>> muhats <- sapply(testablegenes, function(gn) {
>> ...
>> muhat
>> })
> [snip]
>
> names(muhats) <- testablegenes
>
> -steve
> --
> Steve Lianoglou
> Graduate Student: Computational Systems Biology
>  | Memorial Sloan-Kettering Cancer Center
>  | Weill Medical College of Cornell University
> Contact Info: http://cbio.mskcc.org/~lianos/contact



-- 
Steve Lianoglou
Graduate Student: Computational Systems Biology
 | Memorial Sloan-Kettering Cancer Center
 | Weill Medical College of Cornell University
Contact Info: http://cbio.mskcc.org/~lianos/contact



More information about the Bioc-devel mailing list