3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
boundingsphere.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_BOUNDINGSPHERE_HH
25#define DUMUX_GEOMETRY_BOUNDINGSPHERE_HH
26
27#include <algorithm>
28#include <vector>
29
31
32namespace Dumux {
33
40template<class ConvexGeometry>
41static inline Sphere<typename ConvexGeometry::ctype, ConvexGeometry::coorddimension>
42boundingSphere(const ConvexGeometry& geometry)
43{
44 constexpr int dimWorld = ConvexGeometry::coorddimension;
45 assert(geometry.corners() >= 1);
46
47 auto corner = geometry.corner(0);
48 auto xMin = corner, xMax = corner;
49 using std::max; using std::min;
50
51 // Compute the min and max over the remaining vertices
52 for (std::size_t i = 1; i < geometry.corners(); ++i)
53 {
54 corner = geometry.corner(i);
55 for (std::size_t dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
56 {
57 xMin[dimIdx] = min(xMin[dimIdx], corner[dimIdx]);
58 xMax[dimIdx] = max(xMax[dimIdx], corner[dimIdx]);
59 }
60 }
61
62 auto center = 0.5*xMax + xMin;
63 auto radius = (center-xMax).two_norm();
64 return { std::move(center), std::move(radius) };
65}
66
67namespace Detail {
68template<class Scalar, int dim>
70{
71 const std::vector<Dune::FieldVector<Scalar, dim>>& points_;
72 PointsToGeometryWrapper(const std::vector<Dune::FieldVector<Scalar, dim>>& points)
73 : points_(points) {}
74
75 static constexpr int coorddimension = dim;
76 using ctype = Scalar;
77
78 auto corners() const { return points_.size(); }
79 const auto& corner(std::size_t i) const { return points_[i]; }
80};
81} // end namespace Detail
82
89template<class Scalar, int dim>
90static inline Sphere<Scalar, dim>
91boundingSphere(const std::vector<Dune::FieldVector<Scalar, dim>>& points)
92{
94 return boundingSphere(geometry);
95}
96
97} // end namespace Dumux
98
99#endif
A function to compute bounding spheres of points clouds or convex polytopes.
static Sphere< typename ConvexGeometry::ctype, ConvexGeometry::coorddimension > boundingSphere(const ConvexGeometry &geometry)
Computes a bounding sphere of a convex polytope geometry (Dune::Geometry interface)
Definition: boundingsphere.hh:42
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:36
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: boundingsphere.hh:70
static constexpr int coorddimension
Definition: boundingsphere.hh:75
const std::vector< Dune::FieldVector< Scalar, dim > > & points_
Definition: boundingsphere.hh:71
const auto & corner(std::size_t i) const
Definition: boundingsphere.hh:79
Scalar ctype
Definition: boundingsphere.hh:76
auto corners() const
Definition: boundingsphere.hh:78
PointsToGeometryWrapper(const std::vector< Dune::FieldVector< Scalar, dim > > &points)
Definition: boundingsphere.hh:72