]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Added QSorted quadrature rule.
authorheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 7 May 2014 17:17:14 +0000 (17:17 +0000)
committerheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 7 May 2014 17:17:14 +0000 (17:17 +0000)
git-svn-id: https://svn.dealii.org/trunk@32890 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/doc/users/config.sample
deal.II/include/deal.II/base/quadrature_lib.h
deal.II/source/base/quadrature_lib.cc
tests/base/quadrature_sorted_test.cc [new file with mode: 0644]
tests/base/quadrature_sorted_test.output [new file with mode: 0644]

index a6597cf1d4479df2c067cd5b815d16ce0037d760..e52d5cf9ad4dd3ef53cdb3bb756b69ce26a3efbf 100644 (file)
@@ -148,6 +148,14 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li> New: There is now a QSorted<dim> quadrature which takes an 
+  arbitrary quadrature at construction time and reorders the quadrature
+  points according to the weigths, from smaller to bigger. This should
+  improve stability of higher order polynomial integration.
+  <br>
+  (Luca Heltai, 2014/05/07)
+  </li>  
+
   <li> New: The class VectorizedArray<Number> now provides methods
   VectorizedArray::load(ptr) to read from arbitrary pointer addresses and
   VectorizedArray::store(ptr) to write to arbitrary pointer addresses,
index a64c3dcbb49ea72a6f7b5d9575e5ab1a9228a094..b37132c37046809ef14838d0b721c30c272092ec 100644 (file)
@@ -1,37 +1,9 @@
-## ---------------------------------------------------------------------
-## $Id$
-##
-## Copyright (C) 2013 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 at
-## the top level of the deal.II distribution.
-##
-## ---------------------------------------------------------------------
-
-#
-# Example configuration file
-#
-# This file can be used as a script to preseed the initial cache of a build
-# directory. To do so, invoke cmake with
-#
-#   $ cmake -C configuration.cmake [...]
-#
-#
-# Please note that this file is actually intrepeted as a full featured
-# CMake script file, so in order to set variables in the cache a statement
-# of the form
-#
-#   SET(<variable> <value> CACHE <type> "<string>")
-#
-# has to be used. See doc/readme.html and doc/users/cmake.html for further
-# details on how to use the cmake build system of deal.II.
-#
+##                                                                       ##
+#                      Example configuration file                         #
+#                                                                         #
+#        See doc/readme.html and doc/users/cmake.html for further         #
+#        details on how to use the cmake build system of deal.II.         #
+##                                                                       ##
 
 
 ###########################################################################
 #   "Fortran Compiler."
 #   )
 #
-# SET(DEAL_II_CXX_FLAGS "" CACHE STRING
-#   "The user supplied cache variable will be appended _at the end_ of the auto generated DEAL_II_CXX_FLAGS variable"
+# SET(CMAKE_CXX_FLAGS "" CACHE STRING
+#   "The user supplied cache variable will be appended _at the end_ of the auto generated CMAKE_CXX_FLAGS variable"
 #   )
 #
 # SET(DEAL_II_CXX_FLAGS_DEBUG "" CACHE STRING
 #   "The user supplied cache variable will be appended _at the end_ of the auto generated DEAL_II_LINKER_FLAGS_RELEASE variable"
 #   )
 #
-# SET(DEAL_II_DEFINITIONS "" CACHE STRING
-#   "Additional, user supplied compile definitions"
-#   )
-#
-# SET(DEAL_II_DEFINITIONS_DEBUG "" CACHE STRING
-#   "Additional, user supplied compile definitions"
-#   )
-#
-# SET(DEAL_II_DEFINITIONS_RELEASE "" CACHE STRING
-#   "Additional, user supplied compile definitions"
-#   )
-#
 # SET(BUILD_SHARED_LIBS "ON" CACHE BOOL
 #   "Build a shared library"
 #   )
 #   "If set to ON, then use 64-bit data types to represent global degree of freedom indices. The default is to OFF. You only want to set this to ON if you will solve problems with more than 2^31 (approximately 2 billion) unknowns. If set to ON, you also need to ensure that both Trilinos and/or PETSc support 64-bit indices."
 #   )
 #
+#
 
 
 ###########################################################################
 #
 
 
