Decide how many cells will be p-adapted from all cells flagged for adaptation.
* Each cell flagged for h-refinement will also be flagged for p-refinement.
* The same applies to coarsening.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, int spacedim>
void
* Each entry of the parameter @p p_flags needs to correspond to an active
* cell.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, int spacedim>
void
* Each entry of the parameter @p criteria needs to correspond to an active
* cell.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, typename Number, int spacedim>
void
* cell. Parameters @p p_refine_fraction and @p p_coarsen_fraction need to be
* in the interval $[0,1]$.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, typename Number, int spacedim>
void
const ComparisonFunction<typename identity<Number>::type>
&compare_coarsen = std::less_equal<Number>());
+ /**
+ * Adapt which finite element to use on a given fraction of cells.
+ *
+ * Out of all cells flagged for a certain type of adaptation, be it
+ * refinement or coarsening, we will determine a fixed number of cells among
+ * this subset that will be flagged for the corresponding p-adaptive
+ * variant.
+ *
+ * For each of both refinement and coarsening subsets, we will determine a
+ * threshold based on the provided parameter @p criteria containing
+ * indicators for every active cell. In the default case for refinement, all
+ * cells with an indicator larger than or equal to the corresponding
+ * threshold will be considered for p-refinement, while for coarsening all
+ * cells with an indicator less than or equal to the matching threshold are
+ * taken into account. However, different compare function objects can be
+ * supplied via the parameters @p compare_refine and @p compare_coarsen to
+ * impose different decision strategies.
+ *
+ * For refinement, the threshold will be associated with the cell that has
+ * the @p p_refine_fraction times Triangulation::n_active_cells() largest
+ * indicator, while it is the cell with the @p p_refine_coarsen times
+ * Triangulation::n_active_cells() lowest indicator for coarsening.
+ *
+ * Each entry of the parameter @p criteria needs to correspond to an active
+ * cell. Parameters @p p_refine_fraction and @p p_coarsen_fraction need to be
+ * in the interval $[0,1]$.
+ *
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
+ */
+ template <int dim, typename Number, int spacedim>
+ void
+ p_adaptivity_fixed_number(
+ const hp::DoFHandler<dim, spacedim> &dof_handler,
+ const Vector<Number> & criteria,
+ const double p_refine_fraction = 0.5,
+ const double p_coarsen_fraction = 0.5,
+ const ComparisonFunction<typename identity<Number>::type>
+ &compare_refine = std::greater_equal<Number>(),
+ const ComparisonFunction<typename identity<Number>::type>
+ &compare_coarsen = std::less_equal<Number>());
+
/**
* Adapt which finite element to use on cells based on the regularity of the
* (unknown) analytical solution.
* }
* @endcode
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, typename Number, int spacedim>
void
* Each entry of the parameters @p criteria and @p references needs to
* correspond to an active cell.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, typename Number, int spacedim>
void
* Removes all refine and coarsen flags on cells that have a
* @p future_fe_index assigned.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, int spacedim>
void
* the decision that Triangulation::prepare_coarsening_and_refinement()
* would have made later on.
*
- * @note Preceeding calls of Triangulation::prepare_for_coarsening_and_refinement()
- * may change refine and coarsen flags, which will ultimately change the
- * results of this function.
+ * @note Triangulation::prepare_for_coarsening_and_refinement() may change
+ * refine and coarsen flags. Avoid calling it before this particular
+ * function.
*/
template <int dim, int spacedim>
void
#include <deal.II/base/mpi.h>
+#include <deal.II/distributed/grid_refinement.h>
+#include <deal.II/distributed/shared_tria.h>
#include <deal.II/distributed/tria_base.h>
#include <deal.II/grid/grid_refinement.h>
DEAL_II_NAMESPACE_OPEN
+namespace
+{
+ /**
+ * ComparisonFunction returning 'true' or 'false' for any set of parameters.
+ *
+ * These will be used to overwrite user-provided comparison functions whenever
+ * no actual comparison is required in the decision process, i.e. when no or
+ * all cells will be refined or coarsened.
+ */
+ template <typename Number>
+ hp::Refinement::ComparisonFunction<Number> compare_false =
+ [](const Number &, const Number &) { return false; };
+ template <typename Number>
+ hp::Refinement::ComparisonFunction<Number> compare_true =
+ [](const Number &, const Number &) { return true; };
+} // namespace
+
+
+
namespace hp
{
namespace Refinement
}
}
- if (const parallel::TriangulationBase<dim, spacedim> *parallel_tria =
- dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
- &dof_handler.get_triangulation()))
+ const parallel::TriangulationBase<dim, spacedim> *parallel_tria =
+ dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
+ &dof_handler.get_triangulation());
+ if (parallel_tria != nullptr &&
+ dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
+ &dof_handler.get_triangulation()) == nullptr)
{
max_criterion_refine =
Utilities::MPI::max(max_criterion_refine,
+ template <int dim, typename Number, int spacedim>
+ void
+ p_adaptivity_fixed_number(
+ const hp::DoFHandler<dim, spacedim> &dof_handler,
+ const Vector<Number> & criteria,
+ const double p_refine_fraction,
+ const double p_coarsen_fraction,
+ const ComparisonFunction<typename identity<Number>::type> &compare_refine,
+ const ComparisonFunction<typename identity<Number>::type>
+ &compare_coarsen)
+ {
+ AssertDimension(dof_handler.get_triangulation().n_active_cells(),
+ criteria.size());
+ Assert((p_refine_fraction >= 0) && (p_refine_fraction <= 1),
+ dealii::GridRefinement::ExcInvalidParameterValue());
+ Assert((p_coarsen_fraction >= 0) && (p_coarsen_fraction <= 1),
+ dealii::GridRefinement::ExcInvalidParameterValue());
+
+ // 1.) First extract from the vector of indicators the ones that
+ // correspond to cells that we locally own.
+ unsigned int n_flags_refinement = 0;
+ unsigned int n_flags_coarsening = 0;
+ Vector<Number> indicators_refinement(
+ dof_handler.get_triangulation().n_active_cells());
+ Vector<Number> indicators_coarsening(
+ dof_handler.get_triangulation().n_active_cells());
+ for (const auto &cell :
+ dof_handler.get_triangulation().active_cell_iterators())
+ if (!cell->is_artificial() && cell->is_locally_owned())
+ {
+ if (cell->refine_flag_set())
+ indicators_refinement(n_flags_refinement++) =
+ criteria(cell->active_cell_index());
+ else if (cell->coarsen_flag_set())
+ indicators_coarsening(n_flags_coarsening++) =
+ criteria(cell->active_cell_index());
+ }
+ indicators_refinement.grow_or_shrink(n_flags_refinement);
+ indicators_coarsening.grow_or_shrink(n_flags_coarsening);
+
+ // 2.) Determine the number of cells for p-refinement and p-coarsening on
+ // basis of the flagged cells.
+ //
+ // 3.) Find thresholds for p-refinment and p-coarsening on only those
+ // cells flagged for adaptation.
+ //
+ // For cases in which no or all cells flagged for refinement and/or
+ // coarsening are subject to p-adaptation, we usually pick thresholds
+ // that apply to all or none of the cells at once. However here, we
+ // do not know which threshold would suffice for this task because the
+ // user could provide any comparison function. Thus if necessary, we
+ // overwrite the user's choice with suitable functions simplying
+ // returning 'true' and 'false' for any cell with reference wrappers.
+ // Thus, no function object copies are stored.
+ //
+ // 4.) Perform p-adaptation with absolute thresholds.
+ Number threshold_refinement = 0.;
+ Number threshold_coarsening = 0.;
+ auto reference_compare_refine = std::cref(compare_refine);
+ auto reference_compare_coarsen = std::cref(compare_coarsen);
+
+ const parallel::TriangulationBase<dim, spacedim> *parallel_tria =
+ dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
+ &dof_handler.get_triangulation());
+ if (parallel_tria != nullptr &&
+ dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
+ &dof_handler.get_triangulation()) == nullptr)
+ {
+#ifndef DEAL_II_WITH_P4EST
+ Assert(false, ExcInternalError());
+#else
+ //
+ // parallel implementation with distributed memory
+ //
+
+ MPI_Comm mpi_communicator = parallel_tria->get_communicator();
+
+ // 2.) Communicate the number of cells scheduled for p-adaptation
+ // globally.
+ const unsigned int n_global_flags_refinement =
+ Utilities::MPI::sum(n_flags_refinement, mpi_communicator);
+ const unsigned int n_global_flags_coarsening =
+ Utilities::MPI::sum(n_flags_coarsening, mpi_communicator);
+
+ const unsigned int target_index_refinement =
+ static_cast<unsigned int>(
+ std::floor(p_refine_fraction * n_global_flags_refinement));
+ const unsigned int target_index_coarsening =
+ static_cast<unsigned int>(
+ std::ceil((1 - p_coarsen_fraction) * n_global_flags_coarsening));
+
+ // 3.) Figure out the global max and min of the criteria. We don't
+ // need it here, but it's a collective communication call.
+ const std::pair<Number, Number> global_min_max_refinement =
+ dealii::internal::parallel::distributed::GridRefinement::
+ compute_global_min_and_max_at_root(indicators_refinement,
+ mpi_communicator);
+
+ const std::pair<Number, Number> global_min_max_coarsening =
+ dealii::internal::parallel::distributed::GridRefinement::
+ compute_global_min_and_max_at_root(indicators_coarsening,
+ mpi_communicator);
+
+ // 3.) Compute thresholds if necessary.
+ if (target_index_refinement == 0)
+ reference_compare_refine = std::cref(compare_false<Number>);
+ else if (target_index_refinement == n_global_flags_refinement)
+ reference_compare_refine = std::cref(compare_true<Number>);
+ else
+ threshold_refinement = dealii::internal::parallel::distributed::
+ GridRefinement::RefineAndCoarsenFixedNumber::compute_threshold(
+ indicators_refinement,
+ global_min_max_refinement,
+ target_index_refinement,
+ mpi_communicator);
+
+ if (target_index_coarsening == n_global_flags_coarsening)
+ reference_compare_coarsen = std::cref(compare_false<Number>);
+ else if (target_index_coarsening == 0)
+ reference_compare_coarsen = std::cref(compare_true<Number>);
+ else
+ threshold_coarsening = dealii::internal::parallel::distributed::
+ GridRefinement::RefineAndCoarsenFixedNumber::compute_threshold(
+ indicators_coarsening,
+ global_min_max_coarsening,
+ target_index_coarsening,
+ mpi_communicator);
+#endif
+ }
+ else
+ {
+ //
+ // serial implementation (and parallel::shared implementation)
+ //
+
+ // 2.) Determine the number of cells scheduled for p-adaptation.
+ const unsigned int n_p_refine_cells = static_cast<unsigned int>(
+ std::floor(p_refine_fraction * n_flags_refinement));
+ const unsigned int n_p_coarsen_cells = static_cast<unsigned int>(
+ std::floor(p_coarsen_fraction * n_flags_coarsening));
+
+ // 3.) Compute thresholds if necessary.
+ if (n_p_refine_cells == 0)
+ reference_compare_refine = std::cref(compare_false<Number>);
+ else if (n_p_refine_cells == n_flags_refinement)
+ reference_compare_refine = std::cref(compare_true<Number>);
+ else
+ {
+ std::nth_element(indicators_refinement.begin(),
+ indicators_refinement.begin() +
+ n_p_refine_cells - 1,
+ indicators_refinement.end(),
+ std::greater<Number>());
+ threshold_refinement =
+ *(indicators_refinement.begin() + n_p_refine_cells - 1);
+ }
+
+ if (n_p_coarsen_cells == 0)
+ reference_compare_coarsen = std::cref(compare_false<Number>);
+ else if (n_p_coarsen_cells == n_flags_coarsening)
+ reference_compare_coarsen = std::cref(compare_true<Number>);
+ else
+ {
+ std::nth_element(indicators_coarsening.begin(),
+ indicators_coarsening.begin() +
+ n_p_coarsen_cells - 1,
+ indicators_coarsening.end(),
+ std::less<Number>());
+ threshold_coarsening =
+ *(indicators_coarsening.begin() + n_p_coarsen_cells - 1);
+ }
+ }
+
+ // 4.) Finally perform adaptation.
+ p_adaptivity_from_absolute_threshold(dof_handler,
+ criteria,
+ threshold_refinement,
+ threshold_coarsening,
+ std::cref(reference_compare_refine),
+ std::cref(
+ reference_compare_coarsen));
+ }
+
+
+
template <int dim, typename Number, int spacedim>
void
p_adaptivity_from_regularity(
const ComparisonFunction<typename identity<S>::type> &,
const ComparisonFunction<typename identity<S>::type> &);
+ template void
+ p_adaptivity_fixed_number<deal_II_dimension,
+ S,
+ deal_II_space_dimension>(
+ const hp::DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
+ const Vector<S> &,
+ const double,
+ const double,
+ const ComparisonFunction<typename identity<S>::type> &,
+ const ComparisonFunction<typename identity<S>::type> &);
+
template void
p_adaptivity_from_regularity<deal_II_dimension,
S,
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// validate algorithms that will flag cells for p-adaptivity:
+// - hp::Refinement::p_adaptivity_fixed_number
+
+
+#include <deal.II/base/geometry_info.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/hp/dof_handler.h>
+#include <deal.II/hp/fe_collection.h>
+#include <deal.II/hp/refinement.h>
+
+#include <deal.II/lac/vector.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+validate(const hp::DoFHandler<dim> &dh)
+{
+ deallog << " (cellid,feidx):";
+ for (const auto &cell : dh.active_cell_iterators())
+ if (!cell->is_artificial() && cell->is_locally_owned())
+ {
+ const std::string cellid = cell->id().to_string();
+ const unsigned int coarse_cellid =
+ std::stoul(cellid.substr(0, cellid.find("_")));
+
+ deallog << " (" << coarse_cellid << "," << cell->future_fe_index()
+ << ")";
+ }
+ deallog << std::endl;
+}
+
+
+
+template <int dim>
+void
+setup(Triangulation<dim> & tria,
+ hp::DoFHandler<dim> & dh,
+ const hp::FECollection<dim> &fes)
+{
+ // Initialize triangulation and dofhandler.
+ GridGenerator::subdivided_hyper_cube(tria, 4);
+ Assert(tria.n_cells(0) == tria.n_global_active_cells(), ExcInternalError());
+
+ dh.initialize(tria, fes);
+
+ // Set all active fe indices to 1.
+ // Flag first half of cells for refinement, and the other half for coarsening.
+ for (const auto &cell : dh.active_cell_iterators())
+ if (!cell->is_artificial() && cell->is_locally_owned())
+ {
+ cell->set_active_fe_index(1);
+
+ const std::string cellid = cell->id().to_string();
+ const unsigned int coarse_cellid =
+ std::stoul(cellid.substr(0, cellid.find('_')));
+
+ if (coarse_cellid < 0.5 * tria.n_global_active_cells())
+ cell->set_refine_flag();
+ else
+ cell->set_coarsen_flag();
+ }
+}
+
+
+
+template <int dim>
+void
+test()
+{
+ hp::FECollection<dim> fes;
+ for (unsigned int d = 1; d <= 3; ++d)
+ fes.push_back(FE_Q<dim>(d));
+
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ hp::DoFHandler<dim> dh;
+ setup(tria, dh, fes);
+
+ deallog << "starting situation" << std::endl;
+ validate(dh);
+
+
+ // We flag the first half of all cells to be refined and the last half of all
+ // cells to be coarsened for p adapativity. Ultimately, the first quarter of
+ // all cells will be flagged for p refinement, and the last quarter for p
+ // coarsening.
+
+ Vector<double> indicators(tria.n_active_cells());
+ for (const auto &cell : tria.active_cell_iterators())
+ if (!cell->is_artificial() && cell->is_locally_owned())
+ {
+ const std::string cellid = cell->id().to_string();
+ const unsigned int coarse_cellid =
+ std::stoul(cellid.substr(0, cellid.find('_')));
+
+ if (coarse_cellid < .25 * tria.n_global_active_cells())
+ indicators[cell->active_cell_index()] = 2.;
+ else if (coarse_cellid < .75 * tria.n_global_active_cells())
+ indicators[cell->active_cell_index()] = 1.;
+ else
+ indicators[cell->active_cell_index()] = 0.;
+ }
+
+ hp::Refinement::p_adaptivity_fixed_number(dh, indicators);
+
+ deallog << "p-adaptivity fixed number" << std::endl;
+ validate(dh);
+
+
+ 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();
+}
--- /dev/null
+
+DEAL:0:2d::starting situation
+DEAL:0:2d:: (cellid,feidx): (0,1) (1,1) (4,1) (5,1) (2,1) (3,1) (6,1) (7,1) (8,1) (9,1) (12,1) (13,1) (10,1) (11,1) (14,1) (15,1)
+DEAL:0:2d::p-adaptivity fixed number
+DEAL:0:2d:: (cellid,feidx): (0,2) (1,2) (4,1) (5,1) (2,2) (3,2) (6,1) (7,1) (8,1) (9,1) (12,0) (13,0) (10,1) (11,1) (14,0) (15,0)
+DEAL:0:2d::OK
+DEAL:0:3d::starting situation
+DEAL:0:3d:: (cellid,feidx): (0,1) (1,1) (8,1) (9,1) (2,1) (3,1) (10,1) (11,1) (16,1) (17,1) (24,1) (25,1) (18,1) (19,1) (26,1) (27,1) (4,1) (5,1) (12,1) (13,1) (6,1) (7,1) (14,1) (15,1) (20,1) (21,1) (28,1) (29,1) (22,1) (23,1) (30,1) (31,1) (32,1) (33,1) (40,1) (41,1) (34,1) (35,1) (42,1) (43,1) (48,1) (49,1) (56,1) (57,1) (50,1) (51,1) (58,1) (59,1) (36,1) (37,1) (44,1) (45,1) (38,1) (39,1) (46,1) (47,1) (52,1) (53,1) (60,1) (61,1) (54,1) (55,1) (62,1) (63,1)
+DEAL:0:3d::p-adaptivity fixed number
+DEAL:0:3d:: (cellid,feidx): (0,2) (1,2) (8,2) (9,2) (2,2) (3,2) (10,2) (11,2) (16,1) (17,1) (24,1) (25,1) (18,1) (19,1) (26,1) (27,1) (4,2) (5,2) (12,2) (13,2) (6,2) (7,2) (14,2) (15,2) (20,1) (21,1) (28,1) (29,1) (22,1) (23,1) (30,1) (31,1) (32,1) (33,1) (40,1) (41,1) (34,1) (35,1) (42,1) (43,1) (48,0) (49,0) (56,0) (57,0) (50,0) (51,0) (58,0) (59,0) (36,1) (37,1) (44,1) (45,1) (38,1) (39,1) (46,1) (47,1) (52,0) (53,0) (60,0) (61,0) (54,0) (55,0) (62,0) (63,0)
+DEAL:0:3d::OK
--- /dev/null
+
+DEAL:0:2d::starting situation
+DEAL:0:2d:: (cellid,feidx): (0,1) (1,1) (4,1) (5,1) (2,1) (3,1) (6,1) (7,1)
+DEAL:0:2d::p-adaptivity fixed number
+DEAL:0:2d:: (cellid,feidx): (0,2) (1,2) (4,1) (5,1) (2,2) (3,2) (6,1) (7,1)
+DEAL:0:2d::OK
+DEAL:0:3d::starting situation
+DEAL:0:3d:: (cellid,feidx): (0,1) (1,1) (8,1) (9,1) (2,1) (3,1) (10,1) (11,1) (16,1) (17,1) (24,1) (25,1) (18,1) (19,1) (26,1) (27,1) (4,1) (5,1) (12,1) (13,1) (6,1) (7,1) (14,1) (15,1) (20,1) (21,1) (28,1) (29,1) (22,1) (23,1) (30,1) (31,1)
+DEAL:0:3d::p-adaptivity fixed number
+DEAL:0:3d:: (cellid,feidx): (0,2) (1,2) (8,2) (9,2) (2,2) (3,2) (10,2) (11,2) (16,1) (17,1) (24,1) (25,1) (18,1) (19,1) (26,1) (27,1) (4,2) (5,2) (12,2) (13,2) (6,2) (7,2) (14,2) (15,2) (20,1) (21,1) (28,1) (29,1) (22,1) (23,1) (30,1) (31,1)
+DEAL:0:3d::OK
+
+DEAL:1:2d::starting situation
+DEAL:1:2d:: (cellid,feidx): (8,1) (9,1) (12,1) (13,1) (10,1) (11,1) (14,1) (15,1)
+DEAL:1:2d::p-adaptivity fixed number
+DEAL:1:2d:: (cellid,feidx): (8,1) (9,1) (12,0) (13,0) (10,1) (11,1) (14,0) (15,0)
+DEAL:1:2d::OK
+DEAL:1:3d::starting situation
+DEAL:1:3d:: (cellid,feidx): (32,1) (33,1) (40,1) (41,1) (34,1) (35,1) (42,1) (43,1) (48,1) (49,1) (56,1) (57,1) (50,1) (51,1) (58,1) (59,1) (36,1) (37,1) (44,1) (45,1) (38,1) (39,1) (46,1) (47,1) (52,1) (53,1) (60,1) (61,1) (54,1) (55,1) (62,1) (63,1)
+DEAL:1:3d::p-adaptivity fixed number
+DEAL:1:3d:: (cellid,feidx): (32,1) (33,1) (40,1) (41,1) (34,1) (35,1) (42,1) (43,1) (48,0) (49,0) (56,0) (57,0) (50,0) (51,0) (58,0) (59,0) (36,1) (37,1) (44,1) (45,1) (38,1) (39,1) (46,1) (47,1) (52,0) (53,0) (60,0) (61,0) (54,0) (55,0) (62,0) (63,0)
+DEAL:1:3d::OK
+
--- /dev/null
+
+DEAL:0:2d::starting situation
+DEAL:0:2d:: (cellid,feidx): (0,1) (1,1)
+DEAL:0:2d::p-adaptivity fixed number
+DEAL:0:2d:: (cellid,feidx): (0,2) (1,2)
+DEAL:0:2d::OK
+DEAL:0:3d::starting situation
+DEAL:0:3d:: (cellid,feidx): (0,1) (1,1) (2,1) (3,1) (4,1) (5,1) (6,1) (7,1)
+DEAL:0:3d::p-adaptivity fixed number
+DEAL:0:3d:: (cellid,feidx): (0,2) (1,2) (2,2) (3,2) (4,2) (5,2) (6,2) (7,2)
+DEAL:0:3d::OK
+
+DEAL:1:2d::starting situation
+DEAL:1:2d:: (cellid,feidx): (2,1) (3,1)
+DEAL:1:2d::p-adaptivity fixed number
+DEAL:1:2d:: (cellid,feidx): (2,2) (3,2)
+DEAL:1:2d::OK
+DEAL:1:3d::starting situation
+DEAL:1:3d:: (cellid,feidx): (8,1) (9,1) (10,1) (11,1) (12,1) (13,1) (14,1) (15,1)
+DEAL:1:3d::p-adaptivity fixed number
+DEAL:1:3d:: (cellid,feidx): (8,2) (9,2) (10,2) (11,2) (12,2) (13,2) (14,2) (15,2)
+DEAL:1:3d::OK
+
+
+DEAL:2:2d::starting situation
+DEAL:2:2d:: (cellid,feidx): (4,1) (5,1)
+DEAL:2:2d::p-adaptivity fixed number
+DEAL:2:2d:: (cellid,feidx): (4,1) (5,1)
+DEAL:2:2d::OK
+DEAL:2:3d::starting situation
+DEAL:2:3d:: (cellid,feidx): (16,1) (17,1) (18,1) (19,1) (20,1) (21,1) (22,1) (23,1)
+DEAL:2:3d::p-adaptivity fixed number
+DEAL:2:3d:: (cellid,feidx): (16,1) (17,1) (18,1) (19,1) (20,1) (21,1) (22,1) (23,1)
+DEAL:2:3d::OK
+
+
+DEAL:3:2d::starting situation
+DEAL:3:2d:: (cellid,feidx): (6,1) (7,1)
+DEAL:3:2d::p-adaptivity fixed number
+DEAL:3:2d:: (cellid,feidx): (6,1) (7,1)
+DEAL:3:2d::OK
+DEAL:3:3d::starting situation
+DEAL:3:3d:: (cellid,feidx): (24,1) (25,1) (26,1) (27,1) (28,1) (29,1) (30,1) (31,1)
+DEAL:3:3d::p-adaptivity fixed number
+DEAL:3:3d:: (cellid,feidx): (24,1) (25,1) (26,1) (27,1) (28,1) (29,1) (30,1) (31,1)
+DEAL:3:3d::OK
+
+
+DEAL:4:2d::starting situation
+DEAL:4:2d:: (cellid,feidx): (8,1) (9,1)
+DEAL:4:2d::p-adaptivity fixed number
+DEAL:4:2d:: (cellid,feidx): (8,1) (9,1)
+DEAL:4:2d::OK
+DEAL:4:3d::starting situation
+DEAL:4:3d:: (cellid,feidx): (32,1) (33,1) (34,1) (35,1) (36,1) (37,1) (38,1) (39,1)
+DEAL:4:3d::p-adaptivity fixed number
+DEAL:4:3d:: (cellid,feidx): (32,1) (33,1) (34,1) (35,1) (36,1) (37,1) (38,1) (39,1)
+DEAL:4:3d::OK
+
+
+DEAL:5:2d::starting situation
+DEAL:5:2d:: (cellid,feidx): (10,1) (11,1)
+DEAL:5:2d::p-adaptivity fixed number
+DEAL:5:2d:: (cellid,feidx): (10,1) (11,1)
+DEAL:5:2d::OK
+DEAL:5:3d::starting situation
+DEAL:5:3d:: (cellid,feidx): (40,1) (41,1) (42,1) (43,1) (44,1) (45,1) (46,1) (47,1)
+DEAL:5:3d::p-adaptivity fixed number
+DEAL:5:3d:: (cellid,feidx): (40,1) (41,1) (42,1) (43,1) (44,1) (45,1) (46,1) (47,1)
+DEAL:5:3d::OK
+
+
+DEAL:6:2d::starting situation
+DEAL:6:2d:: (cellid,feidx): (12,1) (13,1)
+DEAL:6:2d::p-adaptivity fixed number
+DEAL:6:2d:: (cellid,feidx): (12,0) (13,0)
+DEAL:6:2d::OK
+DEAL:6:3d::starting situation
+DEAL:6:3d:: (cellid,feidx): (48,1) (49,1) (50,1) (51,1) (52,1) (53,1) (54,1) (55,1)
+DEAL:6:3d::p-adaptivity fixed number
+DEAL:6:3d:: (cellid,feidx): (48,0) (49,0) (50,0) (51,0) (52,0) (53,0) (54,0) (55,0)
+DEAL:6:3d::OK
+
+
+DEAL:7:2d::starting situation
+DEAL:7:2d:: (cellid,feidx): (14,1) (15,1)
+DEAL:7:2d::p-adaptivity fixed number
+DEAL:7:2d:: (cellid,feidx): (14,0) (15,0)
+DEAL:7:2d::OK
+DEAL:7:3d::starting situation
+DEAL:7:3d:: (cellid,feidx): (56,1) (57,1) (58,1) (59,1) (60,1) (61,1) (62,1) (63,1)
+DEAL:7:3d::p-adaptivity fixed number
+DEAL:7:3d:: (cellid,feidx): (56,0) (57,0) (58,0) (59,0) (60,0) (61,0) (62,0) (63,0)
+DEAL:7:3d::OK
+