Hello,

I am having the following code to generate a mesh (code at the end). It 
crashes if I use more than 1 MPI processes (-np 1 works fine, but still 
have the WARNING message), with the following error message:

--------------------------------------------------------------------------

WARNING: There is at least non-excluded one OpenFabrics device found,

but there are no active ports detected (or Open MPI was unable to use

them).  This is most certainly not what you wanted.  Check your

cables, subnet manager configuration, etc.  The openib BTL will be

ignored for this job.


  Local host: ganymede

--------------------------------------------------------------------------

[ganymede:21801] 1 more process has sent help message help-mpi-btl-openib.txt 
/ no active ports found


[ganymede:21801] Set MCA parameter "orte_base_help_aggregate" to 0 to see 
all help / error messages


[ganymede:21803] *** Process received signal ***

[ganymede:21803] Signal: Segmentation fault (11)

[ganymede:21803] Signal code: Address not mapped (1)

[ganymede:21803] Failing at address: 0x40171e9b8

[ganymede:21803] [ 0] /lib64/libpthread.so.0(+0xf5e0)[0x7fc11d30d5e0]

[ganymede:21803] [ 1] /data/yiyang.zhang/deal_II_trilinos/deal.II-v8.5.1/lib
/libdeal_II.so.8.5.1(
_ZN6dealii13TriangulationILi2ELi2EE33prepare_coarsening_and_refinementEv+
0x731)[0x7fc1238247f1]

[ganymede:21803] [ 2] /data/yiyang.zhang/deal_II_trilinos/deal.II-v8.5.1/lib
/libdeal_II.so.8.5.1(
_ZN6dealii8parallel11distributed13TriangulationILi2ELi2EE33prepare_coarsening_and_refinementEv
+0xe9)[0x7fc123ae2949]

[ganymede:21803] [ 3] /data/yiyang.zhang/deal_II_trilinos/deal.II-v8.5.1/lib
/libdeal_II.so.8.5.1(
_ZN6dealii8parallel11distributed13TriangulationILi2ELi2EE33execute_coarsening_and_refinementEv
+0x4e)[0x7fc123ae212e]

[ganymede:21803] [ 4] /data/yiyang.zhang/deal_II_trilinos/deal.II-v8.5.1/lib
/libdeal_II.so.8.5.1(_ZN6dealii13TriangulationILi2ELi2EE13refine_globalEj+
0xee)[0x7fc1237a163e]

[ganymede:21803] [ 5] ./string[0x44d864]

[ganymede:21803] [ 6] ./string[0x4427cb]

[ganymede:21803] [ 7] ./string[0x43fc8d]

[ganymede:21803] [ 8] ./string[0x43d07b]

[ganymede:21803] [ 9] ./string[0x42ee46]

[ganymede:21803] [10] /lib64/libc.so.6(__libc_start_main+0xf5)[
0x7fc11cf5cc05]

[ganymede:21803] [11] ./string[0x42dc99]

[ganymede:21803] *** End of error message ***

--------------------------------------------------------------------------

mpirun noticed that process rank 0 with PID 21803 on node ganymede exited 
on signal 11 (Segmentation fault).

With std::cout method, it shows the code crashes at tria.refine_global(
refine);

If this line is commented out, this function will finish and return, but 
the code will still crash in the following refine_mesh() function.

I tried the following test code in dealii test folder. All of them passed.

/mpi/refinement_listener_01.cc : PASSED
/mpi/refinement_listener_02.cc : PASSED

/mpi/refine_and_coarsen_fixed_number_01.cc : PASSED
/mpi/refine_and_coarsen_fixed_number_07.cc : PASSED

/distributed_grids/solution_transfer_01.cc : PASSED (1 process)
/distributed_grids/solution_transfer_02.cc : PASSED (1 process)
/distributed_grids/solution_transfer_03.cc : PASSED (1 process)
/distributed_grids/solution_transfer_04.cc : PASSED (1 process, 2 processes)

Thank you!


template<int dim>

    void MyClass<dim>::make_ball_mesh_2d(const Real radius, const int refine
){

        //generate a hyper_ball tria with center at (0,0), and radius = 
radius

        const Real isCenter = 1e-3; //some small value

        const types::manifold_id MFD_ID_SPH = 1;

        GridGenerator::hyper_ball(tria, Point<dim>(), radius);



        //set up post_refinement signal

        //such that boundary indicators are refreshed after every 
refinement.

        //Ref. parallel::distributed::Triangulation<> class

        tria.signals.post_refinement.connect(std::bind(&MyClass<dim>::
set_boundary_ids,

                                                       std::cref(*this),

                                                       std::ref(tria) ) );



        //set manifold id to get spherical manifold description

        tria.set_all_manifold_ids_on_boundary(MFD_ID_SPH);

        tria.set_manifold(MFD_ID_SPH, spherical_mfd);



        //set boundary ID

        set_boundary_ids(tria);



        const Point<dim> mesh_center;

        const Real core_radius = 1.0 / 100.0 * radius;

        const Real inner_radius = 1.0 / 50.0 * radius;

        // Step 1: Shrink the inner cell


        for (typename Triangulation<dim>::active_cell_iterator cell = tria.
begin_active();

             cell != tria.end(); ++cell) {

            if(cell->is_locally_owned()){

                if (mesh_center.distance(cell->center()) < isCenter)

                    //cell center almost on the mesh_center

                {

                    for (unsigned int v = 0;

                         v < GeometryInfo<dim>::vertices_per_cell;

                         ++v)

                        cell->vertex(v) *= core_radius / mesh_center.
distance(cell->vertex(v));

                    //the vertices of this cell now have a radius of 
core_radius.

                }

            }

        }

        //when cells are moved locally, all the MPI processes need to be 
notified

        //by using void Triangulation< dim, spacedim 
>::communicate_locally_moved_vertices().

        sync_moved_cells();



        // Step 2: Refine all cells except the central one

        for (typename Triangulation<dim>::active_cell_iterator cell = tria.
begin_active();

             cell != tria.end(); ++cell) {

            if(cell->is_locally_owned()){

                if (mesh_center.distance(cell->center()) >= isCenter)

                    cell->set_refine_flag();

            }

        }

        tria.execute_coarsening_and_refinement();



        // Step 3: Resize the inner children of the outer cells

        for (typename Triangulation<dim>::active_cell_iterator cell = tria.
begin_active();

             cell != tria.end(); ++cell) {

            if(cell->is_locally_owned()){

                for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_
per_cell; ++v)

                {

                    const double dist = mesh_center.distance(cell->vertex(v
));

                    if (dist > core_radius*1.0001 && dist < 0.9999*radius)

                        cell->vertex(v) *= inner_radius / dist;

                }

            }

        }

        sync_moved_cells();



        // Step 4: Apply curved manifold description

        for (typename Triangulation<dim>::active_cell_iterator cell = tria.
begin_active();

             cell != tria.end(); ++cell) {

            if(cell->is_locally_owned()){

                bool is_in_inner_circle = false;

                for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_
per_cell; ++v)

                    if (mesh_center.distance(cell->vertex(v)) < inner_radius
)

                    {

                        is_in_inner_circle = true;

                        break;

                    }

                if (is_in_inner_circle == false)

                    cell->set_all_manifold_ids(MFD_ID_SPH);

            }

        }

        tria.refine_global(refine);

        return;

  }

    template<int dim>
  void MyClass<dim>::set_boundary_ids(parallel::distributed::Triangulation
<dim>& tria_) const {
      for (typename Triangulation<dim>::active_cell_iterator
           cell = tria_.begin_active();
           cell != tria_.end();
           ++cell){
          if(cell->is_locally_owned()){
              for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell
; ++f)
                  if (cell->face(f)->at_boundary())
                      cell->face(f)->set_all_boundary_ids(BD_INFTY);
          }
      }
      return;
  }

   template<int dim>
  void MyClass<dim>::sync_moved_cells(){
      std::vector<bool> locally_owned_vertices = GridTools::
get_locally_owned_vertices(tria);
      tria.communicate_locally_moved_vertices(locally_owned_vertices);
      return;
 }


-- 
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