]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Update PolynomialsAdini to derive from ScalarPolynomialsBase and add a test for it 9447/head
authorgrahambenharper <grahambenharper@gmail.com>
Mon, 27 Jan 2020 23:59:57 +0000 (16:59 -0700)
committergrahambenharper <grahambenharper@gmail.com>
Wed, 29 Jan 2020 19:52:53 +0000 (12:52 -0700)
include/deal.II/base/polynomials_adini.h
include/deal.II/base/scalar_polynomials_base.h
source/base/polynomials_adini.cc
tests/base/polynomials_adini_01.cc [new file with mode: 0644]
tests/base/polynomials_adini_01.output [new file with mode: 0644]

index 535be7ee157fffec2e5e1f0d3c37e252bfbb96cd..4b9c0b6a06cdf2c37962210ab942a6bc6d5bff0e 100644 (file)
@@ -20,6 +20,7 @@
 #include <deal.II/base/config.h>
 
 #include <deal.II/base/point.h>
+#include <deal.II/base/scalar_polynomials_base.h>
 #include <deal.II/base/table.h>
 #include <deal.II/base/tensor.h>
 
@@ -34,18 +35,22 @@ DEAL_II_NAMESPACE_OPEN
  * The basis of the space is chosen to match the node functionals of the Adini
  * element.
  *
- * @todo This polynomial space is implemented in 2D only.
+ * @todo This polynomial space is implemented in 2D only and does not compute
+ * derivatives of order 3 or higher.
  *
- * @author Bärbel Janssen, 2007
+ * @ingroup Polynomials
+ * @author Bärbel Janssen
+ * @date 2007
  */
