From c3765195dfa96f6e9e527377c35ded88a3d9002d Mon Sep 17 00:00:00 2001 From: Peter Munch Date: Sun, 7 Nov 2021 09:27:28 +0100 Subject: [PATCH] Enabel FETools::get_projection_matrix() for multiple components --- include/deal.II/fe/fe_tools.templates.h | 34 +- .../fe/get_projection_matrix_fe_system_01.cc | 113 ++++++ .../get_projection_matrix_fe_system_01.output | 347 ++++++++++++++++++ 3 files changed, 475 insertions(+), 19 deletions(-) create mode 100644 tests/fe/get_projection_matrix_fe_system_01.cc create mode 100644 tests/fe/get_projection_matrix_fe_system_01.output diff --git a/include/deal.II/fe/fe_tools.templates.h b/include/deal.II/fe/fe_tools.templates.h index 03fccbee81..c8f7678dc6 100644 --- a/include/deal.II/fe/fe_tools.templates.h +++ b/include/deal.II/fe/fe_tools.templates.h @@ -1454,19 +1454,19 @@ namespace FETools const FiniteElement &fe2, FullMatrix & matrix) { - Assert(fe1.n_components() == 1, ExcNotImplemented()); - Assert(fe1.n_components() == fe2.n_components(), - ExcDimensionMismatch(fe1.n_components(), fe2.n_components())); Assert(matrix.m() == fe2.n_dofs_per_cell() && matrix.n() == fe1.n_dofs_per_cell(), ExcMatrixDimensionMismatch(matrix.m(), matrix.n(), fe2.n_dofs_per_cell(), fe1.n_dofs_per_cell())); + AssertDimension(fe1.n_components(), fe2.n_components()); + matrix = 0; const unsigned int n1 = fe1.n_dofs_per_cell(); const unsigned int n2 = fe2.n_dofs_per_cell(); + const unsigned int nd = fe1.n_components(); const ReferenceCell reference_cell = fe1.reference_cell(); @@ -1486,6 +1486,8 @@ namespace FETools const auto quadrature = reference_cell.get_gauss_type_quadrature(degree + 1); + const unsigned int nq = quadrature.size(); + // Set up FEValues. const UpdateFlags flags = update_values | update_quadrature_points | update_JxW_values; @@ -1498,16 +1500,13 @@ namespace FETools // Integrate and invert mass matrix. This happens in the target space FullMatrix mass(n2, n2); - for (unsigned int k = 0; k < quadrature.size(); ++k) - { - const double dx = val2.JxW(k); - for (unsigned int i = 0; i < n2; ++i) - { - const double v = val2.shape_value(i, k); - for (unsigned int j = 0; j < n2; ++j) - mass(i, j) += v * val2.shape_value(j, k) * dx; - } - } + for (unsigned int i = 0; i < n2; ++i) + for (unsigned int j = 0; j < n2; ++j) + for (unsigned int d = 0; d < nd; ++d) + for (unsigned int k = 0; k < nq; ++k) + mass(i, j) += val2.JxW(k) * val2.shape_value_component(i, k, d) * + val2.shape_value_component(j, k, d); + // Invert the matrix. Gauss-Jordan should be sufficient since we expect // the mass matrix to be well-conditioned mass.gauss_jordan(); @@ -1522,12 +1521,9 @@ namespace FETools b = 0.; for (unsigned int i = 0; i < n2; ++i) for (unsigned int k = 0; k < quadrature.size(); ++k) - { - const double dx = val2.JxW(k); - const double u = val1.shape_value(j, k); - const double v = val2.shape_value(i, k); - b(i) += u * v * dx; - } + for (unsigned int d = 0; d < nd; ++d) + b(i) += val1.shape_value_component(j, k, d) * + val2.shape_value_component(i, k, d) * val2.JxW(k); // Multiply by the inverse mass.vmult(x, b); diff --git a/tests/fe/get_projection_matrix_fe_system_01.cc b/tests/fe/get_projection_matrix_fe_system_01.cc new file mode 100644 index 0000000000..f32e0402dd --- /dev/null +++ b/tests/fe/get_projection_matrix_fe_system_01.cc @@ -0,0 +1,113 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2021 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 FETools::get_projection_matrix() for FESystem, FE_RaviartThomas, +// and FE_Nedelec. + +#include +#include +#include +#include +#include +#include +#include + +#include + +#include "../tests.h" + +template +void +test(const FiniteElement &fe_0, + const FiniteElement &fe_1, + const bool do_scalar_test) +{ + { + FullMatrix matrix(fe_1.n_dofs_per_cell(), fe_0.n_dofs_per_cell()); + FETools::get_projection_matrix(fe_0, fe_1, matrix); + matrix.print_formatted(deallog.get_file_stream(), 8, false, 5, "0"); + deallog << std::endl; + + if (do_scalar_test) + { + FullMatrix matrix_1(fe_1.base_element(0).n_dofs_per_cell(), + fe_0.base_element(0).n_dofs_per_cell()); + FETools::get_projection_matrix(fe_0.base_element(0), + fe_1.base_element(0), + matrix_1); + + for (unsigned i = 0; i < matrix.m(); ++i) + for (unsigned j = 0; j < matrix.n(); ++j) + { + const auto component_index_i = fe_1.system_to_component_index(i); + const auto component_index_j = fe_0.system_to_component_index(j); + + if (component_index_i.first != component_index_j.first) + { + Assert(std::abs(matrix[i][j]) < 1e-8, ExcInternalError()); + } + else + { + Assert(std::abs(matrix[i][j] - + matrix_1[component_index_i.second] + [component_index_j.second]) < 1e-8, + ExcInternalError()); + } + } + } + } +} + + +template +void +test(const unsigned int fe_degree_0, const unsigned int fe_degree_1) +{ + test(FESystem(FE_Q(fe_degree_0), dim), + FESystem(FE_Q(fe_degree_1), dim), + true); + test(FESystem(FE_Q(fe_degree_1), dim), + FESystem(FE_Q(fe_degree_0), dim), + true); + + test(FE_RaviartThomas(fe_degree_0), + FE_RaviartThomas(fe_degree_1), + false); + test(FE_RaviartThomas(fe_degree_1), + FE_RaviartThomas(fe_degree_0), + false); + + test(FE_RaviartThomasNodal(fe_degree_0), + FE_RaviartThomas(fe_degree_1), + false); + test(FE_RaviartThomasNodal(fe_degree_1), + FE_RaviartThomas(fe_degree_0), + false); + + test(FE_Nedelec(fe_degree_0), FE_Nedelec(fe_degree_1), false); + test(FE_Nedelec(fe_degree_1), FE_Nedelec(fe_degree_0), false); +} + + +int +main() +{ + initlog(); + + for (unsigned int i = 1; i <= 3; ++i) + for (unsigned int j = i + i; j <= 3; ++j) + test<2>(i, j); +} diff --git a/tests/fe/get_projection_matrix_fe_system_01.output b/tests/fe/get_projection_matrix_fe_system_01.output new file mode 100644 index 0000000000..19c2f2d36c --- /dev/null +++ b/tests/fe/get_projection_matrix_fe_system_01.output @@ -0,0 +1,347 @@ + +1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 +0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 +0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 +0 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 +0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 0 +0 0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 +0.50000000 0 0.00000000 0 0.50000000 0 0.00000000 0 +0 0.50000000 0 0.00000000 0 0.50000000 0 0.00000000 +0.00000000 0 0.50000000 0 0.00000000 0 0.50000000 0 +0 0.00000000 0 0.50000000 0 0.00000000 0 0.50000000 +0.50000000 0 0.50000000 0 0.00000000 0 0.00000000 0 +0 0.50000000 0 0.50000000 0 0.00000000 0 0.00000000 +0.00000000 0 0.00000000 0 0.50000000 0 0.50000000 0 +0 0.00000000 0 0.00000000 0 0.50000000 0 0.50000000 +0.25000000 0 0.25000000 0 0.25000000 0 0.25000000 0 +0 0.25000000 0 0.25000000 0 0.25000000 0 0.25000000 +DEAL:: +0.44444444 0 -0.22222222 0 -0.22222222 0 0.11111111 0 0.44444444 0 -0.22222222 0 0.44444444 0 -0.22222222 0 0.44444444 0 +0 0.44444444 0 -0.22222222 0 -0.22222222 0 0.11111111 0 0.44444444 0 -0.22222222 0 0.44444444 0 -0.22222222 0 0.44444444 +-0.22222222 0 0.44444444 0 0.11111111 0 -0.22222222 0 -0.22222222 0 0.44444444 0 0.44444444 0 -0.22222222 0 0.44444444 0 +0 -0.22222222 0 0.44444444 0 0.11111111 0 -0.22222222 0 -0.22222222 0 0.44444444 0 0.44444444 0 -0.22222222 0 0.44444444 +-0.22222222 0 0.11111111 0 0.44444444 0 -0.22222222 0 0.44444444 0 -0.22222222 0 -0.22222222 0 0.44444444 0 0.44444444 0 +0 -0.22222222 0 0.11111111 0 0.44444444 0 -0.22222222 0 0.44444444 0 -0.22222222 0 -0.22222222 0 0.44444444 0 0.44444444 +0.11111111 0 -0.22222222 0 -0.22222222 0 0.44444444 0 -0.22222222 0 0.44444444 0 -0.22222222 0 0.44444444 0 0.44444444 0 +0 0.11111111 0 -0.22222222 0 -0.22222222 0 0.44444444 0 -0.22222222 0 0.44444444 0 -0.22222222 0 0.44444444 0 0.44444444 +DEAL:: +1.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 1.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 1.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 1.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 1.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 1.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 1.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 1.00000000 0 0.00000000 +-0.28867513 0.00000000 0.28867513 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 1.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 1.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 -0.28867513 0.00000000 0.28867513 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 -0.28867513 0.00000000 0.28867513 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 -0.28867513 0.00000000 0.28867513 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +DEAL:: +0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0 -1.73205081 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0.00000000 0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0 0 0 0 0 0 0.00000000 0 0.00000000 0 0.00000000 0 -1.73205081 0 0.00000000 0 0.00000000 0 +0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0 1.73205081 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0.00000000 0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0 0 0 0 0 0 0.00000000 0 0.00000000 0 0.00000000 0 1.73205081 0 0.00000000 0 0.00000000 0 +0 0 0 0 0 0 0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 -1.73205081 0 0.00000000 0 0.00000000 +0 0 0 0 0 0 0.00000000 0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 -1.73205081 0 0.00000000 +0 0 0 0 0 0 0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 1.73205081 0 0.00000000 0 0.00000000 +0 0 0 0 0 0 0.00000000 0.50000000 0.00000000 0.00000000 0.50000000 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 1.73205081 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 +DEAL:: +0.50000000 0.50000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +-0.28867513 0.28867513 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.50000000 0.50000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 -0.28867513 0.28867513 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.50000000 0.50000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 -0.28867513 0.28867513 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.50000000 0.50000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 -0.28867513 0.28867513 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.08333333 0.08333333 0.08333333 0.08333333 0 0 0 0 0.33333333 0.33333333 0 0 +0 0 0 0 0.08333333 0.08333333 0.08333333 0.08333333 0 0 0.33333333 0.33333333 +-0.14433757 -0.14433757 0.14433757 0.14433757 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 -0.04811252 0.04811252 -0.04811252 0.04811252 0 0 -0.19245009 0.19245009 +-0.04811252 0.04811252 -0.04811252 0.04811252 0 0 0 0 -0.19245009 0.19245009 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.08333333 -0.08333333 -0.08333333 0.08333333 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 -0.14433757 -0.14433757 0.14433757 0.14433757 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.08333333 -0.08333333 -0.08333333 0.08333333 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +DEAL:: +0.12500000 0.50000000 0.12500000 0.04166667 0.16666667 0.04166667 0 0 0 0 0 0 0.09316950 -0.09316950 0.37267800 -0.37267800 0.09316950 -0.09316950 0 0 0 0 0 0 +-0.21650635 0.00000000 0.21650635 -0.07216878 0.00000000 0.07216878 0 0 0 0 0 0 -0.16137431 0.16137431 0.00000000 0.00000000 0.16137431 -0.16137431 0 0 0 0 0 0 +0.04166667 0.16666667 0.04166667 0.12500000 0.50000000 0.12500000 0 0 0 0 0 0 -0.09316950 0.09316950 -0.37267800 0.37267800 -0.09316950 0.09316950 0 0 0 0 0 0 +-0.07216878 0.00000000 0.07216878 -0.21650635 0.00000000 0.21650635 0 0 0 0 0 0 0.16137431 -0.16137431 0.00000000 0.00000000 -0.16137431 0.16137431 0 0 0 0 0 0 +0 0 0 0 0 0 0.12500000 0.50000000 0.12500000 0.04166667 0.16666667 0.04166667 0 0 0 0 0 0 0.09316950 0.37267800 0.09316950 -0.09316950 -0.37267800 -0.09316950 +0 0 0 0 0 0 -0.21650635 0.00000000 0.21650635 -0.07216878 0.00000000 0.07216878 0 0 0 0 0 0 -0.16137431 0.00000000 0.16137431 0.16137431 0.00000000 -0.16137431 +0 0 0 0 0 0 0.04166667 0.16666667 0.04166667 0.12500000 0.50000000 0.12500000 0 0 0 0 0 0 -0.09316950 -0.37267800 -0.09316950 0.09316950 0.37267800 0.09316950 +0 0 0 0 0 0 -0.07216878 0.00000000 0.07216878 -0.21650635 0.00000000 0.21650635 0 0 0 0 0 0 0.16137431 0.00000000 -0.16137431 -0.16137431 0.00000000 0.16137431 +0.01388889 0.05555556 0.01388889 0.01388889 0.05555556 0.01388889 0 0 0 0 0 0 0.06944444 0.06944444 0.27777778 0.27777778 0.06944444 0.06944444 0 0 0 0 0 0 +0 0 0 0 0 0 0.01388889 0.05555556 0.01388889 0.01388889 0.05555556 0.01388889 0 0 0 0 0 0 0.06944444 0.27777778 0.06944444 0.06944444 0.27777778 0.06944444 +-0.02405626 0.00000000 0.02405626 -0.02405626 0.00000000 0.02405626 0 0 0 0 0 0 -0.12028131 -0.12028131 0.00000000 0.00000000 0.12028131 0.12028131 0 0 0 0 0 0 +0 0 0 0 0 0 -0.02405626 0.00000000 0.02405626 -0.02405626 0.00000000 0.02405626 0 0 0 0 0 0 -0.12028131 0.00000000 0.12028131 -0.12028131 0.00000000 0.12028131 +DEAL:: +1.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 1.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 1.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 1.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 1.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 1.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +DEAL:: +1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.22360680 0.00000000 0.00000000 +0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.22360680 0.00000000 +0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 -0.22360680 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 0.00000000 +0 0 0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.22360680 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 +0 0 0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0.22360680 0.00000000 0.00000000 0 0 0 0 0 0 +0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 +0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 0.00000000 0.00000000 0 0 0 0 0 0 +0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 +0 0 0 0 0 0 0.00000000 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 +0.00000000 0 0.00000000 0.00000000 0 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 +DEAL:: +1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0 1.00000000 0 0.00000000 0 0 0 +0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 +0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 +0 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 +0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 0 +0 0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 +0.72360680 0 0.00000000 0 0.27639320 0 0.00000000 0 +0.27639320 0 0.00000000 0 0.72360680 0 0 0 +0 0.72360680 0 0.00000000 0 0.27639320 0 0.00000000 +0 0.27639320 0 0.00000000 0 0.72360680 0 0 +0.00000000 0 0.72360680 0 0.00000000 0 0.27639320 0 +0 0 0.27639320 0 0.00000000 0 0.72360680 0 +0 0.00000000 0 0.72360680 0 0.00000000 0 0.27639320 +0 0 0 0.27639320 0 0.00000000 0 0.72360680 +0.72360680 0 0.27639320 0 0.00000000 0 0 0 +0.27639320 0 0.72360680 0 0.00000000 0 0.00000000 0 +0 0.72360680 0 0.27639320 0 0.00000000 0 0.00000000 +0 0.27639320 0 0.72360680 0 0.00000000 0 0.00000000 +0.00000000 0 0.00000000 0 0.72360680 0 0.27639320 0 +0.00000000 0 0.00000000 0 0.27639320 0 0.72360680 0 +0 0.00000000 0 0 0 0.72360680 0 0.27639320 +0 0.00000000 0 0.00000000 0 0.27639320 0 0.72360680 +0.52360680 0 0.20000000 0 0.20000000 0 0.07639320 0 +0.20000000 0 0.52360680 0 0.07639320 0 0.20000000 0 +0.20000000 0 0.07639320 0 0.52360680 0 0.20000000 0 +0.07639320 0 0.20000000 0 0.20000000 0 0.52360680 0 +0 0.52360680 0 0.20000000 0 0.20000000 0 0.07639320 +0 0.20000000 0 0.52360680 0 0.07639320 0 0.20000000 +0 0.20000000 0 0.07639320 0 0.52360680 0 0.20000000 +0 0.07639320 0 0.20000000 0 0.20000000 0 0.52360680 +DEAL:: +0.11111111 0 -0.05555556 0 -0.05555556 0 0.02777778 0 0.32522789 -0.04745011 0 0 -0.16261394 0.02372505 0 0 0.32522789 -0.04745011 0 0 -0.16261394 0.02372505 0 0 0.95195861 -0.13888889 -0.13888889 0.02026362 0 0 0 0 +0 0.11111111 0 -0.05555556 0 -0.05555556 0 0.02777778 0 0 0.32522789 -0.04745011 0 0 -0.16261394 0.02372505 0 0 0.32522789 -0.04745011 0 0 -0.16261394 0.02372505 0 0 0 0 0.95195861 -0.13888889 -0.13888889 0.02026362 +-0.05555556 0 0.11111111 0 0.02777778 0 -0.05555556 0 -0.16261394 0.02372505 0 0 0.32522789 -0.04745011 0 0 -0.04745011 0.32522789 0 0 0.02372505 -0.16261394 0 0 -0.13888889 0.95195861 0.02026362 -0.13888889 0 0 0 0 +0 -0.05555556 0 0.11111111 0 0.02777778 0 -0.05555556 0 0 -0.16261394 0.02372505 0 0 0.32522789 -0.04745011 0 0 -0.04745011 0.32522789 0 0 0.02372505 -0.16261394 0 0 0 0 -0.13888889 0.95195861 0.02026362 -0.13888889 +-0.05555556 0 0.02777778 0 0.11111111 0 -0.05555556 0 -0.04745011 0.32522789 0 0 0.02372505 -0.16261394 0 0 -0.16261394 0.02372505 0 0 0.32522789 -0.04745011 0 0 -0.13888889 0.02026362 0.95195861 -0.13888889 0 0 0 0 +0 -0.05555556 0 0.02777778 0 0.11111111 0 -0.05555556 0 0 -0.04745011 0.32522789 0 0 0.02372505 -0.16261394 0 0 -0.16261394 0.02372505 0 0 0.32522789 -0.04745011 0 0 0 0 -0.13888889 0.02026362 0.95195861 -0.13888889 +0.02777778 0 -0.05555556 0 -0.05555556 0 0.11111111 0 0.02372505 -0.16261394 0 0 -0.04745011 0.32522789 0 0 0.02372505 -0.16261394 0 0 -0.04745011 0.32522789 0 0 0.02026362 -0.13888889 -0.13888889 0.95195861 0 0 0 0 +0 0.02777778 0 -0.05555556 0 -0.05555556 0 0.11111111 0 0 0.02372505 -0.16261394 0 0 -0.04745011 0.32522789 0 0 0.02372505 -0.16261394 0 0 -0.04745011 0.32522789 0 0 0 0 0.02026362 -0.13888889 -0.13888889 0.95195861 +DEAL:: +1.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 1.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 1.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 1.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 1.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 1.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 1.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 1.00000000 0 0.00000000 +-0.28867513 0.00000000 0.28867513 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 1.00000000 +0.22360680 0.00000000 0.22360680 0.00000000 0 0 0 0 -0.44721360 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 1.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 -0.28867513 0.00000000 0.28867513 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 -0.28867513 0.00000000 0.28867513 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.22360680 0.00000000 0.22360680 0 0 0 0 0.00000000 0 -0.44721360 0 +0 0 0 0 0.00000000 -0.28867513 0.00000000 0.28867513 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.22360680 0.00000000 0.22360680 0.00000000 0 -0.44721360 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.22360680 0.00000000 0.22360680 0 0.00000000 0 -0.44721360 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0 0.00000000 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 0.00000000 +DEAL:: +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 1.00000000 0 -1.73205081 0 2.23606798 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 0 -1.73205081 0 2.23606798 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 1.00000000 0 1.73205081 0 2.23606798 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 0 1.73205081 0 2.23606798 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 -1.73205081 0 0.00000000 0 0.00000000 0 0.00000000 0 2.23606798 0 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 -1.73205081 0 0.00000000 0 0.00000000 0 0.00000000 0 2.23606798 0 0.00000000 0 0.00000000 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 1.73205081 0 0.00000000 0 0.00000000 0 0.00000000 0 2.23606798 0 0.00000000 0 0.00000000 0 0.00000000 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 1.73205081 0 0.00000000 0 0.00000000 0 0.00000000 0 2.23606798 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0.00000000 0 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0 1.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 0 0.00000000 +DEAL:: +0.50000000 0.50000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +-0.28867513 0.28867513 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.50000000 0.50000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 -0.28867513 0.28867513 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.50000000 0.50000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 -0.28867513 0.28867513 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.50000000 0.50000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 -0.28867513 0.28867513 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.08333333 0.08333333 0.08333333 0.08333333 0 0 0 0 0.33333333 0.33333333 0 0 +0 0 0 0 0.08333333 0.08333333 0.08333333 0.08333333 0 0 0.33333333 0.33333333 +-0.14433757 -0.14433757 0.14433757 0.14433757 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 -0.04811252 0.04811252 -0.04811252 0.04811252 0 0 -0.19245009 0.19245009 +0.07453560 0.07453560 0.07453560 0.07453560 0 0 0 0 -0.14907120 -0.14907120 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +-0.04811252 0.04811252 -0.04811252 0.04811252 0 0 0 0 -0.19245009 0.19245009 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.08333333 -0.08333333 -0.08333333 0.08333333 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 -0.14433757 -0.14433757 0.14433757 0.14433757 0 0 0.00000000 0.00000000 +-0.04303315 0.04303315 -0.04303315 0.04303315 0 0 0 0 0.08606630 -0.08606630 0 0 +0 0 0 0 0.08333333 -0.08333333 -0.08333333 0.08333333 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.07453560 0.07453560 0.07453560 0.07453560 0 0 -0.14907120 -0.14907120 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 -0.04303315 0.04303315 -0.04303315 0.04303315 0 0 0.08606630 -0.08606630 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0.00000000 0.00000000 +DEAL:: +0.03750000 0.18750000 0.18750000 0.03750000 0.01250000 0.06250000 0.06250000 0.01250000 0 0 0 0 0 0 0 0 0.08344171 -0.04444444 -0.00566393 0.41720854 -0.22222222 -0.02831965 0.41720854 -0.22222222 -0.02831965 0.08344171 -0.04444444 -0.00566393 0 0 0 0 0 0 0 0 0 0 0 0 +-0.06495191 -0.14523688 0.14523688 0.06495191 -0.02165064 -0.04841229 0.04841229 0.02165064 0 0 0 0 0 0 0 0 -0.14452528 0.07698004 0.00981022 -0.32316835 0.17213259 0.02193631 0.32316835 -0.17213259 -0.02193631 0.14452528 -0.07698004 -0.00981022 0 0 0 0 0 0 0 0 0 0 0 0 +0.01250000 0.06250000 0.06250000 0.01250000 0.03750000 0.18750000 0.18750000 0.03750000 0 0 0 0 0 0 0 0 -0.00566393 -0.04444444 0.08344171 -0.02831965 -0.22222222 0.41720854 -0.02831965 -0.22222222 0.41720854 -0.00566393 -0.04444444 0.08344171 0 0 0 0 0 0 0 0 0 0 0 0 +-0.02165064 -0.04841229 0.04841229 0.02165064 -0.06495191 -0.14523688 0.14523688 0.06495191 0 0 0 0 0 0 0 0 0.00981022 0.07698004 -0.14452528 0.02193631 0.17213259 -0.32316835 -0.02193631 -0.17213259 0.32316835 -0.00981022 -0.07698004 0.14452528 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.03750000 0.18750000 0.18750000 0.03750000 0.01250000 0.06250000 0.06250000 0.01250000 0 0 0 0 0 0 0 0 0 0 0 0 0.08344171 0.41720854 0.41720854 0.08344171 -0.04444444 -0.22222222 -0.22222222 -0.04444444 -0.00566393 -0.02831965 -0.02831965 -0.00566393 +0 0 0 0 0 0 0 0 -0.06495191 -0.14523688 0.14523688 0.06495191 -0.02165064 -0.04841229 0.04841229 0.02165064 0 0 0 0 0 0 0 0 0 0 0 0 -0.14452528 -0.32316835 0.32316835 0.14452528 0.07698004 0.17213259 -0.17213259 -0.07698004 0.00981022 0.02193631 -0.02193631 -0.00981022 +0 0 0 0 0 0 0 0 0.01250000 0.06250000 0.06250000 0.01250000 0.03750000 0.18750000 0.18750000 0.03750000 0 0 0 0 0 0 0 0 0 0 0 0 -0.00566393 -0.02831965 -0.02831965 -0.00566393 -0.04444444 -0.22222222 -0.22222222 -0.04444444 0.08344171 0.41720854 0.41720854 0.08344171 +0 0 0 0 0 0 0 0 -0.02165064 -0.04841229 0.04841229 0.02165064 -0.06495191 -0.14523688 0.14523688 0.06495191 0 0 0 0 0 0 0 0 0 0 0 0 0.00981022 0.02193631 -0.02193631 -0.00981022 0.07698004 0.17213259 -0.17213259 -0.07698004 -0.14452528 -0.32316835 0.32316835 0.14452528 +0.00416667 0.02083333 0.02083333 0.00416667 0.00416667 0.02083333 0.02083333 0.00416667 0 0 0 0 0 0 0 0 0.02268519 0.02962963 0.02268519 0.11342593 0.14814815 0.11342593 0.11342593 0.14814815 0.11342593 0.02268519 0.02962963 0.02268519 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.00416667 0.02083333 0.02083333 0.00416667 0.00416667 0.02083333 0.02083333 0.00416667 0 0 0 0 0 0 0 0 0 0 0 0 0.02268519 0.11342593 0.11342593 0.02268519 0.02962963 0.14814815 0.14814815 0.02962963 0.02268519 0.11342593 0.11342593 0.02268519 +-0.00721688 -0.01613743 0.01613743 0.00721688 -0.00721688 -0.01613743 0.01613743 0.00721688 0 0 0 0 0 0 0 0 -0.03929189 -0.05132002 -0.03929189 -0.08785934 -0.11475506 -0.08785934 0.08785934 0.11475506 0.08785934 0.03929189 0.05132002 0.03929189 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 -0.00721688 -0.01613743 0.01613743 0.00721688 -0.00721688 -0.01613743 0.01613743 0.00721688 0 0 0 0 0 0 0 0 0 0 0 0 -0.03929189 -0.08785934 0.08785934 0.03929189 -0.05132002 -0.11475506 0.11475506 0.05132002 -0.03929189 -0.08785934 0.08785934 0.03929189 +DEAL:: +1.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 1.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 1.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 1.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 1.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 1.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0.00000000 0.00000000 +DEAL:: +1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.22360680 0.00000000 0.00000000 0.00000000 -0.18898224 0.00000000 0.00000000 0.00000000 +0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.22360680 0.00000000 0.00000000 0.00000000 -0.18898224 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 0.00000000 0.00000000 0.00000000 -0.18898224 0.00000000 0.00000000 0.00000000 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 0.00000000 0.00000000 0.00000000 -0.18898224 0.00000000 0.00000000 +0 0 0 0 0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.22360680 -0.18898224 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.22360680 -0.18898224 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 -0.18898224 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 -0.22360680 -0.18898224 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 -0.65465367 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 +0 0 0 0 0 0 0 0 0.00000000 0 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 1.00000000 0.00000000 -0.65465367 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 +0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 -0.65465367 0.00000000 0.00000000 0.00000000 +0.00000000 0 0.00000000 0.00000000 0.00000000 0 0.00000000 0.00000000 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0.00000000 1.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 0.00000000 -0.65465367 0.00000000 0.00000000 +DEAL:: -- 2.39.5