[Bioc-devel] OpenMP error when using ShortRead::FastqStreamer with parallel::parLapply

Martin Morgan mtmorgan at fhcrc.org
Mon May 27 17:29:06 CEST 2013


On 05/26/2013 11:57 PM, Michael Stadler wrote:
> Dear bioc core developers,
>
> We have recently (a few months back?) started to get error messages of
> the type:
>
> libgomp: Thread creation failed: Resource temporarily unavailable
>
> when using ShortRead::FastqStreamer in combination with
> parallel::parLapply (for example in QuasR::preprocessReads or
> QuasR::qQCReport). These errors are not fully reproducible, are less
> frequent when reducing the number of nodes in the cluster object and
> disappear completely when FastqStreamer is used directly rather than
> within a parLapply call of the sort:
>
> cl <- makeCluster(n)
> results <- parLapply(cl, filenames, function(filename) {
>    fs <- FastqStreamer(filename, n=1e6)
>    while (length(chunks <- yield(fs)) != 0L) {
>      # do work
>      ...
>    }
>    close(fs)
> })
>
> I guess that this could be a conflict of nested parallelization and/or a
> resource limit of our setup. If I am not mistaken, FastqStreamer makes
> use of OpenMP for parsing fastq files.
>
> My questions is how the "libgomp" error message can be avoided in the
> above situations. Is it possible to switch off the use of OpenMP in
> FastqStreamer and friends at run-time, other than switching off OpenMP
> at compile time for R or the whole ShortRead package?

Hi Michael --

Your diagnosis sounds on target. It is possible to set the number of openMP 
threads, which I think will address the problem. You'd like to not interfere 
with others' use of openMP, so reset the thread count when your function exits. 
Along the lines of

     nthreads <- .Call(ShortRead:::.set_omp_threads, 1L)
     on.exit(.Call(ShortRead:::.set_omp_threads, nthreads))

ShortRead's srapply function does this internally when it detects that the 
`parallel` package is in use (reasoning that mclapply and forked processes would 
be the common use case for single-machine parallel evaluation). I'll give some 
thought to more aggressively discovering when it is necessary to invoke 
FastqStreamer with this option set.

Sorry for the trouble,

Martin

>
> With kind regards,
> Michael
>
>
> My sessionInfo():
> R version 3.0.1 (2013-05-16)
> Platform: x86_64-unknown-linux-gnu (64-bit)
>
> locale:
> [1] C
>
> attached base packages:
> [1] parallel  stats     graphics  grDevices utils     datasets  methods
> [8] base
>
> other attached packages:
> [1] ShortRead_1.18.0     latticeExtra_0.6-24  Rsamtools_1.12.3
> [4] lattice_0.20-15      Biostrings_2.28.0    GenomicRanges_1.12.3
> [7] IRanges_1.18.1       BiocGenerics_0.7.0   RColorBrewer_1.0-5
>
> loaded via a namespace (and not attached):
> [1] Biobase_2.20.0 bitops_1.0-5   grid_3.0.1     hwriter_1.3
> stats4_3.0.1
> [6] tcltk_3.0.1    tools_3.0.1    zlibbioc_1.6.0
>
> _______________________________________________
> Bioc-devel at r-project.org mailing list
> https://stat.ethz.ch/mailman/listinfo/bioc-devel
>


-- 
Computational Biology / Fred Hutchinson Cancer Research Center
1100 Fairview Ave. N.
PO Box 19024 Seattle, WA 98109

Location: Arnold Building M1 B861
Phone: (206) 667-2793



More information about the Bioc-devel mailing list