File: blacklist.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 (186 lines) | stat: -rw-r--r-- 6,406 bytes parent folder | download | duplicates (4)
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
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
// -*- 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::BlackList
 */
#ifndef EWOMS_BLACK_LIST_HH
#define EWOMS_BLACK_LIST_HH

#include "overlaptypes.hh"

#if HAVE_MPI
#include <opm/models/parallel/mpibuffer.hh>

#include <dune/grid/common/datahandleif.hh>
#include <dune/grid/common/gridenums.hh>
#endif // HAVE_MPI

#include <iostream>
#include <algorithm>

namespace Opm {
namespace Linear {
/*!
 * \brief Expresses which degrees of freedom are blacklisted for the parallel linear
 *        solvers and which domestic indices they correspond to.
 */
class BlackList
{
public:
    struct PeerBlackListedEntry {
        Index nativeIndexOfPeer;
        Index myOwnNativeIndex;
    };
    using PeerBlackList = std::vector<PeerBlackListedEntry>;
    using PeerBlackLists = std::map<ProcessRank, PeerBlackList>;

    BlackList()
    { }

    BlackList(const BlackList&) = default;

    bool hasIndex(Index nativeIdx) const
    { return nativeBlackListedIndices_.count(nativeIdx) > 0; }

    void addIndex(Index nativeIdx)
    { nativeBlackListedIndices_.insert(nativeIdx); }

    Index nativeToDomestic(Index nativeIdx) const
    {
        auto it = nativeToDomesticMap_.find(nativeIdx);
        if (it == nativeToDomesticMap_.end())
            return -1;
        return it->second;
    }

    void setPeerList(ProcessRank peerRank, const PeerBlackList& peerBlackList)
    { peerBlackLists_[peerRank] = peerBlackList; }

    template <class DomesticOverlap>
    void updateNativeToDomesticMap([[maybe_unused]] const DomesticOverlap& domesticOverlap)
    {
#if HAVE_MPI
        auto peerListIt = peerBlackLists_.begin();
        const auto& peerListEndIt = peerBlackLists_.end();
        for (; peerListIt != peerListEndIt; ++peerListIt) {
            sendGlobalIndices_(peerListIt->first,
                               peerListIt->second,
                               domesticOverlap);
        }

        peerListIt = peerBlackLists_.begin();
        for (; peerListIt != peerListEndIt; ++peerListIt) {
            receiveGlobalIndices_(peerListIt->first, domesticOverlap);
        }

        peerListIt = peerBlackLists_.begin();
        for (; peerListIt != peerListEndIt; ++peerListIt) {
            numGlobalIdxSendBuff_.at(peerListIt->first).wait();
            globalIdxSendBuff_.at(peerListIt->first).wait();
        }
#endif // HAVE_MPI
    }

    void print() const
    {
        std::cout << "my own blacklisted indices:\n";
        auto idxIt = nativeBlackListedIndices_.begin();
        const auto& idxEndIt = nativeBlackListedIndices_.end();
        for (; idxIt != idxEndIt; ++idxIt)
            std::cout << " (native index: " << *idxIt
                      << ", domestic index: " << nativeToDomestic(*idxIt) << ")\n";
        std::cout << "blacklisted indices of the peers in my own domain:\n";
        auto peerListIt = peerBlackLists_.begin();
        const auto& peerListEndIt = peerBlackLists_.end();
        for (; peerListIt != peerListEndIt; ++peerListIt) {
            ProcessRank peerRank = peerListIt->first;
            std::cout << " peer " << peerRank << ":\n";
            auto idx2It = peerListIt->second.begin();
            const auto& idx2EndIt = peerListIt->second.end();
            for (; idx2It != idx2EndIt; ++ idx2It)
                std::cout << "   (native index: " << idx2It->myOwnNativeIndex
                          << ", native peer index: " << idx2It->nativeIndexOfPeer << ")\n";
        }
    }

private:
#if HAVE_MPI
    template <class DomesticOverlap>
    void sendGlobalIndices_(ProcessRank peerRank,
                            const PeerBlackList& peerIndices,
                            const DomesticOverlap& domesticOverlap)
    {
        auto& numIdxBuff = numGlobalIdxSendBuff_[peerRank];
        auto& idxBuff = globalIdxSendBuff_[peerRank];

        numIdxBuff.resize(1);
        numIdxBuff[0] = static_cast<unsigned>(peerIndices.size());
        numIdxBuff.send(peerRank);

        idxBuff.resize(2*peerIndices.size());
        for (size_t i = 0; i < peerIndices.size(); ++i) {
            // global index
            Index myNativeIdx = peerIndices[i].myOwnNativeIndex;
            Index myDomesticIdx = domesticOverlap.nativeToDomestic(myNativeIdx);
            idxBuff[2*i + 0] = domesticOverlap.domesticToGlobal(myDomesticIdx);

            // native peer index
            idxBuff[2*i + 1] = peerIndices[i].nativeIndexOfPeer;
        }
        idxBuff.send(peerRank);
    }

    template <class DomesticOverlap>
    void receiveGlobalIndices_(ProcessRank peerRank,
                               const DomesticOverlap& domesticOverlap)
    {
        MpiBuffer<unsigned> numGlobalIdxBuf(1);
        numGlobalIdxBuf.receive(peerRank);
        unsigned numIndices = numGlobalIdxBuf[0];

        MpiBuffer<Index> globalIdxBuf(2*numIndices);
        globalIdxBuf.receive(peerRank);
        for (unsigned i = 0; i < numIndices; ++i) {
            Index globalIdx = globalIdxBuf[2*i + 0];
            Index nativeIdx = globalIdxBuf[2*i + 1];

            nativeToDomesticMap_[nativeIdx] = domesticOverlap.globalToDomestic(globalIdx);
        }
    }
#endif // HAVE_MPI

    std::set<Index> nativeBlackListedIndices_;
    std::map<Index, Index> nativeToDomesticMap_;
#if HAVE_MPI
    std::map<ProcessRank, MpiBuffer<unsigned>> numGlobalIdxSendBuff_;
    std::map<ProcessRank, MpiBuffer<Index>> globalIdxSendBuff_;
#endif // HAVE_MPI

    PeerBlackLists peerBlackLists_;
};

} // namespace Linear
} // namespace Opm

#endif