]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added tests for AffineConstraints diagnostics. 12623/head
authorMarc Fehling <mafehling.git@gmail.com>
Wed, 4 Aug 2021 23:42:05 +0000 (17:42 -0600)
committerMarc Fehling <marc.fehling@gmx.net>
Thu, 5 Aug 2021 03:35:35 +0000 (21:35 -0600)
tests/hp/identity_constraints_01.cc [new file with mode: 0644]
tests/hp/identity_constraints_01.output [new file with mode: 0644]
tests/hp/identity_constraints_02.cc [new file with mode: 0644]
tests/hp/identity_constraints_02.output [new file with mode: 0644]
tests/lac/inhomogeneous_constraints.cc
tests/lac/inhomogeneous_constraints.output

diff --git a/tests/hp/identity_constraints_01.cc b/tests/hp/identity_constraints_01.cc
new file mode 100644 (file)
index 0000000..28470d9
--- /dev/null
@@ -0,0 +1,85 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Create a 2x1x1 grid with FE_Q elements with degrees 2 & 4 assigned.
+// Verify that we do not unify dofs on lines in 3D.
+// On each of the four lines on the interface, the central dofs are identical
+// and will be treated with constraints.
+//
+// +----+----+
+// | Q2 | Q4 |
+// +----+----+
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include <deal.II/lac/affine_constraints.h>
+
+#include "../tests.h"
+
+#include "../test_grids.h"
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> tria;
+  TestGrids::hyper_line(tria, 2);
+
+  hp::FECollection<dim> fe;
+  fe.push_back(FE_Q<dim>(4));
+  fe.push_back(FE_Q<dim>(2));
+
+  DoFHandler<dim> dh(tria);
+  dh.begin_active()->set_active_fe_index(1);
+  dh.distribute_dofs(fe);
+
+  AffineConstraints<double> constraints;
+  DoFTools::make_hanging_node_constraints(dh, constraints);
+  constraints.close();
+
+  deallog << "Total constraints:          " << constraints.n_constraints()
+          << std::endl
+          << "  Inhomogenous constraints: " << constraints.n_inhomogeneities()
+          << std::endl
+          << "  Identity constraints:     " << constraints.n_identities()
+          << std::endl;
+}
+
+
+int
+main()
+{
+  initlog();
+
+  deallog.push("1d");
+  test<1>();
+  deallog.pop();
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+  deallog.push("3d");
+  test<3>();
+  deallog.pop();
+}
diff --git a/tests/hp/identity_constraints_01.output b/tests/hp/identity_constraints_01.output
new file mode 100644 (file)
index 0000000..163e459
--- /dev/null
@@ -0,0 +1,10 @@
+
+DEAL:1d::Total constraints:          0
+DEAL:1d::  Inhomogenous constraints: 0
+DEAL:1d::  Identity constraints:     0
+DEAL:2d::Total constraints:          2
+DEAL:2d::  Inhomogenous constraints: 0
+DEAL:2d::  Identity constraints:     0
+DEAL:3d::Total constraints:          20
+DEAL:3d::  Inhomogenous constraints: 0
+DEAL:3d::  Identity constraints:     4
diff --git a/tests/hp/identity_constraints_02.cc b/tests/hp/identity_constraints_02.cc
new file mode 100644 (file)
index 0000000..14103e1
--- /dev/null
@@ -0,0 +1,128 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Create a 2x2x1 grid with FE_Q elements with degrees 1 - 4 and 2 - 5 assigned
+// in two respective scenarios.
+// Verify that we do not unify dofs on lines in 3D.
+// On each of the four lines on the interface between the Q2 and Q4 element, the
+// central dofs are identical and will be treated with constraints.
+//
+// If the dominating element in the collection of finite elements on the central
+// line does not have a central dof, the dof duplicates on the elements Q2 and
+// Q4 will not be recognized with constraints.
+// This is the reason why there is one identity constraint less in scenario 1
+// than in scenario 2 -- the dominating Q1 element has no dofs on lines.
+// This is a bug.
+//
+// Scenario 1:    Scenario 2:
+// +----+----+    +----+----+
+// | Q3 | Q4 |    | Q4 | Q5 |
+// +----+----+    +----+----+
+// | Q1 | Q2 |    | Q2 | Q3 |
+// +----+----+    +----+----+
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include <deal.II/lac/affine_constraints.h>
+
+#include <array>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test(std::array<unsigned int, 4> fe_degrees)
+{
+  Triangulation<dim> tria;
+  {
+    std::vector<unsigned int> repetitions(dim);
+    Point<dim>                bottom_left, top_right;
+    for (unsigned int d = 0; d < dim; ++d)
+      if (d < 2)
+        {
+          repetitions[d] = 2;
+          bottom_left[d] = -1.;
+          top_right[d]   = 1.;
+        }
+      else
+        {
+          repetitions[d] = 1;
+          bottom_left[d] = 0.;
+          top_right[d]   = 1.;
+        }
+    GridGenerator::subdivided_hyper_rectangle(tria,
+                                              repetitions,
+                                              bottom_left,
+                                              top_right);
+  }
+
+  hp::FECollection<dim> fe;
+  for (const auto d : fe_degrees)
+    fe.push_back(FE_Q<dim>(d));
+
+  DoFHandler<dim> dh(tria);
+  {
+    unsigned int i = 0;
+    for (const auto &cell : dh.active_cell_iterators())
+      cell->set_active_fe_index(i++);
+    Assert(i == 4, ExcInternalError());
+  }
+  dh.distribute_dofs(fe);
+
+  AffineConstraints<double> constraints;
+  DoFTools::make_hanging_node_constraints(dh, constraints);
+  constraints.close();
+
+  deallog << "Total constraints:          " << constraints.n_constraints()
+          << std::endl
+          << "  Inhomogenous constraints: " << constraints.n_inhomogeneities()
+          << std::endl
+          << "  Identity constraints:     " << constraints.n_identities()
+          << std::endl;
+}
+
+
+int
+main()
+{
+  initlog();
+
+  deallog << "FE degrees: 1 - 4" << std::endl;
+  deallog.push("2d");
+  test<2>({{1, 2, 3, 4}});
+  deallog.pop();
+  deallog.push("3d");
+  test<3>({{1, 2, 3, 4}});
+  deallog.pop();
+
+  deallog << "FE degrees: 2 - 5" << std::endl;
+  deallog.push("2d");
+  test<2>({{2, 3, 4, 5}});
+  deallog.pop();
+  deallog.push("3d");
+  test<3>({{2, 3, 4, 5}});
+  deallog.pop();
+}
diff --git a/tests/hp/identity_constraints_02.output b/tests/hp/identity_constraints_02.output
new file mode 100644 (file)
index 0000000..7eee87e
--- /dev/null
@@ -0,0 +1,15 @@
+
+DEAL::FE degrees: 1 - 4
+DEAL:2d::Total constraints:          8
+DEAL:2d::  Inhomogenous constraints: 0
+DEAL:2d::  Identity constraints:     0
+DEAL:3d::Total constraints:          55
+DEAL:3d::  Inhomogenous constraints: 0
+DEAL:3d::  Identity constraints:     3
+DEAL::FE degrees: 2 - 5
+DEAL:2d::Total constraints:          12
+DEAL:2d::  Inhomogenous constraints: 0
+DEAL:2d::  Identity constraints:     0
+DEAL:3d::Total constraints:          92
+DEAL:3d::  Inhomogenous constraints: 0
+DEAL:3d::  Identity constraints:     4
index 088dcc6c48c8d6c6517622ed9f1632505a332561..72bf5dc3b4772cbd9e4c9c7948dbd0ab234efdd8 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2009 - 2020 by the deal.II authors
+// Copyright (C) 2009 - 2021 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -635,7 +635,11 @@ LaplaceProblem<dim>::run()
               << "   Number of hanging node constraints: "
               << hanging_nodes_only.n_constraints() << std::endl
               << "   Total number of constraints:        "
