SatfuncConsistencyCheckManager.hpp
Go to the documentation of this file.
1/*
2 Copyright 2024 Equinor 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
20#ifndef SATFUNC_CONSISTENCY_CHECK_MANAGER_HPP_INCLUDED
21#define SATFUNC_CONSISTENCY_CHECK_MANAGER_HPP_INCLUDED
22
23#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
24
25#include <opm/input/eclipse/EclipseState/Grid/SatfuncPropertyInitializers.hpp>
26
27#include <opm/material/fluidmatrixinteractions/EclEpsGridProperties.hpp>
28
30
33
34#include <dune/grid/common/partitionset.hh>
35#include <dune/grid/common/rangegenerators.hh>
36
37#include <cstddef>
38#include <functional>
39#include <memory>
40#include <string>
41#include <string_view>
42#include <type_traits>
43#include <vector>
44
46 template <typename Scalar>
47 class UnscaledSatfuncCheckPoint;
48} // namespace Opm::Satfunc::PhaseChecks
49
51
56 template <typename Scalar>
58 {
59 public:
62 using LocalToGlobal = std::function<std::size_t(int)>;
63
66 using ReportRecordOutput = typename
68
70 using ViolationLevel = typename
72
89 explicit SatfuncConsistencyCheckManager(const std::size_t numSamplePoints,
90 const EclipseState& eclipseState,
91 const LocalToGlobal& localToGlobal);
92
100 {
101 this->root_ = root;
102 return *this;
103 }
104
127 template <typename GridView, typename GetCellIndex>
128 void run(const GridView& gv, GetCellIndex&& getCellIndex)
129 {
130 this->isRoot_ = gv.comm().rank() == this->root_;
131
132 this->warnIfDirectionalOrIrreversibleEPS();
133
134 for (const auto& elem : elements(gv, Dune::Partitions::interior)) {
135 this->runCellChecks(getCellIndex(elem));
136 }
137
138 gv.comm().barrier();
139
140 this->collectFailures(gv.comm());
141 }
142
145
148
166 const ReportRecordOutput& emitReportRecord) const;
167
168 private:
172 struct CurveCollection
173 {
196 explicit CurveCollection(std::unique_ptr<SatfuncCheckPointInterface<Scalar>> point,
197 std::string_view pointName,
198 const std::size_t numSamplePoints);
199
205 std::unique_ptr<SatfuncCheckPointInterface<Scalar>> point;
206
209 };
210
216 std::reference_wrapper<const EclipseState> eclipseState_;
217
223 LocalToGlobal localToGlobal_;
224
229 satfunc::RawTableEndPoints rtep_{};
230
235 satfunc::RawFunctionValues rfunc_{};
236
242 std::vector<EclEpsGridProperties> gridProps_{};
243
246 std::vector<CurveCollection> curves_{};
247
249 int root_{0};
250
253 bool isRoot_{false};
254
260 void warnIfDirectionalOrIrreversibleEPS() const;
261
267 void runCellChecks(const int cellIdx);
268
273 void configureCurveChecks(const std::size_t numSamplePoints);
274
291 std::unique_ptr<UnscaledSatfuncCheckPoint<Scalar>>
292 configureUnscaledCurveChecks(const std::string& regionName,
293 const std::size_t numSamplePoints);
294
308 void configureScaledCurveChecks(const UnscaledSatfuncCheckPoint<Scalar>& unscaledChecks,
309 const bool useImbibition,
310 const std::size_t numSamplePoints);
311
313 void addChecks();
314
322 void collectFailures(const Parallel::Communication& comm);
323
334 template <typename Body>
335 void curveLoop(Body&& body);
336
347 template <typename Body>
348 void curveLoop(Body&& body) const;
349 };
350
351} // namespace Opm::Satfunc::PhaseChecks
352
353#endif // SATFUNC_CONSISTENCY_CHECK_MANAGER_HPP_INCLUDED
Definition: SatfuncConsistencyCheckManager.hpp:58
SatfuncConsistencyCheckManager & collectFailuresTo(const int root)
Definition: SatfuncConsistencyCheckManager.hpp:99
bool anyFailedStandardChecks() const
Whether or not any checks failed at the Standard level.
void run(const GridView &gv, GetCellIndex &&getCellIndex)
Definition: SatfuncConsistencyCheckManager.hpp:128
std::function< std::size_t(int)> LocalToGlobal
Definition: SatfuncConsistencyCheckManager.hpp:62
typename SatfuncConsistencyChecks< Scalar >::ViolationLevel ViolationLevel
Severity level for consistency condition violation.
Definition: SatfuncConsistencyCheckManager.hpp:71
SatfuncConsistencyCheckManager(const std::size_t numSamplePoints, const EclipseState &eclipseState, const LocalToGlobal &localToGlobal)
bool anyFailedCriticalChecks() const
Whether or not any checks failed at the Critical level.
void reportFailures(const ViolationLevel level, const ReportRecordOutput &emitReportRecord) const
typename SatfuncConsistencyChecks< Scalar >::ReportRecordOutput ReportRecordOutput
Definition: SatfuncConsistencyCheckManager.hpp:67
Definition: SatfuncConsistencyChecks.hpp:46
ViolationLevel
Severity level for consistency condition violation.
Definition: SatfuncConsistencyChecks.hpp:107
std::function< void(std::string_view)> ReportRecordOutput
Definition: SatfuncConsistencyChecks.hpp:122
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: GasPhaseConsistencyChecks.hpp:29
Definition: SatfuncCheckPointInterface.hpp:40