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
34class Well;
35
38template<class Scalar>
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<Scalar> communicateAbove(Scalar first_value,
87 const Scalar* current,
88 std::size_t size);
89
95 std::vector<Scalar> communicateBelow(Scalar first_value,
96 const Scalar* 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
134template<class Scalar>
136{
137public:
140 using GlobalIndex = typename IndexSet::IndexPair::GlobalIndex;
141
145 const Parallel::Communication comm,
146 int num_local_perfs);
147
153 std::vector<Scalar> createGlobal(const std::vector<Scalar>& local_perf_container,
154 std::size_t num_components) const;
155
160 void copyGlobalToLocal(const std::vector<Scalar>& global, std::vector<Scalar>& local,
161 std::size_t num_components) const;
162
163 int numGlobalPerfs() const;
164
165private:
166 const IndexSet& local_indices_;
168 int num_global_perfs_;
170 std::vector<int> sizes_;
172 std::vector<int> displ_;
174 std::vector<int> map_received_;
178 std::vector<int> perf_ecl_index_;
179};
180
184template<class Scalar>
186{
187public:
188 static constexpr int INVALID_ECL_INDEX = -1;
189
191 ParallelWellInfo(const std::string& name = {""},
192 bool hasLocalCells = true);
193
200 ParallelWellInfo(const std::pair<std::string,bool>& well_info,
202
204 {
205 return *comm_;
206 }
207
209 void communicateFirstPerforation(bool hasFirst);
210
211
215 template<class T>
217
223 std::vector<Scalar> communicateAboveValues(Scalar first_value,
224 const Scalar* current,
225 std::size_t size) const;
226
230 std::vector<Scalar> communicateAboveValues(Scalar first_value,
231 const std::vector<Scalar>& current) const;
232
238 std::vector<Scalar> communicateBelowValues(Scalar last_value,
239 const Scalar* current,
240 std::size_t size) const;
241
245 std::vector<Scalar> communicateBelowValues(Scalar last_value,
246 const std::vector<Scalar>& current) const;
247
255 void pushBackEclIndex(int above, int current);
256
258 const std::string& name() const
259 {
260 return name_;
261 }
262
264 bool hasLocalCells() const
265 {
266 return hasLocalCells_;
267 }
268 bool isOwner() const
269 {
270 return isOwner_;
271 }
272
277
279 void endReset();
280
282 template<typename It>
283 typename It::value_type sumPerfValues(It begin, It end) const;
284
292 template<class RAIterator>
293 void partialSumPerfValues(RAIterator begin, RAIterator end) const
294 {
295 commAboveBelow_->partialSumPerfValues(begin, end);
296 }
297
299 void clear();
300
307
308private:
309
311 struct DestroyComm
312 {
313 void operator()(Parallel::Communication* comm);
314 };
315
316
318 std::string name_;
320 bool hasLocalCells_;
322 bool isOwner_;
324 int rankWithFirstPerf_;
328 std::unique_ptr<Parallel::Communication, DestroyComm> comm_;
329
331 std::unique_ptr<CommunicateAboveBelow<Scalar>> commAboveBelow_;
332
333 std::unique_ptr<GlobalPerfContainerFactory<Scalar>> globalPerfCont_;
334};
335
339template<class Scalar>
341{
342public:
344 const ParallelWellInfo<Scalar>& info);
345
350 void connectionFound(std::size_t index);
351
353
354private:
355 std::vector<std::size_t> foundConnections_;
356 const Well& well_;
357 const ParallelWellInfo<Scalar>& pwinfo_;
358};
359
360template<class Scalar>
362
363template<class Scalar>
365
366template<class Scalar>
368
369template<class Scalar>
370bool operator<(const std::pair<std::string, bool>& pair, const ParallelWellInfo<Scalar>& well);
371
372template<class Scalar>
373bool operator<( const ParallelWellInfo<Scalar>& well, const std::pair<std::string, bool>& pair);
374
375template<class Scalar>
376bool operator==(const std::pair<std::string, bool>& pair, const ParallelWellInfo<Scalar>& well);
377
378template<class Scalar>
379bool operator==(const ParallelWellInfo<Scalar>& well, const std::pair<std::string, bool>& pair);
380
381template<class Scalar>
382bool operator!=(const std::pair<std::string, bool>& pair, const ParallelWellInfo<Scalar>& well);
383
384template<class Scalar>
385bool operator!=(const ParallelWellInfo<Scalar>& well, const std::pair<std::string, bool>& pair);
386
387} // end namespace Opm
388
389#endif // OPM_PARALLELWELLINFO_HEADER_INCLUDED
Class checking that all connections are on active cells.
Definition: ParallelWellInfo.hpp:341
void connectionFound(std::size_t index)
Inidicate that the i-th completion was found.
CheckDistributedWellConnections(const Well &well, const ParallelWellInfo< Scalar > &info)
Class to facilitate getting values associated with the above/below perforation.
Definition: ParallelWellInfo.hpp:40
Dune::RemoteIndices< IndexSet > RI
Definition: ParallelWellInfo.hpp:54
std::vector< Scalar > communicateBelow(Scalar first_value, const Scalar *current, std::size_t size)
Creates an array of values for the perforation below.
CommunicateAboveBelow(const Parallel::Communication &comm)
void clear()
Clear all the parallel information.
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
int endReset()
Indicates that the index information is complete.
const IndexSet & getIndexSet() const
Get index set for the local perforations.
std::vector< Scalar > communicateAbove(Scalar first_value, const Scalar *current, std::size_t size)
Creates an array of values for the perforation above.
void beginReset()
Indicates that we will add the index information.
void partialSumPerfValues(RAIterator begin, RAIterator end) const
Do a (in place) partial sum on values attached to all perforations.
Dune::ParallelIndexSet< int, LocalIndex, 50 > IndexSet
Definition: ParallelWellInfo.hpp:52
Attribute
Definition: ParallelWellInfo.hpp:42
@ ownerAbove
Definition: ParallelWellInfo.hpp:48
@ overlap
Definition: ParallelWellInfo.hpp:44
@ overlapAbove
Definition: ParallelWellInfo.hpp:49
@ owner
Definition: ParallelWellInfo.hpp:43
A factory for creating a global data representation for distributed wells.
Definition: ParallelWellInfo.hpp:136
GlobalPerfContainerFactory(const IndexSet &local_indices, const Parallel::Communication comm, int num_local_perfs)
Constructor.
std::vector< Scalar > createGlobal(const std::vector< Scalar > &local_perf_container, std::size_t num_components) const
Creates a container that holds values for all perforations.
void copyGlobalToLocal(const std::vector< Scalar > &global, std::vector< Scalar > &local, std::size_t num_components) const
Copies the values of the global perforation to the local representation.
typename CommunicateAboveBelow< Scalar >::Attribute Attribute
Definition: ParallelWellInfo.hpp:139
typename IndexSet::IndexPair::GlobalIndex GlobalIndex
Definition: ParallelWellInfo.hpp:140
typename CommunicateAboveBelow< Scalar >::IndexSet IndexSet
Definition: ParallelWellInfo.hpp:138
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
std::vector< Scalar > communicateBelowValues(Scalar last_value, const Scalar *current, std::size_t size) const
Creates an array of values for the perforation below.
std::vector< Scalar > communicateBelowValues(Scalar last_value, const std::vector< Scalar > &current) const
Creates an array of values for the perforation above.
ParallelWellInfo(const std::string &name={""}, bool hasLocalCells=true)
Constructs object using MPI_COMM_SELF.
const std::string & name() const
Name of the well.
Definition: ParallelWellInfo.hpp:258
T broadcastFirstPerforationValue(const T &t) const
ParallelWellInfo(const std::pair< std::string, bool > &well_info, Parallel::Communication allComm)
Constructs object with communication between all rank sharing a well.
bool isOwner() const
Definition: ParallelWellInfo.hpp:268
std::vector< Scalar > communicateAboveValues(Scalar first_value, const Scalar *current, std::size_t size) const
Creates an array of values for the perforation above.
void beginReset()
Inidicate that we will reset the ecl index information.
void pushBackEclIndex(int above, int current)
Adds information about the ecl indices of the perforations.
bool hasLocalCells() const
Whether local cells are perforated somewhen.
Definition: ParallelWellInfo.hpp:264
It::value_type sumPerfValues(It begin, It end) const
Sum all the values of the perforations.
void partialSumPerfValues(RAIterator begin, RAIterator end) const
Do a (in place) partial sum on values attached to all perforations.
Definition: ParallelWellInfo.hpp:293
void clear()
Free data of communication data structures.
const GlobalPerfContainerFactory< Scalar > & getGlobalPerfContainerFactory() const
Get a factor to create a global representation of peforation data.
void communicateFirstPerforation(bool hasFirst)
Collectively decide which rank has first perforation.
static constexpr int INVALID_ECL_INDEX
Definition: ParallelWellInfo.hpp:188
std::vector< Scalar > communicateAboveValues(Scalar first_value, const std::vector< Scalar > &current) const
Creates an array of values for the perforation above.
void endReset()
Inidicate completion of reset of the ecl index information.
const Parallel::Communication & communication() const
Definition: ParallelWellInfo.hpp:203
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:37
bool operator<(const ParallelWellInfo< Scalar > &well1, const ParallelWellInfo< Scalar > &well2)
bool operator==(const aligned_allocator< T1, Alignment > &, const aligned_allocator< T2, Alignment > &) noexcept
Definition: alignedallocator.hh:210
bool operator!=(const aligned_allocator< T1, Alignment > &, const aligned_allocator< T2, Alignment > &) noexcept
Definition: alignedallocator.hh:218