]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
New function DoFTools::extract_dofs_with_support_on_boundary(). New tests. Function...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 7 May 2010 20:29:42 +0000 (20:29 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 7 May 2010 20:29:42 +0000 (20:29 +0000)
git-svn-id: https://svn.dealii.org/trunk@21098 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/dofs/dof_tools.h
deal.II/deal.II/source/dofs/dof_tools.cc
deal.II/doc/news/changes.h
tests/bits/dof_tools_20.cc [new file with mode: 0644]
tests/bits/dof_tools_20/cmp/generic [new file with mode: 0644]

index 4823207c444289bedb3c2c319c5d97d1b8d9e6f5..5aa9a647e36747ed12227a413e85e01c0f6575d3 100644 (file)
@@ -1058,6 +1058,44 @@ class DoFTools
                           const std::vector<bool>    &component_select,
                           std::vector<bool>          &selected_dofs,
                           const std::set<unsigned char> &boundary_indicators = std::set<unsigned char>());
+
+                                    /**
+                                     * This function is similar to
+                                     * the extract_boundary_dofs()
+                                     * function but it extracts those
+                                     * degrees of freedom whose shape
+                                     * functions are nonzero on at
+                                     * least part of the selected
+                                     * boundary. For continuous
+                                     * elements, this is exactly the
+                                     * set of shape functions whose
+                                     * degrees of freedom are defined
+                                     * on boundary faces. On the
+                                     * other hand, if the finite
+                                     * element in used is a
+                                     * discontinuous element, all
+                                     * degrees of freedom are defined
+                                     * in the inside of cells and
+                                     * consequently none would be
+                                     * boundary degrees of
+                                     * freedom. Several of those
+                                     * would have shape functions
+                                     * that are nonzero on the
+                                     * boundary, however. This
+                                     * function therefore extracts
+                                     * all those for which the
+                                     * FiniteElement::has_support_on_face
+                                     * function says that it is
+                                     * nonzero on any face on one of
+                                     * the selected boundary parts.
+                                     */
+    template <class DH>
+    static void
+    extract_dofs_with_support_on_boundary (const DH                   &dof_handler,
+                                          const std::vector<bool>    &component_select,
+                                          std::vector<bool>          &selected_dofs,
+                                          const std::set<unsigned char> &boundary_indicators = std::set<unsigned char>());
+
                                     /**
                                      * @name Hanging Nodes
                                      * @{
index 15fcc3f5eab7916d37c5e3803cfdfdcb456a3d26..1552b8c31b02249783d55a0811aad2231735f802 100644 (file)
@@ -3462,6 +3462,113 @@ DoFTools::extract_boundary_dofs (const DH                      &dof_handler,
 }
 
 
+
+template <class DH>
+void
+DoFTools::extract_dofs_with_support_on_boundary (const DH                      &dof_handler,
+                                                const std::vector<bool>       &component_select,
+                                                std::vector<bool>             &selected_dofs,
+                                                const std::set<unsigned char> &boundary_indicators)
+{
+  Assert (component_select.size() == n_components(dof_handler),
+         ExcWrongSize (component_select.size(),
+                       n_components(dof_handler)));
+  Assert (boundary_indicators.find (255) == boundary_indicators.end(),
+         ExcInvalidBoundaryIndicator());
+
+                                  // let's see whether we have to
+                                  // check for certain boundary
+                                  // indicators or whether we can
+                                  // accept all
+  const bool check_boundary_indicator = (boundary_indicators.size() != 0);
+
+                                   // also see whether we have to
+                                   // check whether a certain vector
+                                   // component is selected, or all
+  const bool check_vector_component
+    = (component_select != std::vector<bool>(component_select.size(),
+                                             true));
+
+                                  // clear and reset array by default
+                                  // values
+  selected_dofs.clear ();
+  selected_dofs.resize (dof_handler.n_dofs(), false);
+  std::vector<unsigned int> cell_dof_indices;
+  cell_dof_indices.reserve (max_dofs_per_cell(dof_handler));
+
+                                  // now loop over all cells and
+                                  // check whether their faces are at
+                                  // the boundary. note that we need
+                                  // not take special care of single
+                                  // lines being at the boundary
+                                  // (using
+                                  // @p{cell->has_boundary_lines}),
+                                  // since we do not support
+                                  // boundaries of dimension dim-2,
+                                  // and so every isolated boundary
+                                  // line is also part of a boundary
+                                  // face which we will be visiting
+                                  // sooner or later
+  for (typename DH::active_cell_iterator cell=dof_handler.begin_active();
+       cell!=dof_handler.end(); ++cell)
+    for (unsigned int face=0;
+        face<GeometryInfo<DH::dimension>::faces_per_cell; ++face)
+      if (cell->at_boundary(face))
+       if (! check_boundary_indicator ||
+           (boundary_indicators.find (cell->face(face)->boundary_indicator())
+            != boundary_indicators.end()))
+         {
+            const FiniteElement<DH::dimension> &fe = cell->get_fe();
+
+            const unsigned int dofs_per_cell = fe.dofs_per_cell;
+            cell_dof_indices.resize (dofs_per_cell);
+           cell->get_dof_indices (cell_dof_indices);
+
+           for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
+             if (fe.has_support_on_face(i,face))
+               {
+                 if (!check_vector_component)
+                   selected_dofs[cell_dof_indices[i]] = true;
+                 else
+                                                    // check for
+                                                    // component is
+                                                    // required. somewhat
+                                                    // tricky as usual
+                                                    // for the case that
+                                                    // the shape function
+                                                    // is non-primitive,
+                                                    // but use usual
+                                                    // convention (see
+                                                    // docs)
+                   {
+                     if (fe.is_primitive (i))
+                       selected_dofs[cell_dof_indices[i]]
+                         = (component_select[fe.system_to_component_index(i).first]
+                            == true);
+                     else // not primitive
+                       {
+                         const unsigned int first_nonzero_comp
+                           = (std::find (fe.get_nonzero_components(i).begin(),
+                                         fe.get_nonzero_components(i).end(),
+                                         true)
+                              -
+                              fe.get_nonzero_components(i).begin());
+                         Assert (first_nonzero_comp < fe.n_components(),
+                                 ExcInternalError());
+
+                         selected_dofs[cell_dof_indices[i]]
+                           = (component_select[first_nonzero_comp]
+                              == true);
+                       }
+                   }
+               }
+         }
+}
+
+
+
+
+
 #else  // 1d
 
 
@@ -3580,9 +3687,130 @@ DoFTools::extract_boundary_dofs (const DH                      &dof_handler,
 }
 
 
+
+template <class DH>
+void
+DoFTools::extract_dofs_with_support_on_boundary (const DH                      &dof_handler,
+                                                const std::vector<bool>       &component_select,
+                                                std::vector<bool>             &selected_dofs,
+                                                const std::set<unsigned char> &boundary_indicators)
+{
+  Assert (component_select.size() == n_components(dof_handler),
+         ExcWrongSize (component_select.size(),
+                       n_components(dof_handler)));
+
+                                  // clear and reset array by default
+                                  // values
+  selected_dofs.clear ();
+  selected_dofs.resize (dof_handler.n_dofs(), false);
+
+                                  // let's see whether we have to
+                                  // check for certain boundary
+                                  // indicators or whether we can
+                                  // accept all
+  const bool check_left_vertex  = ((boundary_indicators.size() == 0) ||
+                                  (boundary_indicators.find(0) !=
+                                   boundary_indicators.end()));
+  const bool check_right_vertex = ((boundary_indicators.size() == 0) ||
+                                  (boundary_indicators.find(1) !=
+                                   boundary_indicators.end()));
+
+                                   // see whether we have to check
+                                   // whether a certain vector
+                                   // component is selected, or all
+  const bool check_vector_component
+    = (component_select != std::vector<bool>(component_select.size(),
+                                             true));
+
+  std::vector<unsigned int> dof_indices;
+
+                                  // loop over cells
+  for (typename DH::active_cell_iterator cell=dof_handler.begin_active(0);
+       cell!=dof_handler.end(0); ++cell)
+    {
+      const FiniteElement<1> &fe = cell->get_fe();
+
+                                      // check left-most vertex
+      if (check_left_vertex)
+       if (cell->neighbor(0) == dof_handler.end())
+          {
+           dof_indices.resize (fe.dofs_per_cell);
+           cell->get_dof_indices(dof_indices);
+            for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
+             if (cell->get_fe().has_support_on_face (i,0) == true)
+               {
+                 if (!check_vector_component)
+                   selected_dofs[dof_indices[i]] = true;
+                 else
+                                                    // check
+                                                    // component. make sure
+                                                    // we don't ask the
+                                                    // wrong question
+                                                    // (leading to an
+                                                    // exception) in case
+                                                    // the shape function
+                                                    // is non-primitive. note
+                                                    // that the face dof
+                                                    // index i is also the
+                                                    // cell dof index of a
+                                                    // corresponding dof in 1d
+                   {
+                     const unsigned int component =
+                       (fe.is_primitive(i) ?
+                        fe.system_to_component_index(i).first :
+                        (std::find (fe.get_nonzero_components(i).begin(),
+                                    fe.get_nonzero_components(i).end(),
+                                    true)
+                         -
+                         fe.get_nonzero_components(i).begin()));
+                     Assert (component < fe.n_components(),
+                             ExcInternalError());
+
+                     if (component_select[component] == true)
+                       selected_dofs[dof_indices[i]] = true;
+                   }
+               }
+          }
+
+                                      // check right-most
+                                      // vertex. same procedure here
+                                      // as above
+      if (check_right_vertex)
+       if (cell->neighbor(1) == dof_handler.end())
+          {
+           dof_indices.resize (fe.dofs_per_cell);
+           cell->get_dof_indices (dof_indices);
+            for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
+             if (cell->get_fe().has_support_on_face (i,1) == true)
+               {
+                 if (!check_vector_component)
+                   selected_dofs[dof_indices[i]] = true;
+                 else
+                   {
+                     const unsigned int component =
+                       (fe.is_primitive(i) ?
+                        fe.face_system_to_component_index(i).first :
+                        (std::find (fe.get_nonzero_components(i).begin(),
+                                    fe.get_nonzero_components(i).end(),
+                                    true)
+                         -
+                         fe.get_nonzero_components(i).begin()));
+                     Assert (component < fe.n_components(),
+                             ExcInternalError());
+
+                     if (component_select[component] == true)
+                       selected_dofs[dof_indices[i]] = true;
+                   }
+               }
+         }
+    }
+}
+
+
 #endif
 
 
+
 namespace internal
 {
   namespace DoFTools
@@ -5643,6 +5871,21 @@ DoFTools::extract_boundary_dofs<hp::DoFHandler<deal_II_dimension> >
  std::vector<bool>                        &,
  const std::set<unsigned char> &);
 
+template
+void
+DoFTools::extract_dofs_with_support_on_boundary<DoFHandler<deal_II_dimension> >
+(const DoFHandler<deal_II_dimension> &,
+ const std::vector<bool>                  &,
+ std::vector<bool>                        &,
+ const std::set<unsigned char> &);
+template
+void
+DoFTools::extract_dofs_with_support_on_boundary<hp::DoFHandler<deal_II_dimension> >
+(const hp::DoFHandler<deal_II_dimension> &,
+ const std::vector<bool>                  &,
+ std::vector<bool>                        &,
+ const std::set<unsigned char> &);
+
 template
 void
 DoFTools::extract_hanging_node_dofs
index f1e62f48280db3492f2a0365fe7a9fa8d8b0ce03..9afaefa2ca45f9b424bf0c1a6dfa684d8285c4a2 100644 (file)
@@ -621,6 +621,12 @@ inconvenience this causes.
 
 <ol>
   <li>
+  <p>New: There is a new function DoFTools::extract_dofs_with_support_on_boundary().
+  fixed.
+  <br>
+  (WB 2010/05/07)
+  </p></li>
+
   <p>Fixed: FE_DGQ::has_support_on_face() returned the wrong value in 1d if the
   polynomial degree of the finite element equals zero (i.e. for piecewise
   constants) where the lone shape function is nonzero on all faces. This is now
diff --git a/tests/bits/dof_tools_20.cc b/tests/bits/dof_tools_20.cc
new file mode 100644 (file)
index 0000000..4055fd8
--- /dev/null
@@ -0,0 +1,62 @@
+//----------------------------  dof_tools_5.cc  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2003, 2004, 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  dof_tools_5.cc  ---------------------------
+
+#include "../tests.h"
+#include "dof_tools_common.cc"
+
+// check
+//   DoFTools::extract_dofs_with_support_on_boundary
+
+
+std::string output_file_name = "dof_tools_20/output";
+
+
+template <int dim>
+void
+check_this (const DoFHandler<dim> &dof_handler)
+{
+  std::vector<bool> component_select (dof_handler.get_fe().n_components(),
+                                      true);
+  std::vector<bool> boundary_dofs (dof_handler.n_dofs());
+
+                                   // first with all components
+  {
+    DoFTools::extract_dofs_with_support_on_boundary (dof_handler,
+                                     component_select,
+                                     boundary_dofs);
+    output_bool_vector (boundary_dofs);
+  }
+
+                                   // next with only every second
+                                   // component
+  for (unsigned int i=1; i<component_select.size(); i+=2)
+    component_select[i] = false;
+  {
+    DoFTools::extract_dofs_with_support_on_boundary (dof_handler,
+                                     component_select,
+                                     boundary_dofs);
+    output_bool_vector (boundary_dofs);
+  }
+
+                                   // third further restrict to
+                                   // boundary indicator 0
+  {
+    std::set<unsigned char> boundary_ids;
+    boundary_ids.insert (0);
+    DoFTools::extract_dofs_with_support_on_boundary (dof_handler,
+                                     component_select,
+                                     boundary_dofs,
+                                     boundary_ids);
+    output_bool_vector (boundary_dofs);
+  }
+}
diff --git a/tests/bits/dof_tools_20/cmp/generic b/tests/bits/dof_tools_20/cmp/generic
new file mode 100644 (file)
index 0000000..04112bb
--- /dev/null
@@ -0,0 +1,257 @@
+
+DEAL::Checking Q1 in 1d:
+DEAL::10001
+DEAL::10001
+DEAL::10000
+DEAL::Checking Q1 in 2d:
+DEAL::101111111010010110
+DEAL::101111111010010110
+DEAL::101000101000000000
+DEAL::Checking Q1 in 3d:
+DEAL::111110111111111111111111111111101110110011000111011111011010
+DEAL::111110111111111111111111111111101110110011000111011111011010
+DEAL::101010100000101010000100101010100000010001000000000000000000
+DEAL::Checking Q2 in 1d:
+DEAL::100000010
+DEAL::100000010
+DEAL::100000000
+DEAL::Checking Q2 in 2d:
+DEAL::101110010111010111010100100100010000001001001111000000100
+DEAL::101110010111010111010100100100010000001001001111000000100
+DEAL::101010000000000101010000000000010000000000000000000000000
+DEAL::Checking Q2 in 3d:
+DEAL::111110111111100110111001100111111110111101100111111101111111101010010111111111111001011111111011011111110111111101111101011101010100111011101010010100110111100101001000110000010011010101110101000001010010000100101000000000001110111010100101001111111110111101001011000001001111011010010010100100011101111000000000001011000
+DEAL::111110111111100110111001100111111110111101100111111101111111101010010111111111111001011111111011011111110111111101111101011101010100111011101010010100110111100101001000110000010011010101110101000001010010000100101000000000001110111010100101001111111110111101001011000001001111011010010010100100011101111000000000001011000
+DEAL::101010101000100010101000000000000000000000000101010100100010101000000000000000000000010100101000000000000101010101000100010101000000000000000000000000010100100101000000000000000001010001010100000000000000000100101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking Q3 in 1d:
+DEAL::1000000000100
+DEAL::1000000000100
+DEAL::1000000000000
+DEAL::Checking Q3 in 2d:
+DEAL::10111100001100001111001100001110110011000000100011000000011000000000000000010001100000011111100000000000000011000000
+DEAL::10111100001100001111001100001110110011000000100011000000011000000000000000010001100000011111100000000000000011000000
+DEAL::10101100000000000000000000001010110000000000000000000000011000000000000000000000000000000000000000000000000000000000
+DEAL::Checking Q3 in 3d:
+DEAL::11111011111111111100001111001111111100000000111111110000000000001111111111110011111111110000111111110000000000001111111100111111111111111100111100001111000000001111000000001111111111111111111111111000000001111000000001111111111111111000011111111000000001111111111111111111000000001111111011111111110011001111110011110000111100001111000000000000111011111100110011000000111100001111000000000000110111111110000110011110000000011110000000000000111100000000000000111100000000000011011001100111111001111000011110000000000000000000110011000000111100000000000000000110000110011110000000000000000000000000000000000000000000000111011111100110011000000111100001111000000000000111111111111110011111111111100001111000000000000101111000000000000001111000000000000111111100111111000011110000000000001000110011000000111100000000000000001111100111111111111000000000000000000000000000000000000000000011001111110000000000000000
+DEAL::11111011111111111100001111001111111100000000111111110000000000001111111111110011111111110000111111110000000000001111111100111111111111111100111100001111000000001111000000001111111111111111111111111000000001111000000001111111111111111000011111111000000001111111111111111111000000001111111011111111110011001111110011110000111100001111000000000000111011111100110011000000111100001111000000000000110111111110000110011110000000011110000000000000111100000000000000111100000000000011011001100111111001111000011110000000000000000000110011000000111100000000000000000110000110011110000000000000000000000000000000000000000000000111011111100110011000000111100001111000000000000111111111111110011111111111100001111000000000000101111000000000000001111000000000000111111100111111000011110000000000001000110011000000111100000000000000001111100111111111111000000000000000000000000000000000000000000011001111110000000000000000
+DEAL::10101010110000001100000011001100111100000000000000000000000000000000000000000000000000000000000000000000000000001010101100001100000011001100111100000000000000000000000000000000000000000000000000000000000000000000000001011000011001111000000000000000000000000000000000000000000000001010101011000000110000001100110011110000000000000000000000000000000000000000000000000000000000000000000000000000010110000110000110011110000000000000000000000000000000000000000000000000000000000001011000000110011001111000000000000000000000000000000000000000000000000000000000000110000110011110000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking DGQ0 in 1d:
+DEAL::0000
+DEAL::0000
+DEAL::0000
+DEAL::Checking DGQ0 in 2d:
+DEAL::0000000000
+DEAL::0000000000
+DEAL::0000000000
+DEAL::Checking DGQ0 in 3d:
+DEAL::0000000000000000000000
+DEAL::0000000000000000000000
+DEAL::0000000000000000000000
+DEAL::Checking DGQ1 in 1d:
+DEAL::00000000
+DEAL::00000000
+DEAL::00000000
+DEAL::Checking DGQ1 in 2d:
+DEAL::0000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000
+DEAL::Checking DGQ1 in 3d:
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking DGQ3 in 1d:
+DEAL::0000000000000000
+DEAL::0000000000000000
+DEAL::0000000000000000
+DEAL::Checking DGQ3 in 2d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking DGQ3 in 3d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking DGP0 in 1d:
+DEAL::0000
+DEAL::0000
+DEAL::0000
+DEAL::Checking DGP0 in 2d:
+DEAL::0000000000
+DEAL::0000000000
+DEAL::0000000000
+DEAL::Checking DGP0 in 3d:
+DEAL::0000000000000000000000
+DEAL::0000000000000000000000
+DEAL::0000000000000000000000
+DEAL::Checking DGP1 in 1d:
+DEAL::00000000
+DEAL::00000000
+DEAL::00000000
+DEAL::Checking DGP1 in 2d:
+DEAL::000000000000000000000000000000
+DEAL::000000000000000000000000000000
+DEAL::000000000000000000000000000000
+DEAL::Checking DGP1 in 3d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking DGP3 in 1d:
+DEAL::0000000000000000
+DEAL::0000000000000000
+DEAL::0000000000000000
+DEAL::Checking DGP3 in 2d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking DGP3 in 3d:
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking Nedelec1 in 2d:
+DEAL::10011011010010100000101100010
+DEAL::10011011010010100000101100010
+DEAL::10000001000000100000000000000
+DEAL::Checking Nedelec1 in 3d:
+DEAL::1111100110111111011110111111110111111111111111111110101110111010101111001011000101011100101010010000111010101111101111000111010101011011000101
+DEAL::1111100110111111011110111111110111111111111111111110101110111010101111001011000101011100101010010000111010101111101111000111010101011011000101
+DEAL::1000100010100000000010010001010000000010010000100010001010000000001001001000000100010100000010010000000000000000000000000000000000000000000000
+DEAL::Checking RaviartThomas0 in 2d:
+DEAL::10011011010010100000101100010
+DEAL::10011011010010100000101100010
+DEAL::10000001000000100000000000000
+DEAL::Checking RaviartThomas1 in 2d:
+DEAL::11000011000011001100001100110000000011000000110000000000000000001100000011110000000000000011000000
+DEAL::11000011000011001100001100110000000011000000110000000000000000001100000011110000000000000011000000
+DEAL::11000000000000000000001100000000000000000000110000000000000000000000000000000000000000000000000000
+DEAL::Checking RaviartThomas2 in 2d:
+DEAL::111000000111000000000000111000111000000000000111000111000000000000000000111000000000000000111000000000000000000000000000000000000000111000000000000000111111000000000000000000000000000000000111000000000000000
+DEAL::111000000111000000000000111000111000000000000111000111000000000000000000111000000000000000111000000000000000000000000000000000000000111000000000000000111111000000000000000000000000000000000111000000000000000
+DEAL::111000000000000000000000000000000000000000000111000000000000000000000000000000000000000000111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking RaviartThomasNodal0 in 2d:
+DEAL::10011011010010100000101100010
+DEAL::10011011010010100000101100010
+DEAL::10000001000000100000000000000
+DEAL::Checking RaviartThomasNodal1 in 2d:
+DEAL::11000011000011001100001100110000000011000000110000000000000000001100000011110000000000000011000000
+DEAL::11000011000011001100001100110000000011000000110000000000000000001100000011110000000000000011000000
+DEAL::11000000000000000000001100000000000000000000110000000000000000000000000000000000000000000000000000
+DEAL::Checking RaviartThomasNodal2 in 2d:
+DEAL::111000000111000000000000111000111000000000000111000111000000000000000000111000000000000000111000000000000000000000000000000000000000111000000000000000111111000000000000000000000000000000000111000000000000000
+DEAL::111000000111000000000000111000111000000000000111000111000000000000000000111000000000000000111000000000000000000000000000000000000000111000000000000000111111000000000000000000000000000000000111000000000000000
+DEAL::111000000000000000000000000000000000000000000111000000000000000000000000000000000000000000111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking RaviartThomasNodal0 in 3d:
+DEAL::1001101011010100111001101111110101001010100100010101000100100000001010110100010101001001100000100
+DEAL::1001101011010100111001101111110101001010100100010101000100100000001010110100010101001001100000100
+DEAL::1000000000010000000000100000010000000000100000000100000000100000000000000000000000000000000000000
+DEAL::Checking RaviartThomasNodal1 in 3d:
+DEAL::1111000000001111111100000000000000001111000011111111000000000000000011110000111100000000111100000000000011111111000000001111000000000000111100001111111100000000000011111111111100000000000011110000111100001111000000000000000000001111000011110000000000000000111100000000111100000000000000000000000011110000000000000000111100001111000000000000000000000000111100000000000000000000111100000000000000000000000000000000000000000000000000001111000011110000000000000000111111110000111100000000000000000000000011110000000000000000111100001111000000000000000000001111000000000000000000001111111100000000000000000000000000000000000000000000111100000000000000000000
+DEAL::1111000000001111111100000000000000001111000011111111000000000000000011110000111100000000111100000000000011111111000000001111000000000000111100001111111100000000000011111111111100000000000011110000111100001111000000000000000000001111000011110000000000000000111100000000111100000000000000000000000011110000000000000000111100001111000000000000000000000000111100000000000000000000111100000000000000000000000000000000000000000000000000001111000011110000000000000000111111110000111100000000000000000000000011110000000000000000111100001111000000000000000000001111000000000000000000001111111100000000000000000000000000000000000000000000111100000000000000000000
+DEAL::1111000000000000000000000000000000000000000000000000000000000000000011110000000000000000000000000000000000000000000000000000000000000000111100000000000000000000000000000000000000000000000011110000000000000000000000000000000000000000000000000000000000000000111100000000000000000000000000000000000000000000000000000000111100000000000000000000000000000000000000000000000000000000111100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking RaviartThomasNodal2 in 3d:
+DEAL::111111111000000000000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000111111111000000000000000000000000000000000000000000000000000000111111111111111111000000000000000000111111111000000000000000000000000000000000000000000000000000000111111111000000000111111111111111111000000000000000000000000000000000000000000000000000000111111111111111111111111111000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::111111111000000000000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000111111111000000000000000000000000000000000000000000000000000000111111111111111111000000000000000000111111111000000000000000000000000000000000000000000000000000000111111111000000000111111111111111111000000000000000000000000000000000000000000000000000000111111111111111111111111111000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000111111111000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_Q<1>(1)3 in 1d:
+DEAL::111000000000111
+DEAL::101000000000101
+DEAL::101000000000000
+DEAL::Checking FE_DGQ<1>(2)2 in 1d:
+DEAL::000000000000000000000000
+DEAL::000000000000000000000000
+DEAL::000000000000000000000000
+DEAL::Checking FE_DGP<1>(3)1 in 1d:
+DEAL::0000000000000000
+DEAL::0000000000000000
+DEAL::0000000000000000
+DEAL::Checking FE_Q<2>(1)3 in 2d:
+DEAL::111000111111111111111111111000111000000111000111111000
+DEAL::101000101101101101101101101000101000000101000101101000
+DEAL::101000101000000000101000101000000000000000000000000000
+DEAL::Checking FE_DGQ<2>(2)2 in 2d:
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGP<2>(3)1 in 2d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_Q<3>(1)3 in 3d:
+DEAL::111111111111111000111111111111111111111111111111111111111111111111111111111111111111111111111000111111111000111111000000111111000000000111111111000111111111111111000111111000111000
+DEAL::101101101101101000101101101101101101101101101101101101101101101101101101101101101101101101101000101101101000101101000000101101000000000101101101000101101101101101000101101000101000
+DEAL::101000101000101000101000000000000000101000101000101000000000000101000000101000101000101000101000000000000000000101000000000101000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<3>(2)2 in 3d:
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGP<3>(3)1 in 3d:
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_Q<1>(1)3FE_DGQ<1>(2)2 in 1d:
+DEAL::111000000000000000000000000000111000000
+DEAL::101000000000000000000000000000101000000
+DEAL::101000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<1>(2)2FE_DGP<1>(3)1 in 1d:
+DEAL::0000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000
+DEAL::Checking FE_DGP<1>(3)1FE_DGQ<1>(2)2 in 1d:
+DEAL::0000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000
+DEAL::Checking FE_Q<2>(1)3FE_DGQ<2>(2)2 in 2d:
+DEAL::111000111111000000000000000000111111000000000000000000111111111000000000000000000000111000000000000000000000000000000000000000000000000000000000000111000000000000000000000111111000000000000000000000000000000000000000000000000000000000
+DEAL::101000101101000000000000000000101101000000000000000000101101101000000000000000000000101000000000000000000000000000000000000000000000000000000000000101000000000000000000000101101000000000000000000000000000000000000000000000000000000000
+DEAL::101000101000000000000000000000000000000000000000000000101000101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<2>(2)2FE_DGP<2>(3)1 in 2d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGP<2>(3)1FE_DGQ<2>(2)2 in 2d:
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_Q<1>(1)3FE_DGP<1>(3)1FE_Q<1>(1)3 in 1d:
+DEAL::1111110000000000000000000000000000001111110000
+DEAL::1011010000000000000000000000000000001011010000
+DEAL::1011010000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<1>(2)2FE_DGQ<1>(2)2FE_Q<1>(3)3 in 1d:
+DEAL::111000000000000000000000000000000000000000000000000000000000000000111000000000000000000
+DEAL::101000000000000000000000000000000000000000000000000000000000000000101000000000000000000
+DEAL::101000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGP<1>(3)1FE_DGP<1>(3)1FE_Q<1>(2)3 in 1d:
+DEAL::11100000000000000000000000000000000000000000011100000000000
+DEAL::10100000000000000000000000000000000000000000010100000000000
+DEAL::10100000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_Q<2>(1)3FE_DGP<2>(3)1FE_Q<2>(1)3 in 2d:
+DEAL::1111110000001111111111110000000000111111111111000000000011111111111111111100000000000000001111110000000000000000000000000000000000000000001111110000000000000000111111111111000000000000000000000000000000000000
+DEAL::1011010000001011011011010000000000101101101101000000000010110110110110110100000000000000001011010000000000000000000000000000000000000000001011010000000000000000101101101101000000000000000000000000000000000000
+DEAL::1011010000001011010000000000000000000000000000000000000010110100000010110100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<2>(2)2FE_DGQ<2>(2)2FE_Q<2>(3)3 in 2d:
+DEAL::111000111111111111000000000000111111000000000000000000000000000000000000000000000000111111111111000000111111000000000000000000000000000000000000000000000000111111111000111111000000111111000000000000000000000000000000000000000000000000000000111000000000111111000000000000000000000000000000000000000000000000000000000111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111000000000111111000000000000000000000000000000000000000000000000000000111111111111111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111000000000000000000000000000000000000000000000000000000
+DEAL::101000101101110011000000000000110011000000000000000000000000000000000000000000000000101101110011000000110011000000000000000000000000000000000000000000000000101101101000110011000000110011000000000000000000000000000000000000000000000000000000101000000000110011000000000000000000000000000000000000000000000000000000000110011000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000101000000000110011000000000000000000000000000000000000000000000000000000101101110011110011000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000110011000000000000000000000000000000000000000000000000000000
+DEAL::101000101000110011000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000101000101000110011000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000110011000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGP<2>(3)1FE_DGP<2>(3)1FE_Q<2>(2)3 in 2d:
+DEAL::11100011111111100000011100000000000000000000000111111111000111000000000000000000000001111111110001110001110000000000000000000000000011100000011100000000000000000000000000000111000000000000000000000000000000000000000000000000000000000011100000011100000000000000000000000000111111111111000000000000000000000000000000000000000000000000000000000011100000000000000000000000000
+DEAL::10100010110110100000010100000000000000000000000101101101000101000000000000000000000001011011010001010001010000000000000000000000000010100000010100000000000000000000000000000101000000000000000000000000000000000000000000000000000000000010100000010100000000000000000000000000101101101101000000000000000000000000000000000000000000000000000000000010100000000000000000000000000
+DEAL::10100010100010100000000000000000000000000000000000000000000000000000000000000000000001010001010001010000000000000000000000000000000000000000000000000000000000000000000000000101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<3>(1)3FE_DGP<3>(3)1FE_Q<3>(1)3 in 3d:
+DEAL::11111111111111100011111100000000000000000000000000000000000000000000111111111111000000000000000000000000000000000000000000001111111111111111110000000000000000000000000000000000000000000011111111100000000000000000000000000000000000000000000111111000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000011111111111111111111100000000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000000001111110000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001111111110000000000000000000000000000000000000000000000011111111111100000000000000000000000000000000000000000000111000000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000011100000000000000000000000000000000000000000000000111000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::10110110110110100010110100000000000000000000000000000000000000000000101101101101000000000000000000000000000000000000000000001011011011011011010000000000000000000000000000000000000000000010110110100000000000000000000000000000000000000000000101101000000000000000000000000000000000000000000001010000000000000000000000000000000000000000000010110110110110110110100000000000000000000000000000000000000000000000101101101000000000000000000000000000000000000000000000001011010000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000101101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001011011010000000000000000000000000000000000000000000000010110110110100000000000000000000000000000000000000000000101000000000000000000000000000000000000000000000001010000000000000000000000000000000000000000000010100000000000000000000000000000000000000000000000101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::10100010100010100010100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001010001010001010000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000010100010100010100010100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001010000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000101000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking (FESystem<2>(FE_Q<2>(1),3))3FE_DGQ<2>(0)1FE_Q<2>(1)3 in 2d:
+DEAL::1111111111110000000000001111111111111111111111110111111111111111111111111011111111111111111111111111111111111100000000000001111111111110000000000000000000000000001111111111110000000000000111111111111111111111111000000000000000
+DEAL::1010101011010000000000001010101011011010101011010101010101101101010101101010101010110110101010110110101010110100000000000001010101011010000000000000000000000000001010101011010000000000000101010101101101010101101000000000000000
+DEAL::1010101011010000000000001010101011010000000000000000000000000000000000000010101010110100000000000010101010110100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<2>(3)1FESystem<2>(FE_DGQ<2>(0),3)1FESystem<2>(FE_Q<2>(2),1, FE_DGQ<2>(0),1)2 in 2d:
+DEAL::110011111100001100000000000000000000000111111001100000000000000000000000111111001100110000000000000000000000000110000110000000000000000000000000001100000000000000000000000000000000000000000000000000000011000011000000000000000000000000011111111000000000000000000000000000000000000000000000000000000110000000000000000000000000
+DEAL::110011111100001100000000000000000000000111111001100000000000000000000000111111001100110000000000000000000000000110000110000000000000000000000000001100000000000000000000000000000000000000000000000000000011000011000000000000000000000000011111111000000000000000000000000000000000000000000000000000000110000000000000000000000000
+DEAL::110011001100000000000000000000000000000000000000000000000000000000000000110011001100000000000000000000000000000000000000000000000000000000000000001100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_DGQ<2>(3)1FE_Nedelec<2>(1)2 in 2d:
+DEAL::11000011000000000000000011001100000000000000001100110000000000000000000011000000000000000000110000000000000000000000000000000000000000001100000000000000000011110000000000000000000000000000000000000011000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::Checking FE_Nedelec<2>(1)1FESystem<2>(FE_DGQ<2>(1),2)1FESystem<2>(FE_Q<2>(2),1, FE_Nedelec<2>(1),2)2 in 2d:
+DEAL::110011111111111000000000000001111111000000000011111111111000000011111110000000000111111001111111000000011111110000000000000000011000000000111111100000000000000000001111111000000000000000000000000000000000000000000000000110000000001111111000000000000000001111111111111111110000000000000000000000000000000000000000000111111100000000000000000
+DEAL::100010101100011000000000000001100011000000000010101100011000000011000110000000000101010001100011000000011000110000000000000000010000000000110001100000000000000000001100011000000000000000000000000000000000000000000000000100000000001100011000000000000000001010110001111000110000000000000000000000000000000000000000000110001100000000000000000
+DEAL::100010001100011000000000000000000000000000000000000000000000000000000000000000000100010001100011000000000000000000000000000000000000000000000000000000000000000000001100011000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000

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.