+#
+# muPaser:
+#
+# SET(DEAL_II_WITH_MUPARSER ON CACHE BOOL
+#   "Build deal.II with support for muparser"
+#
+# Automatic detection:
+#
+# Specify a hint with CMAKE_PREFIX_PATH or by setting
+# SET(MUPARSER_DIR "/.../..." CACHE PATH "")
+#
+# Manual setup:
+#
+# SET(MUPARSER_FOUND TRUE CACHE BOOL "")
+# SET(MUPARSER_LIBRARIES "library;and;semicolon;separated;list;of;link;interface" CACHE STRING "")
+# SET(MUPARSER_INCLUDE_DIRS "semicolon;separated;list;of;include;dirs" CACHE STRING "")
+#
+
+
 #
 # Netcdf:
 #
 # C++11 support is autodetected. You can explicitly disable C+11 support by
 # specifying
 #
-#   SET(DEAL_II_WITH_CXX11 FALSE CACHE BOOL "")
+#   SET(DEAL_II_HAVE_CXX11_FLAG FALSE CACHE BOOL "")
 #
 # A custom C++11 flag can be set by setting
 #
-#   SET(DEAL_II_CXX11_FLAG "-std=c++11" CACHE STRING "")
+#   SET(DEAL_II_HAVE_CXX11_FLAG TRUE CACHE BOOL "")
+#   SET(DEAL_II_CXX11_FLAG "-std=c++0x" CACHE STRING "")
 #
 
 
 #   "Library suffix for the debug library"
 #   )
 #
-# SET(DEAL_II_RELEASE_SUFFIX "" CACHE STRING
+# SET_IF_EMPTY(DEAL_II_RELEASE_SUFFIX "" CACHE STRING
 #   "Library suffix for the release library"
 #   )
 #
index a448a80c9e08d9722f62c62f82603567d6112bae..5bfe461aa0f30197386f9e0c4b9e9af4be0a0e1e 100644 (file)
@@ -435,6 +435,29 @@ private:
 
 
 
+/**
+ * Sorted Quadrature. Given an arbitrary quadrature formula, this
+ * class generates a quadrature formula where the quadrature points
+ * are ordered according the weights, from those with smaller
+ * corresponding weight, to those with higher corresponding weights.
+ * This might be necessary, for example, when integrating high order
+ * polynomials, since in these cases you might sum very big numbers
+ * with very small numbers, and summation is not stable if the numbers
+ * to sum are not close to each other.
+ */
+template <int dim>
+class QSorted : public Quadrature<dim>
+{
+public:
+  /** The constructor takes an arbitrary quadrature formula. */
+  QSorted (const Quadrature<dim>);
+  
+  /** A rule to reorder pairs of points and weights.*/
+  bool operator()(const std::pair<double, Point<dim> > &a,
+                 const std::pair<double, Point<dim> > &b);
+};
+
+
 /*@}*/
 
 /* -------------- declaration of explicit specializations ------------- */
index 7752b9d81eb2e346ef242cfb19c6a67b96492414..8bebce9102752bc098c2bc47f0932e9509c41e98 100644 (file)
@@ -936,6 +936,31 @@ QGaussOneOverR<2>::QGaussOneOverR(const unsigned int n,
 }
 
 
+template <int dim>
+QSorted<dim>::QSorted(Quadrature<dim> quad) :
+  Quadrature<dim>(quad.size())
+{
+  std::vector< std::pair<double, Point<dim> > > wp;
+  for (unsigned int i=0; i<quad.size(); ++i)
+    wp.push_back(std::pair<double, Point<dim> >(quad.weight(i),
+                                               quad.point(i)));
+  sort(wp.begin(), wp.end(), *this);
+  for(unsigned int i=0; i<quad.size(); ++i)
+    {
+      this->weights[i] = wp[i].first;
+      this->quadrature_points[i] = wp[i].second;
+    }
+}
+
+
+template <int dim>
+bool QSorted<dim>::operator()(const std::pair<double, Point<dim> > &a,
+                             const std::pair<double, Point<dim> > &b) 
+{
+  return (a.first < b.first);
+}
+
+
 // construct the quadrature formulae in higher dimensions by
 // tensor product of lower dimensions
 
@@ -1009,4 +1034,8 @@ template class QSimpson<3>;
 template class QMilne<3>;
 template class QWeddle<3>;
 
