DUNE PDELab (2.7)

dunefunctionsgridfunctionspace.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_PDELAB_GRIDFUNCTIONSPACE_DUNEFUNCTIONSGRIDFUNCTIONSPACE_HH
4#define DUNE_PDELAB_GRIDFUNCTIONSPACE_DUNEFUNCTIONSGRIDFUNCTIONSPACE_HH
5
6#include <cstddef>
7#include <map>
8#include <bitset>
9
12
13#include <dune/typetree/leafnode.hh>
14#include <dune/typetree/compositenode.hh>
15
16#include <dune/pdelab/gridfunctionspace/gridfunctionspace.hh>
17#include <dune/pdelab/gridfunctionspace/dunefunctionslocalfunctionspace.hh>
18#include <dune/pdelab/gridfunctionspace/dunefunctionslfsindexcache.hh>
19
20#include <dune/pdelab/backend/istl/dunefunctions.hh>
21#include <dune/pdelab/backend/istl.hh>
22
23namespace Dune {
24 namespace PDELab {
25
29
30 namespace Experimental {
31
32 // The following code recognizes whether the given VectorBackend (VBE) is an ISTL backend.
33 // If this is the case, then we need to replace it by ISTL::SimpleVectorBackend,
34 // because we cannot handle anything more complicated at the moment.
35 template<typename VBE>
36 VBE* registerDuneFunctionsCompatibleVBE(VBE*);
37
38 template<std::size_t block_size>
39 ISTL::SimpleVectorBackend<(block_size > 0 ? block_size : 1)>* registerDuneFunctionsCompatibleVBE(ISTL::VectorBackend<ISTL::Blocking::none,block_size>*);
40
41 template<typename VBE>
42 using DuneFunctionsCompatibleVBE = std::decay_t<decltype(*registerDuneFunctionsCompatibleVBE(std::declval<VBE*>()))>;
43
56 template<typename DFBasis, typename VBE, typename CE>
58 : public TypeTree::LeafNode
60 , public DataHandleProvider<GridFunctionSpace<DFBasis,VBE,CE> >
61 {
62 using GV = typename DFBasis::GridView;
63
64 template<typename,typename>
65 friend class GridFunctionSpaceBase;
66
67 public:
69
70 struct Traits {
71
72 using GridView = Dune::PDELab::impl::GridView<typename DFBasis::GridView>;
73 using GridViewType = GridView; // DiscreteGridFunction wants this
74 using EntitySet = Dune::PDELab::impl::EntitySet<typename DFBasis::GridView>;
75
76 using size_type = std::size_t;
77 using SizeType = size_type;
78 using ConstraintsType = CE;
79
80 using FiniteElementType = typename DFBasis::LocalView::Tree::FiniteElement; // DiscreteGridFunction wants this
81
82 using Basis = DFBasis;
83
84 // The following code recognizes whether the given VectorBackend (VBE) is an ISTL backend.
85 // If this is the case, then we replace it by ISTL::SimpleVectorBackend,
86 // because we cannot handle anything more complicated at the moment.
87 using Backend = DuneFunctionsCompatibleVBE<VBE>;
88
90 struct FEM
91 {
92 struct Traits
93 {
94 using FiniteElement = typename DFBasis::LocalView::Tree::FiniteElement;
95 using FiniteElementType = FiniteElement;
96 };
97
98 FEM(const std::shared_ptr<DFBasis>& basis)
99 : _basis(basis)
100 {}
101
112 const typename Traits::FiniteElementType& find (const typename GridView::template Codim<0>::Entity& element) const
113 {
114 auto type = element.type();
115 auto mapEntry = geometryTypeToLocalView_.find(type);
116 if (mapEntry == geometryTypeToLocalView_.end())
117 {
118 auto newLocalView = std::make_shared<typename DFBasis::LocalView>(_basis->localView());
119 newLocalView->bind(element);
120 auto insertedLocalView = geometryTypeToLocalView_.insert(std::make_pair(type, newLocalView));
121 return insertedLocalView.first->second->tree().finiteElement();
122 }
123 else
124 {
125 return mapEntry->second->tree().finiteElement();
126 }
127 }
128
129 void update()
130 {
131 geometryTypeToLocalView_.clear();
132 }
133
134 private:
135 const std::shared_ptr<DFBasis> _basis;
136
137 mutable std::map<GeometryType, std::shared_ptr<typename DFBasis::LocalView> > geometryTypeToLocalView_;
138 };
139
140 using FiniteElementMap = FEM;
141 using FiniteElementMapType = FEM;
142
143 };
144
145 using Basis = DFBasis;
146
155 : public TypeTree::LeafNode
156 {
157
158 struct Traits {
159
162
170 using size_type = std::size_t;
171 using SizeType = size_type;
172
173 using DOFIndexAccessor = Dune::PDELab::DefaultDOFIndexAccessor;
174 };
175
176 using DOFIndex = typename Traits::DOFIndex;
177 using ContainerIndex = typename Traits::ContainerIndex;
178 using size_type = std::size_t;
179
181 : _gfs(gfs)
182 {
183 update();
184 }
185
186 size_type size() const
187 {
188 return _gfs.basis().size();
189 }
190
192 size_type size(const typename DOFIndex::EntityIndex& entity) const
193 {
194 return _containerIndices[entity[0]][entity[1]].size();
195 }
196
197 size_type maxLocalSize() const
198 {
199 return _gfs.basis().localView().maxSize();
200 }
201
204 bool contains(typename Traits::SizeType codim) const
205 {
206 return _contains[codim];
207 }
208
211 bool fixedSize(typename Traits::SizeType codim) const
212 {
213 return _fixedSize[codim];
214 }
215
216 // child_index: Steffen sagt: unklar, im Zweifel einfach ignorieren
217 template<typename CIOutIterator, typename DIOutIterator = DummyDOFIndexIterator>
218 typename Traits::SizeType
219 extract_entity_indices(const typename Traits::DOFIndex::EntityIndex& entityIndex,
220 typename Traits::SizeType child_index,
221 CIOutIterator ci_out, const CIOutIterator ci_end,
222 DIOutIterator dummy) const
223 {
224 for (size_type i=0; i<_containerIndices[entityIndex[0]][entityIndex[1]].size(); i++)
225 {
226 *ci_out = _containerIndices[entityIndex[0]][entityIndex[1]][i];
227 ci_out++;
228 }
229
230 return _containerIndices[entityIndex[0]][entityIndex[1]].size();
231 }
232
233 ContainerIndex containerIndex(const DOFIndex& i) const
234 {
235 return _containerIndices[i.entityIndex()[0]][i.entityIndex()[1]][i.treeIndex()[0]];
236 }
237
238 void update()
239 {
240 _containerIndices.clear();
241
242 constexpr auto dim = GV::dimension;
243 const auto gridView = _gfs.gridView();
244 const auto& indexSet = gridView.indexSet();
245
246 // Count how many dofs there are for each individual entity
247 std::vector<std::vector<size_type> > dofsPerEntity(GlobalGeometryTypeIndex::size(dim));
248 for (size_type codim=0; codim<=dim; codim++)
249 for (auto&& type : indexSet.types(codim))
250 {
251 dofsPerEntity[GlobalGeometryTypeIndex::index(type)].resize(gridView.size(type));
252 std::fill(dofsPerEntity[GlobalGeometryTypeIndex::index(type)].begin(),
253 dofsPerEntity[GlobalGeometryTypeIndex::index(type)].end(), 0);
254 }
255
256 typename DFBasis::LocalView localView = _gfs.basis().localView();
257 for (auto&& element : elements(gridView))
258 {
259 localView.bind(element);
260 const auto refElement = ReferenceElements<typename GV::ctype,dim>::general(element.type());
261
262 const auto& localFiniteElement = localView.tree().finiteElement();
263
264 for (size_type i=0; i<localFiniteElement.size(); i++)
265 {
266 const auto& localKey = localFiniteElement.localCoefficients().localKey(i);
267
268 // Type of the entity that the current local dof is attached to
269 auto subentityTypeIndex = GlobalGeometryTypeIndex::index(refElement.type(localKey.subEntity(), localKey.codim()));
270
271 // Global index of the entity that the current local dof is attached to
272 auto subentityIndex = indexSet.subIndex(element, localKey.subEntity(), localKey.codim());
273
274 dofsPerEntity[subentityTypeIndex][subentityIndex]
275 = std::max(dofsPerEntity[subentityTypeIndex][subentityIndex], (size_type)localKey.index()+1);
276 }
277 }
278
279 // Set up the nested container for the container indices
280 _containerIndices.resize(GlobalGeometryTypeIndex::size(dim));
281
282 for (size_type codim=0; codim<=dim; codim++)
283 for (auto&& type : indexSet.types(codim))
284 {
285 _containerIndices[GlobalGeometryTypeIndex::index(type)].resize(gridView.size(type));
286 for (size_type i=0; i<_containerIndices[GlobalGeometryTypeIndex::index(type)].size(); i++)
287 _containerIndices[GlobalGeometryTypeIndex::index(type)][i].resize(dofsPerEntity[GlobalGeometryTypeIndex::index(type)][i]);
288 }
289
290 // Actually set the container indices for all dofs on all entities
291 for (auto&& element : elements(gridView))
292 {
293 localView.bind(element);
294 const auto refElement = ReferenceElements<typename GV::ctype,dim>::general(element.type());
295
296 const auto& localFiniteElement = localView.tree().finiteElement();
297
298 for (size_type i=0; i<localFiniteElement.size(); i++)
299 {
300 const auto& localKey = localFiniteElement.localCoefficients().localKey(i);
301
302 // Type of the entity that the current local dof is attached to
303 GeometryType subentityType = refElement.type(localKey.subEntity(), localKey.codim());
304
305 // Global index of the entity that the current local dof is attached to
306 auto subentityIndex = indexSet.subIndex(element, localKey.subEntity(), localKey.codim());
307
308 _containerIndices[GlobalGeometryTypeIndex::index(subentityType)][subentityIndex][localKey.index()].set({localView.index(i)});
309 }
310
311 }
312
313 // Precompute for each codimension whether there is at least one entity with degrees of freedom,
314 // and whether all entities have the same number of dofs
315 for (size_type codim=0; codim<=dim; codim++)
316 {
317 _contains[codim] = false;
318 _fixedSize[codim] = true;
319 for (auto&& type : indexSet.types(codim))
320 {
321 const auto& dofs = dofsPerEntity[GlobalGeometryTypeIndex::index(type)];
322
323 if (dofs[0] > 0)
324 _contains[codim] = true;
325
326 for (size_type i=1; i<dofs.size(); i++)
327 {
328 if (dofs[i] > 0)
329 _contains[codim] = true;
330
331 if (dofs[i-1] != dofs[i])
332 _fixedSize[codim] = false;
333 }
334 }
335 }
336 }
337
338 private:
339
340 const GridFunctionSpace& _gfs;
341
342 // Container that contains the ContainerIndices for all dofs, accessible by entities
343 std::vector<std::vector<std::vector<ContainerIndex> > > _containerIndices;
344
348 std::bitset<GV::dimension+1> _contains;
349
352 std::bitset<GV::dimension+1> _fixedSize;
353 };
354
361 struct Ordering
362 : public TypeTree::CompositeNode<LeafOrdering>
363 {
364 friend class LocalFunctionSpace<GridFunctionSpace>;
365
366 using Traits = typename LeafOrdering::Traits;
367
368 static const bool consume_tree_index = false;
369
370 using DOFIndex = typename Traits::DOFIndex;
371 using ContainerIndex = typename Traits::ContainerIndex;
372 using size_type = std::size_t;
373
374 using CacheTag = DuneFunctionsCacheTag;
376
377 Ordering(const GridFunctionSpace& gfs)
379 {}
380
381 size_type size() const
382 {
383 return this->child(Indices::_0).size();
384 }
385
388 size_type blockCount() const
389 {
390 return this->child(Indices::_0).size();
391 }
392
393 size_type maxLocalSize() const
394 {
395 return this->child(Indices::_0).maxLocalSize();
396 }
397
401 bool contains(typename Traits::SizeType codim) const
402 {
403 return this->child(Indices::_0).contains(codim);
404 }
405
408 bool fixedSize(typename Traits::SizeType codim) const
409 {
410 return this->child(Indices::_0).fixedSize(codim);
411 }
412
413 template<typename CIOutIterator, typename DIOutIterator = DummyDOFIndexIterator>
414 typename Traits::SizeType
415 extract_entity_indices(const typename Traits::DOFIndex::EntityIndex& ei,
416 typename Traits::SizeType child_index,
417 CIOutIterator ci_out, const CIOutIterator ci_end) const
418 {
419 return 0;
420 }
421
422 void update()
423 {
424 this->child(Indices::_0).update();
425 }
426
427 private:
428
429 ContainerIndex containerIndex(const DOFIndex& i) const
430 {
431 return this->child(Indices::_0).containerIndex(i);
432 }
433
434 };
435
436
438 template<typename E>
440 {
441
443 using Type = std::conditional_t<
444 std::is_same<
445 CE,
446 NoConstraints
447 >::value,
448 EmptyTransformation,
450 >;
451
452 private:
454 };
455
456 // ****************************************************************************************************
457 // Construct from a dune-functions basis
458 // ****************************************************************************************************
459
461 GridFunctionSpace (std::shared_ptr<DFBasis> df_basis, std::shared_ptr<CE> ce)
462 : _es(df_basis->gridView(), Traits::EntitySet::allCodims())
463 , _df_basis(std::move(df_basis))
464 , _finiteElementMap(_df_basis)
465 , _pce(std::move(ce))
466 , _ordering(*this)
467 {}
468
469 GridFunctionSpace (std::shared_ptr<DFBasis> df_basis)
470 : _es(df_basis->gridView(), Traits::EntitySet::allCodims())
471 , _df_basis(std::move(df_basis))
472 , _finiteElementMap(_df_basis)
473 , _pce(std::make_shared<CE>())
474 , _ordering(*this)
475 {}
476
478 const typename Traits::GridView& gridView () const
479 {
480 return _es.gridView();
481 }
482
484 const typename Traits::EntitySet& entitySet () const
485 {
486 return _es;
487 }
488
490 const auto& finiteElementMap () const
491 {
492 return _finiteElementMap;
493 }
494
496 const typename Traits::ConstraintsType& constraints () const
497 {
498 return *_pce;
499 }
500
502 std::shared_ptr<const CE> constraintsStorage () const
503 {
504 return _pce;
505 }
506
508 const Ordering& ordering() const
509 {
510 return _ordering;
511 }
512
513 typename Traits::SizeType size() const
514 {
515 return _ordering.size();
516 }
517
518 typename Traits::SizeType blockCount() const
519 {
520 return _ordering.blockCount();
521 }
522
523 typename Traits::SizeType globalSize() const
524 {
525 return _ordering.size();
526 }
527
528 typename Traits::SizeType maxLocalSize () const
529 {
530 return _ordering.maxLocalSize();
531 }
532
538 void update(bool force = false)
539 {
540 _es.update(force);
541 _df_basis->update(_es.gridView());
542 _finiteElementMap.update();
543 // Apparently there is no need to update the constraints assembler '_pce';
544 _ordering.update();
545 }
546
547 const std::string& name() const
548 {
549 return _name;
550 }
551
552 void name(const std::string& name)
553 {
554 _name = name;
555 }
556
557 bool isRootSpace() const
558 {
559 return true;
560 }
561
562 const Basis& basis() const
563 {
564 return *_df_basis;
565 }
566
567 private:
568
569 typename Traits::EntitySet _es;
570 std::shared_ptr<DFBasis> _df_basis;
571 typename Traits::FiniteElementMap _finiteElementMap;
572 std::shared_ptr<CE const> _pce;
573 Ordering _ordering;
574 std::string _name;
575 };
576
577 } // namespace Experimental
578
579 } // namespace PDELab
580} // namespace Dune
581
582#endif // DUNE_PDELAB_GRIDFUNCTIONSPACE_DUNEFUNCTIONSGRIDFUNCTIONSPACE_HH
static constexpr std::size_t index(const GeometryType &gt)
Compute the index for the given geometry type over all dimensions.
Definition: typeindex.hh:133
static constexpr std::size_t size(std::size_t maxdim)
Compute total number of geometry types up to and including the given dimension.
Definition: typeindex.hh:120
a class holding transformation for constrained spaces
Definition: constraintstransformation.hh:20
A multi-index representing a degree of freedom in a GridFunctionSpace.
Definition: dofindex.hh:148
EntityIndex & entityIndex()
Returns the index of the grid entity associated with the DOF.
Definition: dofindex.hh:258
TreeIndex & treeIndex()
Returns the tuple of entity-local indices associated with the DOF.
Definition: dofindex.hh:270
A pdelab grid function space implemented by a dune-functions function space basis.
Definition: dunefunctionsgridfunctionspace.hh:61
Mixin base class for specifying output hints to I/O routines like VTK.
Definition: function.hh:126
A class for representing multi-indices.
Definition: multiindex.hh:29
Base class for composite nodes based on variadic templates.
Definition: compositenode.hh:25
Child< k >::Type & child(index_constant< k >={})
Returns the i-th child.
Definition: compositenode.hh:82
Base class for leaf nodes in a dune-typetree.
Definition: leafnode.hh:25
A few common exception classes.
Traits for type conversions and type information.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:180
constexpr index_constant< 0 > _0
Compile time index with value 0.
Definition: indices.hh:51
size_type blockCount() const
Same as size(), because block size is always 1.
Definition: dunefunctionsgridfunctionspace.hh:388
const Ordering & ordering() const
Direct access to the DOF ordering.
Definition: dunefunctionsgridfunctionspace.hh:508
size_type size(const typename DOFIndex::EntityIndex &entity) const
Number of degrees of freedom per entity.
Definition: dunefunctionsgridfunctionspace.hh:192
const Traits::GridView & gridView() const
get grid view
Definition: dunefunctionsgridfunctionspace.hh:478
bool fixedSize(typename Traits::SizeType codim) const
True if for all entities of the given codim the same number of data items has to be communicated.
Definition: dunefunctionsgridfunctionspace.hh:408
bool fixedSize(typename Traits::SizeType codim) const
True if all entities of the given codimension have the same number of dofs.
Definition: dunefunctionsgridfunctionspace.hh:211
const Traits::FiniteElementType & find(const typename GridView::template Codim< 0 >::Entity &element) const
Get local basis functions for entity.
Definition: dunefunctionsgridfunctionspace.hh:112
std::shared_ptr< const CE > constraintsStorage() const
return storage of constraints engine
Definition: dunefunctionsgridfunctionspace.hh:502
const Traits::EntitySet & entitySet() const
get EntitySet
Definition: dunefunctionsgridfunctionspace.hh:484
bool contains(typename Traits::SizeType codim) const
Returns true if there is at least one entity of the given codim for which data needs to be communicat...
Definition: dunefunctionsgridfunctionspace.hh:401
const auto & finiteElementMap() const
get finite element map
Definition: dunefunctionsgridfunctionspace.hh:490
GridFunctionSpace(std::shared_ptr< DFBasis > df_basis, std::shared_ptr< CE > ce)
constructor
Definition: dunefunctionsgridfunctionspace.hh:461
std::conditional_t< std::is_same< CE, NoConstraints >::value, EmptyTransformation, ConstraintsTransformation< typename Ordering::Traits::DOFIndex, typename Ordering::Traits::ContainerIndex, E > > Type
define Type as the Type of a container of E's
Definition: dunefunctionsgridfunctionspace.hh:450
const Traits::ConstraintsType & constraints() const
return constraints engine
Definition: dunefunctionsgridfunctionspace.hh:496
void update(bool force=false)
Update the indexing information of the GridFunctionSpace.
Definition: dunefunctionsgridfunctionspace.hh:538
bool contains(typename Traits::SizeType codim) const
True if there is at least one entity of the given codim that has a dof.
Definition: dunefunctionsgridfunctionspace.hh:204
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Dune namespace.
Definition: alignedallocator.hh:14
STL namespace.
Static tag representing a codimension.
Definition: dimension.hh:22
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:196
extract type for storing constraints
Definition: dunefunctionsgridfunctionspace.hh:440
The actual Ordering object of the grid function space.
Definition: dunefunctionsgridfunctionspace.hh:156
Root of the ordering tree.
Definition: dunefunctionsgridfunctionspace.hh:363
Rudimentary internal implementation of a FiniteElementMap.
Definition: dunefunctionsgridfunctionspace.hh:91
export Traits class
Definition: dunefunctionsgridfunctionspace.hh:70
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)