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 130 131 132 133 134 135 136 137 138 139 140 141
|
// -*- 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 Opm::Linear::VertexBorderListFromGrid
*/
#ifndef EWOMS_VERTEX_BORDER_LIST_FROM_GRID_HH
#define EWOMS_VERTEX_BORDER_LIST_FROM_GRID_HH
#include "overlaptypes.hh"
#include "blacklist.hh"
#include <dune/grid/common/datahandleif.hh>
#include <dune/grid/common/gridenums.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/scalarproducts.hh>
#include <dune/istl/operators.hh>
#include <dune/common/version.hh>
#include <algorithm>
namespace Opm {
namespace Linear {
/*!
* \brief Uses communication on the grid to find the initial seed list
* of indices.
*
* \todo implement this class generically. For this, it must be
* possible to query the mapper whether it contains entities of
* a given codimension without the need to hand it an actual
* entity.
*/
template <class GridView, class VertexMapper>
class VertexBorderListFromGrid
: public Dune::CommDataHandleIF<VertexBorderListFromGrid<GridView, VertexMapper>,
int>
{
static const int dimWorld = GridView::dimensionworld;
public:
VertexBorderListFromGrid(const GridView& gridView, const VertexMapper& map)
: gridView_(gridView), map_(map)
{
gridView.communicate(*this,
Dune::InteriorBorder_InteriorBorder_Interface,
Dune::ForwardCommunication);
auto vIt = gridView.template begin<dimWorld>();
const auto& vEndIt = gridView.template end<dimWorld >();
for (; vIt != vEndIt; ++vIt) {
if (vIt->partitionType() != Dune::InteriorEntity
&& vIt->partitionType() != Dune::BorderEntity)
{
Index vIdx = static_cast<Index>(map_.index(*vIt));
blackList_.addIndex(vIdx);
}
}
}
// data handle methods
bool contains(int dim, int codim) const
{ return dim == codim; }
#if DUNE_VERSION_LT(DUNE_GRID, 2, 8)
bool fixedsize(int, int) const
#else
bool fixedSize(int, int) const
#endif
{ return true; }
template <class EntityType>
size_t size(const EntityType&) const
{ return 2; }
template <class MessageBufferImp, class EntityType>
void gather(MessageBufferImp& buff, const EntityType& e) const
{
buff.write(static_cast<int>(gridView_.comm().rank()));
buff.write(static_cast<int>(map_.index(e)));
}
template <class MessageBufferImp, class EntityType>
void scatter(MessageBufferImp& buff, const EntityType& e, size_t)
{
BorderIndex bIdx;
bIdx.localIdx = static_cast<Index>(map_.index(e));
{
int tmp;
buff.read(tmp);
bIdx.peerRank = static_cast<ProcessRank>(tmp);
}
{
int tmp;
buff.read(tmp);
bIdx.peerIdx = static_cast<Index>(tmp);
}
bIdx.borderDistance = 0;
borderList_.push_back(bIdx);
}
// Access to the border list.
const BorderList& borderList() const
{ return borderList_; }
// Access to the black-list indices.
const BlackList& blackList() const
{ return blackList_; }
private:
const GridView gridView_;
const VertexMapper& map_;
BorderList borderList_;
BlackList blackList_;
};
} // namespace Linear
} // namespace Opm
#endif
|