Seg fault using Extracting hanging node dofs in MPI

24 views
Skip to first unread message

Phani Motamarri

unread,
Dec 8, 2017, 1:55:50 AM12/8/17
to dea...@googlegroups.com
Dear All,

I am trying to create a vector of type bool to identify which dofs correspond to hanging nodes. I  use  the function

DoFTools::extract_hanging_node_dofs(dofHandler,selectedDofsHanging).

As mentioned in the manual, I set the size of the vector selectedDofsHanging to dofHandler.n_dofs().

This works fine in serial but when running in parallel on 8 procs, the above function gives a seg fault at line 854 of dof_tools.cc

  // dof numbers on the centers of the lines bounding this
   // face
   for (unsigned int line=0; line<4; ++line)
      for (unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
            selected_dofs[face->line(line)->child(0)->vertex_dof_index(1,dof)] = true;

The line where the code seg faults has been marked bold.

My hanging node mesh generation process starts with a base uniform mesh containing 64 elements in a cube of domain 40 units and I kept refining just the layer of elements around the center of cube till the mesh size becomes less than 0.05. Final mesh contains around 512 elements.


I am wondering if I am  missing anything here when using the above function in MPI.

Thanks and Regards
Phani

Wolfgang Bangerth

unread,
Dec 8, 2017, 10:52:12 AM12/8/17
to dea...@googlegroups.com
On 12/07/2017 11:55 PM, Phani Motamarri wrote:
>
> This works fine in serial but when running in parallel on 8 procs, the above
> function gives a seg fault at line 854 of dof_tools.cc
>
> // dof numbers on the centers of the lines bounding this
> // face
> for (unsigned int line=0; line<4; ++line)
> for (unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
> *selected_dofs[face->line(line)->child(0)->vertex_dof_index(1,dof)] = true;*
>
> The line where the code seg faults has been marked bold.
>
> My hanging node mesh generation process starts with a base uniform mesh
> containing 64 elements in a cube of domain 40 units and I kept refining just
> the layer of elements around the center of cube till the mesh size becomes
> less than 0.05. Final mesh contains around 512 elements.

This may be a bug, but it's hard to tell without a testcase that allows
someone to look into it. Can you come up with a small program that
demonstrates the issue? It doesn't have to do anything useful -- just build a
DoFHandler in parallel and call this function.

Best
W.

--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Daniel Arndt

unread,
Dec 9, 2017, 6:26:32 AM12/9/17
to deal.II User Group
Phani,

apparently no one used this function for a parallel::distributed::Triangulation before. The problem is that we don't know about DoFs on artficial cells.
The pull request https://github.com/dealii/dealii/pull/5604 fixes this issue by extracting DoFs belonging to locally owned cells.
Essentially, you only have to guard the cell loops in extract_hanging_node_dofs by    

if (cell->is_locally_owned())

Best,
Daniel

Phani Motamarri

unread,
Dec 9, 2017, 8:45:40 AM12/9/17
to dea...@googlegroups.com
Oh thank you very much for the pull request. For now, I was extracting DoFs corresponding to hanging nodes by using is_constrained(dof).in parallel.



--
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+unsubscribe@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply all
Reply to author
Forward
0 new messages