-              << test_all_constraints.n_constraints() << std::endl;
+              << test_all_constraints.n_constraints() << std::endl
+              << "   Number of inhomogenous constraints: "
+              << test_all_constraints.n_inhomogeneities() << std::endl
+              << "   Number of identity constraints:     "
+              << test_all_constraints.n_identities() << std::endl;
 
       assemble_reference();
       assemble_test_1();
index 9a196421b0056808916ba4181634b086c55187f0..4567d42c839be05253c937d1a8db511ececc7f34 100644 (file)
@@ -5,6 +5,8 @@ DEAL:2d::   Number of active cells:             192
 DEAL:2d::   Number of degrees of freedom:       864
 DEAL:2d::   Number of hanging node constraints: 0
 DEAL:2d::   Total number of constraints:        192
+DEAL:2d::   Number of inhomogenous constraints: 188
+DEAL:2d::   Number of identity constraints:     0
 DEAL:2d::  Reference matrix nonzeros: 12672, actually: 8480
 DEAL:2d::  Test matrix 1 nonzeros: 12672, actually: 8480
 DEAL:2d::  Matrix difference norm: 0
@@ -21,6 +23,8 @@ DEAL:2d::   Number of active cells:             360
 DEAL:2d::   Number of degrees of freedom:       1588
 DEAL:2d::   Number of hanging node constraints: 36
 DEAL:2d::   Total number of constraints:        284
