3#ifndef DUNE_PDELAB_GRIDFUNCTIONSPACE_ENTITYINDEXCACHE_HH
4#define DUNE_PDELAB_GRIDFUNCTIONSPACE_ENTITYINDEXCACHE_HH
8#include <dune/typetree/utility.hh>
10#include <dune/pdelab/gridfunctionspace/lfsindexcache.hh>
16 template<
typename GFS,
bool map_dof_indices = false>
17 class EntityIndexCache
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;
28 typedef std::size_t size_type;
30 typedef std::vector<CI> CIVector;
31 typedef std::vector<DI> DIVector;
39 typedef std::array<size_type,leaf_count + 1> Offsets;
41 EntityIndexCache(
const GFS& gfs)
43 , _container_indices(gfs.maxLocalSize())
44 , _dof_indices(map_dof_indices ? gfs.maxLocalSize() : 0)
46 std::fill(_offsets.begin(),_offsets.end(),0);
49 template<
typename Entity>
50 void update(
const Entity& e)
52 std::fill(_offsets.begin(),_offsets.end(),0);
56 _gfs.dataHandleIndices(e,_container_indices,_dof_indices,_offsets.begin(),std::integral_constant<bool,map_dof_indices>());
59 const DI& dofIndex(size_type i)
const
61 assert(map_dof_indices);
62 return _dof_indices[i];
65 const CI& containerIndex(size_type i)
const
67 return _container_indices[i];
70 const GridFunctionSpace& gridFunctionSpace()
const
75 size_type
size()
const
77 return _offsets[leaf_count];
80 const Offsets& offsets()
const
88 CIVector _container_indices;
89 DIVector _dof_indices;
static constexpr int codimension
Know your own codimension.
Definition: entity.hh:106
Dune namespace.
Definition: alignedallocator.hh:13
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
static const std::size_t leafCount
The number of leaf nodes in the TypeTree.
Definition: utility.hh:136