ParallelWellInfo.hpp
Go to the documentation of this file.
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
29
30#include <memory>
31
32namespace Opm
33{
34
35class Well;
36
40{
41public:
42 enum Attribute {
45 // there is a bug in older versions of DUNE that will skip
46 // entries with matching attributes in RemoteIndices that are local
47 // therefore we add one more version for above.
49 overlapAbove = 4
50 };
51 using LocalIndex = Dune::ParallelLocalIndex<Attribute>;
52 using IndexSet = Dune::ParallelIndexSet<int,LocalIndex,50>;
53#if HAVE_MPI
54 using RI = Dune::RemoteIndices<IndexSet>;
55#endif
56
65 void pushBackEclIndex(int above, int current, bool owner=true);
66
68 void clear();
69
72 void beginReset();
73
79 int endReset();
80
86 std::vector<double> communicateAbove(double first_value,
87 const double* current,
88 std::size_t size);
89
95 std::vector<double> communicateBelow(double first_value,
96 const double* current,
97 std::size_t size);
98
106 template<class RAIterator>
107 void partialSumPerfValues(RAIterator begin, RAIterator end) const;
108
110 const IndexSet& getIndexSet() const;
111
112 int numLocalPerfs() const;
113
114private:
117 IndexSet current_indices_;
118#if HAVE_MPI
120 IndexSet above_indices_;
121 RI remote_indices_;
122 Dune::Interface interface_;
123 Dune::BufferedCommunicator communicator_;
124#endif
125 std::size_t num_local_perfs_{};
126};
127
135{
136public:
139 using GlobalIndex = typename IndexSet::IndexPair::GlobalIndex;
140
144 const Parallel::Communication comm,
145 int num_local_perfs);
146
152 std::vector<double> createGlobal(const std::vector<double>& local_perf_container,
153 std::size_t num_components) const;
154
159 void copyGlobalToLocal(const std::vector<double>& global, std::vector<double>& local,
160 std::size_t num_components) const;
161
162 int numGlobalPerfs() const;
163
164private:
165 const IndexSet& local_indices_;
167 int num_global_perfs_;
169 std::vector<int> sizes_;
171 std::vector<int> displ_;
173 std::vector<int> map_received_;
177 std::vector<int> perf_ecl_index_;
178};
179
184{
185public:
186 static constexpr int INVALID_ECL_INDEX = -1;
187
189 ParallelWellInfo(const std::string& name = {""},
190 bool hasLocalCells = true);
191
198 ParallelWellInfo(const std::pair<std::string,bool>& well_info,
200
202 {
203 return *comm_;
204 }
205
207 void communicateFirstPerforation(bool hasFirst);
208
209
213 template<class T>
215
221 std::vector<double> communicateAboveValues(double first_value,
222 const double* current,
223 std::size_t size) const;
224
228 std::vector<double> communicateAboveValues(double first_value,
229 const std::vector<double>& current) const;
230
236 std::vector<double> communicateBelowValues(double last_value,
237 const double* current,
238 std::size_t size) const;
239
243 std::vector<double> communicateBelowValues(double last_value,
244 const std::vector<double>& current) const;
245
253 void pushBackEclIndex(int above, int current);
254
256 const std::string& name() const
257 {
258 return name_;
259 }
260
262 bool hasLocalCells() const
263 {
264 return hasLocalCells_;
265 }
266 bool isOwner() const
267 {
268 return isOwner_;
269 }
270
275
277 void endReset();
278
280 template<typename It>
281 typename It::value_type sumPerfValues(It begin, It end) const;
282
290 template<class RAIterator>
291 void partialSumPerfValues(RAIterator begin, RAIterator end) const
292 {
293 commAboveBelow_->partialSumPerfValues(begin, end);
294 }
295
297 void clear();
298
305
306private:
307
309 struct DestroyComm
310 {
311 void operator()(Parallel::Communication* comm);
312 };
313
314
316 std::string name_;
318 bool hasLocalCells_;
320 bool isOwner_;
322 int rankWithFirstPerf_;
326 std::unique_ptr<Parallel::Communication, DestroyComm> comm_;
327
329 std::unique_ptr<CommunicateAboveBelow> commAboveBelow_;
330
331 std::unique_ptr<GlobalPerfContainerFactory> globalPerfCont_;
332};
333
338{
339public:
341 const ParallelWellInfo& info);
342
347 void connectionFound(std::size_t index);
348
350
351private:
352 std::vector<std::size_t> foundConnections_;
353 const Well& well_;
354 const ParallelWellInfo& pwinfo_;
355};
356
357bool operator<(const ParallelWellInfo& well1, const ParallelWellInfo& well2);
358
359bool operator==(const ParallelWellInfo& well1, const ParallelWellInfo& well2);
360
361bool operator!=(const ParallelWellInfo& well1, const ParallelWellInfo& well2);
362
363bool operator<(const std::pair<std::string, bool>& pair, const ParallelWellInfo& well);
364
365bool operator<( const ParallelWellInfo& well, const std::pair<std::string, bool>& pair);
366
367bool operator==(const std::pair<std::string, bool>& pair, const ParallelWellInfo& well);
368
369bool operator==(const ParallelWellInfo& well, const std::pair<std::string, bool>& pair);
370
371bool operator!=(const std::pair<std::string, bool>& pair, const ParallelWellInfo& well);
372
373bool operator!=(const ParallelWellInfo& well, const std::pair<std::string, bool>& pair);
374
375} // end namespace Opm
376#endif // OPM_PARALLELWELLINFO_HEADER_INCLUDED
Class checking that all connections are on active cells.
Definition: ParallelWellInfo.hpp:338
CheckDistributedWellConnections(const Well &well, const ParallelWellInfo &info)
void connectionFound(std::size_t index)
Inidicate that the i-th completion was found.
Class to facilitate getting values associated with the above/below perforation.
Definition: ParallelWellInfo.hpp:40
void beginReset()
Indicates that we will add the index information.
const IndexSet & getIndexSet() const
Get index set for the local perforations.
void clear()
Clear all the parallel information.
Dune::ParallelIndexSet< int, LocalIndex, 50 > IndexSet
Definition: ParallelWellInfo.hpp:52
int endReset()
Indicates that the index information is complete.
Dune::RemoteIndices< IndexSet > RI
Definition: ParallelWellInfo.hpp:54
std::vector< double > communicateAbove(double first_value, const double *current, std::size_t size)
Creates an array of values for the perforation above.
void partialSumPerfValues(RAIterator begin, RAIterator end) const
Do a (in place) partial sum on values attached to all perforations.
void pushBackEclIndex(int above, int current, bool owner=true)
Adds information about original index of the perforations in ECL Schedule.
Dune::ParallelLocalIndex< Attribute > LocalIndex
Definition: ParallelWellInfo.hpp:51
Attribute
Definition: ParallelWellInfo.hpp:42
@ owner
Definition: ParallelWellInfo.hpp:43
@ overlap
Definition: ParallelWellInfo.hpp:44
@ overlapAbove
Definition: ParallelWellInfo.hpp:49
@ ownerAbove
Definition: ParallelWellInfo.hpp:48
std::vector< double > communicateBelow(double first_value, const double *current, std::size_t size)
Creates an array of values for the perforation below.
CommunicateAboveBelow(const Parallel::Communication &comm)
A factory for creating a global data representation for distributed wells.
Definition: ParallelWellInfo.hpp:135
void copyGlobalToLocal(const std::vector< double > &global, std::vector< double > &local, std::size_t num_components) const
Copies the values of the global perforation to the local representation.
GlobalPerfContainerFactory(const IndexSet &local_indices, const Parallel::Communication comm, int num_local_perfs)
Constructor.
typename IndexSet::IndexPair::GlobalIndex GlobalIndex
Definition: ParallelWellInfo.hpp:139
std::vector< double > createGlobal(const std::vector< double > &local_perf_container, std::size_t num_components) const
Creates a container that holds values for all perforations.
CommunicateAboveBelow::IndexSet IndexSet
Definition: ParallelWellInfo.hpp:137
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
bool hasLocalCells() const
Whether local cells are perforated somewhen.
Definition: ParallelWellInfo.hpp:262
void pushBackEclIndex(int above, int current)
Adds information about the ecl indices of the perforations.
std::vector< double > communicateBelowValues(double last_value, const std::vector< double > &current) const
Creates an array of values for the perforation above.
void clear()
Free data of communication data structures.
const std::string & name() const
Name of the well.
Definition: ParallelWellInfo.hpp:256
ParallelWellInfo(const std::string &name={""}, bool hasLocalCells=true)
Constructs object using MPI_COMM_SELF.
It::value_type sumPerfValues(It begin, It end) const
Sum all the values of the perforations.
void communicateFirstPerforation(bool hasFirst)
Collectively decide which rank has first perforation.
const GlobalPerfContainerFactory & getGlobalPerfContainerFactory() const
Get a factor to create a global representation of peforation data.
static constexpr int INVALID_ECL_INDEX
Definition: ParallelWellInfo.hpp:186
std::vector< double > communicateBelowValues(double last_value, const double *current, std::size_t size) const
Creates an array of values for the perforation below.
std::vector< double > communicateAboveValues(double first_value, const double *current, std::size_t size) const
Creates an array of values for the perforation above.
std::vector< double > communicateAboveValues(double first_value, const std::vector< double > &current) const
Creates an array of values for the perforation above.
T broadcastFirstPerforationValue(const T &t) const
void beginReset()
Inidicate that we will reset the ecl index information.
bool isOwner() const
Definition: ParallelWellInfo.hpp:266
ParallelWellInfo(const std::pair< std::string, bool > &well_info, Parallel::Communication allComm)
Constructs object with communication between all rank sharing a well.
void partialSumPerfValues(RAIterator begin, RAIterator end) const
Do a (in place) partial sum on values attached to all perforations.
Definition: ParallelWellInfo.hpp:291
const Parallel::Communication & communication() const
Definition: ParallelWellInfo.hpp:201
void endReset()
Inidicate completion of reset of the ecl index information.
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27
bool operator<(const ParallelWellInfo &well1, const ParallelWellInfo &well2)
bool operator==(const ParallelWellInfo &well1, const ParallelWellInfo &well2)
bool operator!=(const ParallelWellInfo &well1, const ParallelWellInfo &well2)