--- /dev/null
+New: A new BoundingBoxDataOut class is available, to output a collection
+of objects that can be converted to BoundingBox objects in graphical format.
+<br>
+(Luca Heltai, 2020/05/19)
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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_bounding_box_data_out_h
+#define dealii_bounding_box_data_out_h
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/data_out_base.h>
+#include <deal.II/base/geometry_info.h>
+#include <deal.II/base/iterator_range.h>
+#include <deal.II/base/point.h>
+
+#include <deal.II/boost_adaptors/bounding_box.h>
+#include <deal.II/boost_adaptors/point.h>
+#include <deal.II/boost_adaptors/segment.h>
+
+DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
+#include <boost/geometry/index/rtree.hpp>
+#include <boost/geometry/strategies/strategies.hpp>
+DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
+
+
+DEAL_II_NAMESPACE_OPEN
+
+/**
+ * This class generates graphical output for BoundingBox objects, starting from
+ * any object that can be converted by boost to a BoundingBox.
+ *
+ * @author Luca Heltai, 2020.
+ */
+template <int dim>
+class BoundingBoxDataOut : public DataOutInterface<dim, dim>
+{
+public:
+ BoundingBoxDataOut() = default;
+
+ ~BoundingBoxDataOut() = default;
+
+ /**
+ * Generate patches from a range of objects that can be converted by boost to
+ * a collection of BoundingBox objects.
+ */
+ template <class ConvertibleToBoundingBoxIterator>
+ void
+ build_patches(const ConvertibleToBoundingBoxIterator &begin,
+ const ConvertibleToBoundingBoxIterator &end);
+
+ /**
+ * Generate patches from a container of objects that can be converted by boost
+ * to a collection of BoundingBox objects.
+ */
+ template <class Container>
+ void
+ build_patches(const Container &boxes);
+
+protected:
+ // Copy doc
+ virtual const std::vector<dealii::DataOutBase::Patch<dim, dim>> &
+ get_patches() const override;
+
+ // Copy doc
+ virtual std::vector<std::string>
+ get_dataset_names() const override;
+
+private:
+ /**
+ * The actual boxes.
+ */
+ std::vector<DataOutBase::Patch<dim, dim>> patches;
+
+ /**
+ * Names of datasets.
+ */
+ std::vector<std::string> dataset_names;
+};
+
+
+// Template and inline functions
+#ifndef DOXYGEN
+template <int dim>
+template <class ConvertibleToBoundingBoxIterator>
+void
+BoundingBoxDataOut<dim>::build_patches(
+ const ConvertibleToBoundingBoxIterator &begin,
+ const ConvertibleToBoundingBoxIterator &end)
+{
+ using Getter = boost::geometry::index::indexable<
+ typename ConvertibleToBoundingBoxIterator::value_type>;
+ Getter getter;
+ constexpr unsigned int boxdim =
+ boost::geometry::dimension<typename Getter::result_type>::value;
+ const unsigned int N = std::distance(begin, end);
+ static_assert(boxdim == dim, "Bounding boxes are of the wrong dimension!");
+
+ dataset_names.clear();
+ patches.resize(N);
+
+ unsigned int i = 0;
+ for (const auto &value :
+ IteratorRange<ConvertibleToBoundingBoxIterator>(begin, end))
+ {
+ BoundingBox<dim> box;
+ boost::geometry::convert(getter(*value), box);
+ for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
+ {
+ patches[i].vertices[v] = box.vertex(v);
+ patches[i].patch_index = i;
+ patches[i].n_subdivisions = 1;
+ }
+ ++i;
+ }
+}
+
+
+
+template <int dim>
+template <class Container>
+void
+BoundingBoxDataOut<dim>::build_patches(const Container &boxes)
+{
+ build_patches(boxes.begin(), boxes.end());
+}
+
+
+
+template <int dim>
+const std::vector<DataOutBase::Patch<dim, dim>> &
+BoundingBoxDataOut<dim>::get_patches() const
+{
+ return patches;
+}
+
+
+
+template <int dim>
+std::vector<std::string>
+BoundingBoxDataOut<dim>::get_dataset_names() const
+{
+ return dataset_names;
+}
+
+#endif
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 BoundingBoxDataOut with vector of bounding boxes.
+
+#include <deal.II/base/bounding_box_data_out.h>
+#include <deal.II/base/geometry_info.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test()
+{
+ const unsigned int N = 2;
+ std::vector<BoundingBox<dim>> boxes(N);
+ auto unit = create_unit_bounding_box<dim>();
+
+ Point<dim> ones;
+ for (unsigned int i = 0; i < dim; ++i)
+ ones[i] = 1;
+
+ for (auto &box : boxes)
+ {
+ const auto c = random_point<dim>();
+ const auto d = random_value();
+ box =
+ BoundingBox<dim>({Point<dim>(c - d * ones), Point<dim>(c + d * ones)});
+ }
+
+ std::string fname = "boxes_" + std::to_string(dim) + ".vtk";
+ {
+ std::ofstream ofile(fname);
+ BoundingBoxDataOut<dim> data_out;
+ DataOutBase::VtkFlags flags;
+ flags.print_date_and_time = false;
+ data_out.set_flags(flags);
+ data_out.build_patches(boxes);
+ data_out.write_vtk(ofile);
+ }
+ cat_file(fname.c_str());
+}
+
+
+int
+main()
+{
+ initlog();
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
--- /dev/null
+
+# vtk DataFile Version 3.0
+#This file was generated by the deal.II library.
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 4 double
+0.445805 0 0
+1.23457 0 0
+-0.0153408 0 0
+1.58154 0 0
+
+CELLS 2 6
+2 0 1
+2 2 3
+
+CELL_TYPES 2
+ 3 3
+POINT_DATA 4
+
+# vtk DataFile Version 3.0
+#This file was generated by the deal.II library.
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 8 double
+0.576425 -0.137671 0
+1.24687 -0.137671 0
+0.576425 0.532774 0
+1.24687 0.532774 0
+0.21426 -0.276195 0
+1.3222 -0.276195 0
+0.21426 0.831745 0
+1.3222 0.831745 0
+
+CELLS 2 10
+4 0 1 3 2
+4 4 5 7 6
+
+CELL_TYPES 2
+ 9 9
+POINT_DATA 8
+
+# vtk DataFile Version 3.0
+#This file was generated by the deal.II library.
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 16 double
+-0.0360039 0.11547 -0.148616
+0.990798 0.11547 -0.148616
+-0.0360039 1.14227 -0.148616
+0.990798 1.14227 -0.148616
+-0.0360039 0.11547 0.878185
+0.990798 0.11547 0.878185
+-0.0360039 1.14227 0.878185
+0.990798 1.14227 0.878185
+0.234933 0.198898 -0.0815852
+1.66953 0.198898 -0.0815852
+0.234933 1.63349 -0.0815852
+1.66953 1.63349 -0.0815852
+0.234933 0.198898 1.35301
+1.66953 0.198898 1.35301
+0.234933 1.63349 1.35301
+1.66953 1.63349 1.35301
+
+CELLS 2 18
+8 0 1 3 2 4 5 7 6
+8 8 9 11 10 12 13 15 14
+
+CELL_TYPES 2
+ 12 12
+POINT_DATA 16
+
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 BoundingBoxDataOut with an rtree of std::pair<BoundingBox<dim>, int>.
+
+#include <deal.II/base/bounding_box_data_out.h>
+#include <deal.II/base/geometry_info.h>
+
+#include <deal.II/numerics/rtree.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test()
+{
+ const unsigned int N = 3;
+ // Use a non-trivial object, convertible to BoundingBox by boost:
+ std::vector<std::pair<BoundingBox<dim>, unsigned int>> boxes(N);
+ auto unit = create_unit_bounding_box<dim>();
+
+ Point<dim> ones;
+ for (unsigned int i = 0; i < dim; ++i)
+ ones[i] = 1;
+
+ unsigned int i = 0;
+ for (auto &box : boxes)
+ {
+ const auto c = random_point<dim>();
+ const auto d = random_value();
+ box = std::make_pair(BoundingBox<dim>({Point<dim>(c - d * ones),
+ Point<dim>(c + d * ones)}),
+ i++);
+ }
+
+ const auto tree = pack_rtree(boxes);
+
+ std::string fname = "boxes_" + std::to_string(dim) + ".vtk";
+ {
+ std::ofstream ofile(fname);
+ BoundingBoxDataOut<dim> data_out;
+ DataOutBase::VtkFlags flags;
+ flags.print_date_and_time = false;
+ data_out.set_flags(flags);
+ data_out.build_patches(tree);
+ data_out.write_vtk(ofile);
+ }
+ cat_file(fname.c_str());
+}
+
+
+int
+main()
+{
+ initlog();
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
--- /dev/null
+
+# vtk DataFile Version 3.0
+#This file was generated by the deal.II library.
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 6 double
+0.445805 0 0
+1.23457 0 0
+-0.0153408 0 0
+1.58154 0 0
+0.714096 0 0
+1.1092 0 0
+
+CELLS 3 9
+2 0 1
+2 2 3
+2 4 5
+
+CELL_TYPES 3
+ 3 3 3
+POINT_DATA 6
+
+# vtk DataFile Version 3.0
+#This file was generated by the deal.II library.
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 12 double
+0.057448 0.490455 0
+0.612997 0.490455 0
+0.057448 1.046 0
+0.612997 1.046 0
+-0.074901 -0.151474 0
+1.18284 -0.151474 0
+-0.074901 1.10627 0
+1.18284 1.10627 0
+-0.587445 -0.438829 0
+1.31701 -0.438829 0
+-0.587445 1.46563 0
+1.31701 1.46563 0
+
+CELLS 3 15
+4 0 1 3 2
+4 4 5 7 6
+4 8 9 11 10
+
+CELL_TYPES 3
+ 9 9 9
+POINT_DATA 12
+
+# vtk DataFile Version 3.0
+#This file was generated by the deal.II library.
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 24 double
+0.774593 0.494109 0.575694
+1.0578 0.494109 0.575694
+0.774593 0.777314 0.575694
+1.0578 0.777314 0.575694
+0.774593 0.494109 0.858899
+1.0578 0.494109 0.858899
+0.774593 0.777314 0.858899
+1.0578 0.777314 0.858899
+0.469737 -0.120931 0.105655
+0.7442 -0.120931 0.105655
+0.469737 0.153532 0.105655
+0.7442 0.153532 0.105655
+0.469737 -0.120931 0.380118
+0.7442 -0.120931 0.380118
+0.469737 0.153532 0.380118
+0.7442 0.153532 0.380118
+0.674386 0.0268886 0.271154
+0.933967 0.0268886 0.271154
+0.674386 0.28647 0.271154
+0.933967 0.28647 0.271154
+0.674386 0.0268886 0.530735
+0.933967 0.0268886 0.530735
+0.674386 0.28647 0.530735
+0.933967 0.28647 0.530735
+
+CELLS 3 27
+8 0 1 3 2 4 5 7 6
+8 8 9 11 10 12 13 15 14
+8 16 17 19 18 20 21 23 22
+
+CELL_TYPES 3
+ 12 12 12
+POINT_DATA 24
+