Jack,

are your solvers using MPI?  This looks similar to this problem 
https://github.com/open-mpi/ompi/issues/1081 Which version of MPI are you 
using? How do you initialize MPI? MPI_InitFinalize set 
MPI_THREAD_SERIALIZED which "tells MPI that 
we might use several threads but never call two MPI functions at the same 
time." This is not want you want if your solvers use MPI. Then you need to 
initialize MPI and TBB yourself and use MPI THREAD MULTIPLE.

Best,

Bruno

On Friday, May 19, 2017 at 3:52:39 AM UTC-4, Jack wrote:
>
> Dear all,
>
>  
>
> I’m trying to solve the thermal diffusion and Stokes flow problem 
> simultaneously, similar to the step-32. I opened two threads to solve 
> thermal diffusion and Stokes equations during solving by linear solvers 
> (the former is solved by CG solver and the latter is solved by GMRES 
> solver) after assembling the matrices and RHS. I encountered a problem when 
> I using parallel computing in combination with distributed and shared 
> memory.
>
>  The code is like following:
>
>
> assemble_thermal_diffusion_system();
>
> assemble_stokes_flow_system ();
>
> Threads::Task<void> task =
>
> Threads::new_task (&Problem<dim>::solve_thermal_diffsuion, *this);
>
> solve_stokes_flow ();
>
> task.join();
>
>
> The program ran successfully with single process and showed a high 
> performance of thread parallel computing, unfortunately, when I ran the 
> program using several processes I got following error. I have googled this 
> error and can not find any solution. Could someone help me?
>
> Thanks in advance!
>
>  
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
> [warn] opal_libevent2021_event_base_loop: reentrant invocation.  Only one 
> event_base_loop can run on each event_base at once.
>
>
> Thank you very much.
>
>  
>
> Best,
>
>  
>
> Jack
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to