19#ifndef OPM_PARALLELWELLINFO_HEADER_INCLUDED
20#define OPM_PARALLELWELLINFO_HEADER_INCLUDED
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>
28#include <opm/simulators/utils/ParallelCommunication.hpp>
31#include <unordered_map>
40class CommunicateAboveBelow
52 using LocalIndex = Dune::ParallelLocalIndex<Attribute>;
53 using IndexSet = Dune::ParallelIndexSet<int,LocalIndex,50>;
55 using RI = Dune::RemoteIndices<IndexSet>;
58 explicit CommunicateAboveBelow(
const Parallel::Communication& comm);
89 const Scalar* current,
98 const Scalar* current,
108 template<
class RAIterator>
114 int numLocalPerfs()
const;
117 Parallel::Communication comm_;
119 IndexSet current_indices_;
122 IndexSet above_indices_;
124 Dune::Interface interface_;
125 Dune::BufferedCommunicator communicator_;
127 std::size_t num_local_perfs_{};
136template<
class Scalar>
140 using IndexSet =
typename CommunicateAboveBelow<Scalar>::IndexSet;
141 using Attribute =
typename CommunicateAboveBelow<Scalar>::Attribute;
142 using GlobalIndex =
typename IndexSet::IndexPair::GlobalIndex;
149 const Parallel::Communication comm,
150 int num_local_perfs);
157 std::vector<Scalar>
createGlobal(
const std::vector<Scalar>& local_perf_container,
158 std::size_t num_quantities)
const;
164 void copyGlobalToLocal(
const std::vector<Scalar>& global, std::vector<Scalar>& local,
165 std::size_t num_quantities)
const;
167 int numGlobalPerfs()
const;
168 int globalToLocal(
const int globalIndex)
const;
169 int localToGlobal(std::size_t localIndex)
const;
172 void buildLocalToGlobalMap()
const;
173 void buildGlobalToLocalMap()
const;
174 mutable std::unordered_map<std::size_t, int> local_to_global_map_;
175 mutable std::unordered_map<int, std::size_t> global_to_local_map_;
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_;
196template<
class Scalar>
200 static constexpr int INVALID_ECL_INDEX = -1;
213 Parallel::Communication allComm);
215 const Parallel::Communication& communication()
const
224 void setActivePerfToLocalPerfMap(
const std::unordered_map<int,int> active_to_local_map)
const;
226 int activePerfToLocalPerf(
const int activeIndex)
const;
228 int localPerfToActivePerf(std::size_t localIndex)
const;
230 int globalPerfToLocalPerf(
const int globalIndex)
const;
232 int localPerfToGlobalPerf(std::size_t localIndex)
const;
246 const Scalar* current,
247 std::size_t size)
const;
253 const std::vector<Scalar>& current)
const;
261 const Scalar* current,
262 std::size_t size)
const;
268 const std::vector<Scalar>& current)
const;
280 const std::string&
name()
const
288 return hasLocalCells_;
304 template<
typename It>
305 typename It::value_type
sumPerfValues(It begin, It end)
const;
314 template<
class RAIterator>
317 commAboveBelow_->partialSumPerfValues(begin, end);
335 void operator()(Parallel::Communication* comm);
346 int rankWithFirstPerf_;
350 std::unique_ptr<Parallel::Communication, DestroyComm> comm_;
353 std::unique_ptr<CommunicateAboveBelow<Scalar>> commAboveBelow_;
355 std::unique_ptr<GlobalPerfContainerFactory<Scalar>> globalPerfCont_;
357 mutable std::unordered_map<int,int> active_to_local_map_;
358 mutable std::unordered_map<int,int> local_to_active_map_;
364template<
class Scalar>
365class CheckDistributedWellConnections
368 CheckDistributedWellConnections(
const Well& well,
377 bool checkAllConnectionsFound();
380 std::vector<std::size_t> foundConnections_;
385template<
class Scalar>
388template<
class Scalar>
391template<
class Scalar>
394template<
class Scalar>
397template<
class Scalar>
400template<
class Scalar>
403template<
class Scalar>
406template<
class Scalar>
409template<
class Scalar>
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