+template class QSorted<1>;
+template class QSorted<2>;
+template class QSorted<3>;
+
 DEAL_II_NAMESPACE_CLOSE
diff --git a/tests/base/quadrature_sorted_test.cc b/tests/base/quadrature_sorted_test.cc
new file mode 100644 (file)
index 0000000..6393b76
--- /dev/null
@@ -0,0 +1,246 @@
+// ---------------------------------------------------------------------
+// $Id: quadrature_test.cc 31349 2013-10-20 19:07:06Z maier $
+//
+// Copyright (C) 1998 - 2013 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// sort quadratures according to weights, and verify that they
+// integrate accurately polynomial of increasing order.
+// Verify that weights are actually sorted
+
+
+#include "../tests.h"
+#include <iomanip>
+#include <fstream>
+
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/qprojector.h>
+#include <cmath>
+
+template <int dim>
+void
+fill_vector (std::vector<Quadrature<dim> *> &quadratures)
+{
+  quadratures.push_back (new QSorted<dim>(QMidpoint<dim>()));
+  quadratures.push_back (new QSorted<dim>(QTrapez<dim>()));
+  quadratures.push_back (new QSorted<dim>(QSimpson<dim>()));
+  quadratures.push_back (new QSorted<dim>(QMilne<dim>()));
+  quadratures.push_back (new QSorted<dim>(QWeddle<dim>()));
+  for (unsigned int i=0; i<9; ++i)
+    {
+      quadratures.push_back (new QSorted<dim>(QGauss<dim>(i)));
+    }
+  QMilne<1> q1d;
+  quadratures.push_back (new QSorted<dim>(Quadrature<dim>(q1d)));
+  for (unsigned int i=2; i<8; ++i)
+    {
+      quadratures.push_back (new QSorted<dim>(QGaussLobatto<dim>(i)));
+    }
+}
+
+template <int dim>
+void
+check_cells (std::vector<Quadrature<dim>*> &quadratures)
+{
+  Quadrature<dim> quadrature;
+  for (unsigned int n=0; n<quadratures.size(); ++n)
+    {
+      quadrature = *quadratures[n];
+      const std::vector<Point<dim> > &points=quadrature.get_points();
+      const std::vector<double> &weights=quadrature.get_weights();
+
+      deallog << "Quadrature no." << n;
+
+      unsigned int i=0;
+      double quadrature_int=0;
+      double exact_int=0;
+      double err = 0;
+
+      do
+        {
+          ++i;
+
+          quadrature_int=0;
+          // Check the polynomial x^i*y^i
+
+          for (unsigned int x=0; x<quadrature.size(); ++x)
+            {
+              double f=1.;
+              switch (dim)
+                {
+                case 3:
+                  f *= std::pow(static_cast<double>(points[x](2)), i*1.0);
+                case 2:
+                  f *= std::pow(static_cast<double>(points[x](1)), i*1.0);
+                case 1:
+                  f *= std::pow(static_cast<double>(points[x](0)), i*1.0);
+                }
+              quadrature_int+=f*weights[x];
+            }
+
+          // the exact integral is 1/(i+1)
+          exact_int=1./std::pow(static_cast<double>(i+1),dim);
+          err = std::fabs(quadrature_int-exact_int);
+        }
+      while (err<1e-14);
+      // Uncomment here for testing
+//      deallog << " (Int " << quadrature_int << ',' << exact_int << ")";
+      deallog << " is exact for polynomials of degree " << i-1 << std::endl;
+    }
+}
+
+
+template <int dim>
+void
+check_faces (const std::vector<Quadrature<dim-1>*>& quadratures, const bool sub)
+{
+  if (sub)
+    deallog.push("subfaces");
+  else
+    deallog.push("faces");
+
+  for (unsigned int n=0; n<quadratures.size(); ++n)
+    {
+      Quadrature<dim> quadrature (sub == false?
+                                  QProjector<dim>::project_to_all_faces(*quadratures[n]) :
+                                  QProjector<dim>::project_to_all_subfaces(*quadratures[n]));
+      const std::vector<Point<dim> > &points=quadrature.get_points();
+      const std::vector<double> &weights=quadrature.get_weights();
+
+      deallog << "Quadrature no." << n;
+
+      unsigned int i=0;
+      long double quadrature_int=0;
+      double exact_int=0;
+      double err = 0;
+
+      do
+        {
+          ++i;
+
+          quadrature_int=0;
+          // Check the polynomial
+          // x^i*y^i*z^i
+
+          for (unsigned int x=0; x<quadrature.size(); ++x)
+            {
+              long double f=1.;
+              switch (dim)
+                {
+                case 3:
+                  f *= std::pow((long double) points[x](2), i*1.0L);
+                case 2:
+                  f *= std::pow((long double) points[x](1), i*1.0L);
+                case 1:
+                  f *= std::pow((long double) points[x](0), i*1.0L);
+                }
+              quadrature_int+=f*weights[x];
+            }
+
+          // the exact integral is
+          // 1/(i+1)^(dim-1)
+          switch (dim)
+            {
+            case 2:
+              exact_int = 2 * (sub ? 2:1) / (double) (i+1);
+              break;
+            case 3:
+              exact_int = 3 * (sub ? (4+2+2):1)*8 / (double) (i+1)/(i+1);
+              break;
+            }
+
+          err = std::fabs(quadrature_int-exact_int);
+        }
+      // for comparison: use factor 8 in case
+      // of dim==3, as we integrate 8 times
+      // over the whole surface (all
+      // combinations of face_orientation,
+      // face_flip and face_rotation)
+      while (err < (dim==3 ? 8 : 1) * 2e-14);
+      // Uncomment here for testing
+      //      deallog << " (Int " << quadrature_int << '-' << exact_int << '=' << err << ")";
+      deallog << " is exact for polynomials of degree " << i-1 << std::endl;
+    }
+  deallog.pop();
+}
+
+template <int dim>
+void
+check_quadratures (const std::vector<Quadrature<dim>*>& quadratures) {
+  for(unsigned int i=0; i<quadratures.size(); ++i) {
+    Quadrature<dim> &quad = *quadratures[i];
+    bool check = true;
+    for(unsigned int q=1; q<quad.size(); ++q)
+      if(quad.weight(q) < quad.weight(q-1)) {
+       check = false;
+       break;
+      }
+    if(check == false) {
+      deallog << "Something went wrong. The qudrature is not properly sorted." << std::endl;
+      for(unsigned int q=1; q<quad.size(); ++q) {
+       deallog << "q(" << q << "): " << quad.point(q)
+               << ", w(" << q << "): " << quad.weight(q) << std::endl;
+      }
+    }
+  }
+}
+
+int main()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  std::vector<Quadrature<1> *> q1;
+  std::vector<Quadrature<2> *> q2;
+  std::vector<Quadrature<3> *> q3;
+  fill_vector (q1);
+  fill_vector (q2);
+  fill_vector (q3);
+
+  deallog.push("1d");
+  check_cells(q1);
+  check_quadratures(q1);
+  deallog.pop();
+
+  deallog.push("2d");
+  check_cells(q2);
+  check_faces<2>(q1,false);
+  check_faces<2>(q1,true);
+  check_quadratures(q2);
+  deallog.pop();
+
+  deallog.push("3d");
+  check_cells(q3);
+  check_faces<3>(q2,false);
+  check_faces<3>(q2,true);
+  check_quadratures(q3);
+  deallog.pop();
+
+  // delete objects again to avoid
+  // messages about memory leaks when
+  // using purify or other memory
+  // checkers
+  for (unsigned int i=0; i<q1.size(); ++i)
+    delete q1[i];
+  for (unsigned int i=0; i<q2.size(); ++i)
+    delete q2[i];
+  for (unsigned int i=0; i<q3.size(); ++i)
+    delete q3[i];
+}
+
+
diff --git a/tests/base/quadrature_sorted_test.output b/tests/base/quadrature_sorted_test.output
new file mode 100644 (file)
index 0000000..c93bc2d
--- /dev/null
@@ -0,0 +1,148 @@
+
+DEAL:1d::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:1d::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:1d::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:1d::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:1d::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:1d::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:1d::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:1d::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:1d::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:1d::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:1d::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:1d::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:1d::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:1d::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:1d::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:1d::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:1d::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:1d::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:1d::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:1d::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:1d::Quadrature no.20 is exact for polynomials of degree 11
+DEAL:2d::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:2d::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:2d::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:2d::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:2d::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:2d::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:2d::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:2d::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:2d::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:2d::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:2d::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:2d::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:2d::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:2d::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:2d::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:2d::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:2d::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:2d::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:2d::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:2d::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:2d::Quadrature no.20 is exact for polynomials of degree 11
+DEAL:2d:faces::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:2d:faces::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:2d:faces::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:2d:faces::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:2d:faces::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:2d:faces::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:2d:faces::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:2d:faces::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:2d:faces::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:2d:faces::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:2d:faces::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:2d:faces::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:2d:faces::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:2d:faces::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:2d:faces::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:2d:faces::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:2d:faces::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:2d:faces::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:2d:faces::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:2d:faces::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:2d:faces::Quadrature no.20 is exact for polynomials of degree 11
+DEAL:2d:subfaces::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:2d:subfaces::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:2d:subfaces::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:2d:subfaces::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:2d:subfaces::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:2d:subfaces::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:2d:subfaces::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:2d:subfaces::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:2d:subfaces::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:2d:subfaces::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:2d:subfaces::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:2d:subfaces::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:2d:subfaces::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:2d:subfaces::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:2d:subfaces::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:2d:subfaces::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:2d:subfaces::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:2d:subfaces::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:2d:subfaces::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:2d:subfaces::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:2d:subfaces::Quadrature no.20 is exact for polynomials of degree 11
+DEAL:3d::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:3d::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:3d::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:3d::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:3d::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:3d::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:3d::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:3d::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:3d::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:3d::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:3d::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:3d::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:3d::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:3d::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:3d::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:3d::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:3d::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:3d::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:3d::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:3d::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:3d::Quadrature no.20 is exact for polynomials of degree 11
+DEAL:3d:faces::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:3d:faces::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:3d:faces::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:3d:faces::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:3d:faces::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:3d:faces::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:3d:faces::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:3d:faces::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:3d:faces::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:3d:faces::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:3d:faces::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:3d:faces::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:3d:faces::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:3d:faces::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:3d:faces::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:3d:faces::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:3d:faces::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:3d:faces::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:3d:faces::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:3d:faces::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:3d:faces::Quadrature no.20 is exact for polynomials of degree 11
+DEAL:3d:subfaces::Quadrature no.0 is exact for polynomials of degree 1
+DEAL:3d:subfaces::Quadrature no.1 is exact for polynomials of degree 1
+DEAL:3d:subfaces::Quadrature no.2 is exact for polynomials of degree 3
+DEAL:3d:subfaces::Quadrature no.3 is exact for polynomials of degree 5
+DEAL:3d:subfaces::Quadrature no.4 is exact for polynomials of degree 7
+DEAL:3d:subfaces::Quadrature no.5 is exact for polynomials of degree 0
+DEAL:3d:subfaces::Quadrature no.6 is exact for polynomials of degree 1
+DEAL:3d:subfaces::Quadrature no.7 is exact for polynomials of degree 3
+DEAL:3d:subfaces::Quadrature no.8 is exact for polynomials of degree 5
+DEAL:3d:subfaces::Quadrature no.9 is exact for polynomials of degree 7
+DEAL:3d:subfaces::Quadrature no.10 is exact for polynomials of degree 9
+DEAL:3d:subfaces::Quadrature no.11 is exact for polynomials of degree 11
+DEAL:3d:subfaces::Quadrature no.12 is exact for polynomials of degree 13
+DEAL:3d:subfaces::Quadrature no.13 is exact for polynomials of degree 15
+DEAL:3d:subfaces::Quadrature no.14 is exact for polynomials of degree 5
+DEAL:3d:subfaces::Quadrature no.15 is exact for polynomials of degree 1
+DEAL:3d:subfaces::Quadrature no.16 is exact for polynomials of degree 3
+DEAL:3d:subfaces::Quadrature no.17 is exact for polynomials of degree 5
+DEAL:3d:subfaces::Quadrature no.18 is exact for polynomials of degree 7
+DEAL:3d:subfaces::Quadrature no.19 is exact for polynomials of degree 9
+DEAL:3d:subfaces::Quadrature no.20 is exact for polynomials of degree 11

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.