51 using PeerBlackList = BlackList::PeerBlackList;
52 using PeerBlackLists = BlackList::PeerBlackLists;
54 using Element =
typename GridView::template Codim<0>::Entity;
56 class BorderListHandle_
57 :
public Dune::CommDataHandleIF<BorderListHandle_, ProcessRank>
60 BorderListHandle_(
const GridView& gridView,
61 const ElementMapper&
map,
66 , blackList_(blackList)
67 , borderList_(borderList)
70 if (
elem.partitionType() != Dune::InteriorEntity) {
72 blackList_.addIndex(elemIdx);
78 bool contains(
int,
int codim)
const
79 {
return codim == 0; }
81 bool fixedSize(
int,
int)
const
84 template <
class EntityType>
88 template <
class MessageBufferImp,
class EntityType>
91 unsigned myIdx =
static_cast<unsigned>(map_.index(
e));
92 buff.write(
static_cast<unsigned>(gridView_.comm().rank()));
96 template <
class MessageBufferImp>
103 auto isIt = gridView_.ibegin(
e);
104 const auto&
isEndIt = gridView_.iend(
e);
106 if (!
isIt->neighbor())
108 else if (
isIt->outside().partitionType() == Dune::InteriorEntity) {
123 peerSet_.insert(tmp);
130 bIdx.borderDistance = 1;
132 borderList_.push_back(
bIdx);
138 template <
class MessageBufferImp,
class EntityType>
144 const std::set<ProcessRank>& peerSet()
const
149 const ElementMapper& map_;
150 std::set<ProcessRank> peerSet_;
155 class PeerBlackListHandle_
156 :
public Dune::CommDataHandleIF<PeerBlackListHandle_, int>
159 PeerBlackListHandle_(
const GridView& gridView,
160 const ElementMapper&
map,
162 : gridView_(gridView)
168 bool contains(
int,
int codim)
const
169 {
return codim == 0; }
171 bool fixedSize(
int,
int)
const
174 template <
class EntityType>
178 template <
class MessageBufferImp,
class EntityType>
181 buff.write(
static_cast<int>(gridView_.comm().rank()));
182 buff.write(
static_cast<int>(map_.index(
e)));
185 template <
class MessageBufferImp,
class EntityType>
194 localIdx =
static_cast<Index>(map_.index(
e));
197 pIdx.nativeIndexOfPeer =
static_cast<Index>(peerIdx);
198 pIdx.myOwnNativeIndex =
static_cast<Index>(localIdx);
203 const PeerBlackList& peerBlackList(
ProcessRank peerRank)
const
204 {
return peerBlackLists_.at(peerRank); }
208 const ElementMapper& map_;
209 PeerBlackLists peerBlackLists_;
214 : gridView_(gridView)
217 BorderListHandle_
blh(gridView,
map, blackList_, borderList_);
218 gridView.communicate(
blh,
219 Dune::InteriorBorder_All_Interface,
220 Dune::BackwardCommunication);
222 PeerBlackListHandle_
pblh(gridView,
map, peerBlackLists_);
223 gridView.communicate(
pblh,
224 Dune::InteriorBorder_All_Interface,
225 Dune::BackwardCommunication);
236 {
return borderList_; }
240 {
return blackList_; }
243 const GridView gridView_;
244 const ElementMapper& map_;
249 PeerBlackLists peerBlackLists_;
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
std::list< BorderIndex > BorderList
This class managages a list of indices which are on the border of a process' partition of the grid.
Definition overlaptypes.hh:120
Definition blacklist.hh:51