3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
14#include <dune/geometry/referenceelements.hh>
16#include <dune/localfunctions/common/localbasis.hh>
17#include <dune/localfunctions/common/localfiniteelementtraits.hh>
18#include <dune/localfunctions/common/localinterpolation.hh>
20namespace Dune::Functions::Impl
36 struct ContravariantPiolaTransformator
42 template<
typename Values,
typename LocalCoordinate,
typename Geometry>
43 static auto apply(Values& values,
44 const LocalCoordinate& xi,
45 const Geometry& geometry)
47 auto jacobianTransposed = geometry.jacobianTransposed(xi);
48 auto integrationElement = geometry.integrationElement(xi);
50 for (
auto& value : values)
53 jacobianTransposed.mtv(tmp, value);
54 value /= integrationElement;
67 template<
typename Gradients,
typename LocalCoordinate,
typename Geometry>
68 static auto applyJacobian(Gradients& gradients,
69 const LocalCoordinate& xi,
70 const Geometry& geometry)
72 auto jacobianTransposed = geometry.jacobianTransposed(xi);
73 auto integrationElement = geometry.integrationElement(xi);
74 for (
auto& gradient : gradients)
78 for (
size_t k=0; k<gradient.M(); k++)
79 for (
size_t l=0; l<tmp.N(); l++)
81 for(
auto&& [jacobianTransposed_l_j, j] :
sparseRange(jacobianTransposed[l]))
82 gradient[j][k] += jacobianTransposed_l_j * tmp[l][k];
83 gradient /= integrationElement;
94 template<
class Function,
class LocalCoordinate,
class Element>
95 class LocalValuedFunction
98 const Element& element_;
102 LocalValuedFunction(
const Function& f,
const Element& element)
103 : f_(f), element_(element)
106 auto operator()(
const LocalCoordinate& xi)
const
108 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
109 auto globalValue = f(xi);
112 auto jacobianInverseTransposed = element_.geometry().jacobianInverseTransposed(xi);
113 auto integrationElement = element_.geometry().integrationElement(xi);
115 auto localValue = globalValue;
116 jacobianInverseTransposed.mtv(globalValue, localValue);
117 localValue *= integrationElement;
137 struct CovariantPiolaTransformator
143 template<
typename Values,
typename LocalCoordinate,
typename Geometry>
144 static auto apply(Values& values,
145 const LocalCoordinate& xi,
146 const Geometry& geometry)
148 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
150 for (
auto& value : values)
153 jacobianInverseTransposed.mv(tmp, value);
166 template<
typename Gradients,
typename LocalCoordinate,
typename Geometry>
167 static auto applyJacobian(Gradients& gradients,
168 const LocalCoordinate& xi,
169 const Geometry& geometry)
171 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
173 for (
auto& gradient : gradients)
177 for (
size_t j=0; j<gradient.N(); j++)
178 for (
size_t k=0; k<gradient.M(); k++)
180 for(
auto&& [jacobianInverseTransposed_j_l, l] :
sparseRange(jacobianInverseTransposed[j]))
181 gradient[j][k] += jacobianInverseTransposed_j_l * tmp[l][k];
192 template<
class Function,
class LocalCoordinate,
class Element>
193 class LocalValuedFunction
196 const Element& element_;
200 LocalValuedFunction(
const Function& f,
const Element& element)
201 : f_(f), element_(element)
204 auto operator()(
const LocalCoordinate& xi)
const
206 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
207 auto globalValue = f(xi);
210 auto jacobianTransposed = element_.geometry().jacobianTransposed(xi);
212 auto localValue = globalValue;
213 jacobianTransposed.mv(globalValue, localValue);
226 template<
class Transformator,
class LocalValuedLocalBasis,
class Element>
227 class GlobalValuedLocalBasis
230 using Traits =
typename LocalValuedLocalBasis::Traits;
234 void bind(
const LocalValuedLocalBasis& localValuedLocalBasis,
const Element& element)
236 localValuedLocalBasis_ = &localValuedLocalBasis;
244 return localValuedLocalBasis_->size();
248 void evaluateFunction(
const typename Traits::DomainType& x,
249 std::vector<typename Traits::RangeType>& out)
const
251 localValuedLocalBasis_->evaluateFunction(x,out);
253 Transformator::apply(out, x, element_->geometry());
261 void evaluateJacobian(
const typename Traits::DomainType& x,
262 std::vector<typename Traits::JacobianType>& out)
const
264 localValuedLocalBasis_->evaluateJacobian(x,out);
266 Transformator::applyJacobian(out, x, element_->geometry());
275 void partial(
const std::array<unsigned int,2>& order,
276 const typename Traits::DomainType& x,
277 std::vector<typename Traits::RangeType>& out)
const
280 if (totalOrder == 0) {
281 evaluateFunction(x, out);
282 }
else if (totalOrder == 1) {
283 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
290 std::vector<typename Traits::JacobianType> fullJacobian;
291 localValuedLocalBasis_->evaluateJacobian(x,fullJacobian);
293 Transformator::applyJacobian(fullJacobian, x, element_->geometry());
295 for (std::size_t i=0; i<out.size(); i++)
296 for (std::size_t j=0; j<out[i].size(); j++)
297 out[i][j] = fullJacobian[i][j][direction];
300 DUNE_THROW(NotImplemented,
"Partial derivatives of order 2 or higher");
306 return localValuedLocalBasis_->order();
309 const LocalValuedLocalBasis* localValuedLocalBasis_;
310 const Element* element_;
321 template<
class Transformator,
class LocalValuedLocalInterpolation,
class Element>
322 class GlobalValuedLocalInterpolation
327 void bind(
const LocalValuedLocalInterpolation& localValuedLocalInterpolation,
const Element& element)
329 localValuedLocalInterpolation_ = &localValuedLocalInterpolation;
333 template<
typename F,
typename C>
334 void interpolate (
const F& f, std::vector<C>& out)
const
336 using LocalCoordinate =
typename Element::Geometry::LocalCoordinate;
337 typename Transformator::template LocalValuedFunction<F,LocalCoordinate,Element> localValuedFunction(f, *element_);
338 localValuedLocalInterpolation_->interpolate(localValuedFunction, out);
342 const LocalValuedLocalInterpolation* localValuedLocalInterpolation_;
343 const Element* element_;
353 template<
class Transformator,
class LocalValuedLFE,
class Element>
354 class GlobalValuedLocalFiniteElement
356 using LocalBasis = GlobalValuedLocalBasis<Transformator,
357 typename LocalValuedLFE::Traits::LocalBasisType,
359 using LocalInterpolation = GlobalValuedLocalInterpolation<Transformator,
360 typename LocalValuedLFE::Traits::LocalInterpolationType,
366 using Traits = LocalFiniteElementTraits<LocalBasis,
367 typename LocalValuedLFE::Traits::LocalCoefficientsType,
370 GlobalValuedLocalFiniteElement() {}
372 void bind(
const LocalValuedLFE& localValuedLFE,
const Element& element)
374 globalValuedLocalBasis_.bind(localValuedLFE.localBasis(), element);
375 globalValuedLocalInterpolation_.bind(localValuedLFE.localInterpolation(), element);
376 localValuedLFE_ = &localValuedLFE;
383 return globalValuedLocalBasis_;
390 return localValuedLFE_->localCoefficients();
397 return globalValuedLocalInterpolation_;
401 std::size_t size()
const
403 return localValuedLFE_->size();
410 return localValuedLFE_->type();
417 const LocalValuedLFE* localValuedLFE_;
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
auto sparseRange(Range &&range)
Allow structured-binding for-loops for sparse iterators.
Definition: rangeutilities.hh:821
Some useful basic math stuff.
Utilities for reduction like operations on ranges.
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24