+DEAL:2d::   Number of inhomogenous constraints: 254
+DEAL:2d::   Number of identity constraints:     12
 DEAL:2d::  Reference matrix nonzeros: 24092, actually: 17520
 DEAL:2d::  Test matrix 1 nonzeros: 24092, actually: 17520
 DEAL:2d::  Matrix difference norm: 0
@@ -37,6 +41,8 @@ DEAL:2d::   Number of active cells:             690
 DEAL:2d::   Number of degrees of freedom:       2982
 DEAL:2d::   Number of hanging node constraints: 108
 DEAL:2d::   Total number of constraints:        408
+DEAL:2d::   Number of inhomogenous constraints: 314
+DEAL:2d::   Number of identity constraints:     36
 DEAL:2d::  Reference matrix nonzeros: 46420, actually: 36248
 DEAL:2d::  Test matrix 1 nonzeros: 46420, actually: 36248
 DEAL:2d::  Matrix difference norm: 0
@@ -53,6 +59,8 @@ DEAL:3d::   Number of active cells:             8
 DEAL:3d::   Number of degrees of freedom:       27
 DEAL:3d::   Number of hanging node constraints: 0
 DEAL:3d::   Total number of constraints:        26
+DEAL:3d::   Number of inhomogenous constraints: 25
+DEAL:3d::   Number of identity constraints:     0
 DEAL:3d::  Reference matrix nonzeros: 343, actually: 27
 DEAL:3d::  Test matrix 1 nonzeros: 343, actually: 27
 DEAL:3d::  Matrix difference norm: 0
@@ -69,6 +77,8 @@ DEAL:3d::   Number of active cells:             22
 DEAL:3d::   Number of degrees of freedom:       60
 DEAL:3d::   Number of hanging node constraints: 16
 DEAL:3d::   Total number of constraints:        56
+DEAL:3d::   Number of inhomogenous constraints: 55
+DEAL:3d::   Number of identity constraints:     0
 DEAL:3d::  Reference matrix nonzeros: 1006, actually: 66
 DEAL:3d::  Test matrix 1 nonzeros: 1006, actually: 66
 DEAL:3d::  Matrix difference norm: 0
@@ -85,6 +95,8 @@ DEAL:3d::   Number of active cells:             92
 DEAL:3d::   Number of degrees of freedom:       189
 DEAL:3d::   Number of hanging node constraints: 55
 DEAL:3d::   Total number of constraints:        158
+DEAL:3d::   Number of inhomogenous constraints: 148
+DEAL:3d::   Number of identity constraints:     0
 DEAL:3d::  Reference matrix nonzeros: 3767, actually: 415
 DEAL:3d::  Test matrix 1 nonzeros: 3767, actually: 415
 DEAL:3d::  Matrix difference norm: 0

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.