27#ifndef EWOMS_PFF_GRID_VECTOR_HH
28#define EWOMS_PFF_GRID_VECTOR_HH
30#include <dune/grid/common/mcmgmapper.hh>
31#include <dune/common/version.hh>
48template <
class Gr
idView,
class Stencil,
class Data,
class DofMapper>
51 using Element =
typename GridView::template Codim<0>::Entity;
53 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
58 , elementMapper_(gridView_,
Dune::mcmgElementLayout())
59 , dofMapper_(dofMapper)
62 template <
class DistFn>
65 const unsigned numElements = gridView_.size(0);
66 const unsigned numLocalDofs = computeNumLocalDofs_();
68 elemData_.resize(numElements);
69 data_.resize(numLocalDofs);
73 Data* curElemDataPtr = &data_[0];
74 Stencil stencil(gridView_, dofMapper_);
75 for (
const auto& elem : elements(gridView_)) {
77 const unsigned elemIdx = elementMapper_.index(elem);
78 elemData_[elemIdx] = curElemDataPtr;
81 const unsigned numDof = stencil.numDof();
82 for (
unsigned localDofIdx = 0; localDofIdx < numDof; ++ localDofIdx) {
83 distFn(curElemDataPtr[localDofIdx], stencil, localDofIdx);
88 curElemDataPtr += numDof;
94 const unsigned elemIdx = elementMapper_.index(elem);
101 const Data&
get(
const Element& elem,
unsigned localDofIdx)
const
103 const unsigned elemIdx = elementMapper_.index(elem);
104 return elemData_[elemIdx][localDofIdx];
108 unsigned computeNumLocalDofs_()
const
113 Stencil stencil(gridView_, dofMapper_);
114 for (
const auto& elem : elements(gridView_)) {
115 stencil.update(elem);
116 result += stencil.numDof();
123 ElementMapper elementMapper_;
124 const DofMapper& dofMapper_;
125 std::vector<Data> data_;
126 std::vector<Data*> elemData_;
A random-access container which stores data attached to a grid's degrees of freedom in a prefetch fri...
Definition: pffgridvector.hh:50
void prefetch(const Element &elem) const
Definition: pffgridvector.hh:92
const Data & get(const Element &elem, unsigned localDofIdx) const
Definition: pffgridvector.hh:101
PffGridVector(const GridView &gridView, const DofMapper &dofMapper)
Definition: pffgridvector.hh:56
void update(const DistFn &distFn)
Definition: pffgridvector.hh:63
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:39
void prefetch(const T &val, unsigned n=1)
Template function which emits prefetch instructions for a range of memory.
Definition: prefetch.hh:39