DUNE-FUNCTIONS (2.8)

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
9#include <dune/common/deprecated.hh>
10#include <dune/common/fmatrix.hh>
11#include <dune/common/fvector.hh>
12#include <dune/common/math.hh>
13#include <dune/common/rangeutilities.hh>
14
15#include <dune/geometry/referenceelements.hh>
16
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19#include <dune/localfunctions/common/localinterpolation.hh>
20
21namespace Dune::Functions::Impl
22{
23
37 struct ContravariantPiolaTransformator
38 {
43 template<typename Values, typename LocalCoordinate, typename Geometry>
44 static auto apply(Values& values,
45 const LocalCoordinate& xi,
46 const Geometry& geometry)
47 {
48 auto jacobianTransposed = geometry.jacobianTransposed(xi);
49 auto integrationElement = geometry.integrationElement(xi);
50
51 for (auto& value : values)
52 {
53 auto tmp = value;
54 jacobianTransposed.mtv(tmp, value);
55 value /= integrationElement;
56 }
57 }
58
68 template<typename Gradients, typename LocalCoordinate, typename Geometry>
69 static auto applyJacobian(Gradients& gradients,
70 const LocalCoordinate& xi,
71 const Geometry& geometry)
72 {
73 auto jacobianTransposed = geometry.jacobianTransposed(xi);
74 auto integrationElement = geometry.integrationElement(xi);
75 for (auto& gradient : gradients)
76 {
77 auto tmp = gradient;
78 gradient = 0;
79 for (size_t k=0; k<gradient.M(); k++)
80 for (size_t l=0; l<tmp.N(); l++)
81 // Use sparseRange because jacobianTransposed may be a sparse DiagonalMatrix
82 for(auto&& [jacobianTransposed_l_j, j] : sparseRange(jacobianTransposed[l]))
83 gradient[j][k] += jacobianTransposed_l_j * tmp[l][k];
84 gradient /= integrationElement;
85 }
86 }
87
95 template<class Function, class LocalCoordinate, class Element>
96 class LocalValuedFunction
97 {
98 const Function& f_;
99 const Element& element_;
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&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
110 auto globalValue = f(xi);
111
112 // Apply the inverse Piola transform
113 auto jacobianInverseTransposed = element_.geometry().jacobianInverseTransposed(xi);
114 auto integrationElement = element_.geometry().integrationElement(xi);
115
116 auto localValue = globalValue;
117 jacobianInverseTransposed.mtv(globalValue, localValue);
118 localValue *= integrationElement;
119
120 return localValue;
121 }
122 };
123 };
124
138 struct CovariantPiolaTransformator
139 {
144 template<typename Values, typename LocalCoordinate, typename Geometry>
145 static auto apply(Values& values,
146 const LocalCoordinate& xi,
147 const Geometry& geometry)
148 {
149 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
150
151 for (auto& value : values)
152 {
153 auto tmp = value;
154 jacobianInverseTransposed.mv(tmp, value);
155 }
156 }
157
167 template<typename Gradients, typename LocalCoordinate, typename Geometry>
168 static auto applyJacobian(Gradients& gradients,
169 const LocalCoordinate& xi,
170 const Geometry& geometry)
171 {
172 auto jacobianInverseTransposed = geometry.jacobianInverseTransposed(xi);
173
174 for (auto& gradient : gradients)
175 {
176 auto tmp = gradient;
177 gradient = 0;
178 for (size_t j=0; j<gradient.N(); j++)
179 for (size_t k=0; k<gradient.M(); k++)
180 // Use sparseRange because jacobianTransposed may be a sparse DiagonalMatrix
181 for(auto&& [jacobianInverseTransposed_j_l, l] : sparseRange(jacobianInverseTransposed[j]))
182 gradient[j][k] += jacobianInverseTransposed_j_l * tmp[l][k];
183 }
184 }
185
193 template<class Function, class LocalCoordinate, class Element>
194 class LocalValuedFunction
195 {
196 const Function& f_;
197 const Element& element_;
198
199 public:
200
201 LocalValuedFunction(const Function& f, const Element& element)
202 : f_(f), element_(element)
203 {}
204
205 auto operator()(const LocalCoordinate& xi) const
206 {
207 auto&& f = Dune::Impl::makeFunctionWithCallOperator<LocalCoordinate>(f_);
208 auto globalValue = f(xi);
209
210 // Apply the inverse Piola transform
211 auto jacobianTransposed = element_.geometry().jacobianTransposed(xi);
212
213 auto localValue = globalValue;
214 jacobianTransposed.mv(globalValue, localValue);
215
216 return localValue;
217 }
218 };
219 };
220
227 template<class Transformator, class LocalValuedLocalBasis, class Element>
228 class GlobalValuedLocalBasis
229 {
230 public:
231 using Traits = typename LocalValuedLocalBasis::Traits;
232
235 void bind(const LocalValuedLocalBasis& localValuedLocalBasis, const Element& element)
236 {
237 localValuedLocalBasis_ = &localValuedLocalBasis;
238 element_ = &element;
239 }
240
243 auto size() const
244 {
245 return localValuedLocalBasis_->size();
246 }
247
249 void evaluateFunction(const typename Traits::DomainType& x,
250 std::vector<typename Traits::RangeType>& out) const
251 {
252 localValuedLocalBasis_->evaluateFunction(x,out);
253
254 Transformator::apply(out, x, element_->geometry());
255 }
256
262 void evaluateJacobian(const typename Traits::DomainType& x,
263 std::vector<typename Traits::JacobianType>& out) const
264 {
265 localValuedLocalBasis_->evaluateJacobian(x,out);
266
267 Transformator::applyJacobian(out, x, element_->geometry());
268 }
269
276 void partial(const std::array<unsigned int,2>& order,
277 const typename Traits::DomainType& x,
278 std::vector<typename Traits::RangeType>& out) const
279 {
280 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
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));
285 out.resize(size());
286
287 // TODO: The following is wasteful: We compute the full Jacobian and then return
288 // only a part of it. While we need the full Jacobian of the underlying local-valued LFE,
289 // it should be possible to compute only a partial Piola transform for the requested
290 // partial derivatives.
291 std::vector<typename Traits::JacobianType> fullJacobian;
292 localValuedLocalBasis_->evaluateJacobian(x,fullJacobian);
293
294 Transformator::applyJacobian(fullJacobian, x, element_->geometry());
295
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];
299
300 } else
301 DUNE_THROW(NotImplemented, "Partial derivatives of order 2 or higher");
302 }
303
305 auto order() const
306 {
307 return localValuedLocalBasis_->order();
308 }
309
310 const LocalValuedLocalBasis* localValuedLocalBasis_;
311 const Element* element_;
312 };
313
322 template<class Transformator, class LocalValuedLocalInterpolation, class Element>
323 class GlobalValuedLocalInterpolation
324 {
325 public:
328 void bind(const LocalValuedLocalInterpolation& localValuedLocalInterpolation, const Element& element)
329 {
330 localValuedLocalInterpolation_ = &localValuedLocalInterpolation;
331 element_ = &element;
332 }
333
334 template<typename F, typename C>
335 void interpolate (const F& f, std::vector<C>& out) const
336 {
337 using LocalCoordinate = typename Element::Geometry::LocalCoordinate;
338 typename Transformator::template LocalValuedFunction<F,LocalCoordinate,Element> localValuedFunction(f, *element_);
339 localValuedLocalInterpolation_->interpolate(localValuedFunction, out);
340 }
341
342 private:
343 const LocalValuedLocalInterpolation* localValuedLocalInterpolation_;
344 const Element* element_;
345 };
346
347
354 template<class Transformator, class LocalValuedLFE, class Element>
355 class GlobalValuedLocalFiniteElement
356 {
357 using LocalBasis = GlobalValuedLocalBasis<Transformator,
358 typename LocalValuedLFE::Traits::LocalBasisType,
359 Element>;
360 using LocalInterpolation = GlobalValuedLocalInterpolation<Transformator,
361 typename LocalValuedLFE::Traits::LocalInterpolationType,
362 Element>;
363
364 public:
367 using Traits = LocalFiniteElementTraits<LocalBasis,
368 typename LocalValuedLFE::Traits::LocalCoefficientsType,
369 LocalInterpolation>;
370
371 GlobalValuedLocalFiniteElement() {}
372
373 void bind(const LocalValuedLFE& localValuedLFE, const Element& element)
374 {
375 globalValuedLocalBasis_.bind(localValuedLFE.localBasis(), element);
376 globalValuedLocalInterpolation_.bind(localValuedLFE.localInterpolation(), element);
377 localValuedLFE_ = &localValuedLFE;
378 }
379
382 const typename Traits::LocalBasisType& localBasis() const
383 {
384 return globalValuedLocalBasis_;
385 }
386
389 const typename Traits::LocalCoefficientsType& localCoefficients() const
390 {
391 return localValuedLFE_->localCoefficients();
392 }
393
396 const typename Traits::LocalInterpolationType& localInterpolation() const
397 {
398 return globalValuedLocalInterpolation_;
399 }
400
402 std::size_t size() const
403 {
404 return localValuedLFE_->size();
405 }
406
409 GeometryType type() const
410 {
411 return localValuedLFE_->type();
412 }
413
414 private:
415
416 typename Traits::LocalBasisType globalValuedLocalBasis_;
417 typename Traits::LocalInterpolationType globalValuedLocalInterpolation_;
418 const LocalValuedLFE* localValuedLFE_;
419 };
420
421} // namespace Dune::Functions::Impl
422
423#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_GLOBALVALUEDLOCALFINITEELEMENT_HH
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)