Hello,
I want to use the deal.II hp::DoFHandler in a principally similar way as in 
step-46 - so using it with a coupled problem. I'm at the beginning of 
implementing this and want to use deal.II with Trilinos, so that the 
program can run on distributed machines.
Unfortunately I have a problem when it comes to creating the 
TrilinosWrappers::BlockSparsityPattern: the code runs serial (mpirun -np 1 
./problem_fe_nothing) but crashes with np >= 2.
I could isolate the error a little bit and have also created a "minimal 
working example" (which is more or less an adapted step-40 and attached at 
this post).

To explain my problem in a nutshell:
I need to create a FESystem which has one FE_Nothing block:
FESystem<dim> fe_problem(FE_Q<dim>(2)^dim, FE_Nothing<dim>()^1);
I add this FESystem to a hp::FECollection<dim> and distribute the dofs via 
the hp::DoFHandler<dim>.
Afterwards I reorder the DoFs component wise and extract the locally owned 
and locally relevant DoFs.
I then reinit my TrilinosWrappers::BlockSparsityPattern:
sp.reinit(locally_owned_dofs,
          locally_owned_dofs,
          locally_relevant_dofs,
          mpi_communicator);
create the sparsity pattern afterwards with 
DoFTools::make_sparsity_pattern() and then do a:
sp.compress()
and this is where the program crashes (in parallel) with an uncaught 
exception in main() - so there is no exception thrown in a underlying 
deal.II or Trilinos routine (using debug mode).


I have found out, that the program does not crash, if I don't use 
FE_Nothing, so if I use:
FESystem<dim> fe_works(FE_Q<dim>(2)^dim, FE_W<dim>()^1);
instead of fe_problem the code runs fine - even in parallel.
Also if I use fe_problem in conjunction with:
sp.reinit(locally_owned_dofs,
          locally_owned_dofs,
          mpi_communicator);
the code runs fine (but will according to the documentation result in a not 
multi thread writable matrix).

I can debug the program in serial with gdb, but unfortunately this doesn't 
help me here since the problem does only occur in parallel.
As proposed in debug mpi 
<https://www.open-mpi.org/faq/?category=debugging#serial-debuggers>I tried 
to debug it with:
mpirun -np 2 urxvt -e gdb ./problem_fe_nothing
After crashing here "backtrace" just leads to "no stack" and is of no help 
either (everything compiled in debug mode).

For the sake of completeness:
I use the brand new deal.II 9.1.0 with Trilinos from the current master 
branch, gcc 8.3 and mpich 3.3.

It would be great if anybody could help me finding the problem here or 
debugging the problem further it would be great.

Thanks and greetings,

Mathias

-- 
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 [email protected].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/dcbdfa64-358f-495d-9768-50da65433ebc%40googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

<<attachment: problem_fe_nothing.zip>>

Reply via email to