]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added functions merge_with and volume to BoundingBox
authorGiovanni Alzetta <giovannialzetta@hotmail.it>
Fri, 22 Sep 2017 11:31:47 +0000 (11:31 +0000)
committerGiovanni Alzetta <giovannialzetta@hotmail.it>
Sun, 24 Sep 2017 10:48:32 +0000 (10:48 +0000)
doc/news/changes/minor/20170922GiovanniAlzetta [new file with mode: 0644]
include/deal.II/base/bounding_box.h
source/base/bounding_box.cc
tests/base/bounding_box_2.cc [new file with mode: 0644]
tests/base/bounding_box_2.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20170922GiovanniAlzetta b/doc/news/changes/minor/20170922GiovanniAlzetta
new file mode 100644 (file)
index 0000000..38e9cf1
--- /dev/null
@@ -0,0 +1,3 @@
+New: Added two methods to the class BoundingBox. The first, merge_with(BoundingBox) merges the current BoundingBox with BBox. The second, volume(), returns the volume (dim-dimensional measure) of the current BoundingBox
+<br>
+(Giovanni Alzetta, 2017/09/22)
index bec903d3f12851c276964b508486e7491f40b3b5..6afc94203ab4721bd6dfde24910ab8cd29cfb411 100644 (file)
@@ -72,20 +72,32 @@ public:
   BoundingBox (const std::pair<Point<spacedim,Number>,Point<spacedim,Number>> &boundary_points);
 
   /**
-   * Returns the boundary_points
+   * Return the boundary_points
    */
   const std::pair<Point<spacedim,Number>,Point<spacedim,Number>> &get_boundary_points () const;
 
   /**
-   * Returns true if the point is inside the Bounding Box, false otherwise
+   * Enlarge the current object so that it contains @p other_bbox .
+   * If the current object already contains @p other_bbox then it is not changed
+   * by this function.
+   */
+  void merge_with(const BoundingBox<spacedim,Number> &other_bbox);
+
+  /**
+   * Return true if the point is inside the Bounding Box, false otherwise
    */
   bool point_inside (const Point<spacedim, Number> &p) const;
 
+  /**
+   * Compute the volume (i.e. the dim-dimensional measure) of the BoundingBox.
+   */
+  double volume() const;
+
 private:
   std::pair<Point<spacedim, Number>,Point<spacedim,Number>> boundary_points;
 };
 
-/*------------------------------- Inline functions: Point ---------------------------*/
+/*------------------------------- Inline functions: BoundingBox ---------------------------*/
 
 #ifndef DOXYGEN
 
index f2ef5b478a3643f8c7b1deba884a4ce42b77b46f..9390d4b8ccb1db7a396822fa397e66606026bdcc 100644 (file)
@@ -31,11 +31,30 @@ bool BoundingBox<spacedim,Number>::point_inside (const Point<spacedim, Number> &
   return true;
 }
 
+template <int spacedim, typename Number>
+void BoundingBox<spacedim,Number>::merge_with(const BoundingBox<spacedim,Number> &other_bbox)
+{
+  for (unsigned int i=0; i < spacedim; ++i)
+    {
+      this->boundary_points.first[i] = std::min( this->boundary_points.first[i], other_bbox.boundary_points.first[i] );
+      this->boundary_points.second[i] = std::max( this->boundary_points.second[i], other_bbox.boundary_points.second[i] );
+    }
+}
+
 template <int spacedim, typename Number>
 const std::pair<Point<spacedim,Number>,Point<spacedim,Number>> &BoundingBox<spacedim,Number>::get_boundary_points () const
 {
   return this->boundary_points;
 }
 
+template <int spacedim, typename Number>
+double BoundingBox<spacedim,Number>::volume() const
+{
+  double vol = 1.0;
+  for (unsigned int i=0; i < spacedim; ++i)
+    vol *= ( this->boundary_points.second[i] - this->boundary_points.first[i] );
+  return vol;
+}
+
 #include "bounding_box.inst"
 DEAL_II_NAMESPACE_CLOSE
