3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
15#include <dune/geometry/referenceelements.hh>
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19#include <dune/localfunctions/common/localinterpolation.hh>
21namespace Dune::Functions::Impl
37 struct ContravariantPiolaTransformator
43 template<
typename Values,
typename LocalCoordinate,
typename Geometry>
44 static auto apply(Values& values,
45 const LocalCoordinate& xi,
46 const Geometry& geometry)
48 auto jacobianTransposed = geometry.jacobianTransposed(xi);
49 auto integrationElement = geometry.integrationElement(xi);
51 for (
auto& value : values)
54 jacobianTransposed.mtv(tmp, value);
55 value /= integrationElement;
68 template<
typename Gradients,
typename LocalCoordinate,
typename Geometry>
69 static auto applyJacobian(Gradients& gradients,
70 const LocalCoordinate& xi,
71 const Geometry& geometry)
73 auto jacobianTransposed = geometry.jacobianTransposed(xi);
74 auto integrationElement = geometry.integrationElement(xi);
75 for (
auto& gradient : gradients)
79 for (
size_t k=0; k<gradient.M(); k++)
80 for (
size_t l=0; l<tmp.N(); l++)
82 for(
auto&& [jacobianTransposed_l_j, j] :
sparseRange(jacobianTransposed[l]))
83 gradient[j][k] += jacobianTransposed_l_j * tmp[l][k];
84 gradient /= integrationElement;
95 template<
class Function,
class LocalCoordinate,
class Element>
96 class LocalValuedFunction
99 const Element& element_;
103 LocalValuedFunction(
const Function& f,
const Element& element)
104 : f_(f), element_(element)
107 auto operator()(
const LocalCoordinate& xi)
const
109 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
110 auto globalValue = f(xi);
113 auto jacobianInverseTransposed = element_.geometry().jacobianInverseTransposed(xi);
114 auto integrationElement = element_.geometry().integrationElement(xi);
116 auto localValue = globalValue;
117 jacobianInverseTransposed.mtv(globalValue, localValue);
118 localValue *= integrationElement;
138 struct CovariantPiolaTransformator
144 template<
typename Values,
typename LocalCoordinate,
typename Geometry>
145 static auto apply(Values& values,
146 const LocalCoordinate& xi,
147 const Geometry& geometry)
149 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
151 for (
auto& value : values)
154 jacobianInverseTransposed.mv(tmp, value);
167 template<
typename Gradients,
typename LocalCoordinate,
typename Geometry>
168 static auto applyJacobian(Gradients& gradients,
169 const LocalCoordinate& xi,
170 const Geometry& geometry)
172 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
174 for (
auto& gradient : gradients)
178 for (
size_t j=0; j<gradient.N(); j++)
179 for (
size_t k=0; k<gradient.M(); k++)
181 for(
auto&& [jacobianInverseTransposed_j_l, l] :
sparseRange(jacobianInverseTransposed[j]))
182 gradient[j][k] += jacobianInverseTransposed_j_l * tmp[l][k];
193 template<
class Function,
class LocalCoordinate,
class Element>
194 class LocalValuedFunction
197 const Element& element_;
201 LocalValuedFunction(
const Function& f,
const Element& element)
202 : f_(f), element_(element)
205 auto operator()(
const LocalCoordinate& xi)
const
207 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
208 auto globalValue = f(xi);
211 auto jacobianTransposed = element_.geometry().jacobianTransposed(xi);
213 auto localValue = globalValue;
214 jacobianTransposed.mv(globalValue, localValue);
227 template<
class Transformator,
class LocalValuedLocalBasis,
class Element>
228 class GlobalValuedLocalBasis
231 using Traits =
typename LocalValuedLocalBasis::Traits;
235 void bind(
const LocalValuedLocalBasis& localValuedLocalBasis,
const Element& element)
237 localValuedLocalBasis_ = &localValuedLocalBasis;
245 return localValuedLocalBasis_->size();
249 void evaluateFunction(
const typename Traits::DomainType& x,
250 std::vector<typename Traits::RangeType>& out)
const
252 localValuedLocalBasis_->evaluateFunction(x,out);
254 Transformator::apply(out, x, element_->geometry());
262 void evaluateJacobian(
const typename Traits::DomainType& x,
263 std::vector<typename Traits::JacobianType>& out)
const
265 localValuedLocalBasis_->evaluateJacobian(x,out);
267 Transformator::applyJacobian(out, x, element_->geometry());
276 void partial(
const std::array<unsigned int,2>& order,
277 const typename Traits::DomainType& x,
278 std::vector<typename Traits::RangeType>& out)
const
281 if (totalOrder == 0) {
282 evaluateFunction(x, out);
283 }
else if (totalOrder == 1) {
284 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
291 std::vector<typename Traits::JacobianType> fullJacobian;
292 localValuedLocalBasis_->evaluateJacobian(x,fullJacobian);
294 Transformator::applyJacobian(fullJacobian, x, element_->geometry());
296 for (std::size_t i=0; i<out.size(); i++)
297 for (std::size_t j=0; j<out[i].size(); j++)
298 out[i][j] = fullJacobian[i][j][direction];
301 DUNE_THROW(NotImplemented,
"Partial derivatives of order 2 or higher");
307 return localValuedLocalBasis_->order();
310 const LocalValuedLocalBasis* localValuedLocalBasis_;
311 const Element* element_;
322 template<
class Transformator,
class LocalValuedLocalInterpolation,
class Element>
323 class GlobalValuedLocalInterpolation
328 void bind(
const LocalValuedLocalInterpolation& localValuedLocalInterpolation,
const Element& element)
330 localValuedLocalInterpolation_ = &localValuedLocalInterpolation;
334 template<
typename F,
typename C>
335 void interpolate (
const F& f, std::vector<C>& out)
const
337 using LocalCoordinate =
typename Element::Geometry::LocalCoordinate;
338 typename Transformator::template LocalValuedFunction<F,LocalCoordinate,Element> localValuedFunction(f, *element_);
339 localValuedLocalInterpolation_->interpolate(localValuedFunction, out);
343 const LocalValuedLocalInterpolation* localValuedLocalInterpolation_;
344 const Element* element_;
354 template<
class Transformator,
class LocalValuedLFE,
class Element>
355 class GlobalValuedLocalFiniteElement
357 using LocalBasis = GlobalValuedLocalBasis<Transformator,
358 typename LocalValuedLFE::Traits::LocalBasisType,
360 using LocalInterpolation = GlobalValuedLocalInterpolation<Transformator,
361 typename LocalValuedLFE::Traits::LocalInterpolationType,
367 using Traits = LocalFiniteElementTraits<LocalBasis,
368 typename LocalValuedLFE::Traits::LocalCoefficientsType,
371 GlobalValuedLocalFiniteElement() {}
373 void bind(
const LocalValuedLFE& localValuedLFE,
const Element& element)
375 globalValuedLocalBasis_.bind(localValuedLFE.localBasis(), element);
376 globalValuedLocalInterpolation_.bind(localValuedLFE.localInterpolation(), element);
377 localValuedLFE_ = &localValuedLFE;
384 return globalValuedLocalBasis_;
391 return localValuedLFE_->localCoefficients();
398 return globalValuedLocalInterpolation_;
402 std::size_t size()
const
404 return localValuedLFE_->size();
411 return localValuedLFE_->type();
418 const LocalValuedLFE* localValuedLFE_;
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:130
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:216
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:289
auto sparseRange(Range &&range)
Allow structured-binding for-loops for sparse iterators.
Definition: rangeutilities.hh:819
Some useful basic math stuff.
Utilities for reduction like operations on ranges.
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22