]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
new FE_TraceQ by Angela
authormaier <maier@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 13 Feb 2014 15:46:01 +0000 (15:46 +0000)
committermaier <maier@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 13 Feb 2014 15:46:01 +0000 (15:46 +0000)
git-svn-id: https://svn.dealii.org/trunk@32482 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/fe/fe_poly_face.templates.h
deal.II/include/deal.II/fe/fe_tools.h
deal.II/include/deal.II/fe/fe_trace.h [new file with mode: 0644]
deal.II/source/fe/CMakeLists.txt
deal.II/source/fe/fe_tools.cc
deal.II/source/fe/fe_tools.inst.in
deal.II/source/fe/fe_trace.cc [new file with mode: 0644]
deal.II/source/fe/fe_trace.inst.in [new file with mode: 0644]
tests/fe/shapes_traceq.cc [new file with mode: 0644]
tests/fe/shapes_traceq.output [new file with mode: 0644]

index c1067d752180f38fe9ffa76d8d2406c396d99570..e84a3959e266638ba644c60b6d20d4e0824ac1c3 100644 (file)
@@ -207,14 +207,45 @@ FE_PolyFace<POLY,dim,spacedim>::fill_fe_face_values (
 
   const UpdateFlags flags(fe_data.update_once | fe_data.update_each);
 
-  const unsigned int foffset = fe_data.shape_values.size() * face;
   if (flags & update_values)
     for (unsigned int i=0; i<quadrature.size(); ++i)
       {
-        for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+       for (unsigned int k=0; k<this->dofs_per_cell; ++k)
           data.shape_values(k,i) = 0.;
-        for (unsigned int k=0; k<fe_data.shape_values.size(); ++k)
-          data.shape_values(foffset+k,i) = fe_data.shape_values[k][i];
+       switch (dim)
+         {
+           case 3:
+             {
+               // Fill data for quad shape functions
+               if (this->dofs_per_quad !=0)
+                 {
+                   const unsigned int foffset = this->first_quad_index + this->dofs_per_quad * face;
+                   for (unsigned int k=0; k<this->dofs_per_quad; ++k)
+                     data.shape_values(foffset+k,i) = fe_data.shape_values[k+this->first_face_quad_index][i];
+                 }
+             }
+           case 2:
+             {
+               // Fill data for line shape functions
+               if (this->dofs_per_line != 0)
+                 {
+                   const unsigned int foffset = this->first_line_index;
+                   for (unsigned int line=0; line<GeometryInfo<dim>::lines_per_face; ++line)
+                     {
+                       for (unsigned int k=0; k<this->dofs_per_line; ++k)
+                         data.shape_values(foffset+GeometryInfo<dim>::face_to_cell_lines(face, line)*this->dofs_per_line+k,i) = fe_data.shape_values[k+(line*this->dofs_per_line)+this->first_face_line_index][i];
+                     }
+                 }
+             }
+           case 1:
+             {
+               // Fill data for vertex shape functions
+               if (this->dofs_per_vertex != 0)
+                 for (unsigned int lvertex=0; lvertex<GeometryInfo<dim>::vertices_per_face; ++lvertex)
+                   data.shape_values(GeometryInfo<dim>::face_to_cell_vertices(face, lvertex),i) = fe_data.shape_values[lvertex][i];
+               break;
+             }
+         }
       }
 }
 
index 685cacda705cb75cff3b614fdd5d3284b94df1a8..5b9e220f743149d92486748c3f5e9c1b2e9ef798 100644 (file)
@@ -761,6 +761,12 @@ namespace FETools
    * correct size, which is equal to the number of degrees of freedom in the
    * finite element.
    */
+
+  template <int dim>
+  void
+  hierarchic_to_lexicographic_numbering (unsigned int degree,
+                                        std::vector<unsigned int> &h2l);
+
   template <int dim>
   void
   hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe_data,
