]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Enable MGTwoLevelTransfer::reinit_geometric_transfer for FECollections 11578/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 18 Jan 2021 18:17:00 +0000 (19:17 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 9 Feb 2021 09:42:17 +0000 (10:42 +0100)
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
tests/multigrid-global-coarsening/mg_transfer_a_03.cc [new file with mode: 0644]
tests/multigrid-global-coarsening/mg_transfer_a_03.with_p4est=true.mpirun=2.output [new file with mode: 0644]

index 4985e280006d9c21ea249b5263837e9c18d14e72..04f1f7e534bbdd51d0528af4a231a13d53851526 100644 (file)
@@ -868,15 +868,52 @@ 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);
+      // gather ranges for active FE indices on both fine and coarse dofhandlers
+      std::array<unsigned int, 2> min_active_fe_indices = {
+        {std::numeric_limits<unsigned int>::max(),
+         std::numeric_limits<unsigned int>::max()}};
+      std::array<unsigned int, 2> max_active_fe_indices = {{0, 0}};
+
+      for (const auto &cell : dof_handler_fine.active_cell_iterators())
+        if (cell->is_locally_owned())
+          {
+            min_active_fe_indices[0] =
+              std::min(min_active_fe_indices[0], cell->active_fe_index());
+            max_active_fe_indices[0] =
+              std::max(max_active_fe_indices[0], cell->active_fe_index());
+          }
+
+      for (const auto &cell : dof_handler_coarse.active_cell_iterators())
+        if (cell->is_locally_owned())
+          {
+            min_active_fe_indices[1] =
+              std::min(min_active_fe_indices[1], cell->active_fe_index());
+            max_active_fe_indices[1] =
+              std::max(max_active_fe_indices[1], cell->active_fe_index());
+          }
+
+      const auto comm = get_mpi_comm(dof_handler_fine);
+
+      Assert(comm == get_mpi_comm(dof_handler_coarse), ExcNotImplemented());
+
+      ArrayView<unsigned int> temp_min(min_active_fe_indices);
+      ArrayView<unsigned int> temp_max(max_active_fe_indices);
+      Utilities::MPI::min(temp_min, comm, temp_min);
+      Utilities::MPI::max(temp_max, comm, temp_max);
+
+      // make sure that hp is used neither on the coarse nor on the fine
+      // dofhandler
+      AssertDimension(min_active_fe_indices[0], max_active_fe_indices[0]);
+      AssertDimension(min_active_fe_indices[1], max_active_fe_indices[1]);
+
+      const auto &fe_fine = dof_handler_fine.get_fe(min_active_fe_indices[0]);
+      const auto &fe_coarse =
+        dof_handler_coarse.get_fe(min_active_fe_indices[1]);
 
       // extract number of components
-      AssertDimension(dof_handler_fine.get_fe().n_components(),
-                      dof_handler_coarse.get_fe().n_components());
+      AssertDimension(fe_fine.n_components(), fe_coarse.n_components());
 
-      transfer.n_components = dof_handler_fine.get_fe().n_components();
+      transfer.n_components = fe_fine.n_components();
 
       // create partitioners and vectors for internal purposes
       {
@@ -939,28 +976,24 @@ namespace internal
       transfer.schemes.resize(2);
 
       // check if FE is the same; TODO: better way?
-      AssertDimension(dof_handler_coarse.get_fe(0).dofs_per_cell,
-                      dof_handler_fine.get_fe(0).dofs_per_cell);
+      AssertDimension(fe_coarse.dofs_per_cell, fe_fine.dofs_per_cell);
 
       // number of dofs on coarse and fine cells
       transfer.schemes[0].dofs_per_cell_coarse =
         transfer.schemes[0].dofs_per_cell_fine =
-          transfer.schemes[1].dofs_per_cell_coarse =
-            dof_handler_coarse.get_fe(0).dofs_per_cell;
+          transfer.schemes[1].dofs_per_cell_coarse = fe_coarse.dofs_per_cell;
       transfer.schemes[1].dofs_per_cell_fine =
-        dof_handler_coarse.get_fe(0).dofs_per_cell *
-        GeometryInfo<dim>::max_children_per_cell;
+        fe_coarse.dofs_per_cell * GeometryInfo<dim>::max_children_per_cell;
 
       // degree of FE on coarse and fine cell
       transfer.schemes[0].degree_coarse   = transfer.schemes[0].degree_fine =
-        transfer.schemes[1].degree_coarse = dof_handler_coarse.get_fe(0).degree;
-      transfer.schemes[1].degree_fine =
-        dof_handler_coarse.get_fe(0).degree * 2 + 1;
+        transfer.schemes[1].degree_coarse = fe_coarse.degree;
+      transfer.schemes[1].degree_fine     = fe_coarse.degree * 2 + 1;
 
       // continuous or discontinuous
       transfer.schemes[0].fine_element_is_continuous =
         transfer.schemes[1].fine_element_is_continuous =
-          dof_handler_fine.get_fe(0).dofs_per_vertex > 0;
+          fe_fine.dofs_per_vertex > 0;
 
       // count coarse cells for each scheme (0, 1)
       {
@@ -981,8 +1014,8 @@ namespace internal
 
       const auto cell_local_chilren_indices =
         get_child_offsets<dim>(transfer.schemes[0].dofs_per_cell_coarse,
-                               dof_handler_fine.get_fe(0).degree + 1,
-                               dof_handler_fine.get_fe(0).degree);
+                               fe_fine.degree + 1,
+                               fe_fine.degree);
 
 
       // indices
@@ -1010,7 +1043,7 @@ namespace internal
           const Quadrature<1> dummy_quadrature(
             std::vector<Point<1>>(1, Point<1>()));
           internal::MatrixFreeFunctions::ShapeInfo<Number> shape_info;
-          shape_info.reinit(dummy_quadrature, dof_handler_fine.get_fe(0), 0);
+          shape_info.reinit(dummy_quadrature, fe_fine, 0);
           lexicographic_numbering = shape_info.lexicographic_numbering;
         }
 
