File: pffgridvector.hh

package info (click to toggle)
opm-models 2022.10%2Bds-4
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 6,216 kB
  • sloc: cpp: 37,910; ansic: 1,897; sh: 277; xml: 182; makefile: 10
file content (129 lines) | stat: -rw-r--r-- 4,330 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
  This file is part of the Open Porous Media project (OPM).

  OPM is free software: you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation, either version 2 of the License, or
  (at your option) any later version.

  OPM is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with OPM.  If not, see <http://www.gnu.org/licenses/>.

  Consult the COPYING file in the top-level source directory of this
  module for the precise wording of the license and the list of
  copyright holders.
*/
/*!
 * \file
 * \copydoc PffGridVector
 */
#ifndef EWOMS_PFF_GRID_VECTOR_HH
#define EWOMS_PFF_GRID_VECTOR_HH

#include <opm/models/utils/prefetch.hh>

#include <dune/grid/common/mcmgmapper.hh>
#include <dune/common/version.hh>

#include <vector>

namespace Opm {
/*!
 * \brief A random-access container which stores data attached to a grid's degrees of
 *        freedom in a prefetch friendly manner.
 *
 * This container often reduces the number of cache faults considerably, thus improving
 * performance. On the flipside data cannot be written to on an individual basis and it
 * requires significantly more memory than a plain array. PffVector stands for "PreFetch
 * Friendly Grid Vector".
 */
template <class GridView, class Stencil, class Data, class DofMapper>
class PffGridVector
{
    using Element = typename GridView::template Codim<0>::Entity;

    using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;

public:
    PffGridVector(const GridView& gridView, const DofMapper& dofMapper)
        : gridView_(gridView)
        , elementMapper_(gridView_, Dune::mcmgElementLayout())
        , dofMapper_(dofMapper)
    { }

    template <class DistFn>
    void update(const DistFn& distFn)
    {
        unsigned numElements = gridView_.size(/*codim=*/0);
        unsigned numLocalDofs = computeNumLocalDofs_();

        elemData_.resize(numElements);
        data_.resize(numLocalDofs);

        // update the pointers for the element data: for this, we need to loop over the
        // whole grid and update a stencil for each element
        Data *curElemDataPtr = &data_[0];
        Stencil stencil(gridView_, dofMapper_);
        for (const auto& elem : elements(gridView_)) {
            // set the DOF data pointer for the current element
            unsigned elemIdx = elementMapper_.index(elem);
            elemData_[elemIdx] = curElemDataPtr;

            stencil.update(elem);
            unsigned numDof = stencil.numDof();
            for (unsigned localDofIdx = 0; localDofIdx < numDof; ++ localDofIdx)
                distFn(curElemDataPtr[localDofIdx], stencil, localDofIdx);

            // update the element data pointer to make it point to the beginning of the
            // data for DOFs of the next element
            curElemDataPtr += numDof;
        }
    }

    void prefetch(const Element& elem) const
    {
        unsigned elemIdx = elementMapper_.index(elem);

        // we use 0 as the temporal locality, because it is reasonable to assume that an
        // entry will only be accessed once.
        ::Opm::prefetch</*temporalLocality=*/0>(elemData_[elemIdx]);
    }

    const Data& get(const Element& elem, unsigned localDofIdx) const
    {
        unsigned elemIdx = elementMapper_.index(elem);
        return elemData_[elemIdx][localDofIdx];
    }

private:
    unsigned computeNumLocalDofs_() const
    {
        unsigned result = 0;

        // loop over the whole grid and sum up the number of local DOFs of all Stencils
        Stencil stencil(gridView_, dofMapper_);
        for (const auto& elem : elements(gridView_)) {
            stencil.update(elem);
            result += stencil.numDof();
        }

        return result;
    }

    GridView gridView_;
    ElementMapper elementMapper_;
    const DofMapper& dofMapper_;
    std::vector<Data> data_;
    std::vector<Data*> elemData_;
};

} // namespace Opm

#endif