DUNE PDELab (git)

globalvaluedlocalfiniteelement.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
5
6#include <array>
7#include <numeric>
8
11#include <dune/common/math.hh>
13
14#include <dune/geometry/referenceelements.hh>
15
16#include <dune/localfunctions/common/localbasis.hh>
17#include <dune/localfunctions/common/localfiniteelementtraits.hh>
18
19namespace Dune::Functions::Impl
20{
21
35 struct ContravariantPiolaTransformator
36 {
41 template<typename Values, typename LocalCoordinate, typename Geometry>
42 static auto apply(Values& values,
43 const LocalCoordinate& xi,
44 const Geometry& geometry)
45 {
46 auto jacobianTransposed = geometry.jacobianTransposed(xi);
47 auto integrationElement = geometry.integrationElement(xi);
48
49 for (auto& value : values)
50 {
51 auto tmp = value;
52 jacobianTransposed.mtv(tmp, value);
53 value /= integrationElement;
54 }
55 }
56
66 template<typename Gradients, typename LocalCoordinate, typename Geometry>
67 static auto applyJacobian(Gradients& gradients,
68 const LocalCoordinate& xi,
69 const Geometry& geometry)
70 {
71 auto jacobianTransposed = geometry.jacobianTransposed(xi);
72 auto integrationElement = geometry.integrationElement(xi);
73 for (auto& gradient : gradients)
74 {
75 auto tmp = gradient;
76 gradient = 0;
77 for (size_t k=0; k<gradient.M(); k++)
78 for (size_t l=0; l<tmp.N(); l++)
79 // Use sparseRange because jacobianTransposed may be a sparse DiagonalMatrix
80 for(auto&& [jacobianTransposed_l_j, j] : sparseRange(jacobianTransposed[l]))
81 gradient[j][k] += jacobianTransposed_l_j * tmp[l][k];
82 gradient /= integrationElement;
83 }
84 }
85
93 template<class Function, class LocalCoordinate, class Element>
94 class LocalValuedFunction
95 {
96 const Function& f_;
97 const Element& element_;
98
99 using LocalValue = LocalCoordinate;
100
101 public:
102
103 LocalValuedFunction(const Function& f, const Element& element)
104 : f_(f), element_(element)
105 {}
106
107 auto operator()(const LocalCoordinate& xi) const
108 {
109 auto globalValue = f_(xi);
110
111 // Apply the inverse Piola transform
112 auto jacobianInverseTransposed = element_.geometry().jacobianInverseTransposed(xi);
113 auto integrationElement = element_.geometry().integrationElement(xi);
114
115 auto localValue = LocalValue{};
116 jacobianInverseTransposed.mtv(globalValue, localValue);
117 localValue *= integrationElement;
118
119 return localValue;
120 }
121 };
122 };
123
137 struct CovariantPiolaTransformator
138 {
143 template<typename Values, typename LocalCoordinate, typename Geometry>
144 static auto apply(Values& values,
145 const LocalCoordinate& xi,
146 const Geometry& geometry)
147 {
148 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
149
150 for (auto& value : values)
151 {
152 auto tmp = value;
153 jacobianInverseTransposed.mv(tmp, value);
154 }
155 }
156
166 template<typename Gradients, typename LocalCoordinate, typename Geometry>
167 static auto applyJacobian(Gradients& gradients,
168 const LocalCoordinate& xi,
169 const Geometry& geometry)
170 {
171 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
172
173 for (auto& gradient : gradients)
174 {
175 auto tmp = gradient;
176 gradient = 0;
177 for (size_t j=0; j<gradient.N(); j++)
178 for (size_t k=0; k<gradient.M(); k++)
179 // Use sparseRange because jacobianTransposed may be a sparse DiagonalMatrix
180 for(auto&& [jacobianInverseTransposed_j_l, l] : sparseRange(jacobianInverseTransposed[j]))
181 gradient[j][k] += jacobianInverseTransposed_j_l * tmp[l][k];
182 }
183 }
184
192 template<class Function, class LocalCoordinate, class Element>
193 class LocalValuedFunction
194 {
195 const Function& f_;
196 const Element& element_;
197
198 public:
199
200 LocalValuedFunction(const Function& f, const Element& element)
201 : f_(f), element_(element)
202 {}
203
204 auto operator()(const LocalCoordinate& xi) const
205 {
206 auto globalValue = f_(xi);
207
208 // Apply the inverse Piola transform
209 auto jacobianTransposed = element_.geometry().jacobianTransposed(xi);
210
211 auto localValue = globalValue;
212 jacobianTransposed.mv(globalValue, localValue);
213
214 return localValue;
215 }
216 };
217 };
218
225 template<class Transformator, class LocalValuedLocalBasis, class Element>
226 class GlobalValuedLocalBasis
227 {
228 public:
229 using Traits = typename LocalValuedLocalBasis::Traits;
230
233 void bind(const LocalValuedLocalBasis& localValuedLocalBasis, const Element& element)
234 {
235 localValuedLocalBasis_ = &localValuedLocalBasis;
236 element_ = &element;
237 }
238
241 auto size() const
242 {
243 return localValuedLocalBasis_->size();
244 }
245
247 void evaluateFunction(const typename Traits::DomainType& x,
248 std::vector<typename Traits::RangeType>& out) const
249 {
250 localValuedLocalBasis_->evaluateFunction(x,out);
251
252 Transformator::apply(out, x, element_->geometry());
253 }
254
260 void evaluateJacobian(const typename Traits::DomainType& x,
261 std::vector<typename Traits::JacobianType>& out) const
262 {
263 localValuedLocalBasis_->evaluateJacobian(x,out);
264
265 Transformator::applyJacobian(out, x, element_->geometry());
266 }
267
274 void partial(const std::array<unsigned int,2>& order,
275 const typename Traits::DomainType& x,
276 std::vector<typename Traits::RangeType>& out) const
277 {
278 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
279 if (totalOrder == 0) {
280 evaluateFunction(x, out);
281 } else if (totalOrder == 1) {
282 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
283 out.resize(size());
284
285 // TODO: The following is wasteful: We compute the full Jacobian and then return
286 // only a part of it. While we need the full Jacobian of the underlying local-valued LFE,
287 // it should be possible to compute only a partial Piola transform for the requested
288 // partial derivatives.
289 std::vector<typename Traits::JacobianType> fullJacobian;
290 localValuedLocalBasis_->evaluateJacobian(x,fullJacobian);
291
292 Transformator::applyJacobian(fullJacobian, x, element_->geometry());
293
294 for (std::size_t i=0; i<out.size(); i++)
295 for (std::size_t j=0; j<out[i].size(); j++)
296 out[i][j] = fullJacobian[i][j][direction];
297
298 } else
299 DUNE_THROW(NotImplemented, "Partial derivatives of order 2 or higher");
300 }
301
303 auto order() const
304 {
305 return localValuedLocalBasis_->order();
306 }
307
308 const LocalValuedLocalBasis* localValuedLocalBasis_;
309 const Element* element_;
310 };
311
320 template<class Transformator, class LocalValuedLocalInterpolation, class Element>
321 class GlobalValuedLocalInterpolation
322 {
323 public:
326 void bind(const LocalValuedLocalInterpolation& localValuedLocalInterpolation, const Element& element)
327 {
328 localValuedLocalInterpolation_ = &localValuedLocalInterpolation;
329 element_ = &element;
330 }
331
332 template<typename F, typename C>
333 void interpolate (const F& f, std::vector<C>& out) const
334 {
335 using LocalCoordinate = typename Element::Geometry::LocalCoordinate;
336 typename Transformator::template LocalValuedFunction<F,LocalCoordinate,Element> localValuedFunction(f, *element_);
337 localValuedLocalInterpolation_->interpolate(localValuedFunction, out);
338 }
339
340 private:
341 const LocalValuedLocalInterpolation* localValuedLocalInterpolation_;
342 const Element* element_;
343 };
344
345
352 template<class Transformator, class LocalValuedLFE, class Element>
353 class GlobalValuedLocalFiniteElement
354 {
355 using LocalBasis = GlobalValuedLocalBasis<Transformator,
356 typename LocalValuedLFE::Traits::LocalBasisType,
357 Element>;
358 using LocalInterpolation = GlobalValuedLocalInterpolation<Transformator,
359 typename LocalValuedLFE::Traits::LocalInterpolationType,
360 Element>;
361
362 public:
365 using Traits = LocalFiniteElementTraits<LocalBasis,
366 typename LocalValuedLFE::Traits::LocalCoefficientsType,
367 LocalInterpolation>;
368
369 GlobalValuedLocalFiniteElement() {}
370
371 void bind(const LocalValuedLFE& localValuedLFE, const Element& element)
372 {
373 globalValuedLocalBasis_.bind(localValuedLFE.localBasis(), element);
374 globalValuedLocalInterpolation_.bind(localValuedLFE.localInterpolation(), element);
375 localValuedLFE_ = &localValuedLFE;
376 }
377
380 const typename Traits::LocalBasisType& localBasis() const
381 {
382 return globalValuedLocalBasis_;
383 }
384
387 const typename Traits::LocalCoefficientsType& localCoefficients() const
388 {
389 return localValuedLFE_->localCoefficients();
390 }
391
394 const typename Traits::LocalInterpolationType& localInterpolation() const
395 {
396 return globalValuedLocalInterpolation_;
397 }
398
400 std::size_t size() const
401 {
402 return localValuedLFE_->size();
403 }
404
407 GeometryType type() const
408 {
409 return localValuedLFE_->type();
410 }
411
412 private:
413
414 typename Traits::LocalBasisType globalValuedLocalBasis_;
415 typename Traits::LocalInterpolationType globalValuedLocalInterpolation_;
416 const LocalValuedLFE* localValuedLFE_;
417 };
418
419} // namespace Dune::Functions::Impl
420
421#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
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
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:279
auto sparseRange(Range &&range)
Allow structured-binding for-loops for sparse iterators.
Definition: rangeutilities.hh:722
Some useful basic math stuff.
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)