pffgridvector.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_PFF_GRID_VECTOR_HH
28#define EWOMS_PFF_GRID_VECTOR_HH
29
30#include <dune/grid/common/mcmgmapper.hh>
31#include <dune/common/version.hh>
32
34
35#include <vector>
36
37namespace Opm {
38
48template <class GridView, class Stencil, class Data, class DofMapper>
50{
51 using Element = typename GridView::template Codim<0>::Entity;
52
53 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
54
55public:
56 PffGridVector(const GridView& gridView, const DofMapper& dofMapper)
57 : gridView_(gridView)
58 , elementMapper_(gridView_, Dune::mcmgElementLayout())
59 , dofMapper_(dofMapper)
60 {}
61
62 template <class DistFn>
63 void update(const DistFn& distFn)
64 {
65 const unsigned numElements = gridView_.size(/*codim=*/0);
66 const unsigned numLocalDofs = computeNumLocalDofs_();
67
68 elemData_.resize(numElements);
69 data_.resize(numLocalDofs);
70
71 // update the pointers for the element data: for this, we need to loop over the
72 // whole grid and update a stencil for each element
73 Data* curElemDataPtr = &data_[0];
74 Stencil stencil(gridView_, dofMapper_);
75 for (const auto& elem : elements(gridView_)) {
76 // set the DOF data pointer for the current element
77 const unsigned elemIdx = elementMapper_.index(elem);
78 elemData_[elemIdx] = curElemDataPtr;
79
80 stencil.update(elem);
81 const unsigned numDof = stencil.numDof();
82 for (unsigned localDofIdx = 0; localDofIdx < numDof; ++ localDofIdx) {
83 distFn(curElemDataPtr[localDofIdx], stencil, localDofIdx);
84 }
85
86 // update the element data pointer to make it point to the beginning of the
87 // data for DOFs of the next element
88 curElemDataPtr += numDof;
89 }
90 }
91
92 void prefetch(const Element& elem) const
93 {
94 const unsigned elemIdx = elementMapper_.index(elem);
95
96 // we use 0 as the temporal locality, because it is reasonable to assume that an
97 // entry will only be accessed once.
98 ::Opm::prefetch</*temporalLocality=*/0>(elemData_[elemIdx]);
99 }
100
101 const Data& get(const Element& elem, unsigned localDofIdx) const
102 {
103 const unsigned elemIdx = elementMapper_.index(elem);
104 return elemData_[elemIdx][localDofIdx];
105 }
106
107private:
108 unsigned computeNumLocalDofs_() const
109 {
110 unsigned result = 0;
111
112 // loop over the whole grid and sum up the number of local DOFs of all Stencils
113 Stencil stencil(gridView_, dofMapper_);
114 for (const auto& elem : elements(gridView_)) {
115 stencil.update(elem);
116 result += stencil.numDof();
117 }
118
119 return result;
120 }
121
122 GridView gridView_;
123 ElementMapper elementMapper_;
124 const DofMapper& dofMapper_;
125 std::vector<Data> data_;
126 std::vector<Data*> elemData_;
127};
128
129} // namespace Opm
130
131#endif
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