ParallelWBPCalculation.hpp
Go to the documentation of this file.
1/*
2 Copyright 2023 Equinor ASA.
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
20#ifndef PARALLEL_WBP_CALCULATION_HPP
21#define PARALLEL_WBP_CALCULATION_HPP
22
24
26
28
29#include <opm/input/eclipse/Schedule/Well/PAvgCalculator.hpp>
30#include <opm/input/eclipse/Schedule/Well/PAvgCalculatorCollection.hpp>
31
32#include <cstddef>
33#include <functional>
34#include <memory>
35#include <vector>
36
37namespace Opm {
38 class GridDims;
39 template<class Scalar> class ParallelWellInfo;
40 class PAvg;
41 class Well;
42}
43
44namespace Opm {
45
49template<class Scalar>
51{
52public:
56
59
63 using EvaluatorFactory = std::function<Evaluator()>;
64
70 explicit ParallelWBPCalculation(const GridDims& cellIndexMap,
71 const Parallel::Communication& gridComm);
72
81
89
114 std::size_t
115 createCalculator(const Well& well,
116 const ParallelWellInfo<Scalar>& parallelWellInfo,
117 const std::vector<int>& localConnIdx,
118 EvaluatorFactory makeWellSourceEvaluator);
119
127
141
154 void inferBlockAveragePressures(const std::size_t calcIndex,
155 const PAvg& controls,
156 const Scalar gravity,
157 const Scalar refDepth);
158
167 const typename PAvgCalculator<Scalar>::Result&
168 averagePressures(const std::size_t calcIndex) const;
169
170private:
172 class LocalConnSet
173 {
174 public:
185 explicit LocalConnSet(const std::vector<int>& localConnIdx);
186
195 int localIndex(const std::size_t connIdx) const;
196
197 private:
206 std::vector<int> localConnIdx_{};
207 };
208
210 class SourceData
211 {
212 public:
217 operator const ParallelPAvgDynamicSourceData<Scalar>&() const
218 {
219 return *this->srcData_;
220 }
221
225 explicit SourceData(const Parallel::Communication& comm);
226
234 SourceData& localIndex(GlobalToLocal localIdx);
235
242 SourceData& evaluator(Evaluator eval);
243
252 SourceData& evaluatorFactory(EvaluatorFactory evalFactory);
253
263 void buildStructure(const std::vector<std::size_t>& sourceLocations);
264
275
279 const Parallel::Communication& comm() const
280 {
281 return this->comm_.get();
282 }
283
293 std::vector<int>
294 getLocalIndex(const std::vector<std::size_t>& globalIndex) const;
295
296 private:
298 using DataPtr = std::unique_ptr<ParallelPAvgDynamicSourceData<Scalar>>;
299
301 std::reference_wrapper<const Parallel::Communication> comm_;
302
304 GlobalToLocal localIdx_{};
305
308 Evaluator eval_{};
309
313 EvaluatorFactory evalFactory_{};
314
316 DataPtr srcData_{};
317 };
318
320 using WellID = typename std::vector<SourceData>::size_type;
321
323 std::reference_wrapper<const GridDims> cellIndexMap_;
324
326 SourceData reservoirSrc_;
327
330 PAvgCalculatorCollection<Scalar> calculators_{};
331
333 std::vector<SourceData> wellConnSrc_{};
334
336 std::vector<LocalConnSet> localConnSet_{};
337
340 void pruneInactiveWBPCells();
341
343 void pruneInactiveWBPCellsSerial();
344
346 void pruneInactiveWBPCellsParallel();
347
349 void defineReservoirCommunication();
350
355 void defineWellCommunication(const std::size_t well);
356
364 typename PAvgCalculator<Scalar>::Sources
365 makeEvaluationSources(const WellID well) const;
366};
367
368} // namespace Opm
369
370#endif // PARALLEL_WBP_CALCULATION_HPP
std::function< void(int, SourceDataSpan< Scalar >)> Evaluator
Definition: ParallelPAvgDynamicSourceData.hpp:54
std::function< int(const std::size_t)> GlobalToLocal
Definition: ParallelPAvgDynamicSourceData.hpp:42
Definition: ParallelWBPCalculation.hpp:51
void inferBlockAveragePressures(const std::size_t calcIndex, const PAvg &controls, const Scalar gravity, const Scalar refDepth)
const PAvgCalculator< Scalar >::Result & averagePressures(const std::size_t calcIndex) const
std::size_t createCalculator(const Well &well, const ParallelWellInfo< Scalar > &parallelWellInfo, const std::vector< int > &localConnIdx, EvaluatorFactory makeWellSourceEvaluator)
ParallelWBPCalculation & localCellIndex(GlobalToLocal localCellIdx)
ParallelWBPCalculation(const GridDims &cellIndexMap, const Parallel::Communication &gridComm)
typename ParallelPAvgDynamicSourceData< Scalar >::GlobalToLocal GlobalToLocal
Definition: ParallelWBPCalculation.hpp:55
ParallelWBPCalculation & evalCellSource(Evaluator evalCellSrc)
typename ParallelPAvgDynamicSourceData< Scalar >::Evaluator Evaluator
Callback for evaluating WBPn source terms on the current MPI rank.
Definition: ParallelWBPCalculation.hpp:58
std::function< Evaluator()> EvaluatorFactory
Definition: ParallelWBPCalculation.hpp:63
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:37