]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Enable MGTwoLevelTransfer for FESystem 10967/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 25 Sep 2020 06:40:05 +0000 (08:40 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Sat, 26 Sep 2020 10:24:57 +0000 (12:24 +0200)
include/deal.II/multigrid/mg_transfer_global_coarsening.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
tests/multigrid-global-coarsening/mg_transfer_a_02.cc [new file with mode: 0644]
tests/multigrid-global-coarsening/mg_transfer_a_02.with_p4est=true.mpirun=1.output [new file with mode: 0644]
tests/multigrid-global-coarsening/mg_transfer_p_03.cc [new file with mode: 0644]
tests/multigrid-global-coarsening/mg_transfer_p_03.with_p4est=true.mpirun=1.output [new file with mode: 0644]

index 7cb12faa55452f0f6c281fa162bc888a54b35dc5..ce75df0fd5dde16eef46d2f5595c0af4b63d5f7d 100644 (file)
@@ -227,6 +227,11 @@ private:
    */
   AffineConstraints<Number> constraint_coarse;
 
+  /**
+   * Number of components.
+   */
+  unsigned int n_components;
+
   friend class internal::MGTwoLevelTransferImplementation;
 };
 
index 23b63274e652062154de0ca4feb8143324ca6231..5d4cab9ddb9c46245f33848ce848fd7e9b354e88 100644 (file)
@@ -874,6 +874,16 @@ namespace internal
       // copy constrain matrix; TODO: why only for the coarse level?
       transfer.constraint_coarse.copy_from(constraint_coarse);
 
