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);
75 for (
auto& gradient : gradients)
79 for (
size_t j=0; j<gradient.N(); j++)
80 for (
size_t k=0; k<gradient.M(); k++)
83 for (
size_t l=0; l<tmp.N(); l++)
84 gradient[j][k] += jacobianTransposed[l][j] * tmp[l][k];
85 gradient[j][k] /= integrationElement;
98 template<
class Function,
class LocalCoordinate,
class Element>
99 class LocalValuedFunction
102 const Element& element_;
106 LocalValuedFunction(
const Function& f,
const Element& element)
107 : f_(f), element_(element)
110 auto operator()(
const LocalCoordinate& xi)
const
112 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
113 auto globalValue = f(xi);
116 auto jacobianInverseTransposed = element_.geometry().jacobianInverseTransposed(xi);
117 auto integrationElement = element_.geometry().integrationElement(xi);
119 auto localValue = globalValue;
120 jacobianInverseTransposed.mtv(globalValue, localValue);
121 localValue *= integrationElement;
141 struct CovariantPiolaTransformator
147 template<
typename Values,
typename LocalCoordinate,
typename Geometry>
148 static auto apply(Values& values,
149 const LocalCoordinate& xi,
150 const Geometry& geometry)
152 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
154 for (
auto& value : values)
157 jacobianInverseTransposed.mv(tmp, value);
170 template<
typename Gradients,
typename LocalCoordinate,
typename Geometry>
171 static auto applyJacobian(Gradients& gradients,
172 const LocalCoordinate& xi,
173 const Geometry& geometry)
175 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
177 for (
auto& gradient : gradients)
181 for (
size_t j=0; j<gradient.N(); j++)
182 for (
size_t k=0; k<gradient.M(); k++)
185 for (
size_t l=0; l<tmp.N(); l++)
186 gradient[j][k] += jacobianInverseTransposed[j][l] * tmp[l][k];
199 template<
class Function,
class LocalCoordinate,
class Element>
200 class LocalValuedFunction
203 const Element& element_;
207 LocalValuedFunction(
const Function& f,
const Element& element)
208 : f_(f), element_(element)
211 auto operator()(
const LocalCoordinate& xi)
const
213 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
214 auto globalValue = f(xi);
217 auto jacobianTransposed = element_.geometry().jacobianTransposed(xi);
219 auto localValue = globalValue;
220 jacobianTransposed.mv(globalValue, localValue);
233 template<
class Transformator,
class LocalValuedLocalBasis,
class Element>
234 class GlobalValuedLocalBasis
237 using Traits =
typename LocalValuedLocalBasis::Traits;
241 void bind(
const LocalValuedLocalBasis& localValuedLocalBasis,
const Element& element)
243 localValuedLocalBasis_ = &localValuedLocalBasis;
251 return localValuedLocalBasis_->size();
255 void evaluateFunction(
const typename Traits::DomainType& x,
256 std::vector<typename Traits::RangeType>& out)
const
258 localValuedLocalBasis_->evaluateFunction(x,out);
260 Transformator::apply(out, x, element_->geometry());
268 void evaluateJacobian(
const typename Traits::DomainType& x,
269 std::vector<typename Traits::JacobianType>& out)
const
271 localValuedLocalBasis_->evaluateJacobian(x,out);
273 Transformator::applyJacobian(out, x, element_->geometry());
282 void partial(
const std::array<unsigned int,2>& order,
283 const typename Traits::DomainType& x,
284 std::vector<typename Traits::RangeType>& out)
const
287 if (totalOrder == 0) {
288 evaluateFunction(x, out);
289 }
else if (totalOrder == 1) {
290 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
297 std::vector<typename Traits::JacobianType> fullJacobian;
298 localValuedLocalBasis_->evaluateJacobian(x,fullJacobian);
300 Transformator::applyJacobian(fullJacobian, x, element_->geometry());
302 for (std::size_t i=0; i<out.size(); i++)
303 for (std::size_t j=0; j<out[i].size(); j++)
304 out[i][j] = fullJacobian[i][j][direction];
307 DUNE_THROW(NotImplemented,
"Partial derivatives of order 2 or higher");
313 return localValuedLocalBasis_->order();
316 const LocalValuedLocalBasis* localValuedLocalBasis_;
317 const Element* element_;
328 template<
class Transformator,
class LocalValuedLocalInterpolation,
class Element>
329 class GlobalValuedLocalInterpolation
334 void bind(
const LocalValuedLocalInterpolation& localValuedLocalInterpolation,
const Element& element)
336 localValuedLocalInterpolation_ = &localValuedLocalInterpolation;
340 template<
typename F,
typename C>
341 void interpolate (
const F& f, std::vector<C>& out)
const
343 using LocalCoordinate =
typename Element::Geometry::LocalCoordinate;
344 typename Transformator::template LocalValuedFunction<F,LocalCoordinate,Element> localValuedFunction(f, *element_);
345 localValuedLocalInterpolation_->interpolate(localValuedFunction, out);
349 const LocalValuedLocalInterpolation* localValuedLocalInterpolation_;
350 const Element* element_;
360 template<
class Transformator,
class LocalValuedLFE,
class Element>
361 class GlobalValuedLocalFiniteElement
363 using LocalBasis = GlobalValuedLocalBasis<Transformator,
364 typename LocalValuedLFE::Traits::LocalBasisType,
366 using LocalInterpolation = GlobalValuedLocalInterpolation<Transformator,
367 typename LocalValuedLFE::Traits::LocalInterpolationType,
373 using Traits = LocalFiniteElementTraits<LocalBasis,
374 typename LocalValuedLFE::Traits::LocalCoefficientsType,
377 GlobalValuedLocalFiniteElement() {}
379 void bind(
const LocalValuedLFE& localValuedLFE,
const Element& element)
381 globalValuedLocalBasis_.bind(localValuedLFE.localBasis(), element);
382 globalValuedLocalInterpolation_.bind(localValuedLFE.localInterpolation(), element);
383 localValuedLFE_ = &localValuedLFE;
390 return globalValuedLocalBasis_;
397 return localValuedLFE_->localCoefficients();
404 return globalValuedLocalInterpolation_;
408 std::size_t size()
const
410 return localValuedLFE_->size();
417 return localValuedLFE_->type();
424 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
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:289
Some useful basic math stuff.
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22