diff --git a/tests/base/bounding_box_2.cc b/tests/base/bounding_box_2.cc
new file mode 100644 (file)
index 0000000..f7c281b
--- /dev/null
@@ -0,0 +1,84 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// test for BoundingBox<unsigned int spacedim> which tests the functions
+// volume and merge_with of the class
+
+#include "../tests.h"
+
+#include <deal.II/base/point.h>
+#include <deal.II/base/bounding_box.h>
+
+template <int spacedim>
+void test_bounding_box()
+{
+  BoundingBox<spacedim> a;
+  deallog << "Empty BoundingBox: " << std::endl;
+  deallog << a.volume() << std::endl;
+
+
+  std::pair<Point<spacedim>,Point<spacedim>> unit;
+  for (int i=0; i<spacedim; i++)
+    {
+      unit.first[i] = 0.0;
+      unit.second[i] = 1.0;
+    }
+
+  BoundingBox<spacedim> b(unit);
+  deallog << "Bounding box boundaries: " << std::endl;
+  deallog << b.get_boundary_points().first << std::endl;
+  deallog << b.get_boundary_points().second << std::endl;
+  deallog << " Boundary Box volume: " << b.volume() << std::endl;
+  b.merge_with(a);
+  deallog << "Merging it with 0: volume "  << b.volume() << std::endl;
+  deallog << "Boundary points:" << std::endl;
+  deallog << b.get_boundary_points().first << std::endl;
+  deallog << b.get_boundary_points().second << std::endl;
+
+  std::pair<Point<spacedim>,Point<spacedim>> boundaries;
+  for (int i=0; i<spacedim; i++)
+    {
+      unit.first[i] = 1.0;
+      unit.second[i] = 2.0 + i ;
+    }
+
+  BoundingBox<spacedim> c(unit);
+  deallog << "Bounding box boundaries: " << std::endl;
+  deallog << c.get_boundary_points().first << std::endl;
+  deallog << c.get_boundary_points().second << std::endl;
+  deallog << c.volume() << std::endl;
+  c.merge_with(b);
+  deallog << "Merging it with previous: volume " << c.volume() << std::endl;
+  deallog << "Boundary points:" << std::endl;
+  deallog << c.get_boundary_points().first << std::endl;
+  deallog << c.get_boundary_points().second << std::endl;
+  deallog << "End test for dimension " << spacedim << std::endl;
+}
+
+int main()
+{
+  initlog();
+
+  deallog << "Test: Bounding Box class merge and volume functions" << std::endl;
+  deallog << std::endl << "Test for dimension 1" << std::endl;
+  test_bounding_box<1>();
+
+  deallog << std::endl << "Test for dimension 2" << std::endl;
+  test_bounding_box<2>();
+
+  deallog << std::endl << "Test for dimension 3" << std::endl;
+  test_bounding_box<3>();
+}
diff --git a/tests/base/bounding_box_2.output b/tests/base/bounding_box_2.output
new file mode 100644 (file)
index 0000000..c12cf60
--- /dev/null
@@ -0,0 +1,65 @@
+
+DEAL::Test: Bounding Box class merge and volume functions
+DEAL::
+DEAL::Test for dimension 1
+DEAL::Empty BoundingBox: 
+DEAL::0.00000
+DEAL::Bounding box boundaries: 
+DEAL::0.00000
+DEAL::1.00000
+DEAL:: Boundary Box volume: 1.00000
+DEAL::Merging it with 0: volume 1.00000
+DEAL::Boundary points:
+DEAL::0.00000
+DEAL::1.00000
+DEAL::Bounding box boundaries: 
+DEAL::1.00000
+DEAL::2.00000
+DEAL::1.00000
+DEAL::Merging it with previous: volume 2.00000
+DEAL::Boundary points:
+DEAL::0.00000
+DEAL::2.00000
+DEAL::End test for dimension 1
+DEAL::
+DEAL::Test for dimension 2
+DEAL::Empty BoundingBox: 
+DEAL::0.00000
+DEAL::Bounding box boundaries: 
+DEAL::0.00000 0.00000
+DEAL::1.00000 1.00000
+DEAL:: Boundary Box volume: 1.00000
+DEAL::Merging it with 0: volume 1.00000
+DEAL::Boundary points:
+DEAL::0.00000 0.00000
+DEAL::1.00000 1.00000
+DEAL::Bounding box boundaries: 
+DEAL::1.00000 1.00000
+DEAL::2.00000 3.00000
+DEAL::2.00000
+DEAL::Merging it with previous: volume 6.00000
+DEAL::Boundary points:
+DEAL::0.00000 0.00000
+DEAL::2.00000 3.00000
+DEAL::End test for dimension 2
+DEAL::
+DEAL::Test for dimension 3
+DEAL::Empty BoundingBox: 
+DEAL::0.00000
+DEAL::Bounding box boundaries: 
+DEAL::0.00000 0.00000 0.00000
+DEAL::1.00000 1.00000 1.00000
+DEAL:: Boundary Box volume: 1.00000
+DEAL::Merging it with 0: volume 1.00000
+DEAL::Boundary points:
+DEAL::0.00000 0.00000 0.00000
+DEAL::1.00000 1.00000 1.00000
+DEAL::Bounding box boundaries: 
+DEAL::1.00000 1.00000 1.00000
+DEAL::2.00000 3.00000 4.00000
+DEAL::6.00000
+DEAL::Merging it with previous: volume 24.0000
+DEAL::Boundary points:
+DEAL::0.00000 0.00000 0.00000
+DEAL::2.00000 3.00000 4.00000
+DEAL::End test for dimension 3

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.