AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
DiscreteLocalFunction.inc.hpp
1#pragma once
2
3#include <amdis/Output.hpp>
4#include <amdis/algorithm/ForEach.hpp>
5#include <amdis/common/DerivativeTraits.hpp>
6#include <amdis/common/FieldMatVec.hpp>
7#include <amdis/common/Index.hpp>
8#include <amdis/common/Math.hpp>
9#include <amdis/functions/HierarchicNodeToRangeMap.hpp>
10#include <amdis/functions/NodeIndices.hpp>
11#include <amdis/functions/Order.hpp>
12#include <amdis/typetree/FiniteElementType.hpp>
13#include <amdis/utility/LocalToGlobalAdapter.hpp>
14
15#include <dune/common/ftraits.hh>
16#include <dune/common/typetree/childaccess.hh>
17#include <dune/common/typetree/nodeconcepts.hh>
18
19namespace AMDiS {
20namespace Impl {
21
22// specialization for containers with gather method
23template <class Coeff, class LV, class LC>
24auto gather(Coeff const& coeff, LV const& localView, LC& localCoeff, Dune::PriorityTag<2>)
25 -> std::void_t<decltype(coeff.gather(localView, localCoeff))>
26{
27 coeff.gather(localView, localCoeff);
28}
29
30// implementation for containers with multi-index access
31template <class Coeff, class LV, class LC,
32 class MI = typename LV::MultiIndex>
33auto gather(Coeff const& coeff, LV const& localView, LC& localCoeff, Dune::PriorityTag<1>)
34 -> std::void_t<decltype(coeff[std::declval<MI>()])>
35{
36 localCoeff.resize(localView.size());
37 auto it = localCoeff.begin();
38 for (auto const& idx : nodeIndices(localView))
39 *it++ = coeff[idx];
40}
41
42// fallback implementation
43template <class Coeff, class LV, class LC>
44auto gather(Coeff const& coeff, LV const& localView, LC& localCoeff, Dune::PriorityTag<0>)
45 -> std::void_t<decltype(coeff[index_<0>])>
46{
47 auto backend = Dune::Functions::istlVectorBackend(coeff);
48 localCoeff.resize(localView.size());
49 auto it = localCoeff.begin();
50 for (auto const& idx : nodeIndices(localView))
51 *it++ = backend[idx];
52}
53
54} // end namespace Impl
55
56
57template <class C, class GB, class TP, class R>
58 template <class Type>
59class DiscreteFunction<C const,GB,TP,R>::LocalFunction
60{
61 using DomainType = typename DiscreteFunction::Domain;
62 using RangeGenerator = typename DiscreteFunction::RangeGenerator;
63 using RangeType = typename DiscreteFunction::Range;
64
65 template <class T>
66 using DerivativeRange = typename DerivativeTraits<RangeType(DomainType), T>::Range;
67
68public:
69 using Domain = typename EntitySet::LocalCoordinate;
70 using Range = DerivativeRange<Type>;
71
72 enum { hasDerivative = std::is_same<Type, tag::value>::value };
73
74private:
75 using LocalView = typename GlobalBasis::LocalView;
76 using Element = typename EntitySet::Element;
77 using Geometry = typename Element::Geometry;
78
79public:
81 template <class LV,
82 REQUIRES(Concepts::Similar<LV,LocalView>)>
83 LocalFunction(LV&& localView, TP const& treePath, std::shared_ptr<C const> coefficients, Type type)
84 : localView_(FWD(localView))
85 , treePath_(treePath)
86 , coefficients_(std::move(coefficients))
87 , type_(type)
88 {}
89
91 void bind(Element const& element)
92 {
93 localView_.bind(element);
94 geometry_.emplace(element.geometry());
95 Impl::gather(*coefficients_, localView_, localCoefficients_, Dune::PriorityTag<4>{});
96 }
97
99 void unbind()
100 {
101 localView_.unbind();
102 geometry_.reset();
103 }
104
106 bool bound() const
107 {
108 return !!geometry_;
109 }
110
112 Range operator() (const Domain& local) const
113 {
114 assert(bound());
115 if constexpr (std::is_same_v<Type, tag::value>)
116 return evaluateFunction(local);
117 else if constexpr (std::is_same_v<Type, tag::jacobian>)
118 return evaluateJacobian(local);
119 else if constexpr (std::is_same_v<Type, tag::partial>)
120 return evaluatePartialDerivative(local);
121 else if constexpr (std::is_same_v<Type, tag::divergence>) {
122 static constexpr bool isVectorNode
123 = Dune::TypeTree::Concept::UniformInnerTreeNode<SubTree> && int(SubTree::degree()) == GridView::dimensionworld;
124 static_assert(isVectorNode);
125 if constexpr (isVectorNode)
126 return evaluateDivergence(local);
127 }
128
129 return Range{};
130 }
131
133 template <class DerType>
135 {
136 static_assert(std::is_same_v<Type, tag::value>);
137 LocalFunction<DerType> df {localView_, treePath_, coefficients_, type};
138 if (bound())
139 df.bind(localContext());
140 return df;
141 }
142
144 auto order() const
145 -> decltype(AMDiS::order(std::declval<SubTree>()))
146 {
147 assert(bound());
148 return AMDiS::order(subTree());
149 }
150
152 Element const& localContext() const
153 {
154 assert(bound());
155 return localView_.element();
156 }
157
158private:
159
160 template <class LocalCoordinate>
161 auto evaluateFunction (const LocalCoordinate& local) const;
162
163 template <class LocalCoordinate>
164 auto evaluateJacobian (const LocalCoordinate& local) const;
165
166 template <class LocalCoordinate>
167 auto evaluateDivergence (const LocalCoordinate& local) const;
168
169 template <class LocalCoordinate>
170 auto evaluatePartialDerivative (const LocalCoordinate& local) const;
171
172 // get the geometry of the bound element, that is constructed in the bind() method
173 Geometry const& geometry() const
174 {
175 assert(!!geometry_);
176 return *geometry_;
177 }
178
179 // return the sub-tree of the basis-tree associated to the tree-path
180 SubTree const& subTree() const
181 {
182 return Dune::TypeTree::child(localView_.tree(), treePath_);
183 }
184
185 // return the sub-tree of the basis-tree-cache associated to the tree-path
186 // NOTE: this is only return in case the local-view does provide a treeCache() function
187 template <class LV,
188 class = decltype(std::declval<LV>().treeCache())>
189 auto const& subTreeCache(Dune::PriorityTag<1>) const
190 {
191 return Dune::TypeTree::child(localView_.treeCache(), treePath_);
192 }
193
194 // construct a new tree-cache that stores a reference to the sub-tree
195 template <class>
196 auto subTreeCache(Dune::PriorityTag<0>) const
197 {
198 return makeNodeCache(subTree());
199 }
200
201 decltype(auto) subTreeCache() const
202 {
203 return subTreeCache<LocalView>(Dune::PriorityTag<2>());
204 }
205
206private:
207 LocalView localView_;
208 TP treePath_;
209 std::shared_ptr<Coefficients const> coefficients_;
210 Type type_;
211
212 std::optional<Geometry> geometry_;
213 std::vector<Coefficient> localCoefficients_;
214};
215
216
217template <class C, class GB, class TP, class R>
218 template <class Type>
219 template <class LocalCoordinate>
220auto DiscreteFunction<C const,GB,TP,R>::LocalFunction<Type>
221::evaluateFunction(LocalCoordinate const& local) const
222{
223 Range y = RangeGenerator::create(this->subTree());
224 Recursive::forEach(y, [](auto& v) { v = 0; });
225
226 Traversal::forEachLeafNode(this->subTreeCache(), [&](auto const& node, auto const& tp)
227 {
228 auto const& shapeFunctionValues = node.localBasisValuesAt(local);
229 std::size_t size = node.finiteElement().size();
230
231 // Get range entry associated to this node
232 auto re = Dune::Functions::flatVectorView(hierarchicNodeToRangeMap(tp, y));
233
234 for (std::size_t i = 0; i < size; ++i) {
235 // Get coefficient associated to i-th shape function
236 auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);
237
238 // Get value of i-th shape function
239 auto v = Dune::Functions::flatVectorView(shapeFunctionValues[i]);
240
241 std::size_t dimC = c.size();
242 std::size_t dimV = v.size();
243 assert(dimC*dimV == std::size_t(re.size()));
244 for(std::size_t j = 0; j < dimC; ++j) {
245 auto&& c_j = c[j];
246 for(std::size_t k = 0; k < dimV; ++k)
247 re[j*dimV + k] += c_j*v[k];
248 }
249 }
250 });
251
252 return y;
253}
254
255
256template <class C, class GB, class TP, class R>
257 template <class Type>
258 template <class LocalCoordinate>
259auto DiscreteFunction<C const,GB,TP,R>::LocalFunction<Type>
260::evaluateJacobian(LocalCoordinate const& local) const
261{
262 Range dy;
263 Recursive::forEach(dy, [](auto& v) { v = 0; });
264
265 Traversal::forEachLeafNode(this->subTreeCache(), [&](auto const& node, auto const& tp)
266 {
267 LocalToGlobalBasisAdapter localBasis(node, this->geometry());
268 auto const& gradients = localBasis.gradientsAt(local);
269
270 // Get range entry associated to this node
271 auto re = Dune::Functions::flatVectorView(hierarchicNodeToRangeMap(tp, dy));
272
273 for (std::size_t i = 0; i < localBasis.size(); ++i) {
274 // Get coefficient associated to i-th shape function
275 auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);
276
277 // Get value of i-th transformed reference gradient
278 auto grad = Dune::Functions::flatVectorView(gradients[i]);
279
280 std::size_t dimC = c.size();
281 std::size_t dimV = grad.size();
282 assert(dimC*dimV == std::size_t(re.size()));
283 for(std::size_t j = 0; j < dimC; ++j) {
284 auto&& c_j = c[j];
285 for(std::size_t k = 0; k < dimV; ++k)
286 re[j*dimV + k] += c_j*grad[k];
287 }
288 }
289 });
290
291 return dy;
292}
293
294
295template <class C, class GB, class TP, class R>
296 template <class Type>
297 template <class LocalCoordinate>
298auto DiscreteFunction<C const,GB,TP,R>::LocalFunction<Type>
299::evaluateDivergence(LocalCoordinate const& local) const
300{
301 static_assert(static_size_v<Range> == 1, "Divergence implemented for scalar coefficients only.");
302
303 Range dy;
304 Recursive::forEach(dy, [](auto& v) { v = 0; });
305
306 auto&& node = this->subTreeCache();
307
308 LocalToGlobalBasisAdapter localBasis(node.child(0), this->geometry());
309 auto const& gradients = localBasis.gradientsAt(local);
310
311 auto re = Dune::Functions::flatVectorView(dy);
312 assert(int(re.size()) == 1);
313 for (std::size_t i = 0; i < localBasis.size(); ++i) {
314 auto grad = Dune::Functions::flatVectorView(gradients[i]);
315
316 assert(int(grad.size()) == GridView::dimensionworld);
317 for (std::size_t j = 0; j < GridView::dimensionworld; ++j)
318 re[0] += localCoefficients_[node.child(j).localIndex(i)] * grad[j];
319 }
320
321 return dy;
322}
323
324
325template <class C, class GB, class TP, class R>
326 template <class Type>
327 template <class LocalCoordinate>
328auto DiscreteFunction<C const,GB,TP,R>::LocalFunction<Type>
329::evaluatePartialDerivative(LocalCoordinate const& local) const
330{
331 std::size_t comp = this->type_.comp;
332
333 Range dy;
334 Recursive::forEach(dy, [](auto& v) { v = 0; });
335
336 Traversal::forEachLeafNode(this->subTreeCache(), [&](auto const& node, auto const& tp)
337 {
338 LocalToGlobalBasisAdapter localBasis(node, this->geometry());
339 auto const& partial = localBasis.partialsAt(local, comp);
340
341 // Get range entry associated to this node
342 auto re = Dune::Functions::flatVectorView(hierarchicNodeToRangeMap(tp, dy));
343
344 for (std::size_t i = 0; i < localBasis.size(); ++i) {
345 // Get coefficient associated to i-th shape function
346 auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);
347
348 // Get value of i-th transformed reference partial_derivative
349 auto d_comp = Dune::Functions::flatVectorView(partial[i]);
350
351 std::size_t dimC = c.size();
352 std::size_t dimV = d_comp.size();
353 assert(dimC*dimV == std::size_t(re.size()));
354 for(std::size_t j = 0; j < dimC; ++j) {
355 auto&& c_j = c[j];
356 for(std::size_t k = 0; k < dimV; ++k)
357 re[j*dimV + k] += c_j*d_comp[k];
358 }
359 }
360 });
361
362 return dy;
363}
364
365
366} // end namespace AMDiS
Definition DiscreteLocalFunction.inc.hpp:60
bool bound() const
Check whether this localfunction is bound to an element.
Definition DiscreteLocalFunction.inc.hpp:106
auto order() const -> decltype(AMDiS::order(std::declval< SubTree >()))
The polynomial degree of the LocalFunctions basis-functions.
Definition DiscreteLocalFunction.inc.hpp:144
LocalFunction(LV &&localView, TP const &treePath, std::shared_ptr< C const > coefficients, Type type)
Constructor. Stores a copy of the DiscreteFunction.
Definition DiscreteLocalFunction.inc.hpp:83
void bind(Element const &element)
Bind the LocalView to the element.
Definition DiscreteLocalFunction.inc.hpp:91
Element const & localContext() const
Return the bound element.
Definition DiscreteLocalFunction.inc.hpp:152
void unbind()
Unbind the LocalView from the element.
Definition DiscreteLocalFunction.inc.hpp:99
LocalFunction< DerType > makeDerivative(DerType type) const
Create a LocalFunction representing the derivative.
Definition DiscreteLocalFunction.inc.hpp:134
A mutable view on the subspace of a DOFVector,.
Definition DiscreteFunction.hpp:45
Coefficients & coefficients()
Return the mutable DOFVector.
Definition DiscreteFunction.hpp:113
AMDiS::LocalView< Self > LocalView
Type of the local view on the restriction of the basis to a single element.
Definition GlobalBasis.hpp:60
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition Concepts.hpp:181
Definition DerivativeTraits.hpp:29