+      // make sure that no hp is used
+      AssertDimension(dof_handler_coarse.get_fe_collection().size(), 1);
+      AssertDimension(dof_handler_fine.get_fe_collection().size(), 1);
+
+      // extract number of components
+      AssertDimension(dof_handler_fine.get_fe().n_components(),
+                      dof_handler_coarse.get_fe().n_components());
+
+      transfer.n_components = dof_handler_fine.get_fe().n_components();
+
       // create partitioners and vectors for internal purposes
       {
         // ... for fine mesh
@@ -1264,6 +1274,13 @@ namespace internal
       AssertDimension(dof_handler_fine.get_triangulation().n_active_cells(),
                       dof_handler_coarse.get_triangulation().n_active_cells());
 
+      // extract number of components
+      AssertDimension(dof_handler_fine.get_fe_collection().n_components(),
+                      dof_handler_coarse.get_fe_collection().n_components());
+
+      transfer.n_components =
+        dof_handler_fine.get_fe_collection().n_components();
+
       const auto process_cells = [&](const auto &fu) {
         typename MeshType::active_cell_iterator            //
           cell_coarse = dof_handler_coarse.begin_active(), //
@@ -1615,7 +1632,9 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::prolongate(
   LinearAlgebra::distributed::Vector<Number> &      dst,
   const LinearAlgebra::distributed::Vector<Number> &src) const
 {
-  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
+  using VectorizedArrayType = VectorizedArray<Number>;
+
+  const unsigned int n_lanes = VectorizedArrayType::size();
 
   this->vec_coarse.copy_locally_owned_data_from(src);
   this->vec_coarse.update_ghost_values();
@@ -1625,33 +1644,34 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::prolongate(
 
   this->vec_fine = Number(0.);
 
+  AlignedVector<VectorizedArrayType> evaluation_data_fine;
+  AlignedVector<VectorizedArrayType> evaluation_data_coarse;
+
   for (const auto &scheme : schemes)
     {
-      AlignedVector<VectorizedArray<Number>> evaluation_data_fine(
-        scheme.dofs_per_cell_fine);
-      AlignedVector<VectorizedArray<Number>> evaluation_data_coarse(
-        scheme.dofs_per_cell_fine);
-
-      CellProlongator<dim, VectorizedArray<Number>> cell_prolongator(
-        scheme.prolongation_matrix_1d,
-        evaluation_data_coarse.begin(),
-        evaluation_data_fine.begin());
+      evaluation_data_fine.resize(scheme.dofs_per_cell_fine);
+      evaluation_data_coarse.resize(scheme.dofs_per_cell_fine);
+
       CellTransfer cell_transfer(scheme.degree_fine, scheme.degree_coarse);
 
-      for (unsigned int cell = 0; cell < scheme.n_coarse_cells;
-           cell += vec_size)
+      const unsigned int n_scalar_dofs_fine =
+        Utilities::pow(scheme.degree_fine + 1, dim);
+      const unsigned int n_scalar_dofs_coarse =
+        Utilities::pow(scheme.degree_coarse + 1, dim);
+
+      for (unsigned int cell = 0; cell < scheme.n_coarse_cells; cell += n_lanes)
         {
-          const unsigned int n_lanes =
-            (cell + vec_size > scheme.n_coarse_cells) ?
+          const unsigned int n_lanes_filled =
+            (cell + n_lanes > scheme.n_coarse_cells) ?
               (scheme.n_coarse_cells - cell) :
-              vec_size;
+              n_lanes;
 
           // read from source vector
           {
             const unsigned int *indices =
               &scheme
                  .level_dof_indices_coarse[cell * scheme.dofs_per_cell_coarse];
-            for (unsigned int v = 0; v < n_lanes; ++v)
+            for (unsigned int v = 0; v < n_lanes_filled; ++v)
               {
                 for (unsigned int i = 0; i < scheme.dofs_per_cell_coarse; ++i)
                   evaluation_data_coarse[i][v] =
@@ -1661,14 +1681,22 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::prolongate(
           }
 
           // ---------------------------- coarse -----------------------------
-          cell_transfer.run(cell_prolongator);
+          for (int c = n_components - 1; c >= 0; --c)
+            {
+              CellProlongator<dim, VectorizedArrayType> cell_prolongator(
+                scheme.prolongation_matrix_1d,
+                evaluation_data_coarse.begin() + c * n_scalar_dofs_coarse,
+                evaluation_data_fine.begin() + c * n_scalar_dofs_fine);
+
+              cell_transfer.run(cell_prolongator);
+            }
           // ------------------------------ fine -----------------------------
 
           if (scheme.fine_element_is_continuous)
             {
               const Number *w =
                 &scheme.weights[cell * scheme.dofs_per_cell_fine];
-              for (unsigned int v = 0; v < n_lanes; ++v)
+              for (unsigned int v = 0; v < n_lanes_filled; ++v)
                 {
                   for (unsigned int i = 0; i < scheme.dofs_per_cell_fine; ++i)
                     evaluation_data_fine[i][v] *= w[i];
@@ -1681,7 +1709,7 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::prolongate(
           {
             const unsigned int *indices =
               &scheme.level_dof_indices_fine[cell * scheme.dofs_per_cell_fine];
-            for (unsigned int v = 0; v < n_lanes; ++v)
+            for (unsigned int v = 0; v < n_lanes_filled; ++v)
               {
                 for (unsigned int i = 0; i < scheme.dofs_per_cell_fine; ++i)
                   this->vec_fine.local_element(indices[i]) +=
@@ -1709,39 +1737,42 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
   restrict_and_add(LinearAlgebra::distributed::Vector<Number> &      dst,
                    const LinearAlgebra::distributed::Vector<Number> &src) const
 {
-  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
+  using VectorizedArrayType = VectorizedArray<Number>;
+
+  const unsigned int n_lanes = VectorizedArrayType::size();
 
   this->vec_fine.copy_locally_owned_data_from(src);
   this->vec_fine.update_ghost_values();
 
   this->vec_coarse.copy_locally_owned_data_from(dst);
 
+  AlignedVector<VectorizedArrayType> evaluation_data_fine;
+  AlignedVector<VectorizedArrayType> evaluation_data_coarse;
+
   for (const auto &scheme : schemes)
     {
-      AlignedVector<VectorizedArray<Number>> evaluation_data_fine(
-        scheme.dofs_per_cell_fine);
-      AlignedVector<VectorizedArray<Number>> evaluation_data_coarse(
-        scheme.dofs_per_cell_fine);
-
-      CellRestrictor<dim, VectorizedArray<Number>> cell_restrictor(
-        scheme.prolongation_matrix_1d,
-        evaluation_data_fine.begin(),
-        evaluation_data_coarse.begin());
+      evaluation_data_fine.resize(scheme.dofs_per_cell_fine);
+      evaluation_data_coarse.resize(scheme.dofs_per_cell_fine);
+
       CellTransfer cell_transfer(scheme.degree_fine, scheme.degree_coarse);
 
-      for (unsigned int cell = 0; cell < scheme.n_coarse_cells;
-           cell += vec_size)
+      const unsigned int n_scalar_dofs_fine =
+        Utilities::pow(scheme.degree_fine + 1, dim);
+      const unsigned int n_scalar_dofs_coarse =
+        Utilities::pow(scheme.degree_coarse + 1, dim);
+
+      for (unsigned int cell = 0; cell < scheme.n_coarse_cells; cell += n_lanes)
         {
-          const unsigned int n_lanes =
-            (cell + vec_size > scheme.n_coarse_cells) ?
+          const unsigned int n_lanes_filled =
+            (cell + n_lanes > scheme.n_coarse_cells) ?
               (scheme.n_coarse_cells - cell) :
-              vec_size;
+              n_lanes;
 
           // read from source vector
           {
             const unsigned int *indices =
               &scheme.level_dof_indices_fine[cell * scheme.dofs_per_cell_fine];
-            for (unsigned int v = 0; v < n_lanes; ++v)
+            for (unsigned int v = 0; v < n_lanes_filled; ++v)
               {
                 for (unsigned int i = 0; i < scheme.dofs_per_cell_fine; ++i)
                   evaluation_data_fine[i][v] =
@@ -1754,7 +1785,7 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
             {
               const Number *w =
                 &scheme.weights[cell * scheme.dofs_per_cell_fine];
-              for (unsigned int v = 0; v < n_lanes; ++v)
+              for (unsigned int v = 0; v < n_lanes_filled; ++v)
                 {
                   for (unsigned int i = 0; i < scheme.dofs_per_cell_fine; ++i)
                     evaluation_data_fine[i][v] *= w[i];
@@ -1763,7 +1794,15 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
             }
 
           // ------------------------------ fine -----------------------------
-          cell_transfer.run(cell_restrictor);
+          for (int c = n_components - 1; c >= 0; --c)
+            {
+              CellRestrictor<dim, VectorizedArrayType> cell_restrictor(
+                scheme.prolongation_matrix_1d,
+                evaluation_data_fine.begin() + c * n_scalar_dofs_fine,
+                evaluation_data_coarse.begin() + c * n_scalar_dofs_coarse);
+
+              cell_transfer.run(cell_restrictor);
+            }
           // ----------------------------- coarse ----------------------------
 
 
@@ -1772,7 +1811,7 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
             const unsigned int *indices =
               &scheme
                  .level_dof_indices_coarse[cell * scheme.dofs_per_cell_coarse];
-            for (unsigned int v = 0; v < n_lanes; ++v)
+            for (unsigned int v = 0; v < n_lanes_filled; ++v)
               {
                 for (unsigned int i = 0; i < scheme.dofs_per_cell_coarse; ++i)
                   constraint_coarse.distribute_local_to_global(
diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_02.cc b/tests/multigrid-global-coarsening/mg_transfer_a_02.cc
new file mode 100644 (file)
index 0000000..00ffb5c
--- /dev/null
@@ -0,0 +1,140 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+/**
+ * Like mg_transfer_a_01.cc but for a system of finite elements.
+ */
+
+#include <deal.II/base/conditional_ostream.h>
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/mpi.h>
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.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/fe_tools.h>
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+
+#include <deal.II/matrix_free/fe_evaluation.h>
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include <deal.II/multigrid/mg_constrained_dofs.h>
+#include <deal.II/multigrid/mg_transfer_global_coarsening.h>
+
+#include "test.h"
+
+using namespace dealii;
+
+template <int dim, typename Number>
+void
+do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
+{
+  auto create_fine_grid = [](Triangulation<dim> &tria) {
+    GridGenerator::hyper_cube(tria);
+    tria.refine_global();
+
+    for (auto &cell : tria.active_cell_iterators())
+      if (cell->active() && cell->center()[0] < 0.5)
+        cell->set_refine_flag();
+    tria.execute_coarsening_and_refinement();
+  };
+
+  auto execute_global_coarsening = [](Triangulation<dim> &tria) {
+    for (auto &cell : tria.active_cell_iterators())
+      cell->set_coarsen_flag();
+    tria.execute_coarsening_and_refinement();
+  };
+
+  // create coarse grid
+  parallel::distributed::Triangulation<dim> tria_coarse(MPI_COMM_WORLD);
+  create_fine_grid(tria_coarse);
+  execute_global_coarsening(tria_coarse);
+
+  // create fine grid
+  parallel::distributed::Triangulation<dim> tria_fine(MPI_COMM_WORLD);
+  create_fine_grid(tria_fine);
+
+  // setup dof-handlers
+  DoFHandler<dim> dof_handler_fine(tria_fine);
+  dof_handler_fine.distribute_dofs(FESystem<dim>(fe_fine, dim));
+
+  DoFHandler<dim> dof_handler_coarse(tria_coarse);
+  dof_handler_coarse.distribute_dofs(FESystem<dim>(fe_coarse, dim));
+
+  // setup constraint matrix
+  AffineConstraints<Number> constraint_coarse;
+  DoFTools::make_hanging_node_constraints(dof_handler_coarse,
+                                          constraint_coarse);
+  constraint_coarse.close();
+
+  AffineConstraints<Number> constraint_fine;
+  DoFTools::make_hanging_node_constraints(dof_handler_fine, constraint_fine);
+  constraint_coarse.close();
+
+  // setup transfer operator
+  MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
+  transfer.reinit_geometric_transfer(dof_handler_fine,
+                                     dof_handler_coarse,
+                                     constraint_fine,
+                                     constraint_coarse);
+
+  test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
+}
+
+template <int dim, typename Number>
+void
+test(int fe_degree)
+{
+  const auto str_fine   = std::to_string(fe_degree);
+  const auto str_coarse = std::to_string(fe_degree);
+
+  {
+    deallog.push("CG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")");
+    do_test<dim, double>(FE_Q<dim>(fe_degree), FE_Q<dim>(fe_degree));
+    deallog.pop();
+  }
+
+  {
+    deallog.push("DG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")");
+    do_test<dim, double>(FE_DGQ<dim>(fe_degree), FE_Q<dim>(fe_degree));
+    deallog.pop();
+  }
+
+  {
+    deallog.push("DG<2>(" + str_fine + ")<->DG<2>(" + str_coarse + ")");
+    do_test<dim, double>(FE_DGQ<dim>(fe_degree), FE_DGQ<dim>(fe_degree));
+    deallog.pop();
+  }
+}
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  deallog.precision(8);
+
+  test<2, double>(2 /*=degree*/);
+}
diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_02.with_p4est=true.mpirun=1.output b/tests/multigrid-global-coarsening/mg_transfer_a_02.with_p4est=true.mpirun=1.output
new file mode 100644 (file)
index 0000000..a02a589
--- /dev/null
@@ -0,0 +1,67 @@
+
+DEAL:0:CG<2>(2)<->CG<2>(2)::weights:
+DEAL:0:CG<2>(2)<->CG<2>(2)::0.50000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.25000000 0.50000000 0.50000000 0.50000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.25000000 0.50000000 0.50000000 0.25000000 0.50000000 0.50000000 1.0000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 0.25000000 0.50000000 0.50000000 1.0000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:0:CG<2>(2)<->CG<2>(2)::0 12 2 8 16 10 4 14 6 1 13 3 9 17 11 5 15 7 4 14 6 22 28 24 18 26 20 5 15 7 23 29 25 19 27 21 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:CG<2>(2)<->CG<2>(2)::2 24 18 10 28 22 6 26 20 3 25 19 11 29 23 7 27 21 6 26 20 36 48 44 32 46 42 7 27 21 37 49 45 33 47 43 
+DEAL:0:CG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.0000000 0.0000000 0.0000000 0.0000000 1.0000000 0.0000000 0.0000000 0.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::weights:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.50000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.50000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.25000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 1.0000000 1.0000000 0.50000000 0.50000000 1.0000000 0.50000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:0:CG<2>(2)<->CG<2>(2)::30 42 32 32 52 0 38 46 40 40 56 50 34 44 36 36 54 48 34 44 36 36 54 48 62 68 64 64 74 70 58 66 60 60 72 4 31 43 33 33 53 1 39 47 41 41 57 51 35 45 37 37 55 49 35 45 37 37 55 49 63 69 65 65 75 71 59 67 61 61 73 5 58 66 60 60 72 4 80 86 82 82 94 90 76 84 78 78 92 88 76 84 78 78 92 88 100 106 102 102 112 108 96 104 98 98 110 18 59 67 61 61 73 5 81 87 83 83 95 91 77 85 79 79 93 89 77 85 79 79 93 89 101 107 103 103 113 109 97 105 99 99 111 19 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:CG<2>(2)<->CG<2>(2)::0 12 2 8 16 10 4 14 6 1 13 3 9 17 11 5 15 7 4 14 6 34 40 36 30 38 32 5 15 7 35 41 37 31 39 33 
+DEAL:0:CG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.0000000 0.37500000 0.0000000 0.0000000 -0.12500000 0.0000000 0.0000000 0.75000000 1.0000000 1.0000000 0.75000000 0.0000000 0.0000000 -0.12500000 0.0000000 0.0000000 0.37500000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::10.677078
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::17.903387
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.5625000 1.5625000 1.3125000 1.3125000 1.8750000 1.8750000 1.3750000 1.3750000 3.1250000 3.1250000 1.6250000 1.6250000 3.1250000 3.1250000 3.7500000 3.7500000 6.2500000 6.2500000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.5625000 1.5625000 1.3125000 1.3125000 3.1250000 3.1250000 1.6250000 1.6250000 3.1250000 3.1250000 6.2500000 6.2500000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::weights:
+DEAL:0:DG<2>(2)<->CG<2>(2)::
+DEAL:0:DG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:0:DG<2>(2)<->CG<2>(2)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 
+DEAL:0:DG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->CG<2>(2)::2 3 18 19 6 7 20 21 10 11 22 23 24 25 26 27 28 29 6 7 20 21 32 33 42 43 36 37 44 45 26 27 46 47 48 49 
+DEAL:0:DG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.0000000 0.0000000 0.0000000 0.0000000 1.0000000 0.0000000 0.0000000 0.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::weights:
+DEAL:0:DG<2>(2)<->CG<2>(2)::
+DEAL:0:DG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:0:DG<2>(2)<->CG<2>(2)::36 37 38 54 55 56 39 40 41 57 58 59 42 43 44 60 61 62 72 73 74 90 91 92 75 76 77 93 94 95 78 79 80 96 97 98 45 46 47 63 64 65 48 49 50 66 67 68 51 52 53 69 70 71 81 82 83 99 100 101 84 85 86 102 103 104 87 88 89 105 106 107 108 109 110 126 127 128 111 112 113 129 130 131 114 115 116 132 133 134 144 145 146 162 163 164 147 148 149 165 166 167 150 151 152 168 169 170 117 118 119 135 136 137 120 121 122 138 139 140 123 124 125 141 142 143 153 154 155 171 172 173 156 157 158 174 175 176 159 160 161 177 178 179 
+DEAL:0:DG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->CG<2>(2)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 4 5 6 7 30 31 32 33 34 35 36 37 14 15 38 39 40 41 
+DEAL:0:DG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.0000000 0.37500000 0.0000000 0.0000000 -0.12500000 0.0000000 0.0000000 0.75000000 1.0000000 1.0000000 0.75000000 0.0000000 0.0000000 -0.12500000 0.0000000 0.0000000 0.37500000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::13.416408
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::35.844826
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.5625000 4.3750000 2.5625000 5.3750000 13.812500 8.7500000 5.1250000 10.750000 1.5625000 1.5625000 5.3750000 2.5625000 4.3750000 12.250000 8.7500000 13.812500 4.3750000 1.5625000 1.0000000 1.0000000 2.0000000 2.0000000 1.0000000 1.0000000 1.0000000 1.0000000 2.0000000 2.0000000 1.0000000 1.0000000 12.250000 4.3750000 2.5625000 5.3750000 1.5625000 1.5625000 5.3750000 2.5625000 4.3750000 1.5625000 4.3750000 1.5625000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::weights:
+DEAL:0:DG<2>(2)<->DG<2>(2)::
+DEAL:0:DG<2>(2)<->DG<2>(2)::level_dof_indices_fine:
+DEAL:0:DG<2>(2)<->DG<2>(2)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 
+DEAL:0:DG<2>(2)<->DG<2>(2)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->DG<2>(2)::18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 
+DEAL:0:DG<2>(2)<->DG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.0000000 0.0000000 0.0000000 0.0000000 1.0000000 0.0000000 0.0000000 0.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::weights:
+DEAL:0:DG<2>(2)<->DG<2>(2)::
+DEAL:0:DG<2>(2)<->DG<2>(2)::level_dof_indices_fine:
+DEAL:0:DG<2>(2)<->DG<2>(2)::36 37 38 54 55 56 39 40 41 57 58 59 42 43 44 60 61 62 72 73 74 90 91 92 75 76 77 93 94 95 78 79 80 96 97 98 45 46 47 63 64 65 48 49 50 66 67 68 51 52 53 69 70 71 81 82 83 99 100 101 84 85 86 102 103 104 87 88 89 105 106 107 108 109 110 126 127 128 111 112 113 129 130 131 114 115 116 132 133 134 144 145 146 162 163 164 147 148 149 165 166 167 150 151 152 168 169 170 117 118 119 135 136 137 120 121 122 138 139 140 123 124 125 141 142 143 153 154 155 171 172 173 156 157 158 174 175 176 159 160 161 177 178 179 
+DEAL:0:DG<2>(2)<->DG<2>(2)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->DG<2>(2)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 
+DEAL:0:DG<2>(2)<->DG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.0000000 0.37500000 0.0000000 0.0000000 -0.12500000 0.0000000 0.0000000 0.75000000 1.0000000 1.0000000 0.75000000 0.0000000 0.0000000 -0.12500000 0.0000000 0.0000000 0.37500000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::13.416408
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::31.329898
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.5625000 4.3750000 1.5625000 4.3750000 12.250000 4.3750000 1.5625000 4.3750000 1.5625000 1.5625000 4.3750000 1.5625000 4.3750000 12.250000 4.3750000 1.5625000 4.3750000 1.5625000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.5625000 4.3750000 1.5625000 4.3750000 12.250000 4.3750000 1.5625000 4.3750000 1.5625000 1.5625000 4.3750000 1.5625000 4.3750000 12.250000 4.3750000 1.5625000 4.3750000 1.5625000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_03.cc b/tests/multigrid-global-coarsening/mg_transfer_p_03.cc
new file mode 100644 (file)
index 0000000..2e0b039
--- /dev/null
@@ -0,0 +1,128 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.
+//
+// ---------------------------------------------------------------------
+
+
+/**
+ * Like mg_transfer_p_01.cc but for a system of finite elements.
+ */
+
+#include <deal.II/base/conditional_ostream.h>
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/mpi.h>
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.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/fe_tools.h>
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+
+#include <deal.II/matrix_free/fe_evaluation.h>
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include <deal.II/multigrid/mg_constrained_dofs.h>
+#include <deal.II/multigrid/mg_transfer_global_coarsening.h>
+
+#include "test.h"
+
+using namespace dealii;
+
+template <int dim, typename Number>
+void
+do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
+{
+  parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+
+  // create grid
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global();
+
+  for (auto &cell : tria.active_cell_iterators())
+    if (cell->active() && cell->center()[0] < 0.5)
+      cell->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+
+  // setup dof-handlers
+  DoFHandler<dim> dof_handler_fine(tria);
+  dof_handler_fine.distribute_dofs(FESystem<dim>(fe_fine, dim));
+
+  DoFHandler<dim> dof_handler_coarse(tria);
+  dof_handler_coarse.distribute_dofs(FESystem<dim>(fe_coarse, dim));
+
+  AffineConstraints<Number> constraint_coarse;
+  DoFTools::make_hanging_node_constraints(dof_handler_coarse,
+                                          constraint_coarse);
+  constraint_coarse.close();
+
+  AffineConstraints<Number> constraint_fine;
+  DoFTools::make_hanging_node_constraints(dof_handler_fine, constraint_fine);
+  constraint_fine.close();
+
+  // setup transfer operator
+  MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
+  transfer.reinit_polynomial_transfer(dof_handler_fine,
+                                      dof_handler_coarse,
+                                      constraint_fine,
+                                      constraint_coarse);
+
+  test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
+}
+
+template <int dim, typename Number>
+void
+test(int fe_degree_fine, int fe_degree_coarse)
+{
+  const auto str_fine   = std::to_string(fe_degree_fine);
+  const auto str_coarse = std::to_string(fe_degree_coarse);
+
+  {
+    deallog.push("CG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")");
+    do_test<dim, Number>(FE_Q<dim>(fe_degree_fine),
+                         FE_Q<dim>(fe_degree_coarse));
+    deallog.pop();
+  }
+
+  {
+    deallog.push("DG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")");
+    do_test<dim, Number>(FE_DGQ<dim>(fe_degree_fine),
+                         FE_Q<dim>(fe_degree_coarse));
+    deallog.pop();
+  }
+
+  {
+    deallog.push("DG<2>(" + str_fine + ")<->DG<2>(" + str_coarse + ")");
+    do_test<dim, Number>(FE_DGQ<dim>(fe_degree_fine),
+                         FE_DGQ<dim>(fe_degree_coarse));
+    deallog.pop();
+  }
+}
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  deallog.precision(8);
+
+  test<2, double>(2, 1);
+}
diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_03.with_p4est=true.mpirun=1.output b/tests/multigrid-global-coarsening/mg_transfer_p_03.with_p4est=true.mpirun=1.output
new file mode 100644 (file)
index 0000000..eb71e25
--- /dev/null
@@ -0,0 +1,43 @@
+
+DEAL:0:CG<2>(2)<->CG<2>(1)::weights:
+DEAL:0:CG<2>(2)<->CG<2>(1)::0.50000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.25000000 0.50000000 0.50000000 0.50000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.25000000 0.50000000 0.50000000 0.25000000 0.50000000 0.50000000 1.0000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 0.25000000 0.50000000 0.50000000 1.0000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 0.50000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.50000000 1.0000000 0.50000000 0.50000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 0.50000000 1.0000000 0.0000000 0.25000000 0.50000000 0.25000000 0.25000000 0.50000000 0.0000000 0.50000000 1.0000000 0.0000000 0.25000000 0.50000000 0.25000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 0.25000000 0.50000000 0.25000000 0.50000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.25000000 0.50000000 0.25000000 0.50000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 0.50000000 0.50000000 0.50000000 0.25000000 1.0000000 1.0000000 0.50000000 1.0000000 1.0000000 0.50000000 0.25000000 0.50000000 0.0000000 0.50000000 1.0000000 0.0000000 0.50000000 1.0000000 0.50000000 0.25000000 0.50000000 0.0000000 0.50000000 1.0000000 0.0000000 0.50000000 1.0000000 0.50000000 
+DEAL:0:CG<2>(2)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:0:CG<2>(2)<->CG<2>(1)::0 12 2 8 16 10 4 14 6 1 13 3 9 17 11 5 15 7 4 14 6 22 28 24 18 26 20 5 15 7 23 29 25 19 27 21 30 42 32 38 46 40 34 44 36 31 43 33 39 47 41 35 45 37 32 52 0 40 56 50 36 54 48 33 53 1 41 57 51 37 55 49 34 44 36 62 68 64 58 66 60 35 45 37 63 69 65 59 67 61 36 54 48 64 74 70 60 72 4 37 55 49 65 75 71 61 73 5 58 66 60 80 86 82 76 84 78 59 67 61 81 87 83 77 85 79 60 72 4 82 94 90 78 92 88 61 73 5 83 95 91 79 93 89 76 84 78 100 106 102 96 104 98 77 85 79 101 107 103 97 105 99 78 92 88 102 112 108 98 110 18 79 93 89 103 113 109 99 111 19 
+DEAL:0:CG<2>(2)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:0:CG<2>(2)<->CG<2>(1)::0 2 4 6 1 3 5 7 4 6 8 10 5 7 9 11 12 14 16 18 13 15 17 19 14 0 18 20 15 1 19 21 16 18 22 24 17 19 23 25 18 20 24 4 19 21 25 5 22 24 26 28 23 25 27 29 24 4 28 30 25 5 29 31 26 28 32 34 27 29 33 35 28 30 34 8 29 31 35 9 
+DEAL:0:CG<2>(2)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:0:CG<2>(2)<->CG<2>(1)::1.0000000 0.50000000 1.1102230e-16 1.1102230e-16 0.50000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(1)::10.677078
+DEAL:0:CG<2>(2)<->CG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(1)::18.533753
+DEAL:0:CG<2>(2)<->CG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:CG<2>(2)<->CG<2>(1)::3.5000000 3.5000000 2.2500000 2.2500000 5.0000000 5.0000000 3.0000000 3.0000000 3.5000000 3.5000000 2.2500000 2.2500000 2.2500000 2.2500000 3.0000000 3.0000000 3.0000000 3.0000000 4.0000000 4.0000000 0.0000000 0.0000000 3.0000000 3.0000000 4.0000000 4.0000000 3.0000000 3.0000000 4.0000000 4.0000000 0.0000000 0.0000000 2.2500000 2.2500000 3.0000000 3.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(1)::weights:
+DEAL:0:DG<2>(2)<->CG<2>(1)::
+DEAL:0:DG<2>(2)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:0:DG<2>(2)<->CG<2>(1)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 
+DEAL:0:DG<2>(2)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->CG<2>(1)::0 2 4 6 1 3 5 7 4 6 8 10 5 7 9 11 12 14 16 18 13 15 17 19 14 0 18 20 15 1 19 21 16 18 22 24 17 19 23 25 18 20 24 4 19 21 25 5 22 24 26 28 23 25 27 29 24 4 28 30 25 5 29 31 26 28 32 34 27 29 33 35 28 30 34 8 29 31 35 9 
+DEAL:0:DG<2>(2)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->CG<2>(1)::1.0000000 0.50000000 1.1102230e-16 1.1102230e-16 0.50000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(1)::13.416408
+DEAL:0:DG<2>(2)<->CG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(1)::36.280160
+DEAL:0:DG<2>(2)<->CG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->CG<2>(1)::6.7500000 6.7500000 2.2500000 2.2500000 13.500000 13.500000 4.5000000 4.5000000 6.7500000 6.7500000 2.2500000 2.2500000 2.2500000 2.2500000 4.5000000 4.5000000 4.5000000 4.5000000 9.0000000 9.0000000 0.0000000 0.0000000 4.5000000 4.5000000 9.0000000 9.0000000 4.5000000 4.5000000 9.0000000 9.0000000 0.0000000 0.0000000 2.2500000 2.2500000 4.5000000 4.5000000 
+DEAL:0:DG<2>(2)<->DG<2>(1)::weights:
+DEAL:0:DG<2>(2)<->DG<2>(1)::
+DEAL:0:DG<2>(2)<->DG<2>(1)::level_dof_indices_fine:
+DEAL:0:DG<2>(2)<->DG<2>(1)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 
+DEAL:0:DG<2>(2)<->DG<2>(1)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->DG<2>(1)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 
+DEAL:0:DG<2>(2)<->DG<2>(1)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->DG<2>(1)::1.0000000 0.50000000 1.1102230e-16 1.1102230e-16 0.50000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(1)::13.416408
+DEAL:0:DG<2>(2)<->DG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(1)::20.124612
+DEAL:0:DG<2>(2)<->DG<2>(1)::1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 
+DEAL:0:DG<2>(2)<->DG<2>(1)::2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 2.2500000 

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.