]> https://gitweb.dealii.org/ - dealii.git/commitdiff
New test case
authorMartin Kronbichler <martin.kronbichler@rub.de>
Mon, 22 Jul 2024 14:15:52 +0000 (16:15 +0200)
committerMartin Kronbichler <martin.kronbichler@rub.de>
Mon, 22 Jul 2024 14:15:52 +0000 (16:15 +0200)
tests/multigrid-global-coarsening/mg_transfer_p_11.cc [new file with mode: 0644]
tests/multigrid-global-coarsening/mg_transfer_p_11.output [new file with mode: 0644]

diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_11.cc b/tests/multigrid-global-coarsening/mg_transfer_p_11.cc
new file mode 100644 (file)
index 0000000..6ed90fa
--- /dev/null
@@ -0,0 +1,191 @@
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2024 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+
+/**
+ * Test MGTwoLevelTransfer set up with MatrixFree objects using an FESystem,
+ * otherwise the same as mg_transfer_p_08.cc
+ */
+
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include <deal.II/multigrid/mg_transfer_global_coarsening.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+
+template <int dim, typename Number>
+void
+do_test(const FiniteElement<dim> &fe_fine,
+        const FiniteElement<dim> &fe_coarse,
+        const Triangulation<dim> &tria)
+{
+  deallog << "Testing " << fe_fine.get_name() << " <-> " << fe_coarse.get_name()
+          << std::endl;
+  // setup dof-handlers
+  DoFHandler<dim> dof_handler_fine(tria);
+  dof_handler_fine.distribute_dofs(fe_fine);
+
+  DoFHandler<dim> dof_handler_coarse(tria);
+  dof_handler_coarse.distribute_dofs(fe_coarse);
+
+  FE_DGQ<dim>     fe_dummy(0);
+  DoFHandler<dim> dof_handler_dummy(tria);
+  dof_handler_dummy.distribute_dofs(fe_dummy);
+  AffineConstraints<Number> constraints_dummy;
+  constraints_dummy.close();
+
+  AffineConstraints<Number> constraints_coarse(
+    dof_handler_coarse.locally_owned_dofs(),
+    DoFTools::extract_locally_relevant_dofs(dof_handler_coarse));
+  DoFTools::make_hanging_node_constraints(dof_handler_coarse,
+                                          constraints_coarse);
+  constraints_coarse.close();
+
+  AffineConstraints<Number> constraints_fine(
+    dof_handler_fine.locally_owned_dofs(),
+    DoFTools::extract_locally_relevant_dofs(dof_handler_fine));
+  DoFTools::make_hanging_node_constraints(dof_handler_fine, constraints_fine);
+  constraints_fine.close();
+
+  // setup transfer operator
+  MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
+    transfer_ref;
+  transfer_ref.reinit(dof_handler_fine,
+                      dof_handler_coarse,
+                      constraints_fine,
+                      constraints_coarse);
+
+  MappingQ1<dim>          mapping;
+  MatrixFree<dim, Number> mf_coarse, mf_fine;
+  mf_coarse.reinit(mapping,
+                   std::vector<const DoFHandler<dim> *>{
+                     {&dof_handler_dummy, &dof_handler_coarse}},
+                   std::vector<const AffineConstraints<double> *>{
+                     {&constraints_dummy, &constraints_coarse}},
+                   std::vector<Quadrature<1>>{{QGauss<1>(1)}},
+                   typename MatrixFree<dim, Number>::AdditionalData());
+  mf_fine.reinit(mapping,
+                 std::vector<const DoFHandler<dim> *>{{&dof_handler_coarse,
+                                                       &dof_handler_dummy,
+                                                       &dof_handler_fine,
+                                                       &dof_handler_dummy}},
+                 std::vector<const AffineConstraints<double> *>{
+                   {&constraints_coarse,
+                    &constraints_dummy,
+                    &constraints_fine,
+                    &constraints_dummy}},
+                 std::vector<Quadrature<1>>{{QGauss<1>(1)}},
+                 typename MatrixFree<dim, Number>::AdditionalData());
+  MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
+  transfer.reinit(mf_fine, 2, mf_coarse, 1);
+
+  LinearAlgebra::distributed::Vector<Number> vec_fine_0, vec_fine_1,
+    vec_coarse_0, vec_coarse_1;
+  mf_fine.initialize_dof_vector(vec_fine_0, 2);
+  mf_fine.initialize_dof_vector(vec_fine_1, 2);
+  mf_coarse.initialize_dof_vector(vec_coarse_0, 1);
+  mf_coarse.initialize_dof_vector(vec_coarse_1, 1);
+
+  // prolongate
+  for (auto &d : vec_coarse_0)
+    d = random_value<double>();
+  transfer_ref.prolongate_and_add(vec_fine_0, vec_coarse_0);
+  transfer.prolongate_and_add(vec_fine_1, vec_coarse_0);
+  deallog << "Prolongate: " << vec_fine_0.l2_norm() << " "
+          << vec_fine_1.l2_norm() << " ";
+  vec_fine_1 -= vec_fine_0;
+  deallog << vec_fine_1.l2_norm() << std::endl;
+
+  vec_coarse_0 = 0;
+  transfer_ref.restrict_and_add(vec_coarse_0, vec_fine_0);
+  vec_coarse_1 = 0;
+  transfer.restrict_and_add(vec_coarse_1, vec_fine_0);
+  deallog << "Restrict: " << vec_coarse_0.l2_norm() << " "
+          << vec_coarse_1.l2_norm() << " ";
+  vec_coarse_1 -= vec_coarse_0;
+  deallog << vec_coarse_1.l2_norm() << std::endl << std::endl;
+}
+
+
+
+template <int dim, typename Number>
+void
+test(int fe_degree_fine, int fe_degree_coarse)
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria, -1., 1.);
+  tria.refine_global(5 - dim);
+
+  deallog << "Uniform grid" << std::endl << std::endl;
+  do_test<dim, Number>(FESystem<dim>(FE_Q<dim>(fe_degree_fine), dim),
+                       FESystem<dim>(FE_Q<dim>(fe_degree_coarse), dim),
+                       tria);
+  do_test<dim, Number>(FESystem<dim>(FE_DGQ<dim>(fe_degree_fine), dim),
+                       FESystem<dim>(FE_Q<dim>(fe_degree_coarse), dim),
+                       tria);
+  do_test<dim, Number>(FESystem<dim>(FE_DGQ<dim>(fe_degree_fine), dim),
+                       FESystem<dim>(FE_DGQ<dim>(fe_degree_coarse), dim),
+                       tria);
+
+  for (const auto &cell : tria.active_cell_iterators())
+    if (cell->center()[0] < 0 || cell->center()[1] < 0 ||
+        cell->center()[dim - 1] < 0)
+      cell->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+
+  deallog << "Adaptively refined grid" << std::endl << std::endl;
+  do_test<dim, Number>(FESystem<dim>(FE_Q<dim>(fe_degree_fine), dim),
+                       FESystem<dim>(FE_Q<dim>(fe_degree_coarse), dim),
+                       tria);
+  do_test<dim, Number>(FESystem<dim>(FE_DGQ<dim>(fe_degree_fine), dim),
+                       FESystem<dim>(FE_Q<dim>(fe_degree_coarse), dim),
+                       tria);
+  do_test<dim, Number>(FESystem<dim>(FE_DGQHermite<dim>(fe_degree_fine), dim),
+                       FESystem<dim>(FE_DGQ<dim>(fe_degree_coarse), dim),
+                       tria);
+}
+
+
+
+int
+main()
+{
+  initlog();
+
+  MultithreadInfo::set_thread_limit(1);
+
+  deallog.precision(8);
+
+  test<2, double>(2, 1);
+  test<2, double>(3, 1);
+  test<2, double>(1, 1);
+  test<3, double>(2, 1);
+}
diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_11.output b/tests/multigrid-global-coarsening/mg_transfer_p_11.output
new file mode 100644 (file)
index 0000000..ea9587f
--- /dev/null
@@ -0,0 +1,113 @@
+
+DEAL::Uniform grid
+DEAL::
+DEAL::Testing FESystem<2>[FE_Q<2>(2)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 13.663649 13.663649 2.7962968e-16
+DEAL::Restrict: 25.420254 25.420254 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(2)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 20.269046 20.269046 0.0000000
+DEAL::Restrict: 56.395388 56.395388 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(2)^2] <-> FESystem<2>[FE_DGQ<2>(1)^2]
+DEAL::Prolongate: 18.775923 18.775923 0.0000000
+DEAL::Restrict: 27.500564 27.500564 0.0000000
+DEAL::
+DEAL::Adaptively refined grid
+DEAL::
+DEAL::Testing FESystem<2>[FE_Q<2>(2)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 22.524324 22.524324 5.3797688e-16
+DEAL::Restrict: 42.637382 42.637382 1.1749496e-15
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(2)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 33.899239 33.899239 0.0000000
+DEAL::Restrict: 96.747920 96.747920 3.2023728e-15
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQHermite<2>(2)^2] <-> FESystem<2>[FE_DGQ<2>(1)^2]
+DEAL::Prolongate: 33.672479 33.672479 0.0000000
+DEAL::Restrict: 49.253827 49.253827 0.0000000
+DEAL::
+DEAL::Uniform grid
+DEAL::
+DEAL::Testing FESystem<2>[FE_Q<2>(3)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 19.757343 19.757343 3.2671802e-16
+DEAL::Restrict: 54.074714 54.074714 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(3)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 25.082979 25.082979 0.0000000
+DEAL::Restrict: 92.116340 92.116340 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(3)^2] <-> FESystem<2>[FE_DGQ<2>(1)^2]
+DEAL::Prolongate: 24.715226 24.715226 0.0000000
+DEAL::Restrict: 48.264625 48.264625 0.0000000
+DEAL::
+DEAL::Adaptively refined grid
+DEAL::
+DEAL::Testing FESystem<2>[FE_Q<2>(3)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 34.600965 34.600965 7.6179449e-16
+DEAL::Restrict: 98.563611 98.563611 2.0350724e-15
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(3)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 43.956683 43.956683 0.0000000
+DEAL::Restrict: 166.90548 166.90548 3.3528000e-15
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQHermite<2>(3)^2] <-> FESystem<2>[FE_DGQ<2>(1)^2]
+DEAL::Prolongate: 45.004243 45.004243 0.0000000
+DEAL::Restrict: 87.866557 87.866557 0.0000000
+DEAL::
+DEAL::Uniform grid
+DEAL::
+DEAL::Testing FESystem<2>[FE_Q<2>(1)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 7.3364299 7.3364299 0.0000000
+DEAL::Restrict: 7.3364299 7.3364299 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(1)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 13.244635 13.244635 0.0000000
+DEAL::Restrict: 24.977650 24.977650 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(1)^2] <-> FESystem<2>[FE_DGQ<2>(1)^2]
+DEAL::Prolongate: 12.972939 12.972939 0.0000000
+DEAL::Restrict: 12.972939 12.972939 0.0000000
+DEAL::
+DEAL::Adaptively refined grid
+DEAL::
+DEAL::Testing FESystem<2>[FE_Q<2>(1)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 12.653242 12.653242 0.0000000
+DEAL::Restrict: 12.653242 12.653242 0.0000000
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQ<2>(1)^2] <-> FESystem<2>[FE_Q<2>(1)^2]
+DEAL::Prolongate: 22.803418 22.803418 0.0000000
+DEAL::Restrict: 44.733710 44.733710 1.4054301e-15
+DEAL::
+DEAL::Testing FESystem<2>[FE_DGQHermite<2>(1)^2] <-> FESystem<2>[FE_DGQ<2>(1)^2]
+DEAL::Prolongate: 23.131741 23.131741 0.0000000
+DEAL::Restrict: 23.131741 23.131741 0.0000000
+DEAL::
+DEAL::Uniform grid
+DEAL::
+DEAL::Testing FESystem<3>[FE_Q<3>(2)^3] <-> FESystem<3>[FE_Q<3>(1)^3]
+DEAL::Prolongate: 25.768428 25.768428 6.3656552e-16
+DEAL::Restrict: 61.986736 61.986736 0.0000000
+DEAL::
+DEAL::Testing FESystem<3>[FE_DGQ<3>(2)^3] <-> FESystem<3>[FE_Q<3>(1)^3]
+DEAL::Prolongate: 39.653205 39.653205 0.0000000
+DEAL::Restrict: 162.26875 162.26875 0.0000000
+DEAL::
+DEAL::Testing FESystem<3>[FE_DGQ<3>(2)^3] <-> FESystem<3>[FE_DGQ<3>(1)^3]
+DEAL::Prolongate: 38.581375 38.581375 0.0000000
+DEAL::Restrict: 68.771113 68.771113 0.0000000
+DEAL::
+DEAL::Adaptively refined grid
+DEAL::
+DEAL::Testing FESystem<3>[FE_Q<3>(2)^3] <-> FESystem<3>[FE_Q<3>(1)^3]
+DEAL::Prolongate: 61.402632 61.402632 1.7178432e-15
+DEAL::Restrict: 159.07644 159.07644 5.5288661e-15
+DEAL::
+DEAL::Testing FESystem<3>[FE_DGQ<3>(2)^3] <-> FESystem<3>[FE_Q<3>(1)^3]
+DEAL::Prolongate: 101.67905 101.67905 0.0000000
+DEAL::Restrict: 482.16372 482.16372 2.7818565e-14
+DEAL::
+DEAL::Testing FESystem<3>[FE_DGQHermite<3>(2)^3] <-> FESystem<3>[FE_DGQ<3>(1)^3]
+DEAL::Prolongate: 104.98305 104.98305 0.0000000
+DEAL::Restrict: 187.36024 187.36024 0.0000000
+DEAL::

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.