-
-class PolynomialsAdini
+template <int dim>
+class PolynomialsAdini : public ScalarPolynomialsBase<dim>
 {
 public:
   /**
    * Constructor for the polynomials of the described space
    */
   PolynomialsAdini();
+
   /**
    * Compute the value and the first and second derivatives of each
    * polynomial at <tt>unit_point</tt>.
@@ -59,21 +64,21 @@ public:
    * function, rather than using any of the compute_value(), compute_grad() or
    * compute_grad_grad() functions, see below, in a loop over all polynomials.
    */
-
   void
-  evaluate(const Point<2> &           unit_point,
-           std::vector<double> &      values,
-           std::vector<Tensor<1, 2>> &grads,
-           std::vector<Tensor<2, 2>> &grad_grads) const;
+  evaluate(const Point<dim> &           unit_point,
+           std::vector<double> &        values,
+           std::vector<Tensor<1, dim>> &grads,
+           std::vector<Tensor<2, dim>> &grad_grads,
+           std::vector<Tensor<3, dim>> &third_derivatives,
+           std::vector<Tensor<4, dim>> &fourth_derivatives) const override;
 
   /**
    * Compute the value of the <tt>i</tt>th polynomial at <tt>unit_point</tt>.
    *
    * Consider using evaluate() instead.
    */
-
   double
-  compute_value(const unsigned int i, const Point<2> &p) const;
+  compute_value(const unsigned int i, const Point<dim> &p) const;
 
   /**
    * Compute the gradient of the <tt>i</tt>th polynomial at
@@ -81,20 +86,29 @@ public:
    *
    * Consider using evaluate() instead.
    */
+  Tensor<1, dim>
+  compute_grad(const unsigned int i, const Point<dim> &p) const;
 
-  Tensor<1, 2>
-  compute_grad(const unsigned int i, const Point<2> &p) const;
   /**
    * Compute the second derivative (grad_grad) of the <tt>i</tt>th polynomial
    * at <tt>unit_point</tt>.
    *
    * Consider using evaluate() instead.
    */
+  Tensor<2, dim>
+  compute_grad_grad(const unsigned int i, const Point<dim> &p) const;
+
+  /**
+   * Return the name of the space, which is <tt>PolynomialsAdini</tt>.
+   */
+  std::string
+  name() const override;
 
-  Tensor<2, 2>
-  compute_grad_grad(const unsigned int i, const Point<2> &p) const;
-  Tensor<2, 2>
-  compute_grad_grad_2(const unsigned int i, const Point<2> &p) const;
+  /**
+   * @copydoc ScalarPolynomialsBase<dim>::clone()
+   */
+  virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
+  clone() const override;
 
 private:
   /**
@@ -107,24 +121,26 @@ private:
    * Store the coefficients of the x-derivative of the polynomials in the
    * order $1,x,y,x^2,y^2,xy,x^3,y^3,xy^2,x^2y,x^3y,xy^3$
    */
-
   Table<2, double> dx;
+
   /**
    * Store the coefficients of the y-derivative of the polynomials in the
    * order $1,x,y,x^2,y^2,xy,x^3,y^3,xy^2,x^2y,x^3y,xy^3$
    */
-
   Table<2, double> dy;
+
   /**
    * Store the coefficients of the second x-derivative of the polynomials in
    * the order $1,x,y,x^2,y^2,xy,x^3,y^3,xy^2,x^2y,x^3y,xy^3$
    */
   Table<2, double> dxx;
+
   /**
    * Store the coefficients of the second y-derivative of the polynomials in
    * the order $1,x,y,x^2,y^2,xy,x^3,y^3,xy^2,x^2y,x^3y,xy^3$
    */
   Table<2, double> dyy;
+
   /**
    * Store the coefficients of the second mixed derivative of the polynomials
    * in the order $1,x,y,x^2,y^2,xy,x^3,y^3,xy^2,x^2y,x^3y,xy^3$
@@ -134,6 +150,15 @@ private:
 
 
 
+template <int dim>
+inline std::string
+PolynomialsAdini<dim>::name() const
+{
+  return "PolynomialsAdini";
+}
+
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index d65ef3451e53534e352bea6a71cc6c4e8b8b61cd..5b192a885cbc3ec80814c00db17cb145038baa45 100644 (file)
@@ -46,6 +46,7 @@ DEAL_II_NAMESPACE_OPEN
  *
  * Some classes that derive from this class include
  * <ul>
+ *   <li> <tt>PolynomialsAdini</tt>
  *   <li> <tt>PolynomialsRannacherTurek</tt>
  *   <li> <tt>PolynomialsP</tt>
  *   <li> <tt>PolynomialSpace</tt>
index d3eefc64341b3defa6c877063119beb2c22b61d2..05714c57523afeb189dd6f7d58501616089fb7dd 100644 (file)
@@ -13,8 +13,9 @@
 //
 // ---------------------------------------------------------------------
 
-
+#include <deal.II/base/exceptions.h>
 #include <deal.II/base/polynomials_adini.h>
+#include <deal.II/base/std_cxx14/memory.h>
 
 #define ENTER_COEFFICIENTS(                                   \
   koefs, z, a0, a1, a2, a3, a4, a5, a6, a7, a8, a9, a10, a11) \
 DEAL_II_NAMESPACE_OPEN
 
 
-PolynomialsAdini::PolynomialsAdini()
-  : coef(12, 12)
+
+template <int dim>
+PolynomialsAdini<dim>::PolynomialsAdini()
+  : ScalarPolynomialsBase<dim>(3, 12)
+  , coef(12, 12)
   , dx(12, 12)
   , dy(12, 12)
   , dxx(12, 12)
   , dyy(12, 12)
   , dxy(12, 12)
 {
-  //                       1  x  y  xx yy xy 3x 3y xyy xxy 3xy x3y
-  //                       0  1  2  3  4  5  6  7  8  9 10 11
+  Assert(dim == 2, ExcNotImplemented());
+
+  //                          1  x  y  xx yy xy 3x 3y xyy xxy 3xy x3y
+  //                          0  1  2  3  4  5  6  7  8   9   10  11
   ENTER_COEFFICIENTS(coef, 0, 1, 0, 0, -3, -3, -1, 2, 2, 3, 3, -2, -2);
   ENTER_COEFFICIENTS(coef, 1, 0, 1, 0, -2, 0, -1, 1, 0, 0, 2, -1, 0);
   ENTER_COEFFICIENTS(coef, 2, 0, 0, 1, 0, -2, -1, 0, 1, 2, 0, 0, -1);
@@ -84,7 +90,6 @@ PolynomialsAdini::PolynomialsAdini()
   ENTER_COEFFICIENTS(dy, 10, 0, 0, 0, -1, 0, 0, 1, 0, 0, 0, 0, 0);
   ENTER_COEFFICIENTS(dy, 11, 0, 0, 0, 0, 0, -2, 0, 0, 3, 0, 0, 0);
 
-  //                       0  1  2  3  4  5  6  7  8  9 10 11
   ENTER_COEFFICIENTS(dxx, 0, -6, 12, 6, 0, 0, -12, 0, 0, 0, 0, 0, 0);
   ENTER_COEFFICIENTS(dxx, 1, -4, 6, 4, 0, 0, -6, 0, 0, 0, 0, 0, 0);
   ENTER_COEFFICIENTS(dxx, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);
@@ -125,12 +130,30 @@ PolynomialsAdini::PolynomialsAdini()
   ENTER_COEFFICIENTS(dxy, 11, 0, 0, -2, 0, 3, 0, 0, 0, 0, 0, 0, 0);
 }
 
+
+
+template <int dim>
 void
-PolynomialsAdini::evaluate(const Point<2> &           unit_point,
-                           std::vector<double> &      values,
-                           std::vector<Tensor<1, 2>> &grads,
-                           std::vector<Tensor<2, 2>> &grad_grads) const
+PolynomialsAdini<dim>::evaluate(
+  const Point<dim> &           unit_point,
+  std::vector<double> &        values,
+  std::vector<Tensor<1, dim>> &grads,
+  std::vector<Tensor<2, dim>> &grad_grads,
+  std::vector<Tensor<3, dim>> &third_derivatives,
+  std::vector<Tensor<4, dim>> &fourth_derivatives) const
 {
+  const unsigned int n_pols = this->n();
+  Assert(values.size() == n_pols || values.size() == 0,
+         ExcDimensionMismatch(values.size(), n_pols));
+  Assert(grads.size() == n_pols || grads.size() == 0,
+         ExcDimensionMismatch(grads.size(), n_pols));
+  Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
+         ExcDimensionMismatch(grad_grads.size(), n_pols));
+  Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
+         ExcDimensionMismatch(third_derivatives.size(), n_pols));
+  Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
+         ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+
   if (values.empty() == false) // do not bother if empty
     {
       for (unsigned int i = 0; i < values.size(); ++i)
@@ -154,11 +177,16 @@ PolynomialsAdini::evaluate(const Point<2> &           unit_point,
           grad_grads[i] = compute_grad_grad(i, unit_point);
         }
     }
+
   return;
 }
 
+
+
+template <int dim>
 double
-PolynomialsAdini::compute_value(const unsigned int i, const Point<2> &p) const
+PolynomialsAdini<dim>::compute_value(const unsigned int i,
+                                     const Point<dim> & p) const
 {
   const double x = p(0);
   const double y = p(1);
@@ -169,12 +197,16 @@ PolynomialsAdini::compute_value(const unsigned int i, const Point<2> &p) const
          coef(11, i) * x * y * y * y;
 }
 
-Tensor<1, 2>
-PolynomialsAdini::compute_grad(const unsigned int i, const Point<2> &p) const
+
+
+template <int dim>
+Tensor<1, dim>
+PolynomialsAdini<dim>::compute_grad(const unsigned int i,
+                                    const Point<dim> & p) const
 {
-  const double x = p(0);
-  const double y = p(1);
-  Tensor<1, 2> tensor;
+  const double   x = p(0);
+  const double   y = p(1);
+  Tensor<1, dim> tensor;
   tensor[0] = dx(0, i) + dx(1, i) * x + dx(2, i) * y + dx(3, i) * x * x +
               dx(4, i) * y * y + dx(5, i) * x * y + dx(6, i) * x * x * x +
               dx(7, i) * y * y * y + dx(8, i) * x * y * y +
@@ -189,13 +221,16 @@ PolynomialsAdini::compute_grad(const unsigned int i, const Point<2> &p) const
   return tensor;
 }
 
-Tensor<2, 2>
-PolynomialsAdini::compute_grad_grad(const unsigned int i,
-                                    const Point<2> &   p) const
+
+
+template <int dim>
+Tensor<2, dim>
+PolynomialsAdini<dim>::compute_grad_grad(const unsigned int i,
+                                         const Point<dim> & p) const
 {
-  const double x = p(0);
-  const double y = p(1);
-  Tensor<2, 2> tensor;
+  const double   x = p(0);
+  const double   y = p(1);
+  Tensor<2, dim> tensor;
   tensor[0][0] = dxx(0, i) + dxx(1, i) * x + dxx(2, i) * y + dxx(3, i) * x * x +
                  dxx(4, i) * y * y + dxx(5, i) * x * y + dxx(6, i) * x * x * x +
                  dxx(7, i) * y * y * y + dxx(8, i) * x * y * y +
@@ -216,4 +251,19 @@ PolynomialsAdini::compute_grad_grad(const unsigned int i,
 }
 
 
+
+template <int dim>
+std::unique_ptr<ScalarPolynomialsBase<dim>>
+PolynomialsAdini<dim>::clone() const
+{
+  return std_cxx14::make_unique<PolynomialsAdini<dim>>(*this);
+}
+
+
+
+template class PolynomialsAdini<0>;
+template class PolynomialsAdini<1>;
+template class PolynomialsAdini<2>;
+template class PolynomialsAdini<3>;
+
 DEAL_II_NAMESPACE_CLOSE
diff --git a/tests/base/polynomials_adini_01.cc b/tests/base/polynomials_adini_01.cc
new file mode 100644 (file)
index 0000000..3829338
--- /dev/null
@@ -0,0 +1,94 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// evaluate polynomials_adini on the reference cell using a 4x4 grid
+
+#include <deal.II/base/geometry_info.h>
+#include <deal.II/base/polynomials_adini.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/tensor.h>
+
+#include <vector>
+
+#include "../tests.h"
+
+template <int dim>
+void
+check_poly_q(const PolynomialsAdini<dim> &poly)
+{
+  std::vector<Point<dim>>     points;
+  std::vector<double>         values;
+  std::vector<Tensor<1, dim>> grads;
+  std::vector<Tensor<2, dim>> grads2;
+  std::vector<Tensor<3, dim>> thirds;
+  std::vector<Tensor<4, dim>> fourths;
+
+  // set up evaluation points - 4x4 points for cubic
+  for (unsigned int i = 0; i < 4; ++i)
+    {
+      for (unsigned int j = 0; j < 4; ++j)
+        {
+          Point<dim> p;
+          p[0] = 1. / 3. * j;
+          p[1] = 1. / 3. * i;
+          points.push_back(p);
+        }
+    }
+
+  // loop over evaluation points
+  for (unsigned int i = 0; i < points.size(); ++i)
+    {
+      values.clear();
+      grads.clear();
+      grads2.clear();
+      thirds.clear();
+      fourths.clear();
+      values.resize(12);
+      grads.resize(12);
+
+      deallog << "Adini<" << dim << "> point " << i << " (" << points[i][0];
+      for (unsigned int d = 1; d < dim; ++d)
+        deallog << ", " << points[i][d];
+      deallog << ")" << std::endl;
+
+      poly.evaluate(points[i], values, grads, grads2, thirds, fourths);
+
+      // loop through shape fxns
+      for (unsigned int j = 0; j < 12; ++j)
+        {
+          deallog << "Adini<" << dim << "> shape fxn " << j << ": ";
+          deallog << '\t' << values[j];
+          deallog << std::endl;
+
+          deallog << "Adini<" << dim << "> shape grad " << j << ": ";
+          for (unsigned int d = 0; d < dim; ++d)
+            deallog << '\t' << grads[j][d];
+          deallog << std::endl;
+        }
+    }
+}
+
+
+
+int
+main()
+{
+  initlog();
+  deallog << std::setprecision(3);
+
+  PolynomialsAdini<2> p_2d;
+  check_poly_q(p_2d);
+}
diff --git a/tests/base/polynomials_adini_01.output b/tests/base/polynomials_adini_01.output
new file mode 100644 (file)
index 0000000..30e9336
--- /dev/null
@@ -0,0 +1,401 @@
+
+DEAL::Adini<2> point 0 (0.00, 0.00)
+DEAL::Adini<2> shape fxn 0:    1.00
+DEAL::Adini<2> shape grad 0:   0.00    0.00
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   1.00    0.00
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   0.00    1.00
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   0.00    0.00
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.00
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.00
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   0.00    0.00
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.00
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.00
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   0.00    0.00
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    0.00
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.00
+DEAL::Adini<2> point 1 (0.333, 0.00)
+DEAL::Adini<2> shape fxn 0:    0.741
+DEAL::Adini<2> shape grad 0:   -1.33   -0.0741
+DEAL::Adini<2> shape fxn 1:    0.148
+DEAL::Adini<2> shape grad 1:   5.55e-17        -0.148
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   0.00    0.667
+DEAL::Adini<2> shape fxn 3:    0.259
+DEAL::Adini<2> shape grad 3:   1.33    0.0741
+DEAL::Adini<2> shape fxn 4:    -0.0741
+DEAL::Adini<2> shape grad 4:   -0.333  0.0741
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.333
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   0.00    0.0741
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.148
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.00
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   0.00    -0.0741
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    -0.0741
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.00
+DEAL::Adini<2> point 2 (0.667, 0.00)
+DEAL::Adini<2> shape fxn 0:    0.259
+DEAL::Adini<2> shape grad 0:   -1.33   0.0741
+DEAL::Adini<2> shape fxn 1:    0.0741
+DEAL::Adini<2> shape grad 1:   -0.333  -0.0741
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   0.00    0.333
+DEAL::Adini<2> shape fxn 3:    0.741
+DEAL::Adini<2> shape grad 3:   1.33    -0.0741
+DEAL::Adini<2> shape fxn 4:    -0.148
+DEAL::Adini<2> shape grad 4:   0.00    0.148
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.667
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   0.00    -0.0741
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.0741
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.00
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   0.00    0.0741
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    -0.148
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.00
+DEAL::Adini<2> point 3 (1.00, 0.00)
+DEAL::Adini<2> shape fxn 0:    0.00
+DEAL::Adini<2> shape grad 0:   0.00    0.00
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   0.00    0.00
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   0.00    0.00
+DEAL::Adini<2> shape fxn 3:    1.00
+DEAL::Adini<2> shape grad 3:   0.00    0.00
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   1.00    0.00
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    1.00
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   0.00    0.00
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.00
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.00
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   0.00    0.00
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    0.00
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.00
+DEAL::Adini<2> point 4 (0.00, 0.333)
+DEAL::Adini<2> shape fxn 0:    0.741
+DEAL::Adini<2> shape grad 0:   -0.0741 -1.33
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   0.667   0.00
+DEAL::Adini<2> shape fxn 2:    0.148
+DEAL::Adini<2> shape grad 2:   -0.148  5.55e-17
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   0.0741  0.00
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.00
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.148   0.00
+DEAL::Adini<2> shape fxn 6:    0.259
+DEAL::Adini<2> shape grad 6:   0.0741  1.33
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.333   0.00
+DEAL::Adini<2> shape fxn 8:    -0.0741
+DEAL::Adini<2> shape grad 8:   0.0741  -0.333
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   -0.0741 0.00
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    0.00
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  -0.0741 0.00
+DEAL::Adini<2> point 5 (0.333, 0.333)
+DEAL::Adini<2> shape fxn 0:    0.543
+DEAL::Adini<2> shape grad 0:   -0.963  -0.963
+DEAL::Adini<2> shape fxn 1:    0.0988
+DEAL::Adini<2> shape grad 1:   1.11e-16        -0.148
+DEAL::Adini<2> shape fxn 2:    0.0988
+DEAL::Adini<2> shape grad 2:   -0.148  1.11e-16
+DEAL::Adini<2> shape fxn 3:    0.198
+DEAL::Adini<2> shape grad 3:   0.963   -0.370
+DEAL::Adini<2> shape fxn 4:    -0.0494
+DEAL::Adini<2> shape grad 4:   -0.222  0.0741
+DEAL::Adini<2> shape fxn 5:    0.0494
+DEAL::Adini<2> shape grad 5:   0.148   0.00
+DEAL::Adini<2> shape fxn 6:    0.198
+DEAL::Adini<2> shape grad 6:   -0.370  0.963
+DEAL::Adini<2> shape fxn 7:    0.0494
+DEAL::Adini<2> shape grad 7:   0.00    0.148
+DEAL::Adini<2> shape fxn 8:    -0.0494
+DEAL::Adini<2> shape grad 8:   0.0741  -0.222
+DEAL::Adini<2> shape fxn 9:    0.0617
+DEAL::Adini<2> shape grad 9:   0.370   0.370
+DEAL::Adini<2> shape fxn 10:   -0.0247
+DEAL::Adini<2> shape grad 10:  -0.111  -0.0741
+DEAL::Adini<2> shape fxn 11:   -0.0247
+DEAL::Adini<2> shape grad 11:  -0.0741 -0.111
+DEAL::Adini<2> point 6 (0.667, 0.333)
+DEAL::Adini<2> shape fxn 0:    0.198
+DEAL::Adini<2> shape grad 0:   -0.963  -0.370
+DEAL::Adini<2> shape fxn 1:    0.0494
+DEAL::Adini<2> shape grad 1:   -0.222  -0.0741
+DEAL::Adini<2> shape fxn 2:    0.0494
+DEAL::Adini<2> shape grad 2:   -0.148  1.11e-16
+DEAL::Adini<2> shape fxn 3:    0.543
+DEAL::Adini<2> shape grad 3:   0.963   -0.963
+DEAL::Adini<2> shape fxn 4:    -0.0988
+DEAL::Adini<2> shape grad 4:   0.00    0.148
+DEAL::Adini<2> shape fxn 5:    0.0988
+DEAL::Adini<2> shape grad 5:   0.148   0.00
+DEAL::Adini<2> shape fxn 6:    0.0617
+DEAL::Adini<2> shape grad 6:   -0.370  0.370
+DEAL::Adini<2> shape fxn 7:    0.0247
+DEAL::Adini<2> shape grad 7:   -0.111  0.0741
+DEAL::Adini<2> shape fxn 8:    -0.0247
+DEAL::Adini<2> shape grad 8:   0.0741  -0.111
+DEAL::Adini<2> shape fxn 9:    0.198
+DEAL::Adini<2> shape grad 9:   0.370   0.963
+DEAL::Adini<2> shape fxn 10:   -0.0494
+DEAL::Adini<2> shape grad 10:  0.00    -0.148
+DEAL::Adini<2> shape fxn 11:   -0.0494
+DEAL::Adini<2> shape grad 11:  -0.0741 -0.222
+DEAL::Adini<2> point 7 (1.00, 0.333)
+DEAL::Adini<2> shape fxn 0:    -2.22e-16
+DEAL::Adini<2> shape grad 0:   -0.0741 -1.11e-16
+DEAL::Adini<2> shape fxn 1:    5.55e-17
+DEAL::Adini<2> shape grad 1:   -2.22e-16       0.00
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   -0.148  -5.55e-17
+DEAL::Adini<2> shape fxn 3:    0.741
+DEAL::Adini<2> shape grad 3:   0.0741  -1.33
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.667   0.00
+DEAL::Adini<2> shape fxn 5:    0.148
+DEAL::Adini<2> shape grad 5:   0.148   5.55e-17
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   0.0741  1.11e-16
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.00
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.0741  0.00
+DEAL::Adini<2> shape fxn 9:    0.259
+DEAL::Adini<2> shape grad 9:   -0.0741 1.33
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.333   0.00
+DEAL::Adini<2> shape fxn 11:   -0.0741
+DEAL::Adini<2> shape grad 11:  -0.0741 -0.333
+DEAL::Adini<2> point 8 (0.00, 0.667)
+DEAL::Adini<2> shape fxn 0:    0.259
+DEAL::Adini<2> shape grad 0:   0.0741  -1.33
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   0.333   0.00
+DEAL::Adini<2> shape fxn 2:    0.0741
+DEAL::Adini<2> shape grad 2:   -0.0741 -0.333
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   -0.0741 0.00
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.00
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.0741  0.00
+DEAL::Adini<2> shape fxn 6:    0.741
+DEAL::Adini<2> shape grad 6:   -0.0741 1.33
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.667   0.00
+DEAL::Adini<2> shape fxn 8:    -0.148
+DEAL::Adini<2> shape grad 8:   0.148   0.00
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   0.0741  0.00
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    0.00
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  -0.148  0.00
+DEAL::Adini<2> point 9 (0.333, 0.667)
+DEAL::Adini<2> shape fxn 0:    0.198
+DEAL::Adini<2> shape grad 0:   -0.370  -0.963
+DEAL::Adini<2> shape fxn 1:    0.0494
+DEAL::Adini<2> shape grad 1:   1.11e-16        -0.148
+DEAL::Adini<2> shape fxn 2:    0.0494
+DEAL::Adini<2> shape grad 2:   -0.0741 -0.222
+DEAL::Adini<2> shape fxn 3:    0.0617
+DEAL::Adini<2> shape grad 3:   0.370   -0.370
+DEAL::Adini<2> shape fxn 4:    -0.0247
+DEAL::Adini<2> shape grad 4:   -0.111  0.0741
+DEAL::Adini<2> shape fxn 5:    0.0247
+DEAL::Adini<2> shape grad 5:   0.0741  -0.111
+DEAL::Adini<2> shape fxn 6:    0.543
+DEAL::Adini<2> shape grad 6:   -0.963  0.963
+DEAL::Adini<2> shape fxn 7:    0.0988
+DEAL::Adini<2> shape grad 7:   0.00    0.148
+DEAL::Adini<2> shape fxn 8:    -0.0988
+DEAL::Adini<2> shape grad 8:   0.148   0.00
+DEAL::Adini<2> shape fxn 9:    0.198
+DEAL::Adini<2> shape grad 9:   0.963   0.370
+DEAL::Adini<2> shape fxn 10:   -0.0494
+DEAL::Adini<2> shape grad 10:  -0.222  -0.0741
+DEAL::Adini<2> shape fxn 11:   -0.0494
+DEAL::Adini<2> shape grad 11:  -0.148  0.00
+DEAL::Adini<2> point 10 (0.667, 0.667)
+DEAL::Adini<2> shape fxn 0:    0.0617
+DEAL::Adini<2> shape grad 0:   -0.370  -0.370
+DEAL::Adini<2> shape fxn 1:    0.0247
+DEAL::Adini<2> shape grad 1:   -0.111  -0.0741
+DEAL::Adini<2> shape fxn 2:    0.0247
+DEAL::Adini<2> shape grad 2:   -0.0741 -0.111
+DEAL::Adini<2> shape fxn 3:    0.198
+DEAL::Adini<2> shape grad 3:   0.370   -0.963
+DEAL::Adini<2> shape fxn 4:    -0.0494
+DEAL::Adini<2> shape grad 4:   0.00    0.148
+DEAL::Adini<2> shape fxn 5:    0.0494
+DEAL::Adini<2> shape grad 5:   0.0741  -0.222
+DEAL::Adini<2> shape fxn 6:    0.198
+DEAL::Adini<2> shape grad 6:   -0.963  0.370
+DEAL::Adini<2> shape fxn 7:    0.0494
+DEAL::Adini<2> shape grad 7:   -0.222  0.0741
+DEAL::Adini<2> shape fxn 8:    -0.0494
+DEAL::Adini<2> shape grad 8:   0.148   0.00
+DEAL::Adini<2> shape fxn 9:    0.543
+DEAL::Adini<2> shape grad 9:   0.963   0.963
+DEAL::Adini<2> shape fxn 10:   -0.0988
+DEAL::Adini<2> shape grad 10:  0.00    -0.148
+DEAL::Adini<2> shape fxn 11:   -0.0988
+DEAL::Adini<2> shape grad 11:  -0.148  0.00
+DEAL::Adini<2> point 11 (1.00, 0.667)
+DEAL::Adini<2> shape fxn 0:    4.44e-16
+DEAL::Adini<2> shape grad 0:   0.0741  -4.44e-16
+DEAL::Adini<2> shape fxn 1:    1.11e-16
+DEAL::Adini<2> shape grad 1:   0.00    0.00
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   -0.0741 0.00
+DEAL::Adini<2> shape fxn 3:    0.259
+DEAL::Adini<2> shape grad 3:   -0.0741 -1.33
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.333   0.00
+DEAL::Adini<2> shape fxn 5:    0.0741
+DEAL::Adini<2> shape grad 5:   0.0741  -0.333
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   -0.0741 4.44e-16
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.00
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.148   0.00
+DEAL::Adini<2> shape fxn 9:    0.741
+DEAL::Adini<2> shape grad 9:   0.0741  1.33
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.667   0.00
+DEAL::Adini<2> shape fxn 11:   -0.148
+DEAL::Adini<2> shape grad 11:  -0.148  0.00
+DEAL::Adini<2> point 12 (0.00, 1.00)
+DEAL::Adini<2> shape fxn 0:    0.00
+DEAL::Adini<2> shape grad 0:   0.00    0.00
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   0.00    0.00
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   0.00    0.00
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   0.00    0.00
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.00
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.00
+DEAL::Adini<2> shape fxn 6:    1.00
+DEAL::Adini<2> shape grad 6:   0.00    0.00
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   1.00    0.00
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    1.00
+DEAL::Adini<2> shape fxn 9:    0.00
+DEAL::Adini<2> shape grad 9:   0.00    0.00
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  0.00    0.00
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.00
+DEAL::Adini<2> point 13 (0.333, 1.00)
+DEAL::Adini<2> shape fxn 0:    0.00
+DEAL::Adini<2> shape grad 0:   -1.11e-16       -0.0741
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   -5.55e-17       -0.148
+DEAL::Adini<2> shape fxn 2:    5.55e-17
+DEAL::Adini<2> shape grad 2:   0.00    2.22e-16
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   1.11e-16        0.0741
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.0741
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.00
+DEAL::Adini<2> shape fxn 6:    0.741
+DEAL::Adini<2> shape grad 6:   -1.33   0.0741
+DEAL::Adini<2> shape fxn 7:    0.148
+DEAL::Adini<2> shape grad 7:   5.55e-17        0.148
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.667
+DEAL::Adini<2> shape fxn 9:    0.259
+DEAL::Adini<2> shape grad 9:   1.33    -0.0741
+DEAL::Adini<2> shape fxn 10:   -0.0741
+DEAL::Adini<2> shape grad 10:  -0.333  -0.0741
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.333
+DEAL::Adini<2> point 14 (0.667, 1.00)
+DEAL::Adini<2> shape fxn 0:    4.44e-16
+DEAL::Adini<2> shape grad 0:   -4.44e-16       0.0741
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   0.00    -0.0741
+DEAL::Adini<2> shape fxn 2:    1.11e-16
+DEAL::Adini<2> shape grad 2:   0.00    0.00
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   4.44e-16        -0.0741
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.148
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.00
+DEAL::Adini<2> shape fxn 6:    0.259
+DEAL::Adini<2> shape grad 6:   -1.33   -0.0741
+DEAL::Adini<2> shape fxn 7:    0.0741
+DEAL::Adini<2> shape grad 7:   -0.333  0.0741
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.333
+DEAL::Adini<2> shape fxn 9:    0.741
+DEAL::Adini<2> shape grad 9:   1.33    0.0741
+DEAL::Adini<2> shape fxn 10:   -0.148
+DEAL::Adini<2> shape grad 10:  0.00    -0.148
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    0.667
+DEAL::Adini<2> point 15 (1.00, 1.00)
+DEAL::Adini<2> shape fxn 0:    0.00
+DEAL::Adini<2> shape grad 0:   0.00    0.00
+DEAL::Adini<2> shape fxn 1:    0.00
+DEAL::Adini<2> shape grad 1:   0.00    0.00
+DEAL::Adini<2> shape fxn 2:    0.00
+DEAL::Adini<2> shape grad 2:   0.00    0.00
+DEAL::Adini<2> shape fxn 3:    0.00
+DEAL::Adini<2> shape grad 3:   0.00    0.00
+DEAL::Adini<2> shape fxn 4:    0.00
+DEAL::Adini<2> shape grad 4:   0.00    0.00
+DEAL::Adini<2> shape fxn 5:    0.00
+DEAL::Adini<2> shape grad 5:   0.00    0.00
+DEAL::Adini<2> shape fxn 6:    0.00
+DEAL::Adini<2> shape grad 6:   0.00    0.00
+DEAL::Adini<2> shape fxn 7:    0.00
+DEAL::Adini<2> shape grad 7:   0.00    0.00
+DEAL::Adini<2> shape fxn 8:    0.00
+DEAL::Adini<2> shape grad 8:   0.00    0.00
+DEAL::Adini<2> shape fxn 9:    1.00
+DEAL::Adini<2> shape grad 9:   0.00    0.00
+DEAL::Adini<2> shape fxn 10:   0.00
+DEAL::Adini<2> shape grad 10:  1.00    0.00
+DEAL::Adini<2> shape fxn 11:   0.00
+DEAL::Adini<2> shape grad 11:  0.00    1.00

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.