DUNE-FUNCTIONS (unstable)

facenormalgridfunction.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
5// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
6
7#ifndef DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
8#define DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
9
10#include <type_traits>
11#include <optional>
12
13#include <dune/common/exceptions.hh>
14#include <dune/common/typeutilities.hh>
15#include <dune/common/rangeutilities.hh>
16#include <dune/geometry/referenceelements.hh>
17
18#include <dune/functions/gridfunctions/gridfunction.hh>
19#include <dune/functions/gridfunctions/gridviewentityset.hh>
20
21
22namespace Dune::Functions {
23
24namespace Impl {
25
26// Compute closest face to point
27template<class ReferenceElement, class Coordinate>
28auto closestFaceIndex(const ReferenceElement& re, const Coordinate& x)
29{
30 auto closestFaceIndex = decltype(re.subEntity(0,1,0,1)){};
31 double closestFaceDistance = std::numeric_limits<double>::max();
32 for(auto&& faceIndex : Dune::range(re.size(1)))
33 {
34 // For a face unit outer normal consider the orthogonal projection
35 // Px = x + <c-x,n>*n into the face. Then the distance to the face
36 // is given by |x-Px| = |<c-x,n>||n| = <c-x,n>.
37 auto normal = re.integrationOuterNormal(faceIndex);
38 normal /= normal.two_norm();
39 auto c = re.position(faceIndex,1);
40 c -= x;
41 auto faceDistance = (c*normal);
42 if (faceDistance<closestFaceDistance)
43 {
44 closestFaceDistance = faceDistance;
45 closestFaceIndex = faceIndex;
46 }
47 }
48 return closestFaceIndex;
49}
50
51} // end namespace Impl
52
53
54
55
68template<class GV>
70{
71public:
72 using GridView = GV;
74 using Element = typename EntitySet::Element;
75
76 using LocalDomain = typename EntitySet::LocalCoordinate;
77 using Domain = typename EntitySet::GlobalCoordinate;
78 using Range = typename EntitySet::GlobalCoordinate;
79
80private:
81
82 using Traits = Imp::GridFunctionTraits<Range(Domain), EntitySet, DefaultDerivativeTraits, 16>;
83
84 class LocalFunction
85 {
86 using Geometry = typename Element::Geometry;
87 static const int dimension = GV::dimension;
88 public:
89
101 void bind(const Element& element)
102 {
103 element_ = element;
104 geometry_.emplace(element_.geometry());
105 }
106
107 void unbind()
108 {
109 geometry_.reset();
110 }
111
114 bool bound() const
115 {
116 return static_cast<bool>(geometry_);
117 }
118
128 Range operator()(const LocalDomain& x) const
129 {
130 auto&& re = Dune::referenceElement(*geometry_);
131 // Compute reference normal of closest face to given point
132 auto face = Impl::closestFaceIndex(re, x);
133 auto localNormal = re.integrationOuterNormal(face);
134
135 // Transform reference normal into global unit outer normal using
136 // covariant Piola transformation
137 auto normal = Range{};
138 geometry_->jacobianInverseTransposed(x).mv(localNormal, normal);
139 normal /= normal.two_norm();
140 return normal;
141 }
142
144 const Element& localContext() const
145 {
146 return element_;
147 }
148
150 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction& t)
151 {
152 DUNE_THROW(NotImplemented,"not implemented");
153 }
154
155 private:
156 std::optional<Geometry> geometry_;
157 Element element_;
158 };
159
160public:
162 FaceNormalGridFunction(const GridView& gridView) :
163 entitySet_(gridView)
164 {}
165
167 Range operator()(const Domain& x) const
168 {
169 DUNE_THROW(NotImplemented,"not implemented");
170 }
171
173 friend typename Traits::DerivativeInterface derivative(const FaceNormalGridFunction& t)
174 {
175 DUNE_THROW(NotImplemented,"not implemented");
176 }
177
179 friend LocalFunction localFunction(const FaceNormalGridFunction& t)
180 {
181 return LocalFunction{};
182 }
183
185 const EntitySet& entitySet() const
186 {
187 return entitySet_;
188 }
189
190private:
191 EntitySet entitySet_;
192};
193
194
195
196} // namespace Dune::Functions
197
198#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
Grid function implementing the piecewise element face normal.
Definition: facenormalgridfunction.hh:70
const EntitySet & entitySet() const
Return the stored GridViewEntitySet.
Definition: facenormalgridfunction.hh:185
friend Traits::DerivativeInterface derivative(const FaceNormalGridFunction &t)
Not implemented.
Definition: facenormalgridfunction.hh:173
friend LocalFunction localFunction(const FaceNormalGridFunction &t)
Return a local-function associated to FaceNormalGridFunction.
Definition: facenormalgridfunction.hh:179
FaceNormalGridFunction(const GridView &gridView)
Construct the FaceNormalGridFunction.
Definition: facenormalgridfunction.hh:162
Range operator()(const Domain &x) const
Not implemented.
Definition: facenormalgridfunction.hh:167
GridView::template Codim< codim >::Entity Element
Type of Elements contained in this EntitySet.
Definition: gridviewentityset.hh:36
Element::Geometry::LocalCoordinate LocalCoordinate
Type of local coordinates with respect to the Element.
Definition: gridviewentityset.hh:39
Default implementation for derivative traits.
Definition: defaultderivativetraits.hh:41
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Aug 13, 22:30, 2024)