version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_CIRCLEPOINTS_HH
14#define DUMUX_MULTIDOMAIN_EMBEDDED_CIRCLEPOINTS_HH
15
16#include <vector>
17#include <cmath>
18
19#include <dune/common/exceptions.hh>
20
21#include <dumux/common/math.hh>
23
25
37template<class GlobalPosition, class Scalar>
38void circlePoints(std::vector<GlobalPosition>& points,
39 const std::vector<Scalar>& sincos,
40 const GlobalPosition& center,
41 const GlobalPosition& normal,
42 const Scalar radius)
43{
44 assert(sincos.size() % 2 == 0 && "Sample angles have to be pairs of sin/cos so size needs to be even.");
45 static_assert(GlobalPosition::dimension == 3, "Only implemented for world dimension 3");
46
47 // resize the points vector
48 const std::size_t numPoints = sincos.size()/2;
49 points.resize(numPoints);
50
51 // make sure n is a unit vector
52 auto n = normal;
53 n /= n.two_norm();
54
55 // calculate a vector u perpendicular to n
56 auto u = unitNormal(n); u*= radius;
57
58 // the circle parameterization is p(t) = r*cos(t)*u + r*sin(t)*(n x u) + c
59 auto tangent = crossProduct(u, n);
60 tangent *= radius/tangent.two_norm();
61
62 // insert the vertices
63 for (std::size_t i = 0; i < numPoints; ++i)
64 {
65 points[i] = GlobalPosition({u[0]*sincos[2*i+1] + tangent[0]*sincos[2*i] + center[0],
66 u[1]*sincos[2*i+1] + tangent[1]*sincos[2*i] + center[1],
67 u[2]*sincos[2*i+1] + tangent[2]*sincos[2*i] + center[2]});
68 }
69}
70
76template<class Scalar = double>
77std::vector<Scalar> circlePointsSinCos(const std::size_t numPoints)
78{
79 std::vector<Scalar> sincos(2*numPoints);
80 Scalar t = 0 + 0.1; // add arbitrary offset
81 for (std::size_t i = 0; i < numPoints; ++i)
82 {
83 using std::sin; using std::cos;
84 sincos[2*i] = sin(t);
85 sincos[2*i + 1] = cos(t);
86 t += 2*M_PI/numPoints;
87 if(t > 2*M_PI) t -= 2*M_PI;
88 }
89 return sincos;
90}
91
100template<class GlobalPosition, class Scalar>
101std::vector<GlobalPosition> circlePoints(const GlobalPosition& center,
102 const GlobalPosition& normal,
103 const Scalar radius,
104 const std::size_t numPoints = 20)
105{
106 // precompute the sin/cos values
107 const auto sincos = circlePointsSinCos<Scalar>(numPoints);
108
109 std::vector<GlobalPosition> points;
110 circlePoints(points, sincos, center, normal, radius);
111 return points;
112}
113
114} // end namespace Dumux::EmbeddedCoupling
115
116#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:642
std::vector< Scalar > circlePointsSinCos(const std::size_t numPoints)
returns a vector of sin and cos values of a circle parametrization
Definition: circlepoints.hh:77
void circlePoints(std::vector< GlobalPosition > &points, const std::vector< Scalar > &sincos, const GlobalPosition &center, const GlobalPosition &normal, const Scalar radius)
Definition: circlepoints.hh:38
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
Vector unitNormal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:58
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Define some often used mathematical functions.
Definition: circlepoints.hh:24
Helper functions for normals.