opm-simulators
Loading...
Searching...
No Matches
ParallelWellInfo.hpp
1/*
2 Copyright 2020 OPM-OP AS
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef OPM_PARALLELWELLINFO_HEADER_INCLUDED
20#define OPM_PARALLELWELLINFO_HEADER_INCLUDED
21
22#include <dune/common/parallel/communicator.hh>
23#include <dune/common/parallel/interface.hh>
24#include <dune/common/parallel/mpihelper.hh>
25#include <dune/common/parallel/plocalindex.hh>
26#include <dune/common/parallel/remoteindices.hh>
27
28#include <opm/simulators/utils/ParallelCommunication.hpp>
29
30#include <memory>
31#include <unordered_map>
32
33namespace Opm {
34
35class Well;
36
39template<class Scalar>
40class CommunicateAboveBelow
41{
42public:
43 enum Attribute {
44 owner=1,
45 overlap=2,
46 // there is a bug in older versions of DUNE that will skip
47 // entries with matching attributes in RemoteIndices that are local
48 // therefore we add one more version for above.
49 ownerAbove = 3,
50 overlapAbove = 4
51 };
52 using LocalIndex = Dune::ParallelLocalIndex<Attribute>;
53 using IndexSet = Dune::ParallelIndexSet<int,LocalIndex,50>;
54#if HAVE_MPI
55 using RI = Dune::RemoteIndices<IndexSet>;
56#endif
57
58 explicit CommunicateAboveBelow(const Parallel::Communication& comm);
67 void pushBackEclIndex(int above, int current, bool owner = true);
68
70 void clear();
71
74 void beginReset();
75
81 int endReset();
82
88 std::vector<Scalar> communicateAbove(Scalar first_value,
89 const Scalar* current,
90 std::size_t size);
91
97 std::vector<Scalar> communicateBelow(Scalar first_value,
98 const Scalar* current,
99 std::size_t size);
100
108 template<class RAIterator>
109 void partialSumPerfValues(RAIterator begin, RAIterator end) const;
110
112 const IndexSet& getIndexSet() const;
113
114 int numLocalPerfs() const;
115
116private:
117 Parallel::Communication comm_;
119 IndexSet current_indices_;
120#if HAVE_MPI
122 IndexSet above_indices_;
123 RI remote_indices_;
124 Dune::Interface interface_;
125 Dune::BufferedCommunicator communicator_;
126#endif
127 std::size_t num_local_perfs_{};
128};
129
136template<class Scalar>
138{
139public:
140 using IndexSet = typename CommunicateAboveBelow<Scalar>::IndexSet;
141 using Attribute = typename CommunicateAboveBelow<Scalar>::Attribute;
142 using GlobalIndex = typename IndexSet::IndexPair::GlobalIndex;
143
148 GlobalPerfContainerFactory(const IndexSet& local_indices,
149 const Parallel::Communication comm,
150 int num_local_perfs);
151
157 std::vector<Scalar> createGlobal(const std::vector<Scalar>& local_perf_container,
158 std::size_t num_quantities) const;
159
164 void copyGlobalToLocal(const std::vector<Scalar>& global, std::vector<Scalar>& local,
165 std::size_t num_quantities) const;
166
167 int numGlobalPerfs() const;
168 int globalToLocal(const int globalIndex) const;
169 int localToGlobal(std::size_t localIndex) const;
170
171private:
172 void buildLocalToGlobalMap() const;
173 void buildGlobalToLocalMap() const;
174 mutable std::unordered_map<std::size_t, int> local_to_global_map_; // Cache for L2G mapping
175 mutable std::unordered_map<int, std::size_t> global_to_local_map_; // Cache for G2L mapping
176 mutable bool l2g_map_built_ = false;
177 mutable bool g2l_map_built_ = false;
178 const IndexSet& local_indices_;
179 Parallel::Communication comm_;
180 int num_global_perfs_;
182 std::vector<int> sizes_;
184 std::vector<int> displ_;
186 std::vector<int> map_received_;
190 std::vector<int> perf_ecl_index_;
191};
192
196template<class Scalar>
198{
199public:
200 static constexpr int INVALID_ECL_INDEX = -1;
201
203 explicit ParallelWellInfo(const std::string& name = {""},
204 bool hasLocalCells = true);
205
212 ParallelWellInfo(const std::pair<std::string,bool>& well_info,
213 Parallel::Communication allComm);
214
215 const Parallel::Communication& communication() const
216 {
217 return *comm_;
218 }
219
221 void communicateFirstPerforation(bool hasFirst);
222
223 // \brief Set the activePerfToLocalPerf-Map for multisegment wells, to be called from WellState::initWellStateMSWell
224 void setActivePerfToLocalPerfMap(const std::unordered_map<int,int> active_to_local_map) const;
225 // \brief Convert a global active perforation index to a local active perforation index
226 int activePerfToLocalPerf(const int activeIndex) const;
227 // \brief Convert a local active perforation index to a global active perforation index
228 int localPerfToActivePerf(std::size_t localIndex) const;
229 // \brief Convert a global perforation index to a local perforation index
230 int globalPerfToLocalPerf(const int globalIndex) const;
231 // \brief Convert a local perforation index to a global perforation index
232 int localPerfToGlobalPerf(std::size_t localIndex) const;
233
237 template<class T>
238 T broadcastFirstPerforationValue(const T& t) const;
239
245 std::vector<Scalar> communicateAboveValues(Scalar first_value,
246 const Scalar* current,
247 std::size_t size) const;
248
252 std::vector<Scalar> communicateAboveValues(Scalar first_value,
253 const std::vector<Scalar>& current) const;
254
260 std::vector<Scalar> communicateBelowValues(Scalar last_value,
261 const Scalar* current,
262 std::size_t size) const;
263
267 std::vector<Scalar> communicateBelowValues(Scalar last_value,
268 const std::vector<Scalar>& current) const;
269
277 void pushBackEclIndex(int above, int current);
278
280 const std::string& name() const
281 {
282 return name_;
283 }
284
286 bool hasLocalCells() const
287 {
288 return hasLocalCells_;
289 }
290 bool isOwner() const
291 {
292 return isOwner_;
293 }
294
298 void beginReset();
299
301 void endReset();
302
304 template<typename It>
305 typename It::value_type sumPerfValues(It begin, It end) const;
306
314 template<class RAIterator>
315 void partialSumPerfValues(RAIterator begin, RAIterator end) const
316 {
317 commAboveBelow_->partialSumPerfValues(begin, end);
318 }
319
321 void clear();
322
329
330private:
331
333 struct DestroyComm
334 {
335 void operator()(Parallel::Communication* comm);
336 };
337
338
340 std::string name_;
342 bool hasLocalCells_;
344 bool isOwner_;
346 int rankWithFirstPerf_;
350 std::unique_ptr<Parallel::Communication, DestroyComm> comm_;
351
353 std::unique_ptr<CommunicateAboveBelow<Scalar>> commAboveBelow_;
354
355 std::unique_ptr<GlobalPerfContainerFactory<Scalar>> globalPerfCont_;
356
357 mutable std::unordered_map<int,int> active_to_local_map_; // Cache for active perforation index to local perforation index mapping
358 mutable std::unordered_map<int,int> local_to_active_map_; // Cache for local perforation index to active perforation index mapping
359};
360
364template<class Scalar>
365class CheckDistributedWellConnections
366{
367public:
368 CheckDistributedWellConnections(const Well& well,
369 const ParallelWellInfo<Scalar>& info);
370
375 void connectionFound(std::size_t index);
376
377 bool checkAllConnectionsFound();
378
379private:
380 std::vector<std::size_t> foundConnections_;
381 const Well& well_;
382 const ParallelWellInfo<Scalar>& pwinfo_;
383};
384
385template<class Scalar>
386bool operator<(const ParallelWellInfo<Scalar>& well1, const ParallelWellInfo<Scalar>& well2);
387
388template<class Scalar>
389bool operator==(const ParallelWellInfo<Scalar>& well1, const ParallelWellInfo<Scalar>& well2);
390
391template<class Scalar>
392bool operator!=(const ParallelWellInfo<Scalar>& well1, const ParallelWellInfo<Scalar>& well2);
393
394template<class Scalar>
395bool operator<(const std::pair<std::string, bool>& pair, const ParallelWellInfo<Scalar>& well);
396
397template<class Scalar>
398bool operator<( const ParallelWellInfo<Scalar>& well, const std::pair<std::string, bool>& pair);
399
400template<class Scalar>
401bool operator==(const std::pair<std::string, bool>& pair, const ParallelWellInfo<Scalar>& well);
402
403template<class Scalar>
404bool operator==(const ParallelWellInfo<Scalar>& well, const std::pair<std::string, bool>& pair);
405
406template<class Scalar>
407bool operator!=(const std::pair<std::string, bool>& pair, const ParallelWellInfo<Scalar>& well);
408
409template<class Scalar>
410bool operator!=(const ParallelWellInfo<Scalar>& well, const std::pair<std::string, bool>& pair);
411
412} // end namespace Opm
413
414#endif // OPM_PARALLELWELLINFO_HEADER_INCLUDED
void connectionFound(std::size_t index)
Inidicate that the i-th completion was found.
Definition ParallelWellInfo.cpp:778
std::vector< Scalar > communicateBelow(Scalar first_value, const Scalar *current, std::size_t size)
Creates an array of values for the perforation below.
Definition ParallelWellInfo.cpp:421
void clear()
Clear all the parallel information.
Definition ParallelWellInfo.cpp:271
void pushBackEclIndex(int above, int current, bool owner=true)
Adds information about original index of the perforations in ECL Schedule.
Definition ParallelWellInfo.cpp:451
const IndexSet & getIndexSet() const
Get index set for the local perforations.
Definition ParallelWellInfo.cpp:492
int endReset()
Indicates that the index information is complete.
Definition ParallelWellInfo.cpp:296
std::vector< Scalar > communicateAbove(Scalar first_value, const Scalar *current, std::size_t size)
Creates an array of values for the perforation above.
Definition ParallelWellInfo.cpp:391
void beginReset()
Indicates that we will add the index information.
Definition ParallelWellInfo.cpp:283
void partialSumPerfValues(RAIterator begin, RAIterator end) const
Do a (in place) partial sum on values attached to all perforations.
Definition ParallelWellInfo.cpp:317
A factory for creating a global data representation for distributed wells.
Definition ParallelWellInfo.hpp:138
GlobalPerfContainerFactory(const IndexSet &local_indices, const Parallel::Communication comm, int num_local_perfs)
Constructor.
Definition ParallelWellInfo.cpp:78
void copyGlobalToLocal(const std::vector< Scalar > &global, std::vector< Scalar > &local, std::size_t num_quantities) const
Copies the values of the global perforation to the local representation.
Definition ParallelWellInfo.cpp:225
std::vector< Scalar > createGlobal(const std::vector< Scalar > &local_perf_container, std::size_t num_quantities) const
Creates a container that holds values for all perforations.
Definition ParallelWellInfo.cpp:174
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
std::vector< Scalar > communicateBelowValues(Scalar last_value, const Scalar *current, std::size_t size) const
Creates an array of values for the perforation below.
Definition ParallelWellInfo.cpp:677
ParallelWellInfo(const std::string &name={""}, bool hasLocalCells=true)
Constructs object using MPI_COMM_SELF.
Definition ParallelWellInfo.cpp:504
const std::string & name() const
Name of the well.
Definition ParallelWellInfo.hpp:280
T broadcastFirstPerforationValue(const T &t) const
If the well does not have any open connections the member rankWithFirstPerf is not initialized,...
Definition ParallelWellInfo.cpp:629
std::vector< Scalar > communicateAboveValues(Scalar first_value, const Scalar *current, std::size_t size) const
Creates an array of values for the perforation above.
Definition ParallelWellInfo.cpp:658
void beginReset()
Inidicate that we will reset the ecl index information.
Definition ParallelWellInfo.cpp:594
void pushBackEclIndex(int above, int current)
Adds information about the ecl indices of the perforations.
Definition ParallelWellInfo.cpp:588
bool hasLocalCells() const
Whether local cells are perforated somewhen.
Definition ParallelWellInfo.hpp:286
It::value_type sumPerfValues(It begin, It end) const
Sum all the values of the perforations.
Definition ParallelWellInfo.cpp:612
void partialSumPerfValues(RAIterator begin, RAIterator end) const
Do a (in place) partial sum on values attached to all perforations.
Definition ParallelWellInfo.hpp:315
void clear()
Free data of communication data structures.
Definition ParallelWellInfo.cpp:621
const GlobalPerfContainerFactory< Scalar > & getGlobalPerfContainerFactory() const
Get a factor to create a global representation of peforation data.
Definition ParallelWellInfo.cpp:696
void communicateFirstPerforation(bool hasFirst)
Collectively decide which rank has first perforation.
Definition ParallelWellInfo.cpp:576
void endReset()
Inidicate completion of reset of the ecl index information.
Definition ParallelWellInfo.cpp:600
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43