diff --git a/deal.II/include/deal.II/fe/fe_trace.h b/deal.II/include/deal.II/fe/fe_trace.h
new file mode 100644 (file)
index 0000000..d098aa1
--- /dev/null
@@ -0,0 +1,103 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2000 - 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef __deal2__fe_trace_h
+#define __deal2__fe_trace_h
+
+#include <deal.II/base/config.h>
+#include <deal.II/base/tensor_product_polynomials.h>
+#include <deal.II/fe/fe_poly_face.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+/**
+ * A finite element, which is the trace of Fe_Q elements, that is
+ * a tensor product of polynomials on the faces,
+ * undefined in the interior of the cells and continuous. The basis functions on
+ * the faces are from Polynomials::LagrangeEquidistant
+ *
+ * This finite element is the trace space of FE_Q on the
+ * faces.
+ *
+ * @note Since these are only finite elements on faces, only
+ * FEFaceValues and FESubfaceValues will be able to extract reasonable
+ * values from any face polynomial. In order to make the use of
+ * FESystem simpler, FEValues objects will not fail using this finite
+ * element space, but all shape function values extracted will equal
+ * to zero.
+ *
+ * @todo Polynomials::LagrangeEquidistant should be and will be
+ * replaced by Polynomials::LagrangeGaussLobatto as soon as such a
+ * polynomial set exists.
+ * @todo so far, hanging nodes are not implemented
+ *
+ */
+
+template <int dim, int spacedim=dim>
+class FE_TraceQ : public FE_PolyFace<TensorProductPolynomials<dim-1>, dim, spacedim>
+{
+public:
+  /**
+   * Constructor for tensor product
+   * polynomials of degree
+   * <tt>p</tt>. The shape
+   * functions created using this
+   * constructor correspond to
+   * Legendre polynomials in each
+   * coordinate direction.
+   */
+  FE_TraceQ(unsigned int p);
+
+  virtual FiniteElement<dim,spacedim> *clone() const;
+
+  /**
+   * Return a string that uniquely
+   * identifies a finite
+   * element. This class returns
+   * <tt>FE_DGQ<dim>(degree)</tt>, with
+   * <tt>dim</tt> and <tt>degree</tt>
+   * replaced by appropriate
+   * values.
+   */
+  virtual std::string get_name () const;
+
+  /**
+   * Check for non-zero values on a face.
+   *
+   * This function returns
+   * @p true, if the shape
+   * function @p shape_index has
+   * non-zero values on the face
+   * @p face_index.
+   *
+   * Implementation of the
+   * interface in
+   * FiniteElement
+   */
+  virtual bool has_support_on_face (const unsigned int shape_index,
+                                    const unsigned int face_index) const;
+
+private:
+  /**
+   * Return vector with dofs per
+   * vertex, line, quad, hex.
+   */
+  static std::vector<unsigned int> get_dpo_vector (const unsigned int deg);
+};
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 57b9e686f87b3bdf49c511bfbb9266023301830d..31fe760161c998c80c9e4b2a80073f8cf11198ae 100644 (file)
@@ -43,6 +43,7 @@ SET(_src
   fe_system.cc
   fe_tools.cc
   fe_tools_interpolate.cc
+  fe_trace.cc
   fe_values.cc
   fe_values_inst2.cc
   mapping_c1.cc
@@ -78,6 +79,7 @@ SET(_inst
   fe_system.inst.in
   fe_tools.inst.in
   fe_tools_interpolate.inst.in
+  fe_trace.inst.in
   fe_values.decl.1.inst.in
   fe_values.decl.2.inst.in
   fe_values.impl.1.inst.in
index b72d5a6671d4b6b07ffc9ad664414ad73c63943e..87a9685e71d07592322acdb27af4267f5ce44ef3 100644 (file)
@@ -1944,29 +1944,21 @@ namespace FETools
 
   template <int dim>
   void
-  hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe,
-                                         std::vector<unsigned int> &h2l)
+  hierarchic_to_lexicographic_numbering (unsigned int degree, std::vector<unsigned int>& h2l)
   {
-    Assert (h2l.size() == fe.dofs_per_cell,
-            ExcDimensionMismatch (h2l.size(), fe.dofs_per_cell));
-    h2l = hierarchic_to_lexicographic_numbering (fe);
-  }
-
-
+    // number of support points in each
+    // direction
+    const unsigned int n = degree+1;
 
-  template <int dim>
-  std::vector<unsigned int>
-  hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe)
-  {
-    Assert (fe.n_components() == 1, ExcInvalidFE());
+    unsigned int dofs_per_cell = n;
+    for (unsigned int i=1;i<dim;++i)
+      dofs_per_cell *= n;
 
-    std::vector<unsigned int> h2l (fe.dofs_per_cell);
+    // Assert size maches degree
+    AssertDimension (h2l.size(), dofs_per_cell);
 
-    const unsigned int dofs_per_cell = fe.dofs_per_cell;
     // polynomial degree
-    const unsigned int degree = fe.dofs_per_line+1;
-    // number of grid points in each direction
-    const unsigned int n = degree+1;
+    const unsigned int dofs_per_line = degree - 1;
 
     // the following lines of code are somewhat odd, due to the way the
     // hierarchic numbering is organized. if someone would really want to
@@ -1996,29 +1988,27 @@ namespace FETools
         h2l[next_index++] = n*n-1;
 
         // left   line
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (1+i)*n;
 
         // right  line
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (2+i)*n-1;
 
         // bottom line
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = 1+i;
 
         // top    line
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = n*(n-1)+i+1;
 
         // inside quad
-        Assert (fe.dofs_per_quad == fe.dofs_per_line*fe.dofs_per_line,
-                ExcInternalError());
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = n*(i+1)+j+1;
 
-        Assert (next_index == fe.dofs_per_cell, ExcInternalError());
+        Assert (next_index == dofs_per_cell, ExcInternalError());
 
         break;
       }
