7#ifndef DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
8#define DUNE_FUNCTIONS_GRIDFUNCTIONS_FACENORMALGRIDFUNCTION_HH
13#include <dune/common/exceptions.hh>
14#include <dune/common/typeutilities.hh>
15#include <dune/common/rangeutilities.hh>
16#include <dune/geometry/referenceelements.hh>
27template<
class ReferenceElement,
class Coordinate>
28auto closestFaceIndex(
const ReferenceElement& re,
const Coordinate& x)
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)))
37 auto normal = re.integrationOuterNormal(faceIndex);
38 normal /= normal.two_norm();
39 auto c = re.position(faceIndex,1);
41 auto faceDistance = (c*normal);
42 if (faceDistance<closestFaceDistance)
44 closestFaceDistance = faceDistance;
45 closestFaceIndex = faceIndex;
48 return closestFaceIndex;
86 using Geometry =
typename Element::Geometry;
87 static const int dimension = GV::dimension;
101 void bind(
const Element& element)
104 geometry_.emplace(element_.geometry());
116 return static_cast<bool>(geometry_);
130 auto&& re = Dune::referenceElement(*geometry_);
132 auto face = Impl::closestFaceIndex(re, x);
133 auto localNormal = re.integrationOuterNormal(face);
137 auto normal =
Range{};
138 geometry_->jacobianInverseTransposed(x).mv(localNormal, normal);
139 normal /= normal.two_norm();
144 const Element& localContext()
const
152 DUNE_THROW(NotImplemented,
"not implemented");
156 std::optional<Geometry> geometry_;
169 DUNE_THROW(NotImplemented,
"not implemented");
175 DUNE_THROW(NotImplemented,
"not implemented");
181 return LocalFunction{};
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition trigonometricfunction.hh:43
Definition polynomial.hh:17
Definition polynomial.hh:18
Default implementation for derivative traits.
Definition defaultderivativetraits.hh:41
Grid function implementing the piecewise element face normal.
Definition facenormalgridfunction.hh:70
GridViewEntitySet< GridView, 0 > EntitySet
Definition facenormalgridfunction.hh:73
typename EntitySet::LocalCoordinate LocalDomain
Definition facenormalgridfunction.hh:76
typename EntitySet::GlobalCoordinate Domain
Definition facenormalgridfunction.hh:77
typename EntitySet::GlobalCoordinate Range
Definition facenormalgridfunction.hh:78
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
GV GridView
Definition facenormalgridfunction.hh:72
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
typename EntitySet::Element Element
Definition facenormalgridfunction.hh:74
Range operator()(const Domain &x) const
Not implemented.
Definition facenormalgridfunction.hh:167
Definition gridfunction.hh:36
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
Element::Geometry::GlobalCoordinate GlobalCoordinate
Definition gridviewentityset.hh:40