@@ -1094,9 +1127,8 @@ namespace internal
 
       // ------------- prolongation matrix (0) -> identity matrix --------------
       {
-        AssertDimension(dof_handler_fine.get_fe(0).n_base_elements(), 1);
-        std::string fe_name =
-          dof_handler_fine.get_fe(0).base_element(0).get_name();
+        AssertDimension(fe_fine.n_base_elements(), 1);
+        std::string fe_name = fe_fine.base_element(0).get_name();
         {
           const std::size_t template_starts = fe_name.find_first_of('<');
           Assert(fe_name[template_starts + 1] ==
@@ -1117,9 +1149,8 @@ namespace internal
 
       // ----------------------- prolongation matrix (1) -----------------------
       {
-        AssertDimension(dof_handler_fine.get_fe(0).n_base_elements(), 1);
-        std::string fe_name =
-          dof_handler_fine.get_fe(0).base_element(0).get_name();
+        AssertDimension(fe_fine.n_base_elements(), 1);
+        std::string fe_name = fe_fine.base_element(0).get_name();
         {
           const std::size_t template_starts = fe_name.find_first_of('<');
           Assert(fe_name[template_starts + 1] ==
diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_03.cc b/tests/multigrid-global-coarsening/mg_transfer_a_03.cc
new file mode 100644 (file)
index 0000000..76980bd
--- /dev/null
@@ -0,0 +1,149 @@
+// ---------------------------------------------------------------------
+//
+// 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 working with FECollections with more than one
+ * FE.
+ */
+
+#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_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 "mg_transfer_util.h"
+
+using namespace dealii;
+
+template <int dim, typename Number>
+void
+do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
+{
+  const 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();
+  };
+
+  const auto execute_global_coarsening = [](Triangulation<dim> &tria) {
+    for (auto &cell : tria.active_cell_iterators())
+      cell->set_coarsen_flag();
+    tria.execute_coarsening_and_refinement();
+  };
+
+  const hp::FECollection<dim> fe(fe_fine, fe_coarse);
+
+  // 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);
+  for (const auto &cell : dof_handler_fine.active_cell_iterators())
+    if (cell->is_locally_owned())
+      cell->set_active_fe_index(0);
+  dof_handler_fine.distribute_dofs(fe);
+
+  DoFHandler<dim> dof_handler_coarse(tria_coarse);
+  for (const auto &cell : dof_handler_coarse.active_cell_iterators())
+    if (cell->is_locally_owned())
+      cell->set_active_fe_index(1);
+  dof_handler_coarse.distribute_dofs(fe);
+
+  // 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);
+
+  for (unsigned int i = 1; i < 5; i++)
+    test<2, double>(i);
+}
diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_03.with_p4est=true.mpirun=2.output b/tests/multigrid-global-coarsening/mg_transfer_a_03.with_p4est=true.mpirun=2.output
new file mode 100644 (file)
index 0000000..adadf22
--- /dev/null
@@ -0,0 +1,531 @@
+
+DEAL:0:CG<2>(1)<->CG<2>(1)::weights:
+DEAL:0:CG<2>(1)<->CG<2>(1)::0.500000 1.000000 0.250000 0.500000 0.250000 0.500000 0.500000 1.000000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:0:CG<2>(1)<->CG<2>(1)::0 1 2 3 2 3 11 12 
+DEAL:0:CG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:0:CG<2>(1)<->CG<2>(1)::1 4 3 5 3 5 7 8 
+DEAL:0:CG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:0:CG<2>(1)<->CG<2>(1)::1.000000 0.000000 0.000000 1.000000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::weights:
+DEAL:0:CG<2>(1)<->CG<2>(1)::1.000000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.000000 0.500000 0.250000 0.250000 0.000000 0.500000 0.250000 0.250000 0.250000 0.500000 0.250000 0.250000 0.250000 0.500000 0.250000 0.250000 0.000000 0.500000 0.250000 0.250000 0.000000 1.000000 0.500000 0.500000 0.500000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:0:CG<2>(1)<->CG<2>(1)::4 5 5 0 6 7 7 8 6 7 7 8 9 10 10 2 9 10 10 2 13 14 14 15 13 14 14 15 16 17 17 11 
+DEAL:0:CG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:0:CG<2>(1)<->CG<2>(1)::0 1 2 3 2 3 6 7 
+DEAL:0:CG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:0:CG<2>(1)<->CG<2>(1)::1.000000 0.500000 0.500000 0.000000 0.000000 0.500000 0.500000 1.000000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::4.242641
+DEAL:0:CG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::5.678908
+DEAL:0:CG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(1)<->CG<2>(1)::2.250000 1.750000 3.000000 2.000000 1.000000 1.000000 2.250000 1.750000 1.000000 
+DEAL:0:DG<2>(1)<->CG<2>(1)::weights:
+DEAL:0:DG<2>(1)<->CG<2>(1)::
+DEAL:0:DG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:0:DG<2>(1)<->CG<2>(1)::0 1 2 3 20 21 22 23 
+DEAL:0:DG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:0:DG<2>(1)<->CG<2>(1)::1 4 3 5 3 5 7 8 
+DEAL:0:DG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:0:DG<2>(1)<->CG<2>(1)::1.000000 0.000000 0.000000 1.000000 
+DEAL:0:DG<2>(1)<->CG<2>(1)::weights:
+DEAL:0:DG<2>(1)<->CG<2>(1)::
+DEAL:0:DG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:0:DG<2>(1)<->CG<2>(1)::4 5 8 9 6 7 10 11 12 13 16 17 14 15 18 19 24 25 28 29 26 27 30 31 32 33 36 37 34 35 38 39 
+DEAL:0:DG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:0:DG<2>(1)<->CG<2>(1)::0 1 2 3 2 3 6 7 
+DEAL:0:DG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:0:DG<2>(1)<->CG<2>(1)::1.000000 0.500000 0.500000 0.000000 0.000000 0.500000 0.500000 1.000000 
+DEAL:0:DG<2>(1)<->CG<2>(1)::4.472136
+DEAL:0:DG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(1)<->CG<2>(1)::15.874508
+DEAL:0:DG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(1)<->CG<2>(1)::4.000000 5.000000 8.000000 10.000000 1.000000 2.000000 4.000000 5.000000 1.000000 
+DEAL:0:DG<2>(1)<->DG<2>(1)::weights:
+DEAL:0:DG<2>(1)<->DG<2>(1)::
+DEAL:0:DG<2>(1)<->DG<2>(1)::level_dof_indices_fine:
+DEAL:0:DG<2>(1)<->DG<2>(1)::0 1 2 3 20 21 22 23 
+DEAL:0:DG<2>(1)<->DG<2>(1)::level_dof_indices_coarse:
+DEAL:0:DG<2>(1)<->DG<2>(1)::4 5 6 7 12 13 14 15 
+DEAL:0:DG<2>(1)<->DG<2>(1)::prolongation_matrix_1d:
+DEAL:0:DG<2>(1)<->DG<2>(1)::1.000000 0.000000 0.000000 1.000000 
+DEAL:0:DG<2>(1)<->DG<2>(1)::weights:
+DEAL:0:DG<2>(1)<->DG<2>(1)::
+DEAL:0:DG<2>(1)<->DG<2>(1)::level_dof_indices_fine:
+DEAL:0:DG<2>(1)<->DG<2>(1)::4 5 8 9 6 7 10 11 12 13 16 17 14 15 18 19 24 25 28 29 26 27 30 31 32 33 36 37 34 35 38 39 
+DEAL:0:DG<2>(1)<->DG<2>(1)::level_dof_indices_coarse:
+DEAL:0:DG<2>(1)<->DG<2>(1)::0 1 2 3 8 9 10 11 
+DEAL:0:DG<2>(1)<->DG<2>(1)::prolongation_matrix_1d:
+DEAL:0:DG<2>(1)<->DG<2>(1)::1.000000 0.500000 0.500000 0.000000 0.000000 0.500000 0.500000 1.000000 
+DEAL:0:DG<2>(1)<->DG<2>(1)::4.472136
+DEAL:0:DG<2>(1)<->DG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(1)<->DG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(1)<->DG<2>(1)::11.661904
+DEAL:0:DG<2>(1)<->DG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(1)<->DG<2>(1)::4.000000 4.000000 4.000000 4.000000 1.000000 1.000000 1.000000 1.000000 4.000000 4.000000 4.000000 4.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::weights:
+DEAL:0:CG<2>(2)<->CG<2>(2)::0.500000 1.000000 1.000000 1.000000 1.000000 1.000000 0.250000 0.500000 0.500000 0.250000 0.500000 0.500000 1.000000 1.000000 1.000000 0.500000 1.000000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:0:CG<2>(2)<->CG<2>(2)::0 6 1 4 8 5 2 7 3 2 7 3 34 37 35 32 36 33 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1 12 9 5 14 11 3 13 10 3 13 10 18 24 22 16 23 21 
+DEAL:0:CG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::weights:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.000000 1.000000 0.500000 0.500000 1.000000 0.500000 1.000000 1.000000 0.500000 0.500000 1.000000 0.000000 0.500000 0.500000 0.250000 0.250000 0.500000 0.000000 0.500000 0.500000 0.250000 0.250000 0.500000 0.000000 1.000000 1.000000 0.500000 0.500000 1.000000 0.000000 0.500000 0.500000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.250000 0.250000 0.500000 0.250000 1.000000 1.000000 0.500000 0.500000 1.000000 0.000000 0.500000 0.500000 0.250000 0.250000 0.500000 0.000000 0.500000 0.500000 0.250000 0.250000 0.500000 0.000000 1.000000 1.000000 0.500000 0.500000 1.000000 0.000000 1.000000 1.000000 0.500000 0.500000 1.000000 0.500000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:0:CG<2>(2)<->CG<2>(2)::9 15 10 10 20 0 13 17 14 14 22 19 11 16 12 12 21 18 11 16 12 12 21 18 25 28 26 26 31 29 23 27 24 24 30 2 23 27 24 24 30 2 40 43 41 41 47 45 38 42 39 39 46 44 38 42 39 39 46 44 50 53 51 51 56 54 48 52 49 49 55 32 
+DEAL:0:CG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:CG<2>(2)<->CG<2>(2)::0 6 1 4 8 5 2 7 3 2 7 3 17 20 18 15 19 16 
+DEAL:0:CG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.000000 0.375000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.750000 1.000000 1.000000 0.750000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.375000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::7.549834
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::12.659606
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(2)<->CG<2>(2)::1.562500 1.312500 1.875000 1.375000 3.125000 1.625000 3.125000 3.750000 6.250000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.562500 1.312500 3.125000 1.625000 3.125000 6.250000 1.000000 1.000000 1.000000 1.000000 
+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 45 46 47 48 49 50 51 52 53 
+DEAL:0:DG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->CG<2>(2)::1 9 3 10 5 11 12 13 14 3 10 16 21 18 22 13 23 24 
+DEAL:0:DG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 1.000000 
+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)::9 10 11 18 19 20 12 13 14 21 22 23 15 16 17 24 25 26 27 28 29 36 37 38 30 31 32 39 40 41 33 34 35 42 43 44 54 55 56 63 64 65 57 58 59 66 67 68 60 61 62 69 70 71 72 73 74 81 82 83 75 76 77 84 85 86 78 79 80 87 88 89 
+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 2 3 15 16 17 18 7 19 20 
+DEAL:0:DG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.000000 0.375000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.750000 1.000000 1.000000 0.750000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.375000 1.000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::6.708204
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::24.972328
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(2)<->CG<2>(2)::1.562500 5.375000 3.125000 10.750000 12.250000 5.375000 1.562500 5.937500 1.562500 1.000000 2.000000 1.000000 1.000000 2.000000 1.000000 1.562500 5.375000 12.250000 5.375000 4.375000 1.562500 1.000000 1.000000 1.000000 1.000000 
+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 45 46 47 48 49 50 51 52 53 
+DEAL:0:DG<2>(2)<->DG<2>(2)::level_dof_indices_coarse:
+DEAL:0:DG<2>(2)<->DG<2>(2)::9 10 11 12 13 14 15 16 17 27 28 29 30 31 32 33 34 35 
+DEAL:0:DG<2>(2)<->DG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 1.000000 
+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)::9 10 11 18 19 20 12 13 14 21 22 23 15 16 17 24 25 26 27 28 29 36 37 38 30 31 32 39 40 41 33 34 35 42 43 44 54 55 56 63 64 65 57 58 59 66 67 68 60 61 62 69 70 71 72 73 74 81 82 83 75 76 77 84 85 86 78 79 80 87 88 89 
+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 18 19 20 21 22 23 24 25 26 
+DEAL:0:DG<2>(2)<->DG<2>(2)::prolongation_matrix_1d:
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.000000 0.375000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.750000 1.000000 1.000000 0.750000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.375000 1.000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::6.708204
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::22.153583
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(2)<->DG<2>(2)::1.562500 4.375000 1.562500 4.375000 12.250000 4.375000 1.562500 4.375000 1.562500 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.562500 4.375000 1.562500 4.375000 12.250000 4.375000 1.562500 4.375000 1.562500 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::weights:
+DEAL:0:CG<2>(3)<->CG<2>(3)::0.500000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 0.250000 0.500000 0.500000 0.500000 0.250000 0.500000 0.500000 0.500000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 0.500000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:0:CG<2>(3)<->CG<2>(3)::0 8 9 1 4 12 13 6 5 14 15 7 2 10 11 3 2 10 11 3 65 71 72 67 66 73 74 68 63 69 70 64 
+DEAL:0:CG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:0:CG<2>(3)<->CG<2>(3)::1 20 21 16 6 24 25 18 7 26 27 19 3 22 23 17 3 22 23 17 32 45 46 41 33 47 48 42 29 43 44 40 
+DEAL:0:CG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::weights:
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.500000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.000000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.250000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.250000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.000000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 0.500000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:0:CG<2>(3)<->CG<2>(3)::16 24 25 17 17 35 36 0 20 28 29 22 22 39 40 33 21 30 31 23 23 41 42 34 18 26 27 19 19 37 38 32 18 26 27 19 19 37 38 32 45 51 52 47 47 59 60 55 46 53 54 48 48 61 62 56 43 49 50 44 44 57 58 2 43 49 50 44 44 57 58 2 77 83 84 79 79 92 93 88 78 85 86 80 80 94 95 89 75 81 82 76 76 90 91 87 75 81 82 76 76 90 91 87 98 104 105 100 100 112 113 108 99 106 107 101 101 114 115 109 96 102 103 97 97 110 111 63 
+DEAL:0:CG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:0:CG<2>(3)<->CG<2>(3)::0 8 9 1 4 12 13 6 5 14 15 7 2 10 11 3 2 10 11 3 30 36 37 32 31 38 39 33 28 34 35 29 
+DEAL:0:CG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.000000 0.348607 -0.098607 -0.125000 -0.125000 -0.055902 0.055902 0.000000 0.000000 0.779508 0.934017 0.625000 0.625000 0.220492 -0.184017 0.000000 0.000000 -0.184017 0.220492 0.625000 0.625000 0.934017 0.779508 0.000000 0.000000 0.055902 -0.055902 -0.125000 -0.125000 -0.098607 0.348607 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::10.770330
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::19.134240
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(3)<->CG<2>(3)::1.265625 1.140625 1.406250 1.156250 2.671875 2.671875 1.296875 1.296875 2.671875 2.671875 2.968750 2.968750 5.640625 5.640625 5.640625 5.640625 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.265625 1.140625 2.671875 2.671875 1.296875 1.296875 2.671875 2.671875 5.640625 5.640625 5.640625 5.640625 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->CG<2>(3)::weights:
+DEAL:0:DG<2>(3)<->CG<2>(3)::
+DEAL:0:DG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:0:DG<2>(3)<->CG<2>(3)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 
+DEAL:0:DG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:0:DG<2>(3)<->CG<2>(3)::1 16 3 17 6 7 18 19 20 21 22 23 24 25 26 27 3 17 29 40 32 33 41 42 22 23 43 44 45 46 47 48 
+DEAL:0:DG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:0:DG<2>(3)<->CG<2>(3)::1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:DG<2>(3)<->CG<2>(3)::weights:
+DEAL:0:DG<2>(3)<->CG<2>(3)::
+DEAL:0:DG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:0:DG<2>(3)<->CG<2>(3)::16 17 18 19 32 33 34 35 20 21 22 23 36 37 38 39 24 25 26 27 40 41 42 43 28 29 30 31 44 45 46 47 48 49 50 51 64 65 66 67 52 53 54 55 68 69 70 71 56 57 58 59 72 73 74 75 60 61 62 63 76 77 78 79 96 97 98 99 112 113 114 115 100 101 102 103 116 117 118 119 104 105 106 107 120 121 122 123 108 109 110 111 124 125 126 127 128 129 130 131 144 145 146 147 132 133 134 135 148 149 150 151 136 137 138 139 152 153 154 155 140 141 142 143 156 157 158 159 
+DEAL:0:DG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:0:DG<2>(3)<->CG<2>(3)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 2 3 28 29 30 31 32 33 10 11 34 35 36 37 38 39 
+DEAL:0:DG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:0:DG<2>(3)<->CG<2>(3)::1.000000 0.348607 -0.098607 -0.125000 -0.125000 -0.055902 0.055902 0.000000 0.000000 0.779508 0.934017 0.625000 0.625000 0.220492 -0.184017 0.000000 0.000000 -0.184017 0.220492 0.625000 0.625000 0.934017 0.779508 0.000000 0.000000 0.055902 -0.055902 -0.125000 -0.125000 -0.098607 0.348607 1.000000 
+DEAL:0:DG<2>(3)<->CG<2>(3)::8.944272
+DEAL:0:DG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->CG<2>(3)::32.124757
+DEAL:0:DG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->CG<2>(3)::1.000000 4.000000 4.000000 6.000000 3.000000 9.000000 10.000000 4.000000 3.000000 9.000000 12.000000 12.000000 1.000000 3.000000 3.000000 1.000000 1.000000 2.000000 1.000000 1.000000 1.000000 1.000000 2.000000 2.000000 1.000000 1.000000 1.000000 1.000000 3.000000 2.000000 3.000000 9.000000 10.000000 4.000000 9.000000 3.000000 1.000000 3.000000 3.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->DG<2>(3)::weights:
+DEAL:0:DG<2>(3)<->DG<2>(3)::
+DEAL:0:DG<2>(3)<->DG<2>(3)::level_dof_indices_fine:
+DEAL:0:DG<2>(3)<->DG<2>(3)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 
+DEAL:0:DG<2>(3)<->DG<2>(3)::level_dof_indices_coarse:
+DEAL:0:DG<2>(3)<->DG<2>(3)::16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 
+DEAL:0:DG<2>(3)<->DG<2>(3)::prolongation_matrix_1d:
+DEAL:0:DG<2>(3)<->DG<2>(3)::1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:DG<2>(3)<->DG<2>(3)::weights:
+DEAL:0:DG<2>(3)<->DG<2>(3)::
+DEAL:0:DG<2>(3)<->DG<2>(3)::level_dof_indices_fine:
+DEAL:0:DG<2>(3)<->DG<2>(3)::16 17 18 19 32 33 34 35 20 21 22 23 36 37 38 39 24 25 26 27 40 41 42 43 28 29 30 31 44 45 46 47 48 49 50 51 64 65 66 67 52 53 54 55 68 69 70 71 56 57 58 59 72 73 74 75 60 61 62 63 76 77 78 79 96 97 98 99 112 113 114 115 100 101 102 103 116 117 118 119 104 105 106 107 120 121 122 123 108 109 110 111 124 125 126 127 128 129 130 131 144 145 146 147 132 133 134 135 148 149 150 151 136 137 138 139 152 153 154 155 140 141 142 143 156 157 158 159 
+DEAL:0:DG<2>(3)<->DG<2>(3)::level_dof_indices_coarse:
+DEAL:0:DG<2>(3)<->DG<2>(3)::0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 
+DEAL:0:DG<2>(3)<->DG<2>(3)::prolongation_matrix_1d:
+DEAL:0:DG<2>(3)<->DG<2>(3)::1.000000 0.348607 -0.098607 -0.125000 -0.125000 -0.055902 0.055902 0.000000 0.000000 0.779508 0.934017 0.625000 0.625000 0.220492 -0.184017 0.000000 0.000000 -0.184017 0.220492 0.625000 0.625000 0.934017 0.779508 0.000000 0.000000 0.055902 -0.055902 -0.125000 -0.125000 -0.098607 0.348607 1.000000 
+DEAL:0:DG<2>(3)<->DG<2>(3)::8.944272
+DEAL:0:DG<2>(3)<->DG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->DG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->DG<2>(3)::28.844410
+DEAL:0:DG<2>(3)<->DG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(3)<->DG<2>(3)::1.000000 3.000000 3.000000 1.000000 3.000000 9.000000 9.000000 3.000000 3.000000 9.000000 9.000000 3.000000 1.000000 3.000000 3.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 3.000000 3.000000 1.000000 3.000000 9.000000 9.000000 3.000000 3.000000 9.000000 9.000000 3.000000 1.000000 3.000000 3.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::weights:
+DEAL:0:CG<2>(4)<->CG<2>(4)::0.500000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 0.250000 0.500000 0.500000 0.500000 0.500000 0.250000 0.500000 0.500000 0.500000 0.500000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 0.500000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:0:CG<2>(4)<->CG<2>(4)::0 10 11 12 1 4 16 17 18 7 5 19 20 21 8 6 22 23 24 9 2 13 14 15 3 2 13 14 15 3 106 115 116 117 109 107 118 119 120 110 108 121 122 123 111 104 112 113 114 105 
+DEAL:0:CG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:0:CG<2>(4)<->CG<2>(4)::1 30 31 32 25 7 36 37 38 27 8 39 40 41 28 9 42 43 44 29 3 33 34 35 26 3 33 34 35 26 50 72 73 74 66 51 75 76 77 67 52 78 79 80 68 46 69 70 71 65 
+DEAL:0:CG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::weights:
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.500000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.500000 0.000000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.500000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.500000 0.250000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.500000 0.250000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.500000 0.000000 0.500000 0.500000 0.500000 0.500000 0.250000 0.250000 0.500000 0.500000 0.500000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.000000 1.000000 1.000000 1.000000 1.000000 0.500000 0.500000 1.000000 1.000000 1.000000 0.500000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:0:CG<2>(4)<->CG<2>(4)::25 35 36 37 26 26 54 55 56 0 29 41 42 43 32 32 60 61 62 51 30 44 45 46 33 33 63 64 65 52 31 47 48 49 34 34 66 67 68 53 27 38 39 40 28 28 57 58 59 50 27 38 39 40 28 28 57 58 59 50 71 80 81 82 74 74 95 96 97 89 72 83 84 85 75 75 98 99 100 90 73 86 87 88 76 76 101 102 103 91 69 77 78 79 70 70 92 93 94 2 69 77 78 79 70 70 92 93 94 2 126 135 136 137 129 129 151 152 153 145 127 138 139 140 130 130 154 155 156 146 128 141 142 143 131 131 157 158 159 147 124 132 133 134 125 125 148 149 150 144 124 132 133 134 125 125 148 149 150 144 162 171 172 173 165 165 186 187 188 180 163 174 175 176 166 166 189 190 191 181 164 177 178 179 167 167 192 193 194 182 160 168 169 170 161 161 183 184 185 104 
+DEAL:0:CG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:0:CG<2>(4)<->CG<2>(4)::0 10 11 12 1 4 16 17 18 7 5 19 20 21 8 6 22 23 24 9 2 13 14 15 3 2 13 14 15 3 47 56 57 58 50 48 59 60 61 51 49 62 63 64 52 45 53 54 55 46 
+DEAL:0:CG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.000000 0.338508 -0.117187 -0.070651 0.000000 0.000000 0.049844 0.039062 -0.031987 0.000000 0.000000 0.789852 0.884032 0.282970 0.000000 0.000000 -0.164852 -0.118407 0.092030 0.000000 0.000000 -0.188402 0.312500 0.902688 1.000000 1.000000 0.902688 0.312500 -0.188402 0.000000 0.000000 0.092030 -0.118407 -0.164852 0.000000 0.000000 0.282970 0.884032 0.789852 0.000000 0.000000 -0.031987 0.039062 0.049844 0.000000 0.000000 -0.070651 -0.117187 0.338508 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::13.964240
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::25.823990
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:CG<2>(4)<->CG<2>(4)::1.458272 1.250683 1.708954 1.293776 2.132150 3.687460 2.132150 1.366525 1.633889 1.366525 2.132150 3.687460 2.132150 2.498675 4.321349 2.498675 3.117432 5.391462 3.117432 5.391462 9.324298 5.391462 3.117432 5.391462 3.117432 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.458272 1.250683 2.132150 3.687460 2.132150 1.366525 1.633889 1.366525 2.132150 3.687460 2.132150 3.117432 5.391462 3.117432 5.391462 9.324298 5.391462 3.117432 5.391462 3.117432 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->CG<2>(4)::weights:
+DEAL:0:DG<2>(4)<->CG<2>(4)::
+DEAL:0:DG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:0:DG<2>(4)<->CG<2>(4)::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 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 
+DEAL:0:DG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:0:DG<2>(4)<->CG<2>(4)::1 25 3 26 7 8 9 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 3 26 46 65 50 51 52 66 67 68 33 34 35 69 70 71 72 73 74 75 76 77 78 79 80 
+DEAL:0:DG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:0:DG<2>(4)<->CG<2>(4)::1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:DG<2>(4)<->CG<2>(4)::weights:
+DEAL:0:DG<2>(4)<->CG<2>(4)::
+DEAL:0:DG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:0:DG<2>(4)<->CG<2>(4)::25 26 27 28 29 50 51 52 53 54 30 31 32 33 34 55 56 57 58 59 35 36 37 38 39 60 61 62 63 64 40 41 42 43 44 65 66 67 68 69 45 46 47 48 49 70 71 72 73 74 75 76 77 78 79 100 101 102 103 104 80 81 82 83 84 105 106 107 108 109 85 86 87 88 89 110 111 112 113 114 90 91 92 93 94 115 116 117 118 119 95 96 97 98 99 120 121 122 123 124 150 151 152 153 154 175 176 177 178 179 155 156 157 158 159 180 181 182 183 184 160 161 162 163 164 185 186 187 188 189 165 166 167 168 169 190 191 192 193 194 170 171 172 173 174 195 196 197 198 199 200 201 202 203 204 225 226 227 228 229 205 206 207 208 209 230 231 232 233 234 210 211 212 213 214 235 236 237 238 239 215 216 217 218 219 240 241 242 243 244 220 221 222 223 224 245 246 247 248 249 
+DEAL:0:DG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:0:DG<2>(4)<->CG<2>(4)::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 2 3 45 46 47 48 49 50 51 52 13 14 15 53 54 55 56 57 58 59 60 61 62 63 64 
+DEAL:0:DG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:0:DG<2>(4)<->CG<2>(4)::1.000000 0.338508 -0.117187 -0.070651 0.000000 0.000000 0.049844 0.039062 -0.031987 0.000000 0.000000 0.789852 0.884032 0.282970 0.000000 0.000000 -0.164852 -0.118407 0.092030 0.000000 0.000000 -0.188402 0.312500 0.902688 1.000000 1.000000 0.902688 0.312500 -0.188402 0.000000 0.000000 0.092030 -0.118407 -0.164852 0.000000 0.000000 0.282970 0.884032 0.789852 0.000000 0.000000 -0.031987 0.039062 0.049844 0.000000 0.000000 -0.070651 -0.117187 0.338508 1.000000 
+DEAL:0:DG<2>(4)<->CG<2>(4)::11.180340
+DEAL:0:DG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->CG<2>(4)::40.960362
+DEAL:0:DG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->CG<2>(4)::1.458272 3.132150 6.353321 6.264300 1.458272 2.132150 3.117432 8.157087 4.117432 3.132150 4.895049 7.157087 16.431441 12.052136 12.052136 18.563591 3.117432 7.157087 3.117432 2.132150 1.458272 2.132150 4.895049 2.132150 1.458272 1.000000 2.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 2.000000 2.000000 2.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 4.895049 3.132150 1.458272 2.132150 3.117432 8.157087 4.117432 3.132150 7.157087 4.895049 2.132150 3.117432 7.157087 3.117432 2.132150 1.458272 2.132150 4.895049 2.132150 1.458272 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->DG<2>(4)::weights:
+DEAL:0:DG<2>(4)<->DG<2>(4)::
+DEAL:0:DG<2>(4)<->DG<2>(4)::level_dof_indices_fine:
+DEAL:0:DG<2>(4)<->DG<2>(4)::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 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 
+DEAL:0:DG<2>(4)<->DG<2>(4)::level_dof_indices_coarse:
+DEAL:0:DG<2>(4)<->DG<2>(4)::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 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 
+DEAL:0:DG<2>(4)<->DG<2>(4)::prolongation_matrix_1d:
+DEAL:0:DG<2>(4)<->DG<2>(4)::1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:0:DG<2>(4)<->DG<2>(4)::weights:
+DEAL:0:DG<2>(4)<->DG<2>(4)::
+DEAL:0:DG<2>(4)<->DG<2>(4)::level_dof_indices_fine:
+DEAL:0:DG<2>(4)<->DG<2>(4)::25 26 27 28 29 50 51 52 53 54 30 31 32 33 34 55 56 57 58 59 35 36 37 38 39 60 61 62 63 64 40 41 42 43 44 65 66 67 68 69 45 46 47 48 49 70 71 72 73 74 75 76 77 78 79 100 101 102 103 104 80 81 82 83 84 105 106 107 108 109 85 86 87 88 89 110 111 112 113 114 90 91 92 93 94 115 116 117 118 119 95 96 97 98 99 120 121 122 123 124 150 151 152 153 154 175 176 177 178 179 155 156 157 158 159 180 181 182 183 184 160 161 162 163 164 185 186 187 188 189 165 166 167 168 169 190 191 192 193 194 170 171 172 173 174 195 196 197 198 199 200 201 202 203 204 225 226 227 228 229 205 206 207 208 209 230 231 232 233 234 210 211 212 213 214 235 236 237 238 239 215 216 217 218 219 240 241 242 243 244 220 221 222 223 224 245 246 247 248 249 
+DEAL:0:DG<2>(4)<->DG<2>(4)::level_dof_indices_coarse:
+DEAL:0:DG<2>(4)<->DG<2>(4)::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 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 
+DEAL:0:DG<2>(4)<->DG<2>(4)::prolongation_matrix_1d:
+DEAL:0:DG<2>(4)<->DG<2>(4)::1.000000 0.338508 -0.117187 -0.070651 0.000000 0.000000 0.049844 0.039062 -0.031987 0.000000 0.000000 0.789852 0.884032 0.282970 0.000000 0.000000 -0.164852 -0.118407 0.092030 0.000000 0.000000 -0.188402 0.312500 0.902688 1.000000 1.000000 0.902688 0.312500 -0.188402 0.000000 0.000000 0.092030 -0.118407 -0.164852 0.000000 0.000000 0.282970 0.884032 0.789852 0.000000 0.000000 -0.031987 0.039062 0.049844 0.000000 0.000000 -0.070651 -0.117187 0.338508 1.000000 
+DEAL:0:DG<2>(4)<->DG<2>(4)::11.180340
+DEAL:0:DG<2>(4)<->DG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->DG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->DG<2>(4)::36.864132
+DEAL:0:DG<2>(4)<->DG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:0:DG<2>(4)<->DG<2>(4)::1.458272 2.132150 4.895049 2.132150 1.458272 2.132150 3.117432 7.157087 3.117432 2.132150 4.895049 7.157087 16.431441 7.157087 4.895049 2.132150 3.117432 7.157087 3.117432 2.132150 1.458272 2.132150 4.895049 2.132150 1.458272 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.458272 2.132150 4.895049 2.132150 1.458272 2.132150 3.117432 7.157087 3.117432 2.132150 4.895049 7.157087 16.431441 7.157087 4.895049 2.132150 3.117432 7.157087 3.117432 2.132150 1.458272 2.132150 4.895049 2.132150 1.458272 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+
+DEAL:1:CG<2>(1)<->CG<2>(1)::weights:
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:1:CG<2>(1)<->CG<2>(1)::1.000000 0.000000 0.000000 1.000000 
+DEAL:1:CG<2>(1)<->CG<2>(1)::weights:
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:1:CG<2>(1)<->CG<2>(1)::1.000000 0.500000 0.500000 0.000000 0.000000 0.500000 0.500000 1.000000 
+DEAL:1:CG<2>(1)<->CG<2>(1)::4.242641
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:CG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(1)<->CG<2>(1)::5.678908
+DEAL:1:CG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::weights:
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:1:DG<2>(1)<->CG<2>(1)::1.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(1)<->CG<2>(1)::weights:
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::level_dof_indices_fine:
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::level_dof_indices_coarse:
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::prolongation_matrix_1d:
+DEAL:1:DG<2>(1)<->CG<2>(1)::1.000000 0.500000 0.500000 0.000000 0.000000 0.500000 0.500000 1.000000 
+DEAL:1:DG<2>(1)<->CG<2>(1)::4.472136
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->CG<2>(1)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(1)<->CG<2>(1)::15.874508
+DEAL:1:DG<2>(1)<->CG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(1)<->CG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::weights:
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::level_dof_indices_fine:
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::level_dof_indices_coarse:
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::prolongation_matrix_1d:
+DEAL:1:DG<2>(1)<->DG<2>(1)::1.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(1)<->DG<2>(1)::weights:
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::level_dof_indices_fine:
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::level_dof_indices_coarse:
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::prolongation_matrix_1d:
+DEAL:1:DG<2>(1)<->DG<2>(1)::1.000000 0.500000 0.500000 0.000000 0.000000 0.500000 0.500000 1.000000 
+DEAL:1:DG<2>(1)<->DG<2>(1)::4.472136
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:DG<2>(1)<->DG<2>(1)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(1)<->DG<2>(1)::11.661904
+DEAL:1:DG<2>(1)<->DG<2>(1)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(1)<->DG<2>(1)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::weights:
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:1:CG<2>(2)<->CG<2>(2)::1.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:CG<2>(2)<->CG<2>(2)::weights:
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:1:CG<2>(2)<->CG<2>(2)::1.000000 0.375000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.750000 1.000000 1.000000 0.750000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.375000 1.000000 
+DEAL:1:CG<2>(2)<->CG<2>(2)::7.549834
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:CG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(2)<->CG<2>(2)::12.659606
+DEAL:1:CG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::weights:
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:1:DG<2>(2)<->CG<2>(2)::1.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(2)<->CG<2>(2)::weights:
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::level_dof_indices_fine:
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::level_dof_indices_coarse:
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::prolongation_matrix_1d:
+DEAL:1:DG<2>(2)<->CG<2>(2)::1.000000 0.375000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.750000 1.000000 1.000000 0.750000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.375000 1.000000 
+DEAL:1:DG<2>(2)<->CG<2>(2)::6.708204
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->CG<2>(2)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(2)<->CG<2>(2)::24.972328
+DEAL:1:DG<2>(2)<->CG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(2)<->CG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::weights:
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::level_dof_indices_fine:
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::level_dof_indices_coarse:
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::prolongation_matrix_1d:
+DEAL:1:DG<2>(2)<->DG<2>(2)::1.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(2)<->DG<2>(2)::weights:
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::level_dof_indices_fine:
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::level_dof_indices_coarse:
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::prolongation_matrix_1d:
+DEAL:1:DG<2>(2)<->DG<2>(2)::1.000000 0.375000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.750000 1.000000 1.000000 0.750000 0.000000 0.000000 -0.125000 0.000000 0.000000 0.375000 1.000000 
+DEAL:1:DG<2>(2)<->DG<2>(2)::6.708204
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:DG<2>(2)<->DG<2>(2)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(2)<->DG<2>(2)::22.153583
+DEAL:1:DG<2>(2)<->DG<2>(2)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(2)<->DG<2>(2)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::weights:
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:1:CG<2>(3)<->CG<2>(3)::1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:CG<2>(3)<->CG<2>(3)::weights:
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:1:CG<2>(3)<->CG<2>(3)::1.000000 0.348607 -0.098607 -0.125000 -0.125000 -0.055902 0.055902 0.000000 0.000000 0.779508 0.934017 0.625000 0.625000 0.220492 -0.184017 0.000000 0.000000 -0.184017 0.220492 0.625000 0.625000 0.934017 0.779508 0.000000 0.000000 0.055902 -0.055902 -0.125000 -0.125000 -0.098607 0.348607 1.000000 
+DEAL:1:CG<2>(3)<->CG<2>(3)::10.770330
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:CG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(3)<->CG<2>(3)::19.134240
+DEAL:1:CG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::weights:
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:1:DG<2>(3)<->CG<2>(3)::1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(3)<->CG<2>(3)::weights:
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::level_dof_indices_fine:
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::level_dof_indices_coarse:
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::prolongation_matrix_1d:
+DEAL:1:DG<2>(3)<->CG<2>(3)::1.000000 0.348607 -0.098607 -0.125000 -0.125000 -0.055902 0.055902 0.000000 0.000000 0.779508 0.934017 0.625000 0.625000 0.220492 -0.184017 0.000000 0.000000 -0.184017 0.220492 0.625000 0.625000 0.934017 0.779508 0.000000 0.000000 0.055902 -0.055902 -0.125000 -0.125000 -0.098607 0.348607 1.000000 
+DEAL:1:DG<2>(3)<->CG<2>(3)::8.944272
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->CG<2>(3)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(3)<->CG<2>(3)::32.124757
+DEAL:1:DG<2>(3)<->CG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(3)<->CG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::weights:
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::level_dof_indices_fine:
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::level_dof_indices_coarse:
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::prolongation_matrix_1d:
+DEAL:1:DG<2>(3)<->DG<2>(3)::1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(3)<->DG<2>(3)::weights:
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::level_dof_indices_fine:
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::level_dof_indices_coarse:
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::prolongation_matrix_1d:
+DEAL:1:DG<2>(3)<->DG<2>(3)::1.000000 0.348607 -0.098607 -0.125000 -0.125000 -0.055902 0.055902 0.000000 0.000000 0.779508 0.934017 0.625000 0.625000 0.220492 -0.184017 0.000000 0.000000 -0.184017 0.220492 0.625000 0.625000 0.934017 0.779508 0.000000 0.000000 0.055902 -0.055902 -0.125000 -0.125000 -0.098607 0.348607 1.000000 
+DEAL:1:DG<2>(3)<->DG<2>(3)::8.944272
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:DG<2>(3)<->DG<2>(3)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(3)<->DG<2>(3)::28.844410
+DEAL:1:DG<2>(3)<->DG<2>(3)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(3)<->DG<2>(3)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::weights:
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:1:CG<2>(4)<->CG<2>(4)::1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:CG<2>(4)<->CG<2>(4)::weights:
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:1:CG<2>(4)<->CG<2>(4)::1.000000 0.338508 -0.117187 -0.070651 0.000000 0.000000 0.049844 0.039062 -0.031987 0.000000 0.000000 0.789852 0.884032 0.282970 0.000000 0.000000 -0.164852 -0.118407 0.092030 0.000000 0.000000 -0.188402 0.312500 0.902688 1.000000 1.000000 0.902688 0.312500 -0.188402 0.000000 0.000000 0.092030 -0.118407 -0.164852 0.000000 0.000000 0.282970 0.884032 0.789852 0.000000 0.000000 -0.031987 0.039062 0.049844 0.000000 0.000000 -0.070651 -0.117187 0.338508 1.000000 
+DEAL:1:CG<2>(4)<->CG<2>(4)::13.964240
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:CG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(4)<->CG<2>(4)::25.823990
+DEAL:1:CG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:CG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::weights:
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:1:DG<2>(4)<->CG<2>(4)::1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(4)<->CG<2>(4)::weights:
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::level_dof_indices_fine:
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::level_dof_indices_coarse:
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::prolongation_matrix_1d:
+DEAL:1:DG<2>(4)<->CG<2>(4)::1.000000 0.338508 -0.117187 -0.070651 0.000000 0.000000 0.049844 0.039062 -0.031987 0.000000 0.000000 0.789852 0.884032 0.282970 0.000000 0.000000 -0.164852 -0.118407 0.092030 0.000000 0.000000 -0.188402 0.312500 0.902688 1.000000 1.000000 0.902688 0.312500 -0.188402 0.000000 0.000000 0.092030 -0.118407 -0.164852 0.000000 0.000000 0.282970 0.884032 0.789852 0.000000 0.000000 -0.031987 0.039062 0.049844 0.000000 0.000000 -0.070651 -0.117187 0.338508 1.000000 
+DEAL:1:DG<2>(4)<->CG<2>(4)::11.180340
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->CG<2>(4)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(4)<->CG<2>(4)::40.960362
+DEAL:1:DG<2>(4)<->CG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(4)<->CG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::weights:
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::level_dof_indices_fine:
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::level_dof_indices_coarse:
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::prolongation_matrix_1d:
+DEAL:1:DG<2>(4)<->DG<2>(4)::1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 0.000000 0.000000 0.000000 0.000000 0.000000 1.000000 
+DEAL:1:DG<2>(4)<->DG<2>(4)::weights:
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::level_dof_indices_fine:
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::level_dof_indices_coarse:
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::prolongation_matrix_1d:
+DEAL:1:DG<2>(4)<->DG<2>(4)::1.000000 0.338508 -0.117187 -0.070651 0.000000 0.000000 0.049844 0.039062 -0.031987 0.000000 0.000000 0.789852 0.884032 0.282970 0.000000 0.000000 -0.164852 -0.118407 0.092030 0.000000 0.000000 -0.188402 0.312500 0.902688 1.000000 1.000000 0.902688 0.312500 -0.188402 0.000000 0.000000 0.092030 -0.118407 -0.164852 0.000000 0.000000 0.282970 0.884032 0.789852 0.000000 0.000000 -0.031987 0.039062 0.049844 0.000000 0.000000 -0.070651 -0.117187 0.338508 1.000000 
+DEAL:1:DG<2>(4)<->DG<2>(4)::11.180340
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+DEAL:1:DG<2>(4)<->DG<2>(4)::0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 0.000000 
+DEAL:1:DG<2>(4)<->DG<2>(4)::36.864132
+DEAL:1:DG<2>(4)<->DG<2>(4)::1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 1.000000 
+DEAL:1:DG<2>(4)<->DG<2>(4)::
+

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.