3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
circlepoints.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 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_CIRCLEPOINTS_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_CIRCLEPOINTS_HH
27
28#include <vector>
29#include <cmath>
30
31#include <dune/common/exceptions.hh>
32
33namespace Dumux {
34namespace EmbeddedCoupling {
35
44template<class GlobalPosition>
45std::vector<GlobalPosition> circlePoints(const GlobalPosition& center,
46 const GlobalPosition& normal,
47 const typename GlobalPosition::value_type radius,
48 const std::size_t numPoints = 20)
49{
50 using std::abs; using std::sin; using std::cos;
51 using ctype = typename GlobalPosition::value_type;
52
53 constexpr ctype eps = 1.5e-7;
54 static_assert(GlobalPosition::dimension == 3, "Only implemented for world dimension 3");
55
56 std::vector<GlobalPosition> points(numPoints);
57
58 // make sure n is a unit vector
59 auto n = normal;
60 n /= n.two_norm();
61
62 // caculate a vector u perpendicular to n
63 GlobalPosition u;
64 if (abs(n[0]) < eps && abs(n[1]) < eps)
65 if (abs(n[2]) < eps)
66 DUNE_THROW(Dune::MathError, "The normal vector has to be non-zero!");
67 else
68 u = {0, 1, 0};
69 else
70 u = {-n[1], n[0], 0};
71
72 u *= radius/u.two_norm();
73
74 // the circle parameterization is p(t) = r*cos(t)*u + r*sin(t)*(n x u) + c
75 auto tangent = crossProduct(u, n);
76 tangent *= radius/tangent.two_norm();
77
78 // the parameter with an offset
79 ctype t = 0 + 0.1;
80 // insert the vertices
81 for (std::size_t i = 0; i < numPoints; ++i)
82 {
83 points[i] = GlobalPosition({u[0]*cos(t) + tangent[0]*sin(t) + center[0],
84 u[1]*cos(t) + tangent[1]*sin(t) + center[1],
85 u[2]*cos(t) + tangent[2]*sin(t) + center[2]});
86
87 t += 2*M_PI/numPoints;
88
89 // periodic t
90 if(t > 2*M_PI)
91 t -= 2*M_PI;
92 }
93
94 return points;
95}
96
97} // end namespace EmbeddedCoupling
98} // end namespace Dumux
99
100#endif
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:631
std::vector< GlobalPosition > circlePoints(const GlobalPosition &center, const GlobalPosition &normal, const typename GlobalPosition::value_type radius, const std::size_t numPoints=20)
returns a vector of points on a circle
Definition: circlepoints.hh:45
Definition: adapt.hh:29