--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_simplex_fe_lib_h
+#define dealii_simplex_fe_lib_h
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/fe/fe_poly.h>
+
+#include <deal.II/simplex/polynomials.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Simplex
+{
+ /**
+ * Base class of FE_P and FE_DGP.
+ *
+ * @note Only implemented for 2D and 3D.
+ *
+ * @ingroup simplex
+ */
+ template <int dim, int spacedim = dim>
+ class FE_Poly : public dealii::FE_Poly<dim, spacedim>
+ {
+ public:
+ /**
+ * Constructor.
+ */
+ FE_Poly(const unsigned int degree,
+ const std::vector<unsigned int> &dpo_vector);
+
+ private:
+ /**
+ * @copydoc dealii::FiniteElement::convert_generalized_support_point_values_to_dof_values()
+ */
+ void
+ convert_generalized_support_point_values_to_dof_values(
+ const std::vector<Vector<double>> &support_point_values,
+ std::vector<double> & nodal_values) const override;
+ };
+
+
+
+ /**
+ * Implementation of a scalar Lagrange finite element Pp that yields
+ * the finite element space of continuous, piecewise polynomials of
+ * degree p.
+ *
+ * @ingroup simplex
+ */
+ template <int dim, int spacedim = dim>
+ class FE_P : public FE_Poly<dim, spacedim>
+ {
+ public:
+ /**
+ * Constructor.
+ */
+ FE_P(const unsigned int degree);
+
+ /**
+ * @copydoc dealii::FiniteElement::clone()
+ */
+ std::unique_ptr<FiniteElement<dim, spacedim>>
+ clone() const override;
+
+ /**
+ * Return a string that uniquely identifies a finite element. This class
+ * returns <tt>Simplex::FE_P<dim>(degree)</tt>, with @p dim and @p degree
+ * replaced by appropriate values.
+ */
+ std::string
+ get_name() const override;
+ };
+
+
+
+ /**
+ * Implementation of a scalar Lagrange finite element Pp that yields
+ * the finite element space of discontinuous, piecewise polynomials of
+ * degree p.
+ *
+ * @ingroup simplex
+ */
+ template <int dim, int spacedim = dim>
+ class FE_DGP : public FE_Poly<dim, spacedim>
+ {
+ public:
+ /**
+ * Constructor.
+ */
+ FE_DGP(const unsigned int degree);
+
+ /**
+ * @copydoc dealii::FiniteElement::clone()
+ */
+ std::unique_ptr<FiniteElement<dim, spacedim>>
+ clone() const override;
+
+ /**
+ * Return a string that uniquely identifies a finite element. This class
+ * returns <tt>Simplex::FE_DGP<dim>(degree)</tt>, with @p dim and @p degree
+ * replaced by appropriate values.
+ */
+ std::string
+ get_name() const override;
+ };
+
+} // namespace Simplex
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
INCLUDE_DIRECTORIES(BEFORE ${CMAKE_CURRENT_BINARY_DIR})
SET(_unity_include_src
+ fe_lib.cc
polynomials.cc
quadrature_lib.cc
)
)
SET(_inst
+ fe_lib.inst.in
)
FILE(GLOB _header
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/simplex/fe_lib.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Simplex
+{
+ namespace
+ {
+ /**
+ * Helper function to set up the dpo vector of FE_P for a given @p dim and
+ * @p degree.
+ */
+ std::vector<unsigned int>
+ get_dpo_vector_fe_p(const unsigned int dim, const unsigned int degree)
+ {
+ std::vector<unsigned int> dpo(dim + 1, 0U);
+
+ if (degree == 1)
+ {
+ // one dof at each vertex
+ dpo[0] = 1;
+ }
+ else if (degree == 2)
+ {
+ // one dof at each vertex and in the middle of each line
+ dpo[0] = 1;
+ dpo[1] = 1;
+ dpo[2] = 0;
+ }
+ else
+ {
+ Assert(false, ExcNotImplemented());
+ }
+
+ return dpo;
+ }
+
+ /**
+ * Helper function to set up the dpo vector of FE_DGP for a given @p dim and
+ * @p degree.
+ */
+ std::vector<unsigned int>
+ get_dpo_vector_fe_dgp(const unsigned int dim, const unsigned int degree)
+ {
+ std::vector<unsigned int> dpo(dim + 1, 0U);
+
+ // all dofs are internal
+ if (dim == 2 && degree == 1)
+ dpo[dim] = 3;
+ else if (dim == 2 && degree == 2)
+ dpo[dim] = 6;
+ else if (dim == 3 && degree == 1)
+ dpo[dim] = 4;
+ else if (dim == 3 && degree == 2)
+ dpo[dim] = 10;
+ else
+ {
+ Assert(false, ExcNotImplemented());
+ }
+
+ return dpo;
+ }
+ } // namespace
+
+
+
+ template <int dim, int spacedim>
+ FE_Poly<dim, spacedim>::FE_Poly(const unsigned int degree,
+ const std::vector<unsigned int> &dpo_vector)
+ : dealii::FE_Poly<dim, spacedim>(
+ Simplex::ScalarPolynomial<dim>(degree),
+ FiniteElementData<dim>(dpo_vector,
+ dim == 2 ? ReferenceCell::Type::Tri :
+ ReferenceCell::Type::Tet,
+ 1,
+ degree,
+ FiniteElementData<dim>::L2),
+ std::vector<bool>(FiniteElementData<dim>(dpo_vector,
+ dim == 2 ?
+ ReferenceCell::Type::Tri :
+ ReferenceCell::Type::Tet,
+ 1,
+ degree)
+ .dofs_per_cell,
+ true),
+ std::vector<ComponentMask>(
+ FiniteElementData<dim>(dpo_vector,
+ dim == 2 ? ReferenceCell::Type::Tri :
+ ReferenceCell::Type::Tet,
+ 1,
+ degree)
+ .dofs_per_cell,
+ std::vector<bool>(1, true)))
+ {
+ this->unit_support_points.clear();
+
+ if (dim == 2)
+ {
+ if (degree == 1)
+ {
+ this->unit_support_points.emplace_back(1.0, 0.0);
+ this->unit_support_points.emplace_back(0.0, 1.0);
+ this->unit_support_points.emplace_back(0.0, 0.0);
+
+ // TODO
+ this->unit_face_support_points.emplace_back(0.0);
+ this->unit_face_support_points.emplace_back(1.0);
+ }
+ else if (degree == 2)
+ {
+ this->unit_support_points.emplace_back(1.0, 0.0);
+ this->unit_support_points.emplace_back(0.0, 1.0);
+ this->unit_support_points.emplace_back(0.0, 0.0);
+ this->unit_support_points.emplace_back(0.5, 0.5);
+ this->unit_support_points.emplace_back(0.0, 0.5);
+ this->unit_support_points.emplace_back(0.5, 0.0);
+
+ // TODO
+ this->unit_face_support_points.emplace_back(0.0);
+ this->unit_face_support_points.emplace_back(1.0);
+ this->unit_face_support_points.emplace_back(0.5);
+ }
+ else
+ {
+ Assert(false, ExcNotImplemented());
+ }
+ }
+ else if (dim == 3)
+ {
+ if (degree == 1)
+ {
+ this->unit_support_points.emplace_back(0.0, 0.0, 0.0);
+ this->unit_support_points.emplace_back(1.0, 0.0, 0.0);
+ this->unit_support_points.emplace_back(0.0, 1.0, 0.0);
+ this->unit_support_points.emplace_back(0.0, 0.0, 1.0);
+
+ // TODO
+ this->unit_face_support_points.emplace_back(1.0, 0.0);
+ this->unit_face_support_points.emplace_back(0.0, 1.0);
+ this->unit_face_support_points.emplace_back(0.0, 0.0);
+ }
+ else if (degree == 2)
+ {
+ this->unit_support_points.emplace_back(0.0, 0.0, 0.0);
+ this->unit_support_points.emplace_back(1.0, 0.0, 0.0);
+ this->unit_support_points.emplace_back(0.0, 1.0, 0.0);
+ this->unit_support_points.emplace_back(0.0, 0.0, 1.0);
+ this->unit_support_points.emplace_back(0.5, 0.0, 0.0);
+ this->unit_support_points.emplace_back(0.5, 0.5, 0.0);
+ this->unit_support_points.emplace_back(0.0, 0.5, 0.0);
+ this->unit_support_points.emplace_back(0.0, 0.0, 0.5);
+ this->unit_support_points.emplace_back(0.5, 0.0, 0.5);
+ this->unit_support_points.emplace_back(0.0, 0.5, 0.5);
+
+ // TODO
+ this->unit_face_support_points.emplace_back(1.0, 0.0);
+ this->unit_face_support_points.emplace_back(0.0, 1.0);
+ this->unit_face_support_points.emplace_back(0.0, 0.0);
+ this->unit_face_support_points.emplace_back(0.5, 0.5);
+ this->unit_face_support_points.emplace_back(0.0, 0.5);
+ this->unit_face_support_points.emplace_back(0.5, 0.0);
+ }
+ else
+ {
+ Assert(false, ExcNotImplemented());
+ }
+ }
+ else
+ {
+ Assert(false, ExcNotImplemented());
+ }
+ }
+
+
+
+ template <int dim, int spacedim>
+ void
+ FE_Poly<dim, spacedim>::
+ convert_generalized_support_point_values_to_dof_values(
+ const std::vector<Vector<double>> &support_point_values,
+ std::vector<double> & nodal_values) const
+ {
+ AssertDimension(support_point_values.size(),
+ this->get_unit_support_points().size());
+ AssertDimension(support_point_values.size(), nodal_values.size());
+ AssertDimension(this->dofs_per_cell, nodal_values.size());
+
+ for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
+ {
+ AssertDimension(support_point_values[i].size(), 1);
+
+ nodal_values[i] = support_point_values[i](0);
+ }
+ }
+
+
+
+ template <int dim, int spacedim>
+ FE_P<dim, spacedim>::FE_P(const unsigned int degree)
+ : FE_Poly<dim, spacedim>(degree, get_dpo_vector_fe_p(dim, degree))
+ {}
+
+
+
+ template <int dim, int spacedim>
+ std::unique_ptr<FiniteElement<dim, spacedim>>
+ FE_P<dim, spacedim>::clone() const
+ {
+ return std::make_unique<FE_P<dim, spacedim>>(*this);
+ }
+
+
+
+ template <int dim, int spacedim>
+ std::string
+ FE_P<dim, spacedim>::get_name() const
+ {
+ std::ostringstream namebuf;
+ namebuf << "FE_P<" << dim << ">(" << this->degree << ")";
+
+ return namebuf.str();
+ }
+
+
+ template <int dim, int spacedim>
+ FE_DGP<dim, spacedim>::FE_DGP(const unsigned int degree)
+ : FE_Poly<dim, spacedim>(degree, get_dpo_vector_fe_dgp(dim, degree))
+ {}
+
+
+
+ template <int dim, int spacedim>
+ std::unique_ptr<FiniteElement<dim, spacedim>>
+ FE_DGP<dim, spacedim>::clone() const
+ {
+ return std::make_unique<FE_DGP<dim, spacedim>>(*this);
+ }
+
+
+
+ template <int dim, int spacedim>
+ std::string
+ FE_DGP<dim, spacedim>::get_name() const
+ {
+ std::ostringstream namebuf;
+ namebuf << "FE_DGP<" << dim << ">(" << this->degree << ")";
+
+ return namebuf.str();
+ }
+
+
+
+} // namespace Simplex
+
+// explicit instantiations
+#include "fe_lib.inst"
+
+DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : DIMENSIONS)
+ {
+#if deal_II_dimension <= deal_II_space_dimension
+ template class Simplex::FE_Poly<deal_II_dimension, deal_II_space_dimension>;
+ template class Simplex::FE_P<deal_II_dimension, deal_II_space_dimension>;
+ template class Simplex::FE_DGP<deal_II_dimension, deal_II_space_dimension>;
+#endif
+ }
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Test n_dofs_per-methods of Simplex::FE_P and Simplex::FE_DGP.
+
+
+#include <deal.II/simplex/fe_lib.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim, int spacedim>
+void
+test(const FiniteElement<dim, spacedim> &fe)
+{
+ deallog << fe.get_name() << ": " << std::endl;
+
+ deallog << " n_dofs_per_vertex(): " << fe.n_dofs_per_vertex() << std::endl;
+ deallog << " n_dofs_per_line(): " << fe.n_dofs_per_line() << std::endl;
+ deallog << " n_dofs_per_quad(): " << fe.n_dofs_per_quad() << std::endl;
+ deallog << " n_dofs_per_hex(): " << fe.n_dofs_per_hex() << std::endl;
+ deallog << " n_dofs_per_face(): " << fe.n_dofs_per_face() << std::endl;
+ deallog << " n_dofs_per_cell(): " << fe.n_dofs_per_cell() << std::endl;
+ deallog << " tensor_degree(): " << fe.tensor_degree() << std::endl;
+
+ deallog << std::endl;
+}
+
+int
+main()
+{
+ initlog();
+
+ test(Simplex::FE_P<2>(1));
+ test(Simplex::FE_P<2>(2));
+ test(Simplex::FE_P<3>(1));
+ test(Simplex::FE_P<3>(2));
+ test(Simplex::FE_DGP<2>(1));
+ test(Simplex::FE_DGP<2>(2));
+ test(Simplex::FE_DGP<3>(1));
+ test(Simplex::FE_DGP<3>(2));
+}
--- /dev/null
+
+DEAL::FE_P<2>(1):
+DEAL:: n_dofs_per_vertex(): 1
+DEAL:: n_dofs_per_line(): 0
+DEAL:: n_dofs_per_quad(): 0
+DEAL:: n_dofs_per_hex(): 0
+DEAL:: n_dofs_per_face(): 2
+DEAL:: n_dofs_per_cell(): 3
+DEAL:: tensor_degree(): 1
+DEAL::
+DEAL::FE_P<2>(2):
+DEAL:: n_dofs_per_vertex(): 1
+DEAL:: n_dofs_per_line(): 1
+DEAL:: n_dofs_per_quad(): 0
+DEAL:: n_dofs_per_hex(): 0
+DEAL:: n_dofs_per_face(): 3
+DEAL:: n_dofs_per_cell(): 6
+DEAL:: tensor_degree(): 2
+DEAL::
+DEAL::FE_P<3>(1):
+DEAL:: n_dofs_per_vertex(): 1
+DEAL:: n_dofs_per_line(): 0
+DEAL:: n_dofs_per_quad(): 0
+DEAL:: n_dofs_per_hex(): 0
+DEAL:: n_dofs_per_face(): 3
+DEAL:: n_dofs_per_cell(): 4
+DEAL:: tensor_degree(): 1
+DEAL::
+DEAL::FE_P<3>(2):
+DEAL:: n_dofs_per_vertex(): 1
+DEAL:: n_dofs_per_line(): 1
+DEAL:: n_dofs_per_quad(): 0
+DEAL:: n_dofs_per_hex(): 0
+DEAL:: n_dofs_per_face(): 6
+DEAL:: n_dofs_per_cell(): 10
+DEAL:: tensor_degree(): 2
+DEAL::
+DEAL::FE_DGP<2>(1):
+DEAL:: n_dofs_per_vertex(): 0
+DEAL:: n_dofs_per_line(): 0
+DEAL:: n_dofs_per_quad(): 3
+DEAL:: n_dofs_per_hex(): 0
+DEAL:: n_dofs_per_face(): 0
+DEAL:: n_dofs_per_cell(): 3
+DEAL:: tensor_degree(): 1
+DEAL::
+DEAL::FE_DGP<2>(2):
+DEAL:: n_dofs_per_vertex(): 0
+DEAL:: n_dofs_per_line(): 0
+DEAL:: n_dofs_per_quad(): 6
+DEAL:: n_dofs_per_hex(): 0
+DEAL:: n_dofs_per_face(): 0
+DEAL:: n_dofs_per_cell(): 6
+DEAL:: tensor_degree(): 2
+DEAL::
+DEAL::FE_DGP<3>(1):
+DEAL:: n_dofs_per_vertex(): 0
+DEAL:: n_dofs_per_line(): 0
+DEAL:: n_dofs_per_quad(): 0
+DEAL:: n_dofs_per_hex(): 4
+DEAL:: n_dofs_per_face(): 0
+DEAL:: n_dofs_per_cell(): 4
+DEAL:: tensor_degree(): 1
+DEAL::
+DEAL::FE_DGP<3>(2):
+DEAL:: n_dofs_per_vertex(): 0
+DEAL:: n_dofs_per_line(): 0
+DEAL:: n_dofs_per_quad(): 0
+DEAL:: n_dofs_per_hex(): 10
+DEAL:: n_dofs_per_face(): 0
+DEAL:: n_dofs_per_cell(): 10
+DEAL:: tensor_degree(): 2
+DEAL::
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Evaluate Simplex::FE_P and Simplex::FE_DGP at quadrature points.
+
+
+#include <deal.II/simplex/fe_lib.h>
+#include <deal.II/simplex/quadrature_lib.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim, int spacedim>
+void
+test(const FiniteElement<dim, spacedim> &fe, const Quadrature<dim> &quad)
+{
+ deallog.push(fe.get_name());
+
+ for (const auto &point : quad.get_points())
+ {
+ deallog << point << " : " << std::endl;
+ for (unsigned int i = 0; i < fe.n_dofs_per_cell(); i++)
+ deallog << fe.shape_value(i, point) << " ";
+ for (unsigned int i = 0; i < fe.n_dofs_per_cell(); i++)
+ deallog << fe.shape_grad(i, point) << " ";
+ deallog << std::endl;
+ }
+ deallog << std::endl;
+
+ deallog.pop();
+}
+
+int
+main()
+{
+ initlog();
+
+ test(Simplex::FE_P<2>(1), Simplex::PGauss<2>(3));
+ test(Simplex::FE_P<2>(2), Simplex::PGauss<2>(7));
+ test(Simplex::FE_P<3>(1), Simplex::PGauss<3>(4));
+ test(Simplex::FE_P<3>(2), Simplex::PGauss<3>(10));
+ test(Simplex::FE_DGP<2>(1), Simplex::PGauss<2>(3));
+ test(Simplex::FE_DGP<2>(2), Simplex::PGauss<2>(7));
+ test(Simplex::FE_DGP<3>(1), Simplex::PGauss<3>(4));
+ test(Simplex::FE_DGP<3>(2), Simplex::PGauss<3>(10));
+}
--- /dev/null
+
+DEAL:FE_P<2>(1)::0.666667 0.166667 :
+DEAL:FE_P<2>(1)::0.166667 0.666667 0.166667 -1.00000 -1.00000 1.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<2>(1)::0.166667 0.666667 :
+DEAL:FE_P<2>(1)::0.166667 0.166667 0.666667 -1.00000 -1.00000 1.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<2>(1)::0.166667 0.166667 :
+DEAL:FE_P<2>(1)::0.666667 0.166667 0.166667 -1.00000 -1.00000 1.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<2>(1)::
+DEAL:FE_P<2>(2)::0.333333 0.333333 :
+DEAL:FE_P<2>(2)::-0.111111 -0.111111 -0.111111 0.444444 0.444444 0.444444 -0.333333 -0.333333 0.333333 0.00000 0.00000 0.333333 3.99969e-12 -1.33333 1.33333 1.33333 -1.33333 3.99991e-12
+DEAL:FE_P<2>(2)::0.797427 0.101287 :
+DEAL:FE_P<2>(2)::-0.0807686 0.474353 -0.0807686 0.323074 0.323074 0.0410358 0.594854 0.594854 2.18971 0.00000 0.00000 -0.594854 -2.78456 -3.18971 0.405146 3.18971 -0.405146 3.99980e-12
+DEAL:FE_P<2>(2)::0.101287 0.797427 :
+DEAL:FE_P<2>(2)::-0.0807686 -0.0807686 0.474353 0.0410358 0.323074 0.323074 0.594854 0.594854 -0.594854 0.00000 0.00000 2.18971 3.99991e-12 -0.405146 3.18971 0.405146 -3.18971 -2.78456
+DEAL:FE_P<2>(2)::0.101287 0.101287 :
+DEAL:FE_P<2>(2)::0.474353 -0.0807686 -0.0807686 0.323074 0.0410358 0.323074 -2.18971 -2.18971 -0.594854 0.00000 0.00000 -0.594854 2.78456 -0.405146 0.405146 0.405146 -0.405146 2.78456
+DEAL:FE_P<2>(2)::0.0597159 0.470142 :
+DEAL:FE_P<2>(2)::-0.0280749 -0.0525839 -0.0280749 0.112300 0.112300 0.884134 -0.880568 -0.880568 -0.761137 0.00000 0.00000 0.880568 1.64170 -0.238863 1.88057 0.238863 -1.88057 8.00249e-13
+DEAL:FE_P<2>(2)::0.470142 0.0597159 :
+DEAL:FE_P<2>(2)::-0.0280749 -0.0280749 -0.0525839 0.884134 0.112300 0.112300 -0.880568 -0.880568 0.880568 0.00000 0.00000 -0.761137 8.00165e-13 -1.88057 0.238863 1.88057 -0.238863 1.64170
+DEAL:FE_P<2>(2)::0.470142 0.470142 :
+DEAL:FE_P<2>(2)::-0.0525839 -0.0280749 -0.0280749 0.112300 0.884134 0.112300 0.761137 0.761137 0.880568 0.00000 0.00000 0.880568 -1.64170 -1.88057 1.88057 1.88057 -1.88057 -1.64170
+DEAL:FE_P<2>(2)::
+DEAL:FE_P<3>(1)::0.138197 0.138197 0.138197 :
+DEAL:FE_P<3>(1)::0.585410 0.138197 0.138197 0.138197 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<3>(1)::0.585410 0.138197 0.138197 :
+DEAL:FE_P<3>(1)::0.138197 0.585410 0.138197 0.138197 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<3>(1)::0.138197 0.585410 0.138197 :
+DEAL:FE_P<3>(1)::0.138197 0.138197 0.585410 0.138197 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<3>(1)::0.138197 0.138197 0.585410 :
+DEAL:FE_P<3>(1)::0.138197 0.138197 0.138197 0.585410 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_P<3>(1)::
+DEAL:FE_P<3>(2)::0.568431 0.143856 0.143856 :
+DEAL:FE_P<3>(2)::-0.102467 0.0777961 -0.102467 -0.102467 0.327090 0.327090 0.0827787 0.0827787 0.327090 0.0827787 0.424574 0.424574 0.424574 1.27372 0.00000 0.00000 0.00000 -0.424574 0.00000 0.00000 0.00000 -0.424574 -1.69830 -2.27372 -2.27372 0.575426 2.27372 0.00000 -0.575426 -2.22045e-16 -0.575426 -0.575426 -0.575426 -2.22045e-16 0.575426 0.00000 2.27372 0.00000 0.575426 0.575426
+DEAL:FE_P<3>(2)::0.143856 0.143856 0.143856 :
+DEAL:FE_P<3>(2)::0.0777961 -0.102467 -0.102467 -0.102467 0.327090 0.0827787 0.327090 0.327090 0.0827787 0.0827787 -1.27372 -1.27372 -1.27372 -0.424574 0.00000 0.00000 0.00000 -0.424574 0.00000 0.00000 0.00000 -0.424574 1.69830 -0.575426 -0.575426 0.575426 0.575426 0.00000 -0.575426 1.69830 -0.575426 -0.575426 -0.575426 1.69830 0.575426 0.00000 0.575426 0.00000 0.575426 0.575426
+DEAL:FE_P<3>(2)::0.143856 0.143856 0.568431 :
+DEAL:FE_P<3>(2)::-0.102467 -0.102467 -0.102467 0.0777961 0.0827787 0.0827787 0.0827787 0.327090 0.327090 0.327090 0.424574 0.424574 0.424574 -0.424574 0.00000 0.00000 0.00000 -0.424574 0.00000 0.00000 0.00000 1.27372 -3.33067e-16 -0.575426 -0.575426 0.575426 0.575426 0.00000 -0.575426 -3.33067e-16 -0.575426 -2.27372 -2.27372 -1.69830 2.27372 0.00000 0.575426 0.00000 2.27372 0.575426
+DEAL:FE_P<3>(2)::0.143856 0.568431 0.143856 :
+DEAL:FE_P<3>(2)::-0.102467 -0.102467 0.0777961 -0.102467 0.0827787 0.327090 0.327090 0.0827787 0.0827787 0.327090 0.424574 0.424574 0.424574 -0.424574 0.00000 0.00000 0.00000 1.27372 0.00000 0.00000 0.00000 -0.424574 -2.22045e-16 -0.575426 -0.575426 2.27372 0.575426 0.00000 -2.27372 -1.69830 -2.27372 -0.575426 -0.575426 -2.22045e-16 0.575426 0.00000 0.575426 0.00000 0.575426 2.27372
+DEAL:FE_P<3>(2)::0.00000 0.500000 0.500000 :
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 1.00000 1.00000 1.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 2.00000 0.00000 0.00000 -2.00000 -2.00000 -2.00000 -2.00000 -2.00000 -2.00000 2.00000 0.00000 0.00000 0.00000 2.00000 2.00000
+DEAL:FE_P<3>(2)::0.500000 0.00000 0.500000 :
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 1.00000 1.00000 1.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 -2.00000 -2.00000 -2.00000 0.00000 2.00000 0.00000 0.00000 0.00000 0.00000 -2.00000 -2.00000 -2.00000 2.00000 0.00000 2.00000 0.00000 2.00000 0.00000
+DEAL:FE_P<3>(2)::0.500000 0.500000 0.00000 :
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 0.00000 1.00000 1.00000 1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 -2.00000 -2.00000 -2.00000 2.00000 2.00000 0.00000 -2.00000 -2.00000 -2.00000 0.00000 0.00000 0.00000 0.00000 0.00000 2.00000 0.00000 0.00000 2.00000
+DEAL:FE_P<3>(2)::0.500000 0.00000 0.00000 :
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 0.00000 0.00000 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 -2.00000 -2.00000 0.00000 2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 0.00000 2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 0.00000
+DEAL:FE_P<3>(2)::0.00000 0.500000 0.00000 :
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 -1.00000 -1.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 -2.00000 0.00000 -2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 0.00000 0.00000 0.00000 2.00000
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.500000 :
+DEAL:FE_P<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 -1.00000 -1.00000 -1.00000 -1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 2.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 2.00000 0.00000 -2.00000 -2.00000 0.00000 2.00000 0.00000 0.00000 0.00000 2.00000 0.00000
+DEAL:FE_P<3>(2)::
+DEAL:FE_DGP<2>(1)::0.666667 0.166667 :
+DEAL:FE_DGP<2>(1)::0.166667 0.666667 0.166667 -1.00000 -1.00000 1.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<2>(1)::0.166667 0.666667 :
+DEAL:FE_DGP<2>(1)::0.166667 0.166667 0.666667 -1.00000 -1.00000 1.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<2>(1)::0.166667 0.166667 :
+DEAL:FE_DGP<2>(1)::0.666667 0.166667 0.166667 -1.00000 -1.00000 1.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<2>(1)::
+DEAL:FE_DGP<2>(2)::0.333333 0.333333 :
+DEAL:FE_DGP<2>(2)::-0.111111 -0.111111 -0.111111 0.444444 0.444444 0.444444 -0.333333 -0.333333 0.333333 0.00000 0.00000 0.333333 3.99969e-12 -1.33333 1.33333 1.33333 -1.33333 3.99991e-12
+DEAL:FE_DGP<2>(2)::0.797427 0.101287 :
+DEAL:FE_DGP<2>(2)::-0.0807686 0.474353 -0.0807686 0.323074 0.323074 0.0410358 0.594854 0.594854 2.18971 0.00000 0.00000 -0.594854 -2.78456 -3.18971 0.405146 3.18971 -0.405146 3.99980e-12
+DEAL:FE_DGP<2>(2)::0.101287 0.797427 :
+DEAL:FE_DGP<2>(2)::-0.0807686 -0.0807686 0.474353 0.0410358 0.323074 0.323074 0.594854 0.594854 -0.594854 0.00000 0.00000 2.18971 3.99991e-12 -0.405146 3.18971 0.405146 -3.18971 -2.78456
+DEAL:FE_DGP<2>(2)::0.101287 0.101287 :
+DEAL:FE_DGP<2>(2)::0.474353 -0.0807686 -0.0807686 0.323074 0.0410358 0.323074 -2.18971 -2.18971 -0.594854 0.00000 0.00000 -0.594854 2.78456 -0.405146 0.405146 0.405146 -0.405146 2.78456
+DEAL:FE_DGP<2>(2)::0.0597159 0.470142 :
+DEAL:FE_DGP<2>(2)::-0.0280749 -0.0525839 -0.0280749 0.112300 0.112300 0.884134 -0.880568 -0.880568 -0.761137 0.00000 0.00000 0.880568 1.64170 -0.238863 1.88057 0.238863 -1.88057 8.00249e-13
+DEAL:FE_DGP<2>(2)::0.470142 0.0597159 :
+DEAL:FE_DGP<2>(2)::-0.0280749 -0.0280749 -0.0525839 0.884134 0.112300 0.112300 -0.880568 -0.880568 0.880568 0.00000 0.00000 -0.761137 8.00165e-13 -1.88057 0.238863 1.88057 -0.238863 1.64170
+DEAL:FE_DGP<2>(2)::0.470142 0.470142 :
+DEAL:FE_DGP<2>(2)::-0.0525839 -0.0280749 -0.0280749 0.112300 0.884134 0.112300 0.761137 0.761137 0.880568 0.00000 0.00000 0.880568 -1.64170 -1.88057 1.88057 1.88057 -1.88057 -1.64170
+DEAL:FE_DGP<2>(2)::
+DEAL:FE_DGP<3>(1)::0.138197 0.138197 0.138197 :
+DEAL:FE_DGP<3>(1)::0.585410 0.138197 0.138197 0.138197 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<3>(1)::0.585410 0.138197 0.138197 :
+DEAL:FE_DGP<3>(1)::0.138197 0.585410 0.138197 0.138197 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<3>(1)::0.138197 0.585410 0.138197 :
+DEAL:FE_DGP<3>(1)::0.138197 0.138197 0.585410 0.138197 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<3>(1)::0.138197 0.138197 0.585410 :
+DEAL:FE_DGP<3>(1)::0.138197 0.138197 0.138197 0.585410 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000
+DEAL:FE_DGP<3>(1)::
+DEAL:FE_DGP<3>(2)::0.568431 0.143856 0.143856 :
+DEAL:FE_DGP<3>(2)::-0.102467 0.0777961 -0.102467 -0.102467 0.327090 0.327090 0.0827787 0.0827787 0.327090 0.0827787 0.424574 0.424574 0.424574 1.27372 0.00000 0.00000 0.00000 -0.424574 0.00000 0.00000 0.00000 -0.424574 -1.69830 -2.27372 -2.27372 0.575426 2.27372 0.00000 -0.575426 -2.22045e-16 -0.575426 -0.575426 -0.575426 -2.22045e-16 0.575426 0.00000 2.27372 0.00000 0.575426 0.575426
+DEAL:FE_DGP<3>(2)::0.143856 0.143856 0.143856 :
+DEAL:FE_DGP<3>(2)::0.0777961 -0.102467 -0.102467 -0.102467 0.327090 0.0827787 0.327090 0.327090 0.0827787 0.0827787 -1.27372 -1.27372 -1.27372 -0.424574 0.00000 0.00000 0.00000 -0.424574 0.00000 0.00000 0.00000 -0.424574 1.69830 -0.575426 -0.575426 0.575426 0.575426 0.00000 -0.575426 1.69830 -0.575426 -0.575426 -0.575426 1.69830 0.575426 0.00000 0.575426 0.00000 0.575426 0.575426
+DEAL:FE_DGP<3>(2)::0.143856 0.143856 0.568431 :
+DEAL:FE_DGP<3>(2)::-0.102467 -0.102467 -0.102467 0.0777961 0.0827787 0.0827787 0.0827787 0.327090 0.327090 0.327090 0.424574 0.424574 0.424574 -0.424574 0.00000 0.00000 0.00000 -0.424574 0.00000 0.00000 0.00000 1.27372 -3.33067e-16 -0.575426 -0.575426 0.575426 0.575426 0.00000 -0.575426 -3.33067e-16 -0.575426 -2.27372 -2.27372 -1.69830 2.27372 0.00000 0.575426 0.00000 2.27372 0.575426
+DEAL:FE_DGP<3>(2)::0.143856 0.568431 0.143856 :
+DEAL:FE_DGP<3>(2)::-0.102467 -0.102467 0.0777961 -0.102467 0.0827787 0.327090 0.327090 0.0827787 0.0827787 0.327090 0.424574 0.424574 0.424574 -0.424574 0.00000 0.00000 0.00000 1.27372 0.00000 0.00000 0.00000 -0.424574 -2.22045e-16 -0.575426 -0.575426 2.27372 0.575426 0.00000 -2.27372 -1.69830 -2.27372 -0.575426 -0.575426 -2.22045e-16 0.575426 0.00000 0.575426 0.00000 0.575426 2.27372
+DEAL:FE_DGP<3>(2)::0.00000 0.500000 0.500000 :
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 1.00000 1.00000 1.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 2.00000 0.00000 0.00000 -2.00000 -2.00000 -2.00000 -2.00000 -2.00000 -2.00000 2.00000 0.00000 0.00000 0.00000 2.00000 2.00000
+DEAL:FE_DGP<3>(2)::0.500000 0.00000 0.500000 :
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 1.00000 1.00000 1.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 -2.00000 -2.00000 -2.00000 0.00000 2.00000 0.00000 0.00000 0.00000 0.00000 -2.00000 -2.00000 -2.00000 2.00000 0.00000 2.00000 0.00000 2.00000 0.00000
+DEAL:FE_DGP<3>(2)::0.500000 0.500000 0.00000 :
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 0.00000 1.00000 1.00000 1.00000 1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 -2.00000 -2.00000 -2.00000 2.00000 2.00000 0.00000 -2.00000 -2.00000 -2.00000 0.00000 0.00000 0.00000 0.00000 0.00000 2.00000 0.00000 0.00000 2.00000
+DEAL:FE_DGP<3>(2)::0.500000 0.00000 0.00000 :
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 0.00000 0.00000 -1.00000 -1.00000 -1.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 -2.00000 -2.00000 0.00000 2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 0.00000 2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 0.00000
+DEAL:FE_DGP<3>(2)::0.00000 0.500000 0.00000 :
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 -1.00000 -1.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 0.00000 -1.00000 2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 -2.00000 0.00000 -2.00000 0.00000 0.00000 2.00000 0.00000 0.00000 0.00000 0.00000 0.00000 2.00000
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.500000 :
+DEAL:FE_DGP<3>(2)::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 1.00000 0.00000 0.00000 -1.00000 -1.00000 -1.00000 -1.00000 0.00000 0.00000 0.00000 -1.00000 0.00000 0.00000 0.00000 1.00000 2.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 2.00000 0.00000 -2.00000 -2.00000 0.00000 2.00000 0.00000 0.00000 0.00000 2.00000 0.00000
+DEAL:FE_DGP<3>(2)::