48class ElementBorderListFromGrid
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)
69 for (
const auto& elem : elements(gridView_)) {
70 if (elem.partitionType() != Dune::InteriorEntity) {
71 Index elemIdx =
static_cast<Index>(map_.index(elem));
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>
85 size_t size(
const EntityType&)
const
88 template <
class MessageBufferImp,
class EntityType>
89 void gather(MessageBufferImp& buff,
const EntityType& e)
const
91 unsigned myIdx =
static_cast<unsigned>(map_.index(e));
92 buff.write(
static_cast<unsigned>(gridView_.comm().rank()));
96 template <
class MessageBufferImp>
97 void scatter(MessageBufferImp& buff,
102 bool isInteriorNeighbor =
false;
103 auto isIt = gridView_.ibegin(e);
104 const auto& isEndIt = gridView_.iend(e);
105 for (; isIt != isEndIt; ++isIt) {
106 if (!isIt->neighbor())
108 else if (isIt->outside().partitionType() == Dune::InteriorEntity) {
109 isInteriorNeighbor =
true;
113 if (!isInteriorNeighbor)
123 peerSet_.insert(tmp);
132 borderList_.push_back(bIdx);
138 template <
class MessageBufferImp,
class EntityType>
139 void scatter(MessageBufferImp&,
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,
161 PeerBlackLists& peerBlackLists)
162 : gridView_(gridView)
164 , peerBlackLists_(peerBlackLists)
168 bool contains(
int,
int codim)
const
169 {
return codim == 0; }
171 bool fixedSize(
int,
int)
const
174 template <
class EntityType>
175 size_t size(
const EntityType&)
const
178 template <
class MessageBufferImp,
class EntityType>
179 void gather(MessageBufferImp& buff,
const EntityType& e)
const
181 buff.write(
static_cast<int>(gridView_.comm().rank()));
182 buff.write(
static_cast<int>(map_.index(e)));
185 template <
class MessageBufferImp,
class EntityType>
186 void scatter(MessageBufferImp& buff,
const EntityType& e,
size_t)
194 localIdx =
static_cast<Index>(map_.index(e));
196 PeerBlackListedEntry pIdx;
197 pIdx.nativeIndexOfPeer =
static_cast<Index>(peerIdx);
198 pIdx.myOwnNativeIndex =
static_cast<Index>(localIdx);
200 peerBlackLists_[
static_cast<ProcessRank>(peerRank)].push_back(pIdx);
203 const PeerBlackList& peerBlackList(
ProcessRank peerRank)
const
204 {
return peerBlackLists_.at(peerRank); }
208 const ElementMapper& map_;
209 PeerBlackLists peerBlackLists_;
213 ElementBorderListFromGrid(
const GridView& gridView,
const ElementMapper& map)
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);
227 auto peerIt = blh.peerSet().begin();
228 const auto& peerEndIt = blh.peerSet().end();
229 for (; peerIt != peerEndIt; ++peerIt) {
230 blackList_.setPeerList(*peerIt, pblh.peerBlackList(*peerIt));
236 {
return borderList_; }
240 {
return blackList_; }
243 const GridView gridView_;
244 const ElementMapper& map_;
249 PeerBlackLists peerBlackLists_;