dune-functions 2.9.0
facenormalgridfunction.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
4#define DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
5
6#include <type_traits>
7#include <optional>
8
9#include <dune/common/exceptions.hh>
10#include <dune/common/typeutilities.hh>
11#include <dune/common/rangeutilities.hh>
12#include <dune/geometry/referenceelements.hh>
13
16
17
18namespace Dune::Functions {
19
20namespace Impl {
21
22// Compute closest face to point
23template<class ReferenceElement, class Coordinate>
24auto closestFaceIndex(const ReferenceElement& re, const Coordinate& x)
25{
26 auto closestFaceIndex = decltype(re.subEntity(0,1,0,1)){};
27 double closestFaceDistance = std::numeric_limits<double>::max();
28 for(auto&& faceIndex : Dune::range(re.size(1)))
29 {
30 // For a face unit outer normal consider the orthogonal projection
31 // Px = x + <c-x,n>*n into the face. Then the distance to the face
32 // is given by |x-Px| = |<c-x,n>||n| = <c-x,n>.
33 auto normal = re.integrationOuterNormal(faceIndex);
34 normal /= normal.two_norm();
35 auto c = re.position(faceIndex,1);
36 c -= x;
37 auto faceDistance = (c*normal);
38 if (faceDistance<closestFaceDistance)
39 {
40 closestFaceDistance = faceDistance;
41 closestFaceIndex = faceIndex;
42 }
43 }
44 return closestFaceIndex;
45}
46
47} // end namespace Impl
48
49
50
51
64template<class GV>
66{
67public:
68 using GridView = GV;
70 using Element = typename EntitySet::Element;
71
75
76private:
77
78 using Traits = Imp::GridFunctionTraits<Range(Domain), EntitySet, DefaultDerivativeTraits, 16>;
79
80 class LocalFunction
81 {
82 using Geometry = typename Element::Geometry;
83 static const int dimension = GV::dimension;
84 public:
85
97 void bind(const Element& element)
98 {
99 element_ = element;
100 geometry_.emplace(element_.geometry());
101 }
102
103 void unbind()
104 {
105 geometry_.reset();
106 }
107
110 bool bound() const
111 {
112 return static_cast<bool>(geometry_);
113 }
114
124 Range operator()(const LocalDomain& x) const
125 {
126 auto&& re = Dune::referenceElement(*geometry_);
127 // Compute reference normal of closest face to given point
128 auto face = Impl::closestFaceIndex(re, x);
129 auto localNormal = re.integrationOuterNormal(face);
130
131 // Transform reference normal into global unit outer normal using
132 // covariant Piola transformation
133 auto normal = Range{};
134 geometry_->jacobianInverseTransposed(x).mv(localNormal, normal);
135 normal /= normal.two_norm();
136 return normal;
137 }
138
140 const Element& localContext() const
141 {
142 return element_;
143 }
144
146 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction& t)
147 {
148 DUNE_THROW(NotImplemented,"not implemented");
149 }
150
151 private:
152 std::optional<Geometry> geometry_;
153 Element element_;
154 };
155
156public:
159 entitySet_(gridView)
160 {}
161
163 Range operator()(const Domain& x) const
164 {
165 DUNE_THROW(NotImplemented,"not implemented");
166 }
167
170 {
171 DUNE_THROW(NotImplemented,"not implemented");
172 }
173
175 friend LocalFunction localFunction(const FaceNormalGridFunction& t)
176 {
177 return LocalFunction{};
178 }
179
181 const EntitySet& entitySet() const
182 {
183 return entitySet_;
184 }
185
186private:
187 EntitySet entitySet_;
188};
189
190
191
192} // namespace Dune::Functions
193
194#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
Definition: polynomial.hh:11
Default implementation for derivative traits.
Definition: defaultderivativetraits.hh:37
Grid function implementing the piecewise element face normal.
Definition: facenormalgridfunction.hh:66
GridViewEntitySet< GridView, 0 > EntitySet
Definition: facenormalgridfunction.hh:69
typename EntitySet::LocalCoordinate LocalDomain
Definition: facenormalgridfunction.hh:72
typename EntitySet::GlobalCoordinate Domain
Definition: facenormalgridfunction.hh:73
typename EntitySet::GlobalCoordinate Range
Definition: facenormalgridfunction.hh:74
const EntitySet & entitySet() const
Return the stored GridViewEntitySet.
Definition: facenormalgridfunction.hh:181
friend Traits::DerivativeInterface derivative(const FaceNormalGridFunction &t)
Not implemented.
Definition: facenormalgridfunction.hh:169
GV GridView
Definition: facenormalgridfunction.hh:68
friend LocalFunction localFunction(const FaceNormalGridFunction &t)
Return a local-function associated to FaceNormalGridFunction.
Definition: facenormalgridfunction.hh:175
FaceNormalGridFunction(const GridView &gridView)
Construct the FaceNormalGridFunction.
Definition: facenormalgridfunction.hh:158
typename EntitySet::Element Element
Definition: facenormalgridfunction.hh:70
Range operator()(const Domain &x) const
Not implemented.
Definition: facenormalgridfunction.hh:163
Definition: gridfunction.hh:32
GridView::template Codim< codim >::Entity Element
Type of Elements contained in this EntitySet.
Definition: gridviewentityset.hh:32
Element::Geometry::LocalCoordinate LocalCoordinate
Type of local coordinates with respect to the Element.
Definition: gridviewentityset.hh:35
Element::Geometry::GlobalCoordinate GlobalCoordinate
Definition: gridviewentityset.hh:36