--- /dev/null
+Fixed: GridTools::minimal_cell_diameter() and GridTools::maximal_cell_diameter()
+return the maximal respectively minimal cell diameter of the non-artifical part
+of a Triangulation object.
+<br>
+(Daniel Arndt, 2018/03/25)
const Mapping<dim,spacedim> &mapping = (StaticMappingQ1<dim,spacedim>::mapping));
/**
- * Return the diameter of the smallest active cell of a triangulation. See
- * step-24 for an example of use of this function.
+ * Return the diameter of the smallest active cell of the non-artificial part
+ * of a triangulation. See step-24 for an example of use of this function.
*/
template <int dim, int spacedim>
double
minimal_cell_diameter (const Triangulation<dim, spacedim> &triangulation);
/**
- * Return the diameter of the largest active cell of a triangulation.
+ * Return the diameter of the largest active cell of the non-artificial part
+ * of a triangulation.
*/
template <int dim, int spacedim>
double
double
minimal_cell_diameter (const Triangulation<dim, spacedim> &triangulation)
{
- double min_diameter = triangulation.begin_active()->diameter();
- for (typename Triangulation<dim, spacedim>::active_cell_iterator
- cell = triangulation.begin_active(); cell != triangulation.end();
- ++cell)
- min_diameter = std::min (min_diameter,
- cell->diameter());
+ double min_diameter = std::numeric_limits<double>::max();
+ for (const auto cell: triangulation.active_cell_iterators())
+ if (!cell->is_artificial())
+ min_diameter = std::min (min_diameter,
+ cell->diameter());
return min_diameter;
}
double
maximal_cell_diameter (const Triangulation<dim, spacedim> &triangulation)
{
- double max_diameter = triangulation.begin_active()->diameter();
- for (typename Triangulation<dim, spacedim>::active_cell_iterator
- cell = triangulation.begin_active(); cell != triangulation.end();
- ++cell)
- max_diameter = std::max (max_diameter,
- cell->diameter());
+ double max_diameter = 0.;
+ for (const auto cell: triangulation.active_cell_iterators())
+ if (!cell->is_artificial())
+ max_diameter = std::max (max_diameter,
+ cell->diameter());
return max_diameter;
}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+#include "../tests.h"
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_out.h>
+
+
+template <int dim>
+void test1 ()
+{
+ // test 1: hypercube
+ if (true)
+ {
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::hyper_cube(tria);
+
+ for (unsigned int i=0; i<2; ++i)
+ {
+ tria.refine_global(2);
+ deallog << dim << "d, "
+ << "max diameter: "
+ << GridTools::maximal_cell_diameter (tria)
+ << std::endl;
+ Assert (GridTools::maximal_cell_diameter (tria)
+ >=
+ GridTools::minimal_cell_diameter (tria),
+ ExcInternalError());
+ };
+ };
+
+ // test 2: hyperball
+ if (dim >= 2)
+ {
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::hyper_ball(tria, Point<dim>(), 1);
+
+ for (unsigned int i=0; i<2; ++i)
+ {
+ tria.refine_global(2);
+ deallog << dim << "d, "
+ << "max diameter: "
+ << GridTools::maximal_cell_diameter (tria)
+ << std::endl;
+ Assert (GridTools::maximal_cell_diameter (tria)
+ >=
+ GridTools::minimal_cell_diameter (tria),
+ ExcInternalError());
+ };
+ };
+}
+
+
+int main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+
+ MPILogInitAll mpi_init_log;
+
+ test1<2> ();
+ test1<3> ();
+
+ return 0;
+}
+
--- /dev/null
+
+DEAL:0::2d, max diameter: 0.353553
+DEAL:0::2d, max diameter: 0.0883883
+DEAL:0::2d, max diameter: 0.418349
+DEAL:0::2d, max diameter: 0.114019
+DEAL:0::3d, max diameter: 0.433013
+DEAL:0::3d, max diameter: 0.108253
+DEAL:0::3d, max diameter: 0.481724
+DEAL:0::3d, max diameter: 0.132362
+
+DEAL:1::2d, max diameter: 0.353553
+DEAL:1::2d, max diameter: 0.0883883
+DEAL:1::2d, max diameter: 0.418349
+DEAL:1::2d, max diameter: 0.114019
+DEAL:1::3d, max diameter: 0.433013
+DEAL:1::3d, max diameter: 0.108253
+DEAL:1::3d, max diameter: 0.481724
+DEAL:1::3d, max diameter: 0.132362
+
+
+DEAL:2::2d, max diameter: 0.353553
+DEAL:2::2d, max diameter: 0.0883883
+DEAL:2::2d, max diameter: 0.418349
+DEAL:2::2d, max diameter: 0.114019
+DEAL:2::3d, max diameter: 0.433013
+DEAL:2::3d, max diameter: 0.108253
+DEAL:2::3d, max diameter: 0.481724
+DEAL:2::3d, max diameter: 0.132362
+
+
+DEAL:3::2d, max diameter: 0.353553
+DEAL:3::2d, max diameter: 0.0883883
+DEAL:3::2d, max diameter: 0.418349
+DEAL:3::2d, max diameter: 0.114019
+DEAL:3::3d, max diameter: 0.433013
+DEAL:3::3d, max diameter: 0.108253
+DEAL:3::3d, max diameter: 0.481724
+DEAL:3::3d, max diameter: 0.132362
+
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+#include "../tests.h"
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_out.h>
+
+
+template <int dim>
+void test1 ()
+{
+ // test 1: hypercube
+ if (true)
+ {
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::hyper_cube(tria);
+
+ for (unsigned int i=0; i<2; ++i)
+ {
+ tria.refine_global(2);
+ deallog << dim << "d, "
+ << "min diameter: "
+ << GridTools::minimal_cell_diameter (tria)
+ << std::endl;
+ };
+ };
+
+ // test 2: hyperball
+ if (dim >= 2)
+ {
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::hyper_ball(tria, Point<dim>(), 1);
+
+ for (unsigned int i=0; i<2; ++i)
+ {
+ tria.refine_global(2);
+ deallog << dim << "d, "
+ << "min diameter: "
+ << GridTools::minimal_cell_diameter (tria)
+ << std::endl;
+ };
+ };
+}
+
+
+int main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+
+ MPILogInitAll mpi_init_log;
+
+ test1<2> ();
+ test1<3> ();
+
+ return 0;
+}
+
--- /dev/null
+
+DEAL:0::2d, min diameter: 0.353553
+DEAL:0::2d, min diameter: 0.0883883
+DEAL:0::2d, min diameter: 0.207107
+DEAL:0::2d, min diameter: 0.0475189
+DEAL:0::3d, min diameter: 0.433013
+DEAL:0::3d, min diameter: 0.108253
+DEAL:0::3d, min diameter: 0.183013
+DEAL:0::3d, min diameter: 0.0457532
+
+DEAL:1::2d, min diameter: 0.353553
+DEAL:1::2d, min diameter: 0.0883883
+DEAL:1::2d, min diameter: 0.207107
+DEAL:1::2d, min diameter: 0.0475189
+DEAL:1::3d, min diameter: 0.433013
+DEAL:1::3d, min diameter: 0.108253
+DEAL:1::3d, min diameter: 0.183013
+DEAL:1::3d, min diameter: 0.0457532
+
+
+DEAL:2::2d, min diameter: 0.353553
+DEAL:2::2d, min diameter: 0.0883883
+DEAL:2::2d, min diameter: 0.207107
+DEAL:2::2d, min diameter: 0.0475189
+DEAL:2::3d, min diameter: 0.433013
+DEAL:2::3d, min diameter: 0.108253
+DEAL:2::3d, min diameter: 0.183013
+DEAL:2::3d, min diameter: 0.0457532
+
+
+DEAL:3::2d, min diameter: 0.353553
+DEAL:3::2d, min diameter: 0.0883883
+DEAL:3::2d, min diameter: 0.207107
+DEAL:3::2d, min diameter: 0.0475189
+DEAL:3::3d, min diameter: 0.433013
+DEAL:3::3d, min diameter: 0.108253
+DEAL:3::3d, min diameter: 0.183013
+DEAL:3::3d, min diameter: 0.0457532
+