DUNE PDELab (2.8)

entityindexcache.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_ENTITYINDEXCACHE_HH
4#define DUNE_PDELAB_GRIDFUNCTIONSPACE_ENTITYINDEXCACHE_HH
5
6#include <array>
7
8#include <dune/typetree/utility.hh>
9
10#include <dune/pdelab/gridfunctionspace/lfsindexcache.hh>
11
12namespace Dune {
13 namespace PDELab {
14
15
16 template<typename GFS, bool map_dof_indices = false>
17 class EntityIndexCache
18 {
19
20 public:
21
22 typedef GFS GridFunctionSpace;
23 typedef typename GFS::Ordering Ordering;
24 typedef typename Ordering::Traits::ContainerIndex ContainerIndex;
25 typedef ContainerIndex CI;
26 typedef typename Ordering::Traits::DOFIndex DOFIndex;
27 typedef DOFIndex DI;
28 typedef std::size_t size_type;
29
30 typedef std::vector<CI> CIVector;
31 typedef std::vector<DI> DIVector;
32
33 private:
34
35 static const size_type leaf_count = TypeTree::TreeInfo<Ordering>::leafCount;
36
37 public:
38
39 typedef std::array<size_type,leaf_count + 1> Offsets;
40
41 EntityIndexCache(const GFS& gfs)
42 : _gfs(gfs)
43 , _container_indices(gfs.maxLocalSize())
44 , _dof_indices(map_dof_indices ? gfs.maxLocalSize() : 0)
45 {
46 std::fill(_offsets.begin(),_offsets.end(),0);
47 }
48
49 template<typename Entity>
50 void update(const Entity& e)
51 {
52 std::fill(_offsets.begin(),_offsets.end(),0);
53 if (!_gfs.dataHandleContains(Entity::codimension))
54 return;
55
56 _gfs.dataHandleIndices(e,_container_indices,_dof_indices,_offsets.begin(),std::integral_constant<bool,map_dof_indices>());
57 }
58
59 const DI& dofIndex(size_type i) const
60 {
61 assert(map_dof_indices);
62 return _dof_indices[i];
63 }
64
65 const CI& containerIndex(size_type i) const
66 {
67 return _container_indices[i];
68 }
69
70 const GridFunctionSpace& gridFunctionSpace() const
71 {
72 return _gfs;
73 }
74
75 size_type size() const
76 {
77 return _offsets[leaf_count];
78 }
79
80 const Offsets& offsets() const
81 {
82 return _offsets;
83 }
84
85 private:
86
87 const GFS& _gfs;
88 CIVector _container_indices;
89 DIVector _dof_indices;
90 Offsets _offsets;
91
92 };
93
94 } // namespace PDELab
95} // namespace Dune
96
97#endif // DUNE_PDELAB_GRIDFUNCTIONSPACE_ENTITYINDEXCACHE_HH
@ codimension
Know your own codimension.
Definition: entity.hh:105
Dune namespace.
Definition: alignedallocator.hh:11
static const std::size_t leafCount
The number of leaf nodes in the TypeTree.
Definition: utility.hh:183
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)