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
31
32#include <dune/grid/common/mcmgmapper.hh>
33#include <dune/common/version.hh>
34
35#include <vector>
36
37namespace Opm {
47template <class GridView, class Stencil, class Data, class DofMapper>
49{
50 using Element = typename GridView::template Codim<0>::Entity;
51
52 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
53
54public:
55 PffGridVector(const GridView& gridView, const DofMapper& dofMapper)
56 : gridView_(gridView)
57 , elementMapper_(gridView_, Dune::mcmgElementLayout())
58 , dofMapper_(dofMapper)
59 { }
60
61 template <class DistFn>
62 void update(const DistFn& distFn)
63 {
64 unsigned numElements = gridView_.size(/*codim=*/0);
65 unsigned numLocalDofs = computeNumLocalDofs_();
66
67 elemData_.resize(numElements);
68 data_.resize(numLocalDofs);
69
70 // update the pointers for the element data: for this, we need to loop over the
71 // whole grid and update a stencil for each element
72 Data *curElemDataPtr = &data_[0];
73 Stencil stencil(gridView_, dofMapper_);
74 for (const auto& elem : elements(gridView_)) {
75 // set the DOF data pointer for the current element
76 unsigned elemIdx = elementMapper_.index(elem);
77 elemData_[elemIdx] = curElemDataPtr;
78
79 stencil.update(elem);
80 unsigned numDof = stencil.numDof();
81 for (unsigned localDofIdx = 0; localDofIdx < numDof; ++ localDofIdx)
82 distFn(curElemDataPtr[localDofIdx], stencil, localDofIdx);
83
84 // update the element data pointer to make it point to the beginning of the
85 // data for DOFs of the next element
86 curElemDataPtr += numDof;
87 }
88 }
89
90 void prefetch(const Element& elem) const
91 {
92 unsigned elemIdx = elementMapper_.index(elem);
93
94 // we use 0 as the temporal locality, because it is reasonable to assume that an
95 // entry will only be accessed once.
96 ::Opm::prefetch</*temporalLocality=*/0>(elemData_[elemIdx]);
97 }
98
99 const Data& get(const Element& elem, unsigned localDofIdx) const
100 {
101 unsigned elemIdx = elementMapper_.index(elem);
102 return elemData_[elemIdx][localDofIdx];
103 }
104
105private:
106 unsigned computeNumLocalDofs_() const
107 {
108 unsigned result = 0;
109
110 // loop over the whole grid and sum up the number of local DOFs of all Stencils
111 Stencil stencil(gridView_, dofMapper_);
112 for (const auto& elem : elements(gridView_)) {
113 stencil.update(elem);
114 result += stencil.numDof();
115 }
116
117 return result;
118 }
119
120 GridView gridView_;
121 ElementMapper elementMapper_;
122 const DofMapper& dofMapper_;
123 std::vector<Data> data_;
124 std::vector<Data*> elemData_;
125};
126
127} // namespace Opm
128
129#endif
A random-access container which stores data attached to a grid's degrees of freedom in a prefetch fri...
Definition: pffgridvector.hh:49
void prefetch(const Element &elem) const
Definition: pffgridvector.hh:90
const Data & get(const Element &elem, unsigned localDofIdx) const
Definition: pffgridvector.hh:99
PffGridVector(const GridView &gridView, const DofMapper &dofMapper)
Definition: pffgridvector.hh:55
void update(const DistFn &distFn)
Definition: pffgridvector.hh:62
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:37
void prefetch(const T &val, unsigned n=1)
Template function which emits prefetch instructions for a range of memory.
Definition: prefetch.hh:38