@@ -2037,82 +2027,78 @@ namespace FETools
         h2l[next_index++] = (n*n+n+1)*degree;  // 7
 
         // line 0
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (i+1)*n;
         // line 1
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = n-1+(i+1)*n;
         // line 2
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = 1+i;
         // line 3
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = 1+i+n*(n-1);
 
         // line 4
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (n-1)*n*n+(i+1)*n;
         // line 5
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (n-1)*(n*n+1)+(i+1)*n;
         // line 6
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = n*n*(n-1)+i+1;
         // line 7
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = n*n*(n-1)+i+1+n*(n-1);
 
         // line 8
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (i+1)*n*n;
         // line 9
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = n-1+(i+1)*n*n;
         // line 10
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = (i+1)*n*n+n*(n-1);
         // line 11
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
           h2l[next_index++] = n-1+(i+1)*n*n+n*(n-1);
 
 
         // inside quads
-        Assert (fe.dofs_per_quad == fe.dofs_per_line*fe.dofs_per_line,
-                ExcInternalError());
         // face 0
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = (i+1)*n*n+n*(j+1);
         // face 1
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = (i+1)*n*n+n-1+n*(j+1);
         // face 2, note the orientation!
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = (j+1)*n*n+i+1;
         // face 3, note the orientation!
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = (j+1)*n*n+n*(n-1)+i+1;
         // face 4
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = n*(i+1)+j+1;
         // face 5
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
             h2l[next_index++] = (n-1)*n*n+n*(i+1)+j+1;
 
         // inside hex
-        Assert (fe.dofs_per_hex == fe.dofs_per_quad*fe.dofs_per_line,
-                ExcInternalError());
-        for (unsigned int i=0; i<fe.dofs_per_line; ++i)
-          for (unsigned int j=0; j<fe.dofs_per_line; ++j)
-            for (unsigned int k=0; k<fe.dofs_per_line; ++k)
+        for (unsigned int i=0; i<dofs_per_line; ++i)
+          for (unsigned int j=0; j<dofs_per_line; ++j)
+            for (unsigned int k=0; k<dofs_per_line; ++k)
               h2l[next_index++]       = n*n*(i+1)+n*(j+1)+k+1;
 
-        Assert (next_index == fe.dofs_per_cell, ExcInternalError());
+        Assert (next_index == dofs_per_cell, ExcInternalError());
 
         break;
       }
@@ -2120,12 +2106,32 @@ namespace FETools
       default:
         Assert (false, ExcNotImplemented());
       }
+  }
+
 
-    return h2l;
+
+  template <int dim>
+  void
+  hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe,
+                                         std::vector<unsigned int> &h2l)
+  {
+    Assert (h2l.size() == fe.dofs_per_cell,
+            ExcDimensionMismatch (h2l.size(), fe.dofs_per_cell));
+    hierarchic_to_lexicographic_numbering<dim> (fe.dofs_per_line+1, h2l);
   }
 
 
 
+  template <int dim>
+  std::vector<unsigned int>
+  hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe)
+  {
+    Assert (fe.n_components() == 1, ExcInvalidFE());
+    std::vector<unsigned int> h2l(fe.dofs_per_cell);
+    hierarchic_to_lexicographic_numbering<dim> (fe.dofs_per_line+1, h2l);
+    return (h2l);
+  }
+
   template <int dim>
   void
   lexicographic_to_hierarchic_numbering (const FiniteElementData<dim> &fe,
index 3809bdd380bed56c88d4c3aed02c757cd0ef915d..3d7e9c59ced962f1257599aa1cc3e3a3db867273 100644 (file)
@@ -153,6 +153,12 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS
                                                               FullMatrix<double>       &);
 #endif
 
+      template
+       void
+       hierarchic_to_lexicographic_numbering<deal_II_dimension>
+       (unsigned int,
+        std::vector<unsigned int> &);
+
       template
        void
        hierarchic_to_lexicographic_numbering<deal_II_dimension>
diff --git a/deal.II/source/fe/fe_trace.cc b/deal.II/source/fe/fe_trace.cc
new file mode 100644 (file)
index 0000000..902ae85
--- /dev/null
@@ -0,0 +1,266 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2000 - 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.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/grid/tria.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/dofs/dof_accessor.h>
+
+#include <deal.II/base/config.h>
+#include <deal.II/base/tensor_product_polynomials.h>
+#include <deal.II/fe/fe_poly_face.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_tools.h>
+#include <deal.II/dofs/dof_tools.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/function.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/numerics/matrix_tools.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/full_matrix.h>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/lac/compressed_sparsity_pattern.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/precondition.h>
+#include <deal.II/numerics/data_out.h>
+#include <deal.II/numerics/data_out_faces.h>
+#include <fstream>
+#include <iostream>
+
+
+#include <deal.II/fe/fe_face.h>
+#include <deal.II/fe/fe_poly_face.templates.h>
+#include <sstream>
+
+DEAL_II_NAMESPACE_OPEN
+
+/** FE-Class for continuous face functions
+* adaption of FE_Q and FE_FaceQ
+*/
+
+template <int dim, int spacedim=dim>
+class FE_TraceQ : public FE_PolyFace<TensorProductPolynomials<dim-1>, dim, spacedim>
+{
+  public:
+                                     /**
+                                      * Constructor for tensor product
+                                      * polynomials of degree
+                                      * <tt>p</tt>. The shape
+                                      * functions created using this
+                                      * constructor correspond to
+                                      * Legendre polynomials in each
+                                      * coordinate direction.
+                                      */
+    FE_TraceQ(unsigned int p);
+
+    virtual FiniteElement<dim,spacedim>* clone() const;
+
+                                     /**
+                                      * Return a string that uniquely
+                                      * identifies a finite
+                                      * element. This class returns
+                                      * <tt>FE_DGQ<dim>(degree)</tt> , with
+                                      * <tt>dim</tt> and <tt>degree</tt>
+                                      * replaced by appropriate
+                                      * values.
+                                      */
+    virtual std::string get_name () const;
+
+                                     /**
+                                      * Check for non-zero values on a face.
+                                      *
+                                      * This function returns
+                                      * @p true, if the shape
+                                      * function @p shape_index has
+                                      * non-zero values on the face
+                                      * @p face_index.
+                                      *
+                                      * Implementation of the
+                                      * interface in
+                                      * FiniteElement
+                                      */
+    virtual bool has_support_on_face (const unsigned int shape_index,
+                                      const unsigned int face_index) const;
+
+  private:
+                                     /**
+                                      * Return vector with dofs per
+                                      * vertex, line, quad, hex.
+                                      */
+    static std::vector<unsigned int> get_dpo_vector (const unsigned int deg);
+};
+
+
+template <int dim, int spacedim>
+FE_TraceQ<dim,spacedim>::FE_TraceQ (const unsigned int degree)
+                :
+                FE_PolyFace<TensorProductPolynomials<dim-1>, dim, spacedim> (
+                  TensorProductPolynomials<dim-1>(Polynomials::LagrangeEquidistant::generate_complete_basis(degree)),
+                  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
+                  std::vector<bool>(1,true))
+{
+  Assert (degree > 0,
+          ExcMessage ("FE_Trace can only be used for polynomial degrees "
+                      "greater than zero"));
+  std::vector<unsigned int> renumber (this->dofs_per_face);
+  FETools::hierarchic_to_lexicographic_numbering<dim-1> (degree, renumber);
+  this->poly_space.set_numbering(renumber);
+}
+
+
+template <int dim, int spacedim>
+FiniteElement<dim,spacedim>*
+FE_TraceQ<dim,spacedim>::clone() const
+{
+  return new FE_TraceQ<dim,spacedim>(this->degree);
+}
+
+
+template <int dim, int spacedim>
+std::string
+FE_TraceQ<dim,spacedim>::get_name () const
+{
+                                   // note that the
+                                   // FETools::get_fe_from_name
+                                   // function depends on the
+                                   // particular format of the string
+                                   // this function returns, so they
+                                   // have to be kept in synch
+
+  std::ostringstream namebuf;
+  namebuf << "FE_TraceQ<" << dim << ">(" << this->degree << ")";
+
+  return namebuf.str();
+}
+
+template <int dim, int spacedim>
+bool
+FE_TraceQ<dim,spacedim>::has_support_on_face (const unsigned int shape_index,
+                                const unsigned int face_index) const
+{
+  Assert (shape_index < this->dofs_per_cell,
+          ExcIndexRange (shape_index, 0, this->dofs_per_cell));
+  Assert (face_index < GeometryInfo<dim>::faces_per_cell,
+          ExcIndexRange (face_index, 0, GeometryInfo<dim>::faces_per_cell));
+
+                                   // let's see whether this is a
+                                   // vertex
+  if (shape_index < this->first_line_index)
+    {
+                                       // for Q elements, there is
+                                       // one dof per vertex, so
+                                       // shape_index==vertex_number. check
+                                       // whether this vertex is
+                                       // on the given face. thus,
+                                       // for each face, give a
+                                       // list of vertices
+      const unsigned int vertex_no = shape_index;
+      Assert (vertex_no < GeometryInfo<dim>::vertices_per_cell,
+              ExcInternalError());
+
+      for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
+        if (GeometryInfo<dim>::face_to_cell_vertices(face_index, v) == vertex_no)
+          return true;
+
+      return false;
+    }
+  else if (shape_index < this->first_quad_index)
+                                     // ok, dof is on a line
+    {
+      const unsigned int line_index
+        = (shape_index - this->first_line_index) / this->dofs_per_line;
+      Assert (line_index < GeometryInfo<dim>::lines_per_cell,
+              ExcInternalError());
+
+                                       // in 2d, the line is the
+                                       // face, so get the line
+                                       // index
+      if (dim == 2)
+        return (line_index == face_index);
+      else if (dim == 3)
+        {
+                                           // see whether the
+                                           // given line is on the
+                                           // given face.
+          for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_face; ++l)
+            if (GeometryInfo<3>::face_to_cell_lines(face_index, l) == line_index)
+              return true;
+
+          return false;
+        }
+      else
+        Assert (false, ExcNotImplemented());
+    }
+  else if (shape_index < this->first_hex_index)
+                                     // dof is on a quad
+    {
+      const unsigned int quad_index
+        = (shape_index - this->first_quad_index) / this->dofs_per_quad;
+      Assert (static_cast<signed int>(quad_index) <
+              static_cast<signed int>(GeometryInfo<dim>::quads_per_cell),
+              ExcInternalError());
+
+                                       // in 2d, cell bubble are
+                                       // zero on all faces. but
+                                       // we have treated this
+                                       // case above already
+      Assert (dim != 2, ExcInternalError());
+
+                                       // in 3d,
+                                       // quad_index=face_index
+      if (dim == 3)
+        return (quad_index == face_index);
+      else
+        Assert (false, ExcNotImplemented());
+    }
+  else
+                                     // dof on hex
+    {
+                                       // can only happen in 3d,
+                                       // but this case has
+                                       // already been covered
+                                       // above
+      Assert (false, ExcNotImplemented());
+      return false;
+    }
+
+                                   // we should not have gotten here
+  Assert (false, ExcInternalError());
+  return false;
+}
+
+template <int dim, int spacedim>
+std::vector<unsigned int>
+FE_TraceQ<dim,spacedim>::get_dpo_vector (const unsigned int deg)
+{
+  AssertThrow(deg>0,ExcMessage("FE_TraceQ needs to be of degree > 0."));
+  std::vector<unsigned int> dpo(dim+1, 1U);
+  dpo[dim]=0;
+  dpo[0]=1;
+  for (unsigned int i=1;i<dim;++i)
+    dpo[i] = dpo[i-1]*(deg-1);
+  return dpo;
+}
+
+// explicit instantiations
+#include "fe_trace.inst"
+
+
+DEAL_II_NAMESPACE_CLOSE
diff --git a/deal.II/source/fe/fe_trace.inst.in b/deal.II/source/fe/fe_trace.inst.in
new file mode 100644 (file)
index 0000000..21748c4
--- /dev/null
@@ -0,0 +1,22 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2000 - 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.
+//
+// ---------------------------------------------------------------------
+
+for (deal_II_dimension : DIMENSIONS)
+  {
+#if deal_II_dimension > 1
+    template class FE_TraceQ<deal_II_dimension>;
+#endif
+  }
diff --git a/tests/fe/shapes_traceq.cc b/tests/fe/shapes_traceq.cc
new file mode 100644 (file)
index 0000000..25bcc69
--- /dev/null
@@ -0,0 +1,71 @@
+// ---------------------------------------------------------------------
+// $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.
+//
+// ---------------------------------------------------------------------
+
+#include "../tests.h"
+#include "shapes.h"
+//#include "../deal.II/tests/tests.h"
+//#include "../deal.II/tests/fe/shapes.h"
+#include <deal.II/fe/fe_trace.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/mapping_q1.h>
+#include <fstream>
+#include <string>
+
+#define PRECISION 2
+
+
+template<int dim>
+void plot_FE_TraceQ_shape_functions()
+{
+  MappingQ1<dim> m;
+
+  FE_TraceQ<dim> tq1(1);
+  FE_TraceQ<dim> tq2(2);
+  FE_TraceQ<dim> tq3(3);
+  FE_TraceQ<dim> tq4(4);
+  FE_Q<dim> q1(1);
+  FE_Q<dim> q2(2);
+  FESystem<dim> sys1(tq1,1);
+  FESystem<dim> sys2(tq2,1);
+  FESystem<dim> sys11(tq1,1,q1,1);
+  FESystem<dim> sys22(tq2,1,q2,1);
+  plot_face_shape_functions(m, tq1, "TraceQ1", update_values);
+  plot_face_shape_functions(m, tq2, "TraceQ2", update_values);
+  plot_face_shape_functions(m, tq3, "TraceQ3", update_values);
+  plot_face_shape_functions(m, tq4, "TraceQ4", update_values);
+  plot_face_shape_functions(m, sys1, "System1", update_values);
+  plot_face_shape_functions(m, sys2, "System2", update_values);
+  plot_face_shape_functions(m, sys11, "System1-1", update_values);
+  plot_face_shape_functions(m, sys22, "System2-2", update_values);
+}
+
+
+int
+main()
+{
+  const std::string logname = "output";
+  std::ofstream logfile(logname.c_str());
+  deallog << std::setprecision(PRECISION) << std::fixed;
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  plot_FE_TraceQ_shape_functions<2>();
+  //plot_FE_TraceQ_shape_functions<3>();
+
+  return 0;
+}
diff --git a/tests/fe/shapes_traceq.output b/tests/fe/shapes_traceq.output
new file mode 100644 (file)
index 0000000..b71e809
--- /dev/null
@@ -0,0 +1,281 @@
+
+DEAL:Face2d-TraceQ1::0.00 0.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ1::0.00 0.12 1.75 1.00 1.25 1.00
+DEAL:Face2d-TraceQ1::0.00 0.25 1.50 1.00 1.50 1.00
+DEAL:Face2d-TraceQ1::0.00 0.38 1.25 1.00 1.75 1.00
+DEAL:Face2d-TraceQ1::0.00 0.50 1.00 1.00 2.00 1.00
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::0.50 0.00 1.00 1.00 2.00 1.00
+DEAL:Face2d-TraceQ1::0.50 0.06 1.00 1.00 1.88 1.12
+DEAL:Face2d-TraceQ1::0.50 0.12 1.00 1.00 1.75 1.25
+DEAL:Face2d-TraceQ1::0.50 0.19 1.00 1.00 1.62 1.38
+DEAL:Face2d-TraceQ1::0.50 0.25 1.00 1.00 1.50 1.50
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::0.50 0.25 1.00 1.00 1.50 1.50
+DEAL:Face2d-TraceQ1::0.50 0.31 1.00 1.00 1.38 1.62
+DEAL:Face2d-TraceQ1::0.50 0.38 1.00 1.00 1.25 1.75
+DEAL:Face2d-TraceQ1::0.50 0.44 1.00 1.00 1.12 1.88
+DEAL:Face2d-TraceQ1::0.50 0.50 1.00 1.00 1.00 2.00
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::0.00 0.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ1::0.12 0.00 1.75 1.25 1.00 1.00
+DEAL:Face2d-TraceQ1::0.25 0.00 1.50 1.50 1.00 1.00
+DEAL:Face2d-TraceQ1::0.38 0.00 1.25 1.75 1.00 1.00
+DEAL:Face2d-TraceQ1::0.50 0.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::0.00 0.50 1.00 1.00 2.00 1.00
+DEAL:Face2d-TraceQ1::0.12 0.50 1.00 1.00 1.75 1.25
+DEAL:Face2d-TraceQ1::0.25 0.50 1.00 1.00 1.50 1.50
+DEAL:Face2d-TraceQ1::0.38 0.50 1.00 1.00 1.25 1.75
+DEAL:Face2d-TraceQ1::0.50 0.50 1.00 1.00 1.00 2.00
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ1::
+DEAL:Face2d-TraceQ2::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.00 0.12 1.38 1.00 0.88 1.00 1.75 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.00 0.25 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.00 0.38 0.88 1.00 1.38 1.00 1.75 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::0.50 0.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.06 1.00 1.00 1.00 1.66 0.91 1.44 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.12 1.00 1.00 1.00 1.38 0.88 1.75 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.19 1.00 1.00 1.00 1.16 0.91 1.94 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.25 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::0.50 0.25 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.31 1.00 1.00 1.00 0.91 1.16 1.94 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.38 1.00 1.00 1.00 0.88 1.38 1.75 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.44 1.00 1.00 1.00 0.91 1.66 1.44 1.00 1.00
+DEAL:Face2d-TraceQ2::0.50 0.50 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.12 0.00 1.38 0.88 1.00 1.00 1.00 1.00 1.75 1.00
+DEAL:Face2d-TraceQ2::0.25 0.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00
+DEAL:Face2d-TraceQ2::0.38 0.00 0.88 1.38 1.00 1.00 1.00 1.00 1.75 1.00
+DEAL:Face2d-TraceQ2::0.50 0.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::0.12 0.50 1.00 1.00 1.38 0.88 1.00 1.00 1.00 1.75
+DEAL:Face2d-TraceQ2::0.25 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00
+DEAL:Face2d-TraceQ2::0.38 0.50 1.00 1.00 0.88 1.38 1.00 1.00 1.00 1.75
+DEAL:Face2d-TraceQ2::0.50 0.50 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ2::
+DEAL:Face2d-TraceQ3::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.00 0.12 1.12 1.00 1.04 1.00 2.05 0.79 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.00 0.25 0.94 1.00 0.94 1.00 1.56 1.56 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.00 0.38 1.04 1.00 1.12 1.00 0.79 2.05 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::0.50 0.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.06 1.00 1.00 1.00 1.00 1.44 1.06 1.80 0.69 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.12 1.00 1.00 1.00 1.00 1.12 1.04 2.05 0.79 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.19 1.00 1.00 1.00 1.00 0.97 0.98 1.92 1.13 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.25 1.00 1.00 1.00 1.00 0.94 0.94 1.56 1.56 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::0.50 0.25 1.00 1.00 1.00 1.00 0.94 0.94 1.56 1.56 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.31 1.00 1.00 1.00 1.00 0.98 0.97 1.13 1.92 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.38 1.00 1.00 1.00 1.00 1.04 1.12 0.79 2.05 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.44 1.00 1.00 1.00 1.00 1.06 1.44 0.69 1.80 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.50 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.12 0.00 1.12 1.04 1.00 1.00 1.00 1.00 1.00 1.00 2.05 0.79 1.00 1.00
+DEAL:Face2d-TraceQ3::0.25 0.00 0.94 0.94 1.00 1.00 1.00 1.00 1.00 1.00 1.56 1.56 1.00 1.00
+DEAL:Face2d-TraceQ3::0.38 0.00 1.04 1.12 1.00 1.00 1.00 1.00 1.00 1.00 0.79 2.05 1.00 1.00
+DEAL:Face2d-TraceQ3::0.50 0.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::0.12 0.50 1.00 1.00 1.12 1.04 1.00 1.00 1.00 1.00 1.00 1.00 2.05 0.79
+DEAL:Face2d-TraceQ3::0.25 0.50 1.00 1.00 0.94 0.94 1.00 1.00 1.00 1.00 1.00 1.00 1.56 1.56
+DEAL:Face2d-TraceQ3::0.38 0.50 1.00 1.00 1.04 1.12 1.00 1.00 1.00 1.00 1.00 1.00 0.79 2.05
+DEAL:Face2d-TraceQ3::0.50 0.50 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ3::
+DEAL:Face2d-TraceQ4::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.00 0.12 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.00 0.25 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.00 0.38 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::0.50 0.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.06 1.00 1.00 1.00 1.00 1.00 1.27 0.96 2.09 0.45 1.22 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.12 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.19 1.00 1.00 1.00 1.00 1.00 0.96 1.02 1.47 1.70 0.84 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.25 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::0.50 0.25 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.31 1.00 1.00 1.00 1.00 1.00 1.02 0.96 0.84 1.70 1.47 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.38 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.44 1.00 1.00 1.00 1.00 1.00 0.96 1.27 1.22 0.45 2.09 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.50 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.12 0.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.25 0.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.38 0.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.50 0.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.12 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-TraceQ4::0.25 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00
+DEAL:Face2d-TraceQ4::0.38 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00
+DEAL:Face2d-TraceQ4::0.50 0.50 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-TraceQ4::
+DEAL:Face2d-System1::0.00 0.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-System1::0.00 0.12 1.75 1.00 1.25 1.00
+DEAL:Face2d-System1::0.00 0.25 1.50 1.00 1.50 1.00
+DEAL:Face2d-System1::0.00 0.38 1.25 1.00 1.75 1.00
+DEAL:Face2d-System1::0.00 0.50 1.00 1.00 2.00 1.00
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::0.50 0.00 1.00 1.00 2.00 1.00
+DEAL:Face2d-System1::0.50 0.06 1.00 1.00 1.88 1.12
+DEAL:Face2d-System1::0.50 0.12 1.00 1.00 1.75 1.25
+DEAL:Face2d-System1::0.50 0.19 1.00 1.00 1.62 1.38
+DEAL:Face2d-System1::0.50 0.25 1.00 1.00 1.50 1.50
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::0.50 0.25 1.00 1.00 1.50 1.50
+DEAL:Face2d-System1::0.50 0.31 1.00 1.00 1.38 1.62
+DEAL:Face2d-System1::0.50 0.38 1.00 1.00 1.25 1.75
+DEAL:Face2d-System1::0.50 0.44 1.00 1.00 1.12 1.88
+DEAL:Face2d-System1::0.50 0.50 1.00 1.00 1.00 2.00
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::0.00 0.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-System1::0.12 0.00 1.75 1.25 1.00 1.00
+DEAL:Face2d-System1::0.25 0.00 1.50 1.50 1.00 1.00
+DEAL:Face2d-System1::0.38 0.00 1.25 1.75 1.00 1.00
+DEAL:Face2d-System1::0.50 0.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::0.00 0.50 1.00 1.00 2.00 1.00
+DEAL:Face2d-System1::0.12 0.50 1.00 1.00 1.75 1.25
+DEAL:Face2d-System1::0.25 0.50 1.00 1.00 1.50 1.50
+DEAL:Face2d-System1::0.38 0.50 1.00 1.00 1.25 1.75
+DEAL:Face2d-System1::0.50 0.50 1.00 1.00 1.00 2.00
+DEAL:Face2d-System1::
+DEAL:Face2d-System1::
+DEAL:Face2d-System2::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::0.00 0.12 1.38 1.00 0.88 1.00 1.75 1.00 1.00 1.00
+DEAL:Face2d-System2::0.00 0.25 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-System2::0.00 0.38 0.88 1.00 1.38 1.00 1.75 1.00 1.00 1.00
+DEAL:Face2d-System2::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::0.50 0.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::0.50 0.06 1.00 1.00 1.00 1.66 0.91 1.44 1.00 1.00
+DEAL:Face2d-System2::0.50 0.12 1.00 1.00 1.00 1.38 0.88 1.75 1.00 1.00
+DEAL:Face2d-System2::0.50 0.19 1.00 1.00 1.00 1.16 0.91 1.94 1.00 1.00
+DEAL:Face2d-System2::0.50 0.25 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::0.50 0.25 1.00 1.00 1.00 1.00 1.00 2.00 1.00 1.00
+DEAL:Face2d-System2::0.50 0.31 1.00 1.00 1.00 0.91 1.16 1.94 1.00 1.00
+DEAL:Face2d-System2::0.50 0.38 1.00 1.00 1.00 0.88 1.38 1.75 1.00 1.00
+DEAL:Face2d-System2::0.50 0.44 1.00 1.00 1.00 0.91 1.66 1.44 1.00 1.00
+DEAL:Face2d-System2::0.50 0.50 1.00 1.00 1.00 1.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::0.00 0.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::0.12 0.00 1.38 0.88 1.00 1.00 1.00 1.00 1.75 1.00
+DEAL:Face2d-System2::0.25 0.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 1.00
+DEAL:Face2d-System2::0.38 0.00 0.88 1.38 1.00 1.00 1.00 1.00 1.75 1.00
+DEAL:Face2d-System2::0.50 0.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::0.00 0.50 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::0.12 0.50 1.00 1.00 1.38 0.88 1.00 1.00 1.00 1.75
+DEAL:Face2d-System2::0.25 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00
+DEAL:Face2d-System2::0.38 0.50 1.00 1.00 0.88 1.38 1.00 1.00 1.00 1.75
+DEAL:Face2d-System2::0.50 0.50 1.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2::
+DEAL:Face2d-System2::
+DEAL:Face2d-System1-1::0.00 0.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::0.00 0.12 1.75 1.75 1.00 1.00 1.25 1.25 1.00 1.00
+DEAL:Face2d-System1-1::0.00 0.25 1.50 1.50 1.00 1.00 1.50 1.50 1.00 1.00
+DEAL:Face2d-System1-1::0.00 0.38 1.25 1.25 1.00 1.00 1.75 1.75 1.00 1.00
+DEAL:Face2d-System1-1::0.00 0.50 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::0.50 0.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::0.50 0.06 1.00 1.00 1.00 1.88 1.88 1.00 1.12 1.12
+DEAL:Face2d-System1-1::0.50 0.12 1.00 1.00 1.00 1.75 1.75 1.00 1.25 1.25
+DEAL:Face2d-System1-1::0.50 0.19 1.00 1.00 1.00 1.62 1.62 1.00 1.38 1.38
+DEAL:Face2d-System1-1::0.50 0.25 1.00 1.00 1.00 1.50 1.50 1.00 1.50 1.50
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::0.50 0.25 1.00 1.00 1.00 1.50 1.50 1.00 1.50 1.50
+DEAL:Face2d-System1-1::0.50 0.31 1.00 1.00 1.00 1.38 1.38 1.00 1.62 1.62
+DEAL:Face2d-System1-1::0.50 0.38 1.00 1.00 1.00 1.25 1.25 1.00 1.75 1.75
+DEAL:Face2d-System1-1::0.50 0.44 1.00 1.00 1.00 1.12 1.12 1.00 1.88 1.88
+DEAL:Face2d-System1-1::0.50 0.50 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::0.00 0.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::0.12 0.00 1.75 1.75 1.25 1.25 1.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::0.25 0.00 1.50 1.50 1.50 1.50 1.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::0.38 0.00 1.25 1.25 1.75 1.75 1.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::0.50 0.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::0.00 0.50 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00
+DEAL:Face2d-System1-1::0.12 0.50 1.00 1.00 1.00 1.00 1.75 1.75 1.25 1.25
+DEAL:Face2d-System1-1::0.25 0.50 1.00 1.00 1.00 1.00 1.50 1.50 1.50 1.50
+DEAL:Face2d-System1-1::0.38 0.50 1.00 1.00 1.00 1.00 1.25 1.25 1.75 1.75
+DEAL:Face2d-System1-1::0.50 0.50 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System1-1::
+DEAL:Face2d-System2-2::0.00 0.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.00 0.12 1.38 1.38 1.00 1.00 0.88 0.88 1.00 1.00 1.75 1.75 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.00 0.25 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.00 0.38 0.88 0.88 1.00 1.00 1.38 1.38 1.00 1.00 1.75 1.75 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.00 0.50 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::0.50 0.00 1.00 1.00 1.00 2.00 1.00 1.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.06 1.00 1.00 1.00 1.66 1.00 1.00 1.66 0.91 0.91 1.00 1.44 1.44 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.12 1.00 1.00 1.00 1.38 1.00 1.00 1.38 0.88 0.88 1.00 1.75 1.75 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.19 1.00 1.00 1.00 1.16 1.00 1.00 1.16 0.91 0.91 1.00 1.94 1.94 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.25 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::0.50 0.25 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.31 1.00 1.00 1.00 0.91 1.00 1.00 0.91 1.16 1.16 1.00 1.94 1.94 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.38 1.00 1.00 1.00 0.88 1.00 1.00 0.88 1.38 1.38 1.00 1.75 1.75 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.44 1.00 1.00 1.00 0.91 1.00 1.00 0.91 1.66 1.66 1.00 1.44 1.44 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::0.00 0.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.12 0.00 1.38 1.38 0.88 0.88 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.75 1.75 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.25 0.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.38 0.00 0.88 0.88 1.38 1.38 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.75 1.75 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.50 0.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::0.00 0.50 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::0.12 0.50 1.00 1.00 1.00 1.00 1.38 1.38 0.88 0.88 1.00 1.00 1.00 1.00 1.00 1.00 1.75 1.75 1.00
+DEAL:Face2d-System2-2::0.25 0.50 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00
+DEAL:Face2d-System2-2::0.38 0.50 1.00 1.00 1.00 1.00 0.88 0.88 1.38 1.38 1.00 1.00 1.00 1.00 1.00 1.00 1.75 1.75 1.00
+DEAL:Face2d-System2-2::0.50 0.50 1.00 1.00 1.00 1.00 1.00 1.00 2.00 2.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00 1.00
+DEAL:Face2d-System2-2::
+DEAL:Face2d-System2-2::

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.