]> https://gitweb.dealii.org/ - dealii.git/commitdiff
New test: p::d::tria + hp::dofh + hanging_node_constraints. 9555/head
authorMarc Fehling <marc.fehling@gmx.net>
Thu, 20 Feb 2020 12:50:05 +0000 (13:50 +0100)
committerMarc Fehling <marc.fehling@gmx.net>
Thu, 20 Feb 2020 13:14:39 +0000 (14:14 +0100)
tests/mpi/hp_hanging_node_constraints_01.cc [moved from tests/mpi/hp_hanging_node_constraints.cc with 100% similarity]
tests/mpi/hp_hanging_node_constraints_01.with_p4est=true.mpirun=2.output [moved from tests/mpi/hp_hanging_node_constraints.with_p4est=true.mpirun=2.output with 100% similarity]
tests/mpi/hp_hanging_node_constraints_02.cc [new file with mode: 0644]
tests/mpi/hp_hanging_node_constraints_02.with_p4est=true.mpirun=2.output [new file with mode: 0644]

diff --git a/tests/mpi/hp_hanging_node_constraints_02.cc b/tests/mpi/hp_hanging_node_constraints_02.cc
new file mode 100644 (file)
index 0000000..e14b5cd
--- /dev/null
@@ -0,0 +1,123 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// this test failed at some stage:
+// DoFTools::make_hanging_node_constraints() has not worked with an
+// hp::DoFHandler on ghost cells that neighbor artificial cells.
+//
+// see: https://github.com/dealii/dealii/issues/9517
+
+
+#include <deal.II/base/index_set.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_nothing.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/hp/dof_handler.h>
+#include <deal.II/hp/fe_collection.h>
+
+#include <deal.II/lac/affine_constraints.h>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test()
+{
+  // setup triangulation
+  parallel::distributed::Triangulation<dim> triangulation(
+    MPI_COMM_WORLD,
+    typename Triangulation<dim>::MeshSmoothing(
+      Triangulation<dim>::smoothing_on_refinement |
+      Triangulation<dim>::smoothing_on_coarsening));
+  GridGenerator::hyper_rectangle(
+    triangulation,
+    (dim == 3 ? Point<dim>(0.0, 0.0, 0.0) : Point<dim>(0.0, 0.0)),
+    (dim == 3 ? Point<dim>(1.0, 1.0, 1.0) : Point<dim>(1.0, 1.0)),
+    true);
+  triangulation.refine_global(3);
+
+  // Refine all the cells with y < 0.5 twice
+  for (unsigned int n_ref = 0; n_ref < 2; ++n_ref)
+    {
+      for (const auto &cell : triangulation.active_cell_iterators())
+        if (cell->is_locally_owned() && cell->center()(1) < 0.5)
+          cell->set_refine_flag();
+
+      triangulation.prepare_coarsening_and_refinement();
+      triangulation.execute_coarsening_and_refinement();
+    }
+
+  // setup finite elements
+  FESystem<dim> void_fe(FE_Nothing<dim>(), dim);
+  FESystem<dim> solid_fe(FE_Q<dim>(1), dim);
+
+  hp::FECollection<dim> fe_collection;
+  fe_collection.push_back(void_fe);
+  fe_collection.push_back(solid_fe);
+
+  hp::DoFHandler<dim> dof_handler(triangulation);
+  dof_handler.set_fe(fe_collection);
+
+  // Assign void_fe to all the cells with x < 0.5
+  for (const auto &cell : dof_handler.active_cell_iterators())
+    if (cell->is_locally_owned())
+      {
+        if (cell->center()(0) < 0.5)
+          cell->set_active_fe_index(0);
+        else
+          cell->set_active_fe_index(1);
+      }
+
+  dof_handler.distribute_dofs(fe_collection);
+
+  // make constraints
+  IndexSet locally_relevant_dofs;
+  DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+
+  AffineConstraints<double> hanging_node_constraints;
+  hanging_node_constraints.reinit(locally_relevant_dofs);
+  DoFTools::make_hanging_node_constraints(dof_handler,
+                                          hanging_node_constraints);
+
+  deallog << "OK" << std::endl;
+}
+
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+  MPILogInitAll all;
+
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+  MPI_Barrier(MPI_COMM_WORLD);
+  deallog.push("3d");
+  test<3>();
+  deallog.pop();
+}
diff --git a/tests/mpi/hp_hanging_node_constraints_02.with_p4est=true.mpirun=2.output b/tests/mpi/hp_hanging_node_constraints_02.with_p4est=true.mpirun=2.output
new file mode 100644 (file)
index 0000000..4aa9062
--- /dev/null
@@ -0,0 +1,7 @@
+
+DEAL:0:2d::OK
+DEAL:0:3d::OK
+
+DEAL:1:2d::OK
+DEAL:1:3d::OK
+

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.