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 class ParallelWellInfo;
40 class PAvg;
41 class Well;
42}
43
44namespace Opm {
45
50{
51public:
55
58
62 using EvaluatorFactory = std::function<Evaluator()>;
63
69 explicit ParallelWBPCalculation(const GridDims& cellIndexMap,
70 const Parallel::Communication& gridComm);
71
80
88
113 std::size_t
114 createCalculator(const Well& well,
115 const ParallelWellInfo& parallelWellInfo,
116 const std::vector<int>& localConnIdx,
117 EvaluatorFactory makeWellSourceEvaluator);
118
126
140
153 void inferBlockAveragePressures(const std::size_t calcIndex,
154 const PAvg& controls,
155 const double gravity,
156 const double refDepth);
157
166 const PAvgCalculator::Result&
167 averagePressures(const std::size_t calcIndex) const;
168
169private:
171 class LocalConnSet
172 {
173 public:
184 explicit LocalConnSet(const std::vector<int>& localConnIdx);
185
194 int localIndex(const std::size_t connIdx) const;
195
196 private:
205 std::vector<int> localConnIdx_{};
206 };
207
209 class SourceData
210 {
211 public:
216 operator const ParallelPAvgDynamicSourceData&() const
217 {
218 return *this->srcData_;
219 }
220
224 explicit SourceData(const Parallel::Communication& comm);
225
233 SourceData& localIndex(GlobalToLocal localIdx);
234
241 SourceData& evaluator(Evaluator eval);
242
251 SourceData& evaluatorFactory(EvaluatorFactory evalFactory);
252
262 void buildStructure(const std::vector<std::size_t>& sourceLocations);
263
274
278 const Parallel::Communication& comm() const
279 {
280 return this->comm_.get();
281 }
282
292 std::vector<int>
293 getLocalIndex(const std::vector<std::size_t>& globalIndex) const;
294
295 private:
297 using DataPtr = std::unique_ptr<ParallelPAvgDynamicSourceData>;
298
300 std::reference_wrapper<const Parallel::Communication> comm_;
301
303 GlobalToLocal localIdx_{};
304
307 Evaluator eval_{};
308
312 EvaluatorFactory evalFactory_{};
313
315 DataPtr srcData_{};
316 };
317
319 using WellID = std::vector<SourceData>::size_type;
320
322 std::reference_wrapper<const GridDims> cellIndexMap_;
323
325 SourceData reservoirSrc_;
326
329 PAvgCalculatorCollection calculators_{};
330
332 std::vector<SourceData> wellConnSrc_{};
333
335 std::vector<LocalConnSet> localConnSet_{};
336
339 void pruneInactiveWBPCells();
340
342 void pruneInactiveWBPCellsSerial();
343
345 void pruneInactiveWBPCellsParallel();
346
348 void defineReservoirCommunication();
349
354 void defineWellCommunication(const std::size_t well);
355
363 PAvgCalculator::Sources makeEvaluationSources(const WellID well) const;
364};
365
366} // namespace Opm
367
368#endif // PARALLEL_WBP_CALCULATION_HPP
std::function< void(int, SourceDataSpan< double >)> Evaluator
Definition: ParallelPAvgDynamicSourceData.hpp:52
std::function< int(const std::size_t)> GlobalToLocal
Definition: ParallelPAvgDynamicSourceData.hpp:41
Definition: ParallelWBPCalculation.hpp:50
ParallelPAvgDynamicSourceData::Evaluator Evaluator
Callback for evaluating WBPn source terms on the current MPI rank.
Definition: ParallelWBPCalculation.hpp:57
ParallelPAvgDynamicSourceData::GlobalToLocal GlobalToLocal
Definition: ParallelWBPCalculation.hpp:54
void inferBlockAveragePressures(const std::size_t calcIndex, const PAvg &controls, const double gravity, const double refDepth)
std::size_t createCalculator(const Well &well, const ParallelWellInfo &parallelWellInfo, const std::vector< int > &localConnIdx, EvaluatorFactory makeWellSourceEvaluator)
ParallelWBPCalculation & evalCellSource(Evaluator evalCellSrc)
std::function< Evaluator()> EvaluatorFactory
Definition: ParallelWBPCalculation.hpp:62
ParallelWBPCalculation & localCellIndex(GlobalToLocal localCellIdx)
const PAvgCalculator::Result & averagePressures(const std::size_t calcIndex) const
ParallelWBPCalculation(const GridDims &cellIndexMap, const Parallel::Communication &gridComm)
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27