From 41eb40569acdab09f2f46a9b7b1e9b508d476d7f Mon Sep 17 00:00:00 2001 From: Marc Fehling Date: Fri, 22 Nov 2019 18:06:14 +0100 Subject: [PATCH] Introduced p::d::ErrorPredictor. Predicts error after adaptation based on future_fe_indices and CellStatus. --- include/deal.II/distributed/error_predictor.h | 254 ++++++++++++ include/deal.II/hp/refinement.h | 4 +- source/distributed/CMakeLists.txt | 2 + source/distributed/error_predictor.cc | 360 ++++++++++++++++++ source/distributed/error_predictor.inst.in | 30 ++ tests/mpi/error_predictor_01.cc | 157 ++++++++ ...edictor_01.with_p4est=true.mpirun=1.output | 43 +++ ...edictor_01.with_p4est=true.mpirun=4.output | 67 ++++ tests/mpi/error_predictor_02.cc | 183 +++++++++ ...edictor_02.with_p4est=true.mpirun=1.output | 63 +++ ...edictor_02.with_p4est=true.mpirun=4.output | 87 +++++ 11 files changed, 1248 insertions(+), 2 deletions(-) create mode 100644 include/deal.II/distributed/error_predictor.h create mode 100644 source/distributed/error_predictor.cc create mode 100644 source/distributed/error_predictor.inst.in create mode 100644 tests/mpi/error_predictor_01.cc create mode 100644 tests/mpi/error_predictor_01.with_p4est=true.mpirun=1.output create mode 100644 tests/mpi/error_predictor_01.with_p4est=true.mpirun=4.output create mode 100644 tests/mpi/error_predictor_02.cc create mode 100644 tests/mpi/error_predictor_02.with_p4est=true.mpirun=1.output create mode 100644 tests/mpi/error_predictor_02.with_p4est=true.mpirun=4.output diff --git a/include/deal.II/distributed/error_predictor.h b/include/deal.II/distributed/error_predictor.h new file mode 100644 index 0000000000..474840e13b --- /dev/null +++ b/include/deal.II/distributed/error_predictor.h @@ -0,0 +1,254 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 - 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + +#ifndef dealii_distributed_error_predictor_h +#define dealii_distributed_error_predictor_h + +#include + +#include + +#include + +#include + + +DEAL_II_NAMESPACE_OPEN + + +// forward declarations +#ifndef DOXYGEN +template +class Vector; + +namespace parallel +{ + namespace distributed + { + template + class Triangulation; + } +} // namespace parallel + +namespace hp +{ + template + class DoFHandler; +} +#endif + + +namespace parallel +{ + namespace distributed + { + /** + * Predict how current error indicators will change after refinement and + * coarsening were to happen on the provided hp::DoFHandler in context of a + * parallel::distributed::Triangulation. + * + * This algorithm follows the same logic as the error prediction algorithm + * presented by @cite melenk2001hp and has been implemented in deal.II via + * the function hp::Refinement::predict_error(). See the documentation of + * this particular function for details on the algorithm and all featured + * formulas. + * + * In theory, the combination of hp::Refinement::predict_error() to predict + * how the error will change with parallel::distributed::CellDataTransfer + * to transfer data across meshes forms the intended way to apply the + * prediction algorithm. However, p4est determines the details of grid + * refinement in the parallel distributed case; consequently, it yields more + * reliable and trustworthy results when we determine the predicted errors + * during the adaptation process. This is achieved with this class, similar + * to the parallel::distributed::SolutionTransfer and + * parallel::distributed::CellDataTransfer classes. + * + * Before adaptation happens, one or multiple vectors containing error + * estimates for every cell have to be registered with this class using the + * prepare_for_coarsening_and_refinement() functions. Performing adaptation + * on the triangulation initiates the calculation of the error prediction + * based on how p4est performed grid adaptation. After the refinement + * process, the predicted errors need to be unpacked on the updated mesh via + * one of the unpack() functions. + * + * The following code snippet demonstrates how to impose hp-adaptivity based + * on refinement history in a parallel distributed application: + * @code + * // [initialization...] + * Vector predicted_error_per_cell(triangulation.n_active_cells()); + * for(unsigned int i = 0; i < triangulation.n_active_cells(); ++i) + * predicted_error_per_cell[i] = std::numeric_limits::infinity(); + * + * // [during each refinement step...] + * // set h-adaptivity flags + * Vector estimated_error_per_cell(triangulation.n_active_cells()); + * KellyErrorEstimator::estimate(...); + * parallel::distributed::GridRefinement:: + * refine_and_coarsen_fixed_{number|fraction}(...); + * + * // set p-adaptivity flags + * hp::Refinement::p_adaptivity_from_reference( + * hp_dof_handler, + * estimated_error_per_cell, + * predicted_error_per_cell + * std::less(), + * std::less()); + * hp::Refinement::{choose|force}_p_over_h(hp_dof_handler); + * + * // perform adaptation and predict error for the subsequent adaptation + * parallel::distributed::ErrorPredictor predictor(hp_dof_handler); + * predictor.prepare_for_coarsening_and_refinement(estimated_error_per_cell); + * triangulation.execute_coarsening_and_refinement(); + * + * // unpack predicted errors + * predicted_error_per_cell.reinit(triangulation.n_active_cells()); + * predictor.unpack(predicted_error_per_cell); + * @endcode + * + * For serialization of predicted errors, use the + * parallel::distributed::CellDataTransfer class. + * + * @ingroup distributed + * @author Marc Fehling, 2019 - 2020 + */ + template + class ErrorPredictor + { + public: + /** + * Constructor. + * + * @param[in] dof The hp::DoFHandler on which all operations will + * happen. At the time when this constructor is called, the + * DoFHandler still points to the triangulation before the + * refinement in question happens. + */ + ErrorPredictor(const hp::DoFHandler &dof); + + /** + * Prepare the current object for coarsening and refinement. + * @p all_in includes all vectors that are to be interpolated onto the + * new (refined and/or coarsened) grid. + * + * As specified in the error prediction algorithm in + * hp::Refinement::predict_error(), the control parameters @p gamma_p, + * @p gamma_h, @p gamma_n can be used to influence the predicted errors in + * case of p-adaptation, h-adaptation, and no adapation, respectively. + * Their default values comply to the studies of @cite melenk2001hp. + */ + void + prepare_for_coarsening_and_refinement( + const std::vector *> &all_in, + const double gamma_p = std::sqrt(0.4), + const double gamma_h = 2., + const double gamma_n = 1.); + + /** + * Same as the previous function but for only one discrete function to be + * interpolated. + */ + void + prepare_for_coarsening_and_refinement( + const Vector &in, + const double gamma_p = std::sqrt(0.4), + const double gamma_h = 2., + const double gamma_n = 1.); + + /** + * Interpolate the data previously stored in this object before the mesh + * was refined or coarsened onto the current set of cells. Do so for + * each of the vectors provided to + * prepare_for_coarsening_and_refinement() and write the result into the + * given set of vectors. + */ + void + unpack(std::vector *> &all_out); + + /** + * Same as the previous function. It interpolates only one function. It + * assumes the vectors having the right sizes (i.e. + * in.size()==n_cells_old, out.size()==n_cells_refined) + * + * Multiple calling of this function is NOT allowed. Interpolating + * several functions can be performed in one step by using the above + * function. + */ + void + unpack(Vector &out); + + private: + /** + * Pointer to the degree of freedom handler to work with. + */ + SmartPointer, + ErrorPredictor> + dof_handler; + + /** + * A vector that stores pointers to all the vectors we are supposed to + * copy over from the old to the new mesh. + */ + std::vector *> error_indicators; + + /** + * The handle that the Triangulation has assigned to this object + * with which we can access our memory offset and our pack function. + */ + unsigned int handle; + + /** + * Control parameters for the prediction algorithm. + */ + double gamma_p, gamma_h, gamma_n; + + /** + * A callback function used to pack the data on the current mesh into + * objects that can later be retrieved after refinement, coarsening and + * repartitioning. + */ + std::vector + pack_callback( + const typename Triangulation::cell_iterator &cell, + const typename Triangulation::CellStatus status); + + /** + * A callback function used to unpack the data on the current mesh that + * has been packed up previously on the mesh before refinement, + * coarsening and repartitioning. + */ + void + unpack_callback( + const typename Triangulation::cell_iterator &cell, + const typename Triangulation::CellStatus status, + const boost::iterator_range::const_iterator> + & data_range, + std::vector *> &all_out); + + + /** + * Registers the pack_callback() function to the + * parallel::distributed::Triangulation that has been assigned to the + * DoFHandler class member and stores the returning handle. + */ + void + register_data_attach(); + }; + } // namespace distributed +} // namespace parallel + + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/include/deal.II/hp/refinement.h b/include/deal.II/hp/refinement.h index ba86df3c0d..5cbf2e3c9c 100644 --- a/include/deal.II/hp/refinement.h +++ b/include/deal.II/hp/refinement.h @@ -368,7 +368,7 @@ namespace hp /** * Predict how the current @p error_indicators will adapt after refinement - * and coarsening has happened on the provided @p dof_handler, and write its + * and coarsening were to happen on the provided @p dof_handler, and write its * results to @p predicted_errors. Each entry of @p error_indicators and * @p predicted_errors corresponds to an active cell on the underlying * Triangulation, thus each container has to be of size @@ -517,7 +517,7 @@ namespace hp * // set h-adaptivity flags * Vector estimated_error_per_cell(triangulation.n_active_cells()); * KellyErrorEstimator::estimate(...); - * GridRefinemet::refine_and_coarsen_fixed_fraction(...); + * GridRefinemet::refine_and_coarsen_fixed_{number|fraction}(...); * * // set p-adaptivity flags * hp::Refinement::p_adaptivity_from_reference( diff --git a/source/distributed/CMakeLists.txt b/source/distributed/CMakeLists.txt index a16318931c..57ed62bda1 100644 --- a/source/distributed/CMakeLists.txt +++ b/source/distributed/CMakeLists.txt @@ -19,6 +19,7 @@ SET(_unity_include_src grid_refinement.cc cell_weights.cc cell_data_transfer.cc + error_predictor.cc fully_distributed_tria.cc solution_transfer.cc tria.cc @@ -43,6 +44,7 @@ SET(_inst grid_refinement.inst.in cell_weights.inst.in cell_data_transfer.inst.in + error_predictor.inst.in fully_distributed_tria.inst.in solution_transfer.inst.in tria.inst.in diff --git a/source/distributed/error_predictor.cc b/source/distributed/error_predictor.cc new file mode 100644 index 0000000000..d38da55a55 --- /dev/null +++ b/source/distributed/error_predictor.cc @@ -0,0 +1,360 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 - 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + +#include + +#ifdef DEAL_II_WITH_P4EST + +# include +# include + +# include +# include + +# include +# include + +# include + +# include +# include +# include +# include +# include +# include +# include +# include + +# include +# include + + +DEAL_II_NAMESPACE_OPEN + + +namespace parallel +{ + namespace distributed + { + template + ErrorPredictor::ErrorPredictor( + const hp::DoFHandler &dof) + : dof_handler(&dof, typeid(*this).name()) + , handle(numbers::invalid_unsigned_int) + , gamma_p(0.) + , gamma_h(0.) + , gamma_n(0.) + { + Assert( + (dynamic_cast + *>(&dof_handler->get_triangulation()) != nullptr), + ExcMessage( + "parallel::distributed::ErrorPredictor requires a parallel::distributed::Triangulation object.")); + } + + + + template + void + ErrorPredictor::prepare_for_coarsening_and_refinement( + const std::vector *> &all_in, + const double gamma_p, + const double gamma_h, + const double gamma_n) + { + error_indicators = all_in; + this->gamma_p = gamma_p; + this->gamma_h = gamma_h; + this->gamma_n = gamma_n; + register_data_attach(); + } + + + + template + void + ErrorPredictor::register_data_attach() + { + // TODO: casting away constness is bad + parallel::distributed::Triangulation *tria = + (dynamic_cast *>( + const_cast *>( + &dof_handler->get_triangulation()))); + Assert(tria != nullptr, ExcInternalError()); + + handle = tria->register_data_attach( + [this](const typename Triangulation::cell_iterator &cell, + const typename Triangulation::CellStatus status) { + return this->pack_callback(cell, status); + }, + /*returns_variable_size_data=*/false); + } + + + + template + void + ErrorPredictor::prepare_for_coarsening_and_refinement( + const Vector &in, + const double gamma_p, + const double gamma_h, + const double gamma_n) + { + std::vector *> all_in(1, &in); + prepare_for_coarsening_and_refinement(all_in, gamma_p, gamma_h, gamma_n); + } + + + + template + void + ErrorPredictor::unpack(std::vector *> &all_out) + { + Assert(error_indicators.size() == all_out.size(), + ExcDimensionMismatch(error_indicators.size(), all_out.size())); + + // TODO: casting away constness is bad + parallel::distributed::Triangulation *tria = + (dynamic_cast *>( + const_cast *>( + &dof_handler->get_triangulation()))); + Assert(tria != nullptr, ExcInternalError()); + + tria->notify_ready_to_unpack( + handle, + [this, &all_out]( + const typename Triangulation::cell_iterator &cell, + const typename Triangulation::CellStatus status, + const boost::iterator_range::const_iterator> + &data_range) { + this->unpack_callback(cell, status, data_range, all_out); + }); + + for (const auto &out : all_out) + out->compress(::dealii::VectorOperation::insert); + + error_indicators.clear(); + } + + + + template + void + ErrorPredictor::unpack(Vector &out) + { + std::vector *> all_out(1, &out); + unpack(all_out); + } + + + + template + std::vector + ErrorPredictor::pack_callback( + const typename Triangulation::cell_iterator &cell_, + const typename Triangulation::CellStatus status) + { + typename hp::DoFHandler::cell_iterator cell(*cell_, + dof_handler); + + // create buffer for each individual input vector + std::vector predicted_errors(error_indicators.size()); + + auto predicted_error_it = predicted_errors.begin(); + auto estimated_error_it = error_indicators.cbegin(); + for (; estimated_error_it != error_indicators.cend(); + ++predicted_error_it, ++estimated_error_it) + switch (status) + { + case parallel::distributed::Triangulation::CELL_PERSIST: + { + if (cell->future_fe_index_set()) + { + const int degree_difference = + dof_handler->get_fe_collection()[cell->future_fe_index()] + .degree - + cell->get_fe().degree; + + *predicted_error_it = + (**estimated_error_it)[cell->active_cell_index()] * + std::pow(gamma_p, degree_difference); + } + else + { + *predicted_error_it = + (**estimated_error_it)[cell->active_cell_index()] * + gamma_n; + } + break; + } + + case parallel::distributed::Triangulation::CELL_REFINE: + { + // Determine the exponent by the finite element degree on the + // adapted mesh. + const unsigned int future_fe_degree = + dof_handler->get_fe_collection()[cell->future_fe_index()] + .degree; + + *predicted_error_it = + (**estimated_error_it)[cell->active_cell_index()] * + (gamma_h * std::pow(.5, future_fe_degree)); + + // If the future fe index differs from the active one, also take + // into account p-adaptation. + if (cell->future_fe_index_set()) + *predicted_error_it *= + std::pow(gamma_p, + static_cast(future_fe_degree - + cell->get_fe().degree)); + + break; + } + + case parallel::distributed::Triangulation::CELL_COARSEN: + { + // First figure out which finite element will be assigned to the + // parent cell after h-adaptation analogously to + // dealii::internal::hp::DoFHandlerImplementation::Implementation:: + // collect_fe_indices_on_cells_to_be_refined() + std::set fe_indices_children; + for (unsigned int child_index = 0; + child_index < cell->n_children(); + ++child_index) + { + const auto child = cell->child(child_index); + Assert(child->is_active() && child->coarsen_flag_set(), + typename dealii::Triangulation< + dim>::ExcInconsistentCoarseningFlags()); + + fe_indices_children.insert(child->future_fe_index()); + } + + const unsigned int future_fe_index = + dof_handler->get_fe_collection().find_dominated_fe_extended( + fe_indices_children, /*codim=*/0); + + const unsigned int future_fe_degree = + dof_handler->get_fe_collection()[future_fe_index].degree; + + // Then determine the actually contirbution to the predicted + // error of every single cells that is about to be coarsened. + float sqrsum_of_predicted_errors = 0.; + float predicted_error = 0.; + int degree_difference = 0; + for (unsigned int child_index = 0; + child_index < cell->n_children(); + ++child_index) + { + const auto child = cell->child(child_index); + + predicted_error = + (**estimated_error_it)[child->active_cell_index()] / + (gamma_h * std::pow(.5, future_fe_degree)); + + degree_difference = + future_fe_degree - child->get_fe().degree; + + if (degree_difference != 0) + predicted_error *= std::pow(gamma_p, degree_difference); + + sqrsum_of_predicted_errors += + predicted_error * predicted_error; + } + *predicted_error_it = std::sqrt(sqrsum_of_predicted_errors); + + break; + } + + default: + Assert(false, ExcInternalError()); + break; + } + + // We don't have to pack the whole container if there is just one entry. + if (error_indicators.size() == 1) + return Utilities::pack(predicted_errors[0], + /*allow_compression=*/false); + else + return Utilities::pack(predicted_errors, /*allow_compression=*/false); + } + + + + template + void + ErrorPredictor::unpack_callback( + const typename Triangulation::cell_iterator &cell, + const typename Triangulation::CellStatus status, + const boost::iterator_range::const_iterator> + & data_range, + std::vector *> &all_out) + { + std::vector predicted_errors; + + if (all_out.size() == 1) + predicted_errors.push_back( + Utilities::unpack(data_range.begin(), + data_range.end(), + /*allow_compression=*/false)); + else + predicted_errors = + Utilities::unpack>(data_range.begin(), + data_range.end(), + /*allow_compression=*/false); + + Assert(predicted_errors.size() == all_out.size(), ExcInternalError()); + + auto it_input = predicted_errors.cbegin(); + auto it_output = all_out.begin(); + for (; it_input != predicted_errors.cend(); ++it_input, ++it_output) + switch (status) + { + case parallel::distributed::Triangulation::CELL_PERSIST: + case parallel::distributed::Triangulation::CELL_COARSEN: + (**it_output)[cell->active_cell_index()] = *it_input; + break; + + + case parallel::distributed::Triangulation::CELL_REFINE: + for (unsigned int child_index = 0; + child_index < cell->n_children(); + ++child_index) + (**it_output)[cell->child(child_index)->active_cell_index()] = + (*it_input) / std::sqrt(cell->n_children()); + break; + + default: + Assert(false, ExcInternalError()); + break; + } + } + } // namespace distributed +} // namespace parallel + + +// explicit instantiations +# include "error_predictor.inst" + +DEAL_II_NAMESPACE_CLOSE + +#endif /* DEAL_II_WITH_P4EST */ diff --git a/source/distributed/error_predictor.inst.in b/source/distributed/error_predictor.inst.in new file mode 100644 index 0000000000..25beb3cadf --- /dev/null +++ b/source/distributed/error_predictor.inst.in @@ -0,0 +1,30 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2010 - 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) + { + namespace parallel + \{ + namespace distributed + \{ +#if deal_II_dimension <= deal_II_space_dimension + template class ErrorPredictor; +#endif + \} + \} + } diff --git a/tests/mpi/error_predictor_01.cc b/tests/mpi/error_predictor_01.cc new file mode 100644 index 0000000000..9782e336c8 --- /dev/null +++ b/tests/mpi/error_predictor_01.cc @@ -0,0 +1,157 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 - 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +// Test to check if ErrorPredictor works in parallel with hp::DoFHandler. +// This tests is based on hp/error_prediction.cc + + +#include +#include + +#include + +#include + +#include + +#include + +#include + +#include "../tests.h" + + +template +void +test() +{ + const unsigned int myid = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + const unsigned int n_cells = 4; + + // ------ setup ------ + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + std::vector rep(dim, 1); + rep[0] = n_cells; + Point p1, p2; + for (unsigned int d = 0; d < dim; ++d) + { + p1[d] = 0; + p2[d] = (d == 0) ? n_cells : 1; + } + GridGenerator::subdivided_hyper_rectangle(tria, rep, p1, p2); + + { + auto first = tria.begin(0); + if (first->id().to_string() == "0_0:") + first->set_refine_flag(); + } + tria.execute_coarsening_and_refinement(); + + hp::FECollection fes; + for (unsigned int d = 1; d <= 3; ++d) + fes.push_back(FE_Q(d)); + + hp::DoFHandler dh(tria); + dh.set_fe(fes); + for (auto cell = dh.begin(0); cell != dh.end(0); ++cell) + if (cell->id().to_string() == "0_0:") + { + // h-coarsening + for (unsigned int i = 0; i < cell->n_children(); ++i) + if (cell->child(i)->is_locally_owned()) + cell->child(i)->set_coarsen_flag(); + } + else if (cell->id().to_string() == "1_0:") + { + // h-refinement + if (cell->is_locally_owned()) + cell->set_refine_flag(); + } + else if (cell->id().to_string() == "2_0:") + { + // p-refinement + if (cell->is_locally_owned()) + cell->set_future_fe_index(2); + } + + // ----- prepare error indicators ----- + Vector error_indicators(tria.n_active_cells()); + for (unsigned int i = 0; i < error_indicators.size(); ++i) + error_indicators(i) = 10.; + + // ----- verify ------ + deallog << "pre_adaptation" << std::endl; + for (const auto &cell : dh.active_cell_iterators()) + if (cell->is_locally_owned()) + { + deallog << " cell:" << cell->id().to_string() + << " fe_deg:" << cell->get_fe().degree + << " error:" << error_indicators[cell->active_cell_index()]; + + if (cell->coarsen_flag_set()) + deallog << " coarsening"; + else if (cell->refine_flag_set()) + deallog << " refining"; + + if (cell->future_fe_index_set()) + deallog << " future_fe_deg:" << fes[cell->future_fe_index()].degree; + + deallog << std::endl; + } + + // ----- execute adaptation ----- + parallel::distributed::ErrorPredictor predictor(dh); + + predictor.prepare_for_coarsening_and_refinement(error_indicators, + /*gamma_p=*/0.5, + /*gamma_h=*/1., + /*gamma_n=*/1.); + tria.execute_coarsening_and_refinement(); + + Vector predicted_errors(tria.n_active_cells()); + predictor.unpack(predicted_errors); + + // ------ verify ------ + deallog << "post_adaptation" << std::endl; + for (const auto &cell : dh.active_cell_iterators()) + if (cell->is_locally_owned()) + deallog << " cell:" << cell->id().to_string() + << " predicted:" << predicted_errors(cell->active_cell_index()) + << std::endl; + + // make sure no processor is hanging + MPI_Barrier(MPI_COMM_WORLD); + + deallog << "OK" << std::endl; +} + + + +int +main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll log; + + deallog.push("2d"); + test<2>(); + deallog.pop(); + deallog.push("3d"); + test<3>(); + deallog.pop(); +} diff --git a/tests/mpi/error_predictor_01.with_p4est=true.mpirun=1.output b/tests/mpi/error_predictor_01.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..6b5b6a819d --- /dev/null +++ b/tests/mpi/error_predictor_01.with_p4est=true.mpirun=1.output @@ -0,0 +1,43 @@ + +DEAL:0:2d::pre_adaptation +DEAL:0:2d:: cell:1_0: fe_deg:1 error:10.0000 refining +DEAL:0:2d:: cell:2_0: fe_deg:1 error:10.0000 future_fe_deg:3 +DEAL:0:2d:: cell:3_0: fe_deg:1 error:10.0000 +DEAL:0:2d:: cell:0_1:0 fe_deg:1 error:10.0000 coarsening +DEAL:0:2d:: cell:0_1:1 fe_deg:1 error:10.0000 coarsening +DEAL:0:2d:: cell:0_1:2 fe_deg:1 error:10.0000 coarsening +DEAL:0:2d:: cell:0_1:3 fe_deg:1 error:10.0000 coarsening +DEAL:0:2d::post_adaptation +DEAL:0:2d:: cell:0_0: predicted:40.0000 +DEAL:0:2d:: cell:2_0: predicted:2.50000 +DEAL:0:2d:: cell:3_0: predicted:10.0000 +DEAL:0:2d:: cell:1_1:0 predicted:2.50000 +DEAL:0:2d:: cell:1_1:1 predicted:2.50000 +DEAL:0:2d:: cell:1_1:2 predicted:2.50000 +DEAL:0:2d:: cell:1_1:3 predicted:2.50000 +DEAL:0:2d::OK +DEAL:0:3d::pre_adaptation +DEAL:0:3d:: cell:1_0: fe_deg:1 error:10.0000 refining +DEAL:0:3d:: cell:2_0: fe_deg:1 error:10.0000 future_fe_deg:3 +DEAL:0:3d:: cell:3_0: fe_deg:1 error:10.0000 +DEAL:0:3d:: cell:0_1:0 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:1 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:2 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:3 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:4 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:5 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:6 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d:: cell:0_1:7 fe_deg:1 error:10.0000 coarsening +DEAL:0:3d::post_adaptation +DEAL:0:3d:: cell:0_0: predicted:56.5685 +DEAL:0:3d:: cell:2_0: predicted:2.50000 +DEAL:0:3d:: cell:3_0: predicted:10.0000 +DEAL:0:3d:: cell:1_1:0 predicted:1.76777 +DEAL:0:3d:: cell:1_1:1 predicted:1.76777 +DEAL:0:3d:: cell:1_1:2 predicted:1.76777 +DEAL:0:3d:: cell:1_1:3 predicted:1.76777 +DEAL:0:3d:: cell:1_1:4 predicted:1.76777 +DEAL:0:3d:: cell:1_1:5 predicted:1.76777 +DEAL:0:3d:: cell:1_1:6 predicted:1.76777 +DEAL:0:3d:: cell:1_1:7 predicted:1.76777 +DEAL:0:3d::OK diff --git a/tests/mpi/error_predictor_01.with_p4est=true.mpirun=4.output b/tests/mpi/error_predictor_01.with_p4est=true.mpirun=4.output new file mode 100644 index 0000000000..92cd9ee286 --- /dev/null +++ b/tests/mpi/error_predictor_01.with_p4est=true.mpirun=4.output @@ -0,0 +1,67 @@ + +DEAL:0:2d::pre_adaptation +DEAL:0:2d::post_adaptation +DEAL:0:2d:: cell:0_0: predicted:40.0000 +DEAL:0:2d::OK +DEAL:0:3d::pre_adaptation +DEAL:0:3d::post_adaptation +DEAL:0:3d:: cell:0_0: predicted:56.5685 +DEAL:0:3d::OK + +DEAL:1:2d::pre_adaptation +DEAL:1:2d:: cell:0_1:0 fe_deg:1 error:10.0000 coarsening +DEAL:1:2d:: cell:0_1:1 fe_deg:1 error:10.0000 coarsening +DEAL:1:2d:: cell:0_1:2 fe_deg:1 error:10.0000 coarsening +DEAL:1:2d:: cell:0_1:3 fe_deg:1 error:10.0000 coarsening +DEAL:1:2d::post_adaptation +DEAL:1:2d:: cell:1_1:0 predicted:2.50000 +DEAL:1:2d:: cell:1_1:1 predicted:2.50000 +DEAL:1:2d:: cell:1_1:2 predicted:2.50000 +DEAL:1:2d:: cell:1_1:3 predicted:2.50000 +DEAL:1:2d::OK +DEAL:1:3d::pre_adaptation +DEAL:1:3d:: cell:0_1:0 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:1 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:2 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:3 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:4 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:5 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:6 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d:: cell:0_1:7 fe_deg:1 error:10.0000 coarsening +DEAL:1:3d::post_adaptation +DEAL:1:3d:: cell:1_1:0 predicted:1.76777 +DEAL:1:3d:: cell:1_1:1 predicted:1.76777 +DEAL:1:3d:: cell:1_1:2 predicted:1.76777 +DEAL:1:3d:: cell:1_1:3 predicted:1.76777 +DEAL:1:3d:: cell:1_1:4 predicted:1.76777 +DEAL:1:3d:: cell:1_1:5 predicted:1.76777 +DEAL:1:3d:: cell:1_1:6 predicted:1.76777 +DEAL:1:3d:: cell:1_1:7 predicted:1.76777 +DEAL:1:3d::OK + + +DEAL:2:2d::pre_adaptation +DEAL:2:2d:: cell:1_0: fe_deg:1 error:10.0000 refining +DEAL:2:2d::post_adaptation +DEAL:2:2d::OK +DEAL:2:3d::pre_adaptation +DEAL:2:3d::post_adaptation +DEAL:2:3d::OK + + +DEAL:3:2d::pre_adaptation +DEAL:3:2d:: cell:2_0: fe_deg:1 error:10.0000 future_fe_deg:3 +DEAL:3:2d:: cell:3_0: fe_deg:1 error:10.0000 +DEAL:3:2d::post_adaptation +DEAL:3:2d:: cell:2_0: predicted:2.50000 +DEAL:3:2d:: cell:3_0: predicted:10.0000 +DEAL:3:2d::OK +DEAL:3:3d::pre_adaptation +DEAL:3:3d:: cell:1_0: fe_deg:1 error:10.0000 refining +DEAL:3:3d:: cell:2_0: fe_deg:1 error:10.0000 future_fe_deg:3 +DEAL:3:3d:: cell:3_0: fe_deg:1 error:10.0000 +DEAL:3:3d::post_adaptation +DEAL:3:3d:: cell:2_0: predicted:2.50000 +DEAL:3:3d:: cell:3_0: predicted:10.0000 +DEAL:3:3d::OK + diff --git a/tests/mpi/error_predictor_02.cc b/tests/mpi/error_predictor_02.cc new file mode 100644 index 0000000000..8b3aa9d647 --- /dev/null +++ b/tests/mpi/error_predictor_02.cc @@ -0,0 +1,183 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 - 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +// Test to check if ErrorPredictor works in parallel with hp::DoFHandler. +// This tests is based on hp/error_prediction.cc + + +#include +#include + +#include + +#include + +#include + +#include + +#include + +#include "../tests.h" + + +template +void +test() +{ + const unsigned int myid = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + const unsigned int n_cells = 4; + + // ------ setup ------ + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + std::vector rep(dim, 1); + rep[0] = n_cells; + Point p1, p2; + for (unsigned int d = 0; d < dim; ++d) + { + p1[d] = 0; + p2[d] = (d == 0) ? n_cells : 1; + } + GridGenerator::subdivided_hyper_rectangle(tria, rep, p1, p2); + + for (auto cell = tria.begin(0); cell != tria.end(0); ++cell) + if (cell->id().to_string() == "0_0:" || cell->id().to_string() == "1_0:") + cell->set_refine_flag(); + tria.execute_coarsening_and_refinement(); + + hp::FECollection fes; + for (unsigned int d = 1; d <= 3; ++d) + fes.push_back(FE_Q(d)); + + hp::DoFHandler dh(tria); + dh.set_fe(fes); + for (const auto &cell : dh.active_cell_iterators()) + { + // set active fe index + if (cell->is_locally_owned()) + cell->set_active_fe_index(1); + } + for (auto cell = dh.begin(0); cell != dh.end(0); ++cell) + { + // set refinement/coarsening flags + if (cell->id().to_string() == "0_0:") + { + // h-coarsening and p-refinement + for (unsigned int i = 0; i < cell->n_children(); ++i) + if (cell->child(i)->is_locally_owned()) + { + cell->child(i)->set_coarsen_flag(); + cell->child(i)->set_future_fe_index(2); + } + } + else if (cell->id().to_string() == "1_0:") + { + // h-coarsening and p-coarsening + for (unsigned int i = 0; i < cell->n_children(); ++i) + if (cell->child(i)->is_locally_owned()) + { + cell->child(i)->set_coarsen_flag(); + cell->child(i)->set_future_fe_index(0); + } + } + else if (cell->id().to_string() == "2_0:") + { + // h-refinement and p-refinement + if (cell->is_locally_owned()) + { + cell->set_refine_flag(); + cell->set_future_fe_index(2); + } + } + else if (cell->id().to_string() == "3_0:") + { + // h-refinement and p-coarsening + if (cell->is_locally_owned()) + { + cell->set_refine_flag(); + cell->set_future_fe_index(0); + } + } + } + + // ----- prepare error indicators ----- + Vector error_indicators(tria.n_active_cells()); + for (unsigned int i = 0; i < error_indicators.size(); ++i) + error_indicators(i) = 10.; + + // ----- verify ------ + deallog << "pre_adaptation" << std::endl; + for (const auto &cell : dh.active_cell_iterators()) + if (cell->is_locally_owned()) + { + deallog << " cell:" << cell->id().to_string() + << " fe_deg:" << cell->get_fe().degree + << " error:" << error_indicators[cell->active_cell_index()]; + + if (cell->coarsen_flag_set()) + deallog << " coarsening"; + else if (cell->refine_flag_set()) + deallog << " refining"; + + if (cell->future_fe_index_set()) + deallog << " future_fe_deg:" << fes[cell->future_fe_index()].degree; + + deallog << std::endl; + } + + // ----- execute adaptation ----- + parallel::distributed::ErrorPredictor predictor(dh); + + predictor.prepare_for_coarsening_and_refinement(error_indicators, + /*gamma_p=*/0.5, + /*gamma_h=*/1., + /*gamma_n=*/1.); + tria.execute_coarsening_and_refinement(); + + Vector predicted_errors(tria.n_active_cells()); + predictor.unpack(predicted_errors); + + // ------ verify ------ + deallog << "post_adaptation" << std::endl; + for (const auto &cell : dh.active_cell_iterators()) + if (cell->is_locally_owned()) + deallog << " cell:" << cell->id().to_string() + << " predicted:" << predicted_errors(cell->active_cell_index()) + << std::endl; + + // make sure no processor is hanging + MPI_Barrier(MPI_COMM_WORLD); + + deallog << "OK" << std::endl; +} + + + +int +main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll log; + + deallog.push("2d"); + test<2>(); + deallog.pop(); + deallog.push("3d"); + test<3>(); + deallog.pop(); +} diff --git a/tests/mpi/error_predictor_02.with_p4est=true.mpirun=1.output b/tests/mpi/error_predictor_02.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..a6fa722643 --- /dev/null +++ b/tests/mpi/error_predictor_02.with_p4est=true.mpirun=1.output @@ -0,0 +1,63 @@ + +DEAL:0:2d::pre_adaptation +DEAL:0:2d:: cell:2_0: fe_deg:2 error:10.0000 refining future_fe_deg:3 +DEAL:0:2d:: cell:3_0: fe_deg:2 error:10.0000 refining future_fe_deg:1 +DEAL:0:2d:: cell:0_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:0_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:0_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:0_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:1_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:2d:: cell:1_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:2d:: cell:1_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:2d:: cell:1_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:2d::post_adaptation +DEAL:0:2d:: cell:0_0: predicted:80.0000 +DEAL:0:2d:: cell:1_0: predicted:80.0000 +DEAL:0:2d:: cell:2_1:0 predicted:0.312500 +DEAL:0:2d:: cell:2_1:1 predicted:0.312500 +DEAL:0:2d:: cell:2_1:2 predicted:0.312500 +DEAL:0:2d:: cell:2_1:3 predicted:0.312500 +DEAL:0:2d:: cell:3_1:0 predicted:5.00000 +DEAL:0:2d:: cell:3_1:1 predicted:5.00000 +DEAL:0:2d:: cell:3_1:2 predicted:5.00000 +DEAL:0:2d:: cell:3_1:3 predicted:5.00000 +DEAL:0:2d::OK +DEAL:0:3d::pre_adaptation +DEAL:0:3d:: cell:2_0: fe_deg:2 error:10.0000 refining future_fe_deg:3 +DEAL:0:3d:: cell:3_0: fe_deg:2 error:10.0000 refining future_fe_deg:1 +DEAL:0:3d:: cell:0_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:4 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:5 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:6 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:7 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:1_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:4 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:5 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:6 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d:: cell:1_1:7 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:0:3d::post_adaptation +DEAL:0:3d:: cell:0_0: predicted:113.137 +DEAL:0:3d:: cell:1_0: predicted:113.137 +DEAL:0:3d:: cell:2_1:0 predicted:0.220971 +DEAL:0:3d:: cell:2_1:1 predicted:0.220971 +DEAL:0:3d:: cell:2_1:2 predicted:0.220971 +DEAL:0:3d:: cell:2_1:3 predicted:0.220971 +DEAL:0:3d:: cell:2_1:4 predicted:0.220971 +DEAL:0:3d:: cell:2_1:5 predicted:0.220971 +DEAL:0:3d:: cell:2_1:6 predicted:0.220971 +DEAL:0:3d:: cell:2_1:7 predicted:0.220971 +DEAL:0:3d:: cell:3_1:0 predicted:3.53553 +DEAL:0:3d:: cell:3_1:1 predicted:3.53553 +DEAL:0:3d:: cell:3_1:2 predicted:3.53553 +DEAL:0:3d:: cell:3_1:3 predicted:3.53553 +DEAL:0:3d:: cell:3_1:4 predicted:3.53553 +DEAL:0:3d:: cell:3_1:5 predicted:3.53553 +DEAL:0:3d:: cell:3_1:6 predicted:3.53553 +DEAL:0:3d:: cell:3_1:7 predicted:3.53553 +DEAL:0:3d::OK diff --git a/tests/mpi/error_predictor_02.with_p4est=true.mpirun=4.output b/tests/mpi/error_predictor_02.with_p4est=true.mpirun=4.output new file mode 100644 index 0000000000..a925fa4dbb --- /dev/null +++ b/tests/mpi/error_predictor_02.with_p4est=true.mpirun=4.output @@ -0,0 +1,87 @@ + +DEAL:0:2d::pre_adaptation +DEAL:0:2d:: cell:0_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:0_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:0_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d:: cell:0_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:2d::post_adaptation +DEAL:0:2d:: cell:0_0: predicted:80.0000 +DEAL:0:2d:: cell:1_0: predicted:80.0000 +DEAL:0:2d::OK +DEAL:0:3d::pre_adaptation +DEAL:0:3d:: cell:0_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:4 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:5 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:6 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d:: cell:0_1:7 fe_deg:2 error:10.0000 coarsening future_fe_deg:3 +DEAL:0:3d::post_adaptation +DEAL:0:3d:: cell:0_0: predicted:113.137 +DEAL:0:3d:: cell:1_0: predicted:113.137 +DEAL:0:3d::OK + +DEAL:1:2d::pre_adaptation +DEAL:1:2d::post_adaptation +DEAL:1:2d:: cell:2_1:0 predicted:0.312500 +DEAL:1:2d:: cell:2_1:1 predicted:0.312500 +DEAL:1:2d:: cell:2_1:2 predicted:0.312500 +DEAL:1:2d:: cell:2_1:3 predicted:0.312500 +DEAL:1:2d::OK +DEAL:1:3d::pre_adaptation +DEAL:1:3d::post_adaptation +DEAL:1:3d:: cell:2_1:0 predicted:0.220971 +DEAL:1:3d:: cell:2_1:1 predicted:0.220971 +DEAL:1:3d:: cell:2_1:2 predicted:0.220971 +DEAL:1:3d:: cell:2_1:3 predicted:0.220971 +DEAL:1:3d:: cell:2_1:4 predicted:0.220971 +DEAL:1:3d:: cell:2_1:5 predicted:0.220971 +DEAL:1:3d:: cell:2_1:6 predicted:0.220971 +DEAL:1:3d:: cell:2_1:7 predicted:0.220971 +DEAL:1:3d::OK + + +DEAL:2:2d::pre_adaptation +DEAL:2:2d:: cell:1_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:2d:: cell:1_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:2d:: cell:1_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:2d:: cell:1_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:2d::post_adaptation +DEAL:2:2d::OK +DEAL:2:3d::pre_adaptation +DEAL:2:3d:: cell:1_1:0 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:1 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:2 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:3 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:4 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:5 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:6 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d:: cell:1_1:7 fe_deg:2 error:10.0000 coarsening future_fe_deg:1 +DEAL:2:3d::post_adaptation +DEAL:2:3d::OK + + +DEAL:3:2d::pre_adaptation +DEAL:3:2d:: cell:2_0: fe_deg:2 error:10.0000 refining future_fe_deg:3 +DEAL:3:2d:: cell:3_0: fe_deg:2 error:10.0000 refining future_fe_deg:1 +DEAL:3:2d::post_adaptation +DEAL:3:2d:: cell:3_1:0 predicted:5.00000 +DEAL:3:2d:: cell:3_1:1 predicted:5.00000 +DEAL:3:2d:: cell:3_1:2 predicted:5.00000 +DEAL:3:2d:: cell:3_1:3 predicted:5.00000 +DEAL:3:2d::OK +DEAL:3:3d::pre_adaptation +DEAL:3:3d:: cell:2_0: fe_deg:2 error:10.0000 refining future_fe_deg:3 +DEAL:3:3d:: cell:3_0: fe_deg:2 error:10.0000 refining future_fe_deg:1 +DEAL:3:3d::post_adaptation +DEAL:3:3d:: cell:3_1:0 predicted:3.53553 +DEAL:3:3d:: cell:3_1:1 predicted:3.53553 +DEAL:3:3d:: cell:3_1:2 predicted:3.53553 +DEAL:3:3d:: cell:3_1:3 predicted:3.53553 +DEAL:3:3d:: cell:3_1:4 predicted:3.53553 +DEAL:3:3d:: cell:3_1:5 predicted:3.53553 +DEAL:3:3d:: cell:3_1:6 predicted:3.53553 +DEAL:3:3d:: cell:3_1:7 predicted:3.53553 +DEAL:3:3d::OK + -- 2.39.5