Skip to content
Snippets Groups Projects
Commit 4685f756 authored by Timo Koch's avatar Timo Koch Committed by Dennis Gläser
Browse files

[md][embedded] Use new normal functions in circle points generator

parent 217db3b4
No related branches found
No related tags found
1 merge request!2192[geometry] Add header for distance helper functions
......@@ -28,9 +28,11 @@
#include <vector>
#include <cmath>
#include <dumux/common/math.hh>
#include <dune/common/exceptions.hh>
#include <dumux/common/math.hh>
#include <dumux/common/geometry/normal.hh>
namespace Dumux::EmbeddedCoupling {
/*!
......@@ -54,11 +56,6 @@ void circlePoints(std::vector<GlobalPosition>& points,
assert(sincos.size() % 2 == 0 && "Sample angles have to be pairs of sin/cos so size needs to be even.");
static_assert(GlobalPosition::dimension == 3, "Only implemented for world dimension 3");
using std::abs;
using ctype = typename GlobalPosition::value_type;
constexpr ctype eps = 1.5e-7;
// resize the points vector
const std::size_t numPoints = sincos.size()/2;
points.resize(numPoints);
......@@ -68,16 +65,7 @@ void circlePoints(std::vector<GlobalPosition>& points,
n /= n.two_norm();
// calculate a vector u perpendicular to n
GlobalPosition u;
if (abs(n[0]) < eps && abs(n[1]) < eps)
if (abs(n[2]) < eps)
DUNE_THROW(Dune::MathError, "The normal vector has to be non-zero!");
else
u = {0, 1, 0};
else
u = {-n[1], n[0], 0};
u *= radius/u.two_norm();
auto u = unitNormal(n); u*= radius;
// the circle parameterization is p(t) = r*cos(t)*u + r*sin(t)*(n x u) + c
auto tangent = crossProduct(u, n);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment