25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_CIRCLEPOINTS_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_CIRCLEPOINTS_HH
31#include <dune/common/exceptions.hh>
34namespace EmbeddedCoupling {
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)
50 using std::abs;
using std::sin;
using std::cos;
51 using ctype =
typename GlobalPosition::value_type;
53 constexpr ctype
eps = 1.5e-7;
54 static_assert(GlobalPosition::dimension == 3,
"Only implemented for world dimension 3");
56 std::vector<GlobalPosition> points(numPoints);
64 if (abs(n[0]) <
eps && abs(n[1]) <
eps)
66 DUNE_THROW(Dune::MathError,
"The normal vector has to be non-zero!");
72 u *= radius/u.two_norm();
76 tangent *= radius/tangent.two_norm();
81 for (std::size_t i = 0; i < numPoints; ++i)
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]});
87 t += 2*M_PI/numPoints;
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 ¢er, 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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44