From: Peter Munch Date: Sun, 28 May 2023 10:23:04 +0000 (+0200) Subject: MGTwoLevelTransferNonNested: enable simplex X-Git-Tag: relicensing~554^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=3463a43b7e0f722e09417dc2fd963f6525582c11;p=dealii.git MGTwoLevelTransferNonNested: enable simplex --- diff --git a/doc/news/changes/minor/20230802Munch b/doc/news/changes/minor/20230802Munch new file mode 100644 index 0000000000..d334e08bd6 --- /dev/null +++ b/doc/news/changes/minor/20230802Munch @@ -0,0 +1,3 @@ +New: MGTwoLevelTransferNonNested now also supports FE_SimplexP. +
+(Peter Munch, Marco Feder, 2023/08/02) diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h index aad5cb6f42..3369e456f6 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h @@ -31,6 +31,7 @@ #include #include +#include #include #include #include @@ -4399,9 +4400,12 @@ namespace internal unit_points.end()); } + typename NonMatching::MappingInfo::AdditionalData ad; + ad.store_cells = true; + auto mapping_info = std::make_shared>( - rpe.get_mapping(), update_values); + rpe.get_mapping(), update_values, ad); mapping_info->reinit_cells(cell_iterators, unit_points_vector); return mapping_info; @@ -4436,18 +4440,26 @@ namespace internal Assert((dynamic_cast *>( &dof_handler.get_fe().base_element(0)) != nullptr) || (dynamic_cast *>( + &dof_handler.get_fe().base_element(0)) != nullptr) || + (dynamic_cast *>( + &dof_handler.get_fe().base_element(0)) != nullptr) || + (dynamic_cast *>( &dof_handler.get_fe().base_element(0)) != nullptr), + ExcMessage("Function expects FE_DGQ, FE_Q, FE_SimplexP, or " + "FE_SimplexDGP in dof_handler.")); + + Assert((dynamic_cast *>( + &dof_handler_support_points.get_fe()) != nullptr || + dynamic_cast *>( + &dof_handler_support_points.get_fe()) != nullptr) || + ((dynamic_cast *>( + &dof_handler_support_points.get_fe()) != nullptr || + dynamic_cast *>( + &dof_handler_support_points.get_fe()) != nullptr) && + dof_handler_support_points.get_fe().degree == 0), ExcMessage( - "Function expects FE_DGQ of FE_Q elements in dof_handler.")); - - Assert( - (dynamic_cast *>( - &dof_handler_support_points.get_fe()) != nullptr) || - ((dynamic_cast *>( - &dof_handler_support_points.get_fe()) != nullptr) && - dof_handler_support_points.get_fe().degree == 0), - ExcMessage( - "Function expects FE_DGQ&°ree==0 or FE_Q in dof_handler_support_points.")); + "Function expects (FE_DGQ||FE_SimplexDGP)&°ree==0 or " + "(FE_Q||FE_SimplexP) in dof_handler_support_points.")); Assert( dof_handler_support_points.get_fe().n_components() == 1, @@ -4859,6 +4871,10 @@ MGTwoLevelTransferNonNested>:: n_components); else if (const auto fe = dynamic_cast *>(&fe_base)) fe_coarse = fe->clone(); + else if (const auto fe = dynamic_cast *>(&fe_base)) + fe_coarse = fe->clone(); + else if (const auto fe = dynamic_cast *>(&fe_base)) + fe_coarse = fe->clone(); else AssertThrow(false, ExcMessage(dof_handler_coarse.get_fe().get_name())); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_05.cc b/tests/multigrid-global-coarsening/mg_transfer_a_05.cc new file mode 100644 index 0000000000..3bba3c378f --- /dev/null +++ b/tests/multigrid-global-coarsening/mg_transfer_a_05.cc @@ -0,0 +1,164 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2020 - 2022 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 prolongation with uniformly refined meshes. + * + * To compare with non_nested_transfer_04. + */ + +#include +#include +#include +#include +#include + +#include + +#include + +#include +#include +#include +#include +#include + +#include +#include + +#include +#include + +#include +#include + +#include "mg_transfer_util.h" + +using namespace dealii; + +template +class RightHandSideFunction : public Function +{ +public: + RightHandSideFunction() + : Function(1) + {} + + virtual double + value(const Point &p, const unsigned int component = 0) const + { + (void)component; + return p[0]; + } +}; + +template +void +do_test(const FiniteElement &fe_fine, + const FiniteElement &fe_coarse, + const Function &function) +{ + // create coarse grid + Triangulation tria_coarse; // TODO + if (fe_coarse.reference_cell().is_simplex()) + GridGenerator::subdivided_hyper_cube_with_simplices(tria_coarse, 1); + else + GridGenerator::hyper_cube(tria_coarse); + tria_coarse.refine_global(); + + // create fine grid + Triangulation tria_fine; // TODO + + if (fe_fine.reference_cell().is_simplex()) + GridGenerator::subdivided_hyper_cube_with_simplices(tria_fine, 1); + else + GridGenerator::hyper_cube(tria_fine); + tria_fine.refine_global(); + tria_fine.refine_global(); + + // setup dof-handlers + DoFHandler dof_handler_fine(tria_fine); + dof_handler_fine.distribute_dofs(fe_fine); + + DoFHandler dof_handler_coarse(tria_coarse); + dof_handler_coarse.distribute_dofs(fe_coarse); + + // setup constraint matrix + AffineConstraints constraint_coarse; + constraint_coarse.close(); + + AffineConstraints constraint_fine; + constraint_fine.close(); + + // setup transfer operator + MGTwoLevelTransfer> transfer; + + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); + + test_non_nested_transfer(transfer, + dof_handler_fine, + dof_handler_coarse, + function); +} + +template +void +test(int fe_degree, + const Function &function, + const bool do_simplex_mesh) +{ + const auto str_fine = std::to_string(fe_degree); + const auto str_coarse = std::to_string(fe_degree); + + if ((fe_degree > 0) && (do_simplex_mesh == false)) + { + deallog.push("CG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")"); + do_test(FE_Q(fe_degree), + FE_Q(fe_degree), + function); + deallog.pop(); + } + + if ((fe_degree > 0) && do_simplex_mesh) + { + deallog.push("CG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")"); + do_test(FE_SimplexP(fe_degree), + FE_SimplexP(fe_degree), + function); + deallog.pop(); + } +} + +int +main(int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll all; + + deallog.precision(8); + + // Functions::ConstantFunction<2, double> fu(1.); + RightHandSideFunction<2> fu; + + for (unsigned int i = 0; i <= 4; ++i) + test<2, double>(i, fu, false); + + for (unsigned int i = 0; i <= 2; ++i) + test<2, double>(i, fu, true); +} diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_05.with_p4est=true.mpirun=1.output b/tests/multigrid-global-coarsening/mg_transfer_a_05.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..5c952a2813 --- /dev/null +++ b/tests/multigrid-global-coarsening/mg_transfer_a_05.with_p4est=true.mpirun=1.output @@ -0,0 +1,37 @@ + +DEAL:0:CG<2>(1)<->CG<2>(1)::3.0618622 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.50000000 0.0000000 0.50000000 1.0000000 1.0000000 0.0000000 0.50000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.25000000 0.50000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 0.75000000 1.0000000 1.0000000 0.75000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 1.0000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::4.9702238 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.25000000 0.50000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 0.75000000 1.0000000 1.0000000 0.75000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 1.0000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.18750000 1.5000000 0.25000000 2.0000000 2.0625000 2.7500000 0.18750000 1.5000000 2.0625000 +DEAL:0:CG<2>(2)<->CG<2>(2)::5.3560713 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.50000000 0.0000000 0.50000000 0.0000000 0.50000000 0.25000000 0.25000000 0.25000000 1.0000000 1.0000000 1.0000000 0.75000000 0.75000000 0.75000000 0.0000000 0.50000000 0.0000000 0.50000000 0.25000000 0.25000000 1.0000000 1.0000000 0.75000000 0.75000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.12500000 0.50000000 0.50000000 0.50000000 0.37500000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.75000000 0.62500000 0.62500000 0.62500000 1.0000000 1.0000000 1.0000000 0.87500000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 +DEAL:0:CG<2>(2)<->CG<2>(2)::10.383092 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.12500000 0.50000000 0.50000000 0.50000000 0.37500000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.75000000 0.62500000 0.62500000 0.62500000 1.0000000 1.0000000 1.0000000 0.87500000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.93750000 0.0000000 1.1250000 0.0000000 1.8750000 0.78125000 0.93750000 1.5625000 1.5625000 1.8750000 3.1250000 2.3437500 2.8125000 4.6875000 0.0000000 0.93750000 0.0000000 1.8750000 0.78125000 1.5625000 1.5625000 3.1250000 2.3437500 4.6875000 +DEAL:0:CG<2>(3)<->CG<2>(3)::7.6697458 +DEAL:0:CG<2>(3)<->CG<2>(3)::0.0000000 0.50000000 0.0000000 0.50000000 0.0000000 0.0000000 0.50000000 0.50000000 0.13819660 0.36180340 0.13819660 0.36180340 0.13819660 0.36180340 0.13819660 0.36180340 1.0000000 1.0000000 1.0000000 1.0000000 0.63819660 0.86180340 0.63819660 0.86180340 0.63819660 0.86180340 0.63819660 0.86180340 0.0000000 0.50000000 0.0000000 0.0000000 0.50000000 0.50000000 0.13819660 0.36180340 0.13819660 0.36180340 0.13819660 0.36180340 1.0000000 1.0000000 1.0000000 0.63819660 0.86180340 0.63819660 0.86180340 0.63819660 0.86180340 +DEAL:0:CG<2>(3)<->CG<2>(3)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 +DEAL:0:CG<2>(3)<->CG<2>(3)::15.445755 +DEAL:0:CG<2>(3)<->CG<2>(3)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 +DEAL:0:CG<2>(3)<->CG<2>(3)::-0.021093750 0.70312500 -0.023437500 0.78125000 -0.044531250 -0.044531250 1.4843750 1.4843750 0.41641110 0.91952640 0.46267900 1.0216960 0.87909011 1.9412224 0.87909011 1.9412224 1.2867187 1.4296875 2.7164063 2.7164063 1.7523486 2.2554639 1.9470540 2.5060710 3.6994026 4.7615349 3.6994026 4.7615349 -0.021093750 0.70312500 -0.044531250 -0.044531250 1.4843750 1.4843750 0.41641110 0.91952640 0.87909011 1.9412224 0.87909011 1.9412224 1.2867188 2.7164063 2.7164063 1.7523486 2.2554639 3.6994026 4.7615349 3.6994026 4.7615349 +DEAL:0:CG<2>(4)<->CG<2>(4)::9.9861511 +DEAL:0:CG<2>(4)<->CG<2>(4)::0.0000000 0.50000000 0.0000000 0.50000000 0.0000000 0.0000000 0.0000000 0.50000000 0.50000000 0.50000000 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.0000000 0.50000000 0.0000000 0.0000000 0.0000000 0.50000000 0.50000000 0.50000000 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 1.0000000 1.0000000 1.0000000 1.0000000 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 +DEAL:0:CG<2>(4)<->CG<2>(4)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 +DEAL:0:CG<2>(4)<->CG<2>(4)::20.810567 +DEAL:0:CG<2>(4)<->CG<2>(4)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 +DEAL:0:CG<2>(4)<->CG<2>(4)::-3.4694470e-18 0.85447724 -1.9081958e-17 1.0013652 -1.3877788e-17 -4.8572257e-17 -1.7347235e-17 1.2493373 2.1606744 1.2493373 0.18408253 0.92186504 0.88199239 0.21572703 1.0803372 1.0336103 0.26914839 1.3478655 1.2895674 0.46548041 2.3310746 2.2302506 0.26914839 1.3478655 1.2895674 1.4582719 1.7089545 2.1321498 3.6874601 2.1321498 1.2501574 2.7655951 1.9480673 1.4650644 3.2410116 2.2829476 1.8278642 4.0435965 2.8482832 3.1612114 6.9932239 4.9259816 1.8278642 4.0435965 2.8482832 -1.0408341e-17 0.85447724 -1.3877788e-17 -4.8572257e-17 -1.7347235e-17 1.2493373 2.1606744 1.2493373 0.18408253 0.92186504 0.88199239 0.26914839 1.3478655 1.2895674 0.46548041 2.3310746 2.2302506 0.26914839 1.3478655 1.2895674 1.4582719 2.1321498 3.6874601 2.1321498 1.2501574 2.7655951 1.9480673 1.8278642 4.0435965 2.8482832 3.1612114 6.9932239 4.9259816 1.8278642 4.0435965 2.8482832 +DEAL:0:CG<2>(1)<->CG<2>(1)::3.0618622 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.50000000 0.0000000 1.0000000 0.50000000 0.0000000 1.0000000 0.50000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.50000000 0.25000000 0.0000000 0.75000000 0.50000000 1.0000000 0.75000000 0.50000000 0.25000000 0.0000000 0.25000000 0.0000000 1.0000000 0.75000000 1.0000000 0.50000000 0.75000000 1.0000000 0.25000000 0.50000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::4.9812147 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.50000000 0.25000000 0.0000000 0.75000000 0.50000000 1.0000000 0.75000000 0.50000000 0.25000000 0.0000000 0.25000000 0.0000000 1.0000000 0.75000000 1.0000000 0.50000000 0.75000000 1.0000000 0.25000000 0.50000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.12500000 1.3750000 0.25000000 2.2500000 2.0000000 0.25000000 1.8750000 1.6250000 2.7500000 +DEAL:0:CG<2>(2)<->CG<2>(2)::5.3560713 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.50000000 0.0000000 0.25000000 0.25000000 0.0000000 1.0000000 0.50000000 0.75000000 0.75000000 0.50000000 0.0000000 0.25000000 0.25000000 0.0000000 1.0000000 0.50000000 1.0000000 0.75000000 0.75000000 1.0000000 0.25000000 0.50000000 0.75000000 1.0000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.50000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 1.0000000 0.75000000 0.87500000 0.87500000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.37500000 0.37500000 0.25000000 1.0000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.50000000 0.75000000 0.62500000 0.62500000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.25000000 0.50000000 0.37500000 0.37500000 0.50000000 0.12500000 0.25000000 0.37500000 0.50000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.62500000 0.75000000 0.87500000 1.0000000 0.62500000 0.62500000 0.75000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::10.451291 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.50000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 1.0000000 0.75000000 0.87500000 0.87500000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.37500000 0.37500000 0.25000000 1.0000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.50000000 0.75000000 0.62500000 0.62500000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.25000000 0.50000000 0.37500000 0.37500000 0.50000000 0.12500000 0.25000000 0.37500000 0.50000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.62500000 0.75000000 0.87500000 1.0000000 0.62500000 0.62500000 0.75000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::-0.046875000 0.71875000 -0.18750000 0.84375000 1.2500000 0.18750000 1.3906250 0.50000000 2.7187500 3.7500000 2.5000000 -0.14062500 1.2500000 1.2500000 0.18750000 1.2968750 0.53125000 1.4375000 2.9062500 3.7500000 3.5625000 1.0312500 2.5000000 3.7500000 3.5625000 diff --git a/tests/multigrid-global-coarsening/mg_transfer_util.h b/tests/multigrid-global-coarsening/mg_transfer_util.h index 2c182d628e..b58df2f12d 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_util.h +++ b/tests/multigrid-global-coarsening/mg_transfer_util.h @@ -177,8 +177,7 @@ test_transfer_operator( template void test_non_nested_transfer( - const MGTwoLevelTransferNonNested> + const MGTwoLevelTransferBase> &transfer, const MeshType &dof_handler_fine, const MeshType &dof_handler_coarse, diff --git a/tests/multigrid-global-coarsening/non_nested_multigrid_01.cc b/tests/multigrid-global-coarsening/non_nested_multigrid_01.cc index 393a429ddf..881e02f9f0 100644 --- a/tests/multigrid-global-coarsening/non_nested_multigrid_01.cc +++ b/tests/multigrid-global-coarsening/non_nested_multigrid_01.cc @@ -149,7 +149,7 @@ main(int argc, char **argv) for (unsigned int degree = 1; degree <= 4; ++degree) test<2>(n_refinements, degree, false /*quadrilateral*/); - // for (unsigned int n_refinements = 2; n_refinements <= 4; ++n_refinements) - // for (unsigned int degree = 1; degree <= 2; ++degree) - // test<2>(n_refinements, degree, true /*triangle*/); + for (unsigned int n_refinements = 2; n_refinements <= 4; ++n_refinements) + for (unsigned int degree = 1; degree <= 2; ++degree) + test<2>(n_refinements, degree, true /*triangle*/); } diff --git a/tests/multigrid-global-coarsening/non_nested_multigrid_01.mpirun=1.with_trilinos=true.output b/tests/multigrid-global-coarsening/non_nested_multigrid_01.mpirun=1.with_trilinos=true.output index 6e180875b3..07875b2dbe 100644 --- a/tests/multigrid-global-coarsening/non_nested_multigrid_01.mpirun=1.with_trilinos=true.output +++ b/tests/multigrid-global-coarsening/non_nested_multigrid_01.mpirun=1.with_trilinos=true.output @@ -11,3 +11,9 @@ DEAL:0::2 1 4 quad 3 DEAL:0::2 2 4 quad 3 DEAL:0::2 3 4 quad 3 DEAL:0::2 4 4 quad 3 +DEAL:0::2 1 2 tri 3 +DEAL:0::2 2 2 tri 3 +DEAL:0::2 1 3 tri 3 +DEAL:0::2 2 3 tri 3 +DEAL:0::2 1 4 tri 3 +DEAL:0::2 2 4 tri 3 diff --git a/tests/multigrid-global-coarsening/non_nested_transfer_05.cc b/tests/multigrid-global-coarsening/non_nested_transfer_05.cc new file mode 100644 index 0000000000..97e0313f48 --- /dev/null +++ b/tests/multigrid-global-coarsening/non_nested_transfer_05.cc @@ -0,0 +1,173 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2020 - 2022 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 prolongation with uniformly refined meshes, using the + * non-nested infrastructure. + * + * To compare with mg_transfer_a_05. + */ + +#include +#include +#include +#include +#include + +#include + +#include + +#include +#include +#include +#include +#include + +#include +#include + +#include +#include + +#include +#include + +#include "mg_transfer_util.h" + +using namespace dealii; + +template +class RightHandSideFunction : public Function +{ +public: + RightHandSideFunction() + : Function(1) + {} + + virtual double + value(const Point &p, const unsigned int component = 0) const + { + (void)component; + return p[0]; + } +}; + +template +void +do_test(const FiniteElement &fe_fine, + const FiniteElement &fe_coarse, + const Function &function) +{ + // create coarse grid + Triangulation tria_coarse; // TODO + if (fe_coarse.reference_cell().is_simplex()) + GridGenerator::subdivided_hyper_cube_with_simplices(tria_coarse, 1); + else + GridGenerator::hyper_cube(tria_coarse); + tria_coarse.refine_global(); + + // create fine grid + Triangulation tria_fine; // TODO + + if (fe_fine.reference_cell().is_simplex()) + GridGenerator::subdivided_hyper_cube_with_simplices(tria_fine, 1); + else + GridGenerator::hyper_cube(tria_fine); + tria_fine.refine_global(); + tria_fine.refine_global(); + + // setup dof-handlers + DoFHandler dof_handler_fine(tria_fine); + dof_handler_fine.distribute_dofs(fe_fine); + + DoFHandler dof_handler_coarse(tria_coarse); + dof_handler_coarse.distribute_dofs(fe_coarse); + + // setup constraint matrix + AffineConstraints constraint_coarse; + constraint_coarse.close(); + + AffineConstraints constraint_fine; + constraint_fine.close(); + + // setup transfer operator + MGTwoLevelTransferNonNested> + transfer; + + const auto mapping_fine = + fe_fine.reference_cell().template get_default_mapping(1); + const auto mapping_coarse = + fe_coarse.reference_cell().template get_default_mapping(1); + + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + *mapping_fine, + *mapping_coarse, + constraint_fine, + constraint_coarse); + + test_non_nested_transfer(transfer, + dof_handler_fine, + dof_handler_coarse, + function); +} + +template +void +test(int fe_degree, + const Function &function, + const bool do_simplex_mesh) +{ + const auto str_fine = std::to_string(fe_degree); + const auto str_coarse = std::to_string(fe_degree); + + if ((fe_degree > 0) && (do_simplex_mesh == false)) + { + deallog.push("CG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")"); + do_test(FE_Q(fe_degree), + FE_Q(fe_degree), + function); + deallog.pop(); + } + + if ((fe_degree > 0) && do_simplex_mesh) + { + deallog.push("CG<2>(" + str_fine + ")<->CG<2>(" + str_coarse + ")"); + do_test(FE_SimplexP(fe_degree), + FE_SimplexP(fe_degree), + function); + deallog.pop(); + } +} + +int +main(int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll all; + + deallog.precision(8); + + // Functions::ConstantFunction<2, double> fu(1.); + RightHandSideFunction<2> fu; + + for (unsigned int i = 0; i <= 4; ++i) + test<2, double>(i, fu, false); + + for (unsigned int i = 0; i <= 2; ++i) + test<2, double>(i, fu, true); +} diff --git a/tests/multigrid-global-coarsening/non_nested_transfer_05.with_p4est=true.mpirun=1.output b/tests/multigrid-global-coarsening/non_nested_transfer_05.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..174b4ef411 --- /dev/null +++ b/tests/multigrid-global-coarsening/non_nested_transfer_05.with_p4est=true.mpirun=1.output @@ -0,0 +1,37 @@ + +DEAL:0:CG<2>(1)<->CG<2>(1)::3.0618622 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.50000000 0.0000000 0.50000000 1.0000000 1.0000000 0.0000000 0.50000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.25000000 0.50000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 0.75000000 1.0000000 1.0000000 0.75000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 1.0000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::4.9702238 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.25000000 0.50000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 0.75000000 1.0000000 1.0000000 0.75000000 1.0000000 0.0000000 0.25000000 0.50000000 0.0000000 0.25000000 0.50000000 0.75000000 1.0000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.18750000 1.5000000 0.25000000 2.0000000 2.0625000 2.7500000 0.18750000 1.5000000 2.0625000 +DEAL:0:CG<2>(2)<->CG<2>(2)::5.3560713 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.50000000 0.0000000 0.50000000 0.0000000 0.50000000 0.25000000 0.25000000 0.25000000 1.0000000 1.0000000 1.0000000 0.75000000 0.75000000 0.75000000 0.0000000 0.50000000 0.0000000 0.50000000 0.25000000 0.25000000 1.0000000 1.0000000 0.75000000 0.75000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.12500000 0.50000000 0.50000000 0.50000000 0.37500000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.75000000 0.62500000 0.62500000 0.62500000 1.0000000 1.0000000 1.0000000 0.87500000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 +DEAL:0:CG<2>(2)<->CG<2>(2)::10.383092 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.12500000 0.50000000 0.50000000 0.50000000 0.37500000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.75000000 0.62500000 0.62500000 0.62500000 1.0000000 1.0000000 1.0000000 0.87500000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.0000000 0.25000000 0.0000000 0.25000000 0.12500000 0.12500000 0.50000000 0.50000000 0.37500000 0.37500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 0.75000000 0.75000000 0.62500000 0.62500000 1.0000000 1.0000000 0.87500000 0.87500000 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.93750000 0.0000000 1.1250000 0.0000000 1.8750000 0.78125000 0.93750000 1.5625000 1.5625000 1.8750000 3.1250000 2.3437500 2.8125000 4.6875000 0.0000000 0.93750000 0.0000000 1.8750000 0.78125000 1.5625000 1.5625000 3.1250000 2.3437500 4.6875000 +DEAL:0:CG<2>(3)<->CG<2>(3)::7.6697458 +DEAL:0:CG<2>(3)<->CG<2>(3)::0.0000000 0.50000000 0.0000000 0.50000000 0.0000000 0.0000000 0.50000000 0.50000000 0.13819660 0.36180340 0.13819660 0.36180340 0.13819660 0.36180340 0.13819660 0.36180340 1.0000000 1.0000000 1.0000000 1.0000000 0.63819660 0.86180340 0.63819660 0.86180340 0.63819660 0.86180340 0.63819660 0.86180340 0.0000000 0.50000000 0.0000000 0.0000000 0.50000000 0.50000000 0.13819660 0.36180340 0.13819660 0.36180340 0.13819660 0.36180340 1.0000000 1.0000000 1.0000000 0.63819660 0.86180340 0.63819660 0.86180340 0.63819660 0.86180340 +DEAL:0:CG<2>(3)<->CG<2>(3)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 +DEAL:0:CG<2>(3)<->CG<2>(3)::15.445755 +DEAL:0:CG<2>(3)<->CG<2>(3)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.0000000 0.25000000 0.0000000 0.0000000 0.25000000 0.25000000 0.069098301 0.18090170 0.069098301 0.18090170 0.069098301 0.18090170 0.50000000 0.50000000 0.50000000 0.31909830 0.43090170 0.31909830 0.43090170 0.31909830 0.43090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 0.75000000 0.75000000 0.75000000 0.56909830 0.68090170 0.56909830 0.68090170 0.56909830 0.68090170 1.0000000 1.0000000 1.0000000 0.81909830 0.93090170 0.81909830 0.93090170 0.81909830 0.93090170 +DEAL:0:CG<2>(3)<->CG<2>(3)::-0.021093750 0.70312500 -0.023437500 0.78125000 -0.044531250 -0.044531250 1.4843750 1.4843750 0.41641110 0.91952640 0.46267900 1.0216960 0.87909011 1.9412224 0.87909011 1.9412224 1.2867188 1.4296875 2.7164063 2.7164063 1.7523486 2.2554639 1.9470540 2.5060710 3.6994026 4.7615349 3.6994026 4.7615349 -0.021093750 0.70312500 -0.044531250 -0.044531250 1.4843750 1.4843750 0.41641110 0.91952640 0.87909011 1.9412224 0.87909011 1.9412224 1.2867188 2.7164063 2.7164063 1.7523486 2.2554639 3.6994026 4.7615349 3.6994026 4.7615349 +DEAL:0:CG<2>(4)<->CG<2>(4)::9.9861511 +DEAL:0:CG<2>(4)<->CG<2>(4)::0.0000000 0.50000000 0.0000000 0.50000000 0.0000000 0.0000000 0.0000000 0.50000000 0.50000000 0.50000000 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.0000000 0.50000000 0.0000000 0.0000000 0.0000000 0.50000000 0.50000000 0.50000000 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 0.086336582 0.25000000 0.41366342 1.0000000 1.0000000 1.0000000 1.0000000 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 0.58633658 0.75000000 0.91366342 +DEAL:0:CG<2>(4)<->CG<2>(4)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 +DEAL:0:CG<2>(4)<->CG<2>(4)::20.810567 +DEAL:0:CG<2>(4)<->CG<2>(4)::0.0000000 0.25000000 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.0000000 0.25000000 0.0000000 0.0000000 0.0000000 0.25000000 0.25000000 0.25000000 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.043168291 0.12500000 0.20683171 0.50000000 0.50000000 0.50000000 0.50000000 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.29316829 0.37500000 0.45683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.75000000 0.75000000 0.75000000 0.75000000 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 0.54316829 0.62500000 0.70683171 1.0000000 1.0000000 1.0000000 1.0000000 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 0.79316829 0.87500000 0.95683171 +DEAL:0:CG<2>(4)<->CG<2>(4)::-1.9515639e-17 0.85447724 -1.7780916e-17 1.0013652 -2.6888214e-17 -5.2909066e-17 -2.3418767e-17 1.2493373 2.1606744 1.2493373 0.18408253 0.92186504 0.88199239 0.21572703 1.0803372 1.0336103 0.26914839 1.3478655 1.2895674 0.46548041 2.3310746 2.2302506 0.26914839 1.3478655 1.2895674 1.4582719 1.7089545 2.1321498 3.6874601 2.1321498 1.2501574 2.7655951 1.9480673 1.4650644 3.2410116 2.2829476 1.8278642 4.0435965 2.8482832 3.1612114 6.9932239 4.9259816 1.8278642 4.0435965 2.8482832 -6.9388939e-18 0.85447724 -2.6020852e-17 -4.9873300e-17 -2.0816682e-17 1.2493373 2.1606744 1.2493373 0.18408253 0.92186504 0.88199239 0.26914839 1.3478655 1.2895674 0.46548041 2.3310746 2.2302506 0.26914839 1.3478655 1.2895674 1.4582719 2.1321498 3.6874601 2.1321498 1.2501574 2.7655951 1.9480673 1.8278642 4.0435965 2.8482832 3.1612114 6.9932239 4.9259816 1.8278642 4.0435965 2.8482832 +DEAL:0:CG<2>(1)<->CG<2>(1)::3.0618622 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.50000000 0.0000000 1.0000000 0.50000000 0.0000000 1.0000000 0.50000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.50000000 0.25000000 0.0000000 0.75000000 0.50000000 1.0000000 0.75000000 0.50000000 0.25000000 0.0000000 0.25000000 0.0000000 1.0000000 0.75000000 1.0000000 0.50000000 0.75000000 1.0000000 0.25000000 0.50000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::4.9812147 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.0000000 0.25000000 0.0000000 0.50000000 0.25000000 0.0000000 0.75000000 0.50000000 1.0000000 0.75000000 0.50000000 0.25000000 0.0000000 0.25000000 0.0000000 1.0000000 0.75000000 1.0000000 0.50000000 0.75000000 1.0000000 0.25000000 0.50000000 0.75000000 1.0000000 +DEAL:0:CG<2>(1)<->CG<2>(1)::0.12500000 1.3750000 0.25000000 2.2500000 2.0000000 0.25000000 1.8750000 1.6250000 2.7500000 +DEAL:0:CG<2>(2)<->CG<2>(2)::5.3560713 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.50000000 0.0000000 0.25000000 0.25000000 0.0000000 1.0000000 0.50000000 0.75000000 0.75000000 0.50000000 0.0000000 0.25000000 0.25000000 0.0000000 1.0000000 0.50000000 1.0000000 0.75000000 0.75000000 1.0000000 0.25000000 0.50000000 0.75000000 1.0000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.50000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 1.0000000 0.75000000 0.87500000 0.87500000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.37500000 0.37500000 0.25000000 1.0000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.50000000 0.75000000 0.62500000 0.62500000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.25000000 0.50000000 0.37500000 0.37500000 0.50000000 0.12500000 0.25000000 0.37500000 0.50000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.62500000 0.75000000 0.87500000 1.0000000 0.62500000 0.62500000 0.75000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::10.451291 +DEAL:0:CG<2>(2)<->CG<2>(2)::0.0000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.50000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 1.0000000 0.75000000 0.87500000 0.87500000 0.75000000 0.50000000 0.62500000 0.62500000 0.50000000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.25000000 0.37500000 0.37500000 0.25000000 0.0000000 0.12500000 0.12500000 0.0000000 0.37500000 0.37500000 0.25000000 1.0000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.50000000 0.75000000 0.62500000 0.62500000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.25000000 0.50000000 0.37500000 0.37500000 0.50000000 0.12500000 0.25000000 0.37500000 0.50000000 0.75000000 1.0000000 0.87500000 0.87500000 1.0000000 0.62500000 0.75000000 0.87500000 1.0000000 0.62500000 0.62500000 0.75000000 +DEAL:0:CG<2>(2)<->CG<2>(2)::-0.046875000 0.71875000 -0.18750000 0.84375000 1.2500000 0.18750000 1.3906250 0.50000000 2.7187500 3.7500000 2.5000000 -0.14062500 1.2500000 1.2500000 0.18750000 1.2968750 0.53125000 1.4375000 2.9062500 3.7500000 3.5625000 1.0312500 2.5000000 3.7500000 3.5625000