// -*- 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 . 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 #include #include #include 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 PffGridVector { using Element = typename GridView::template Codim<0>::Entity; using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper; public: PffGridVector(const GridView& gridView, const DofMapper& dofMapper) : gridView_(gridView) , elementMapper_(gridView_, Dune::mcmgElementLayout()) , dofMapper_(dofMapper) { } template 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(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_; std::vector elemData_; }; } // namespace Opm #endif