20#ifndef SATFUNC_CONSISTENCY_CHECK_MANAGER_HPP_INCLUDED
21#define SATFUNC_CONSISTENCY_CHECK_MANAGER_HPP_INCLUDED
23#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
25#include <opm/input/eclipse/EclipseState/Grid/SatfuncPropertyInitializers.hpp>
27#include <opm/material/fluidmatrixinteractions/EclEpsGridProperties.hpp>
34#include <dune/grid/common/partitionset.hh>
35#include <dune/grid/common/rangegenerators.hh>
46 template <
typename Scalar>
47 class UnscaledSatfuncCheckPoint;
56 template <
typename Scalar>
90 const EclipseState& eclipseState,
127 template <
typename Gr
idView,
typename GetCellIndex>
128 void run(
const GridView& gv, GetCellIndex&& getCellIndex)
130 this->isRoot_ = gv.comm().rank() == this->root_;
132 this->warnIfDirectionalOrIrreversibleEPS();
134 for (
const auto& elem : elements(gv, Dune::Partitions::interior)) {
135 this->runCellChecks(getCellIndex(elem));
140 this->collectFailures(gv.comm());
172 struct CurveCollection
197 std::string_view pointName,
198 const std::size_t numSamplePoints);
205 std::unique_ptr<SatfuncCheckPointInterface<Scalar>> point;
216 std::reference_wrapper<const EclipseState> eclipseState_;
229 satfunc::RawTableEndPoints rtep_{};
235 satfunc::RawFunctionValues rfunc_{};
242 std::vector<EclEpsGridProperties> gridProps_{};
246 std::vector<CurveCollection> curves_{};
260 void warnIfDirectionalOrIrreversibleEPS()
const;
267 void runCellChecks(
const int cellIdx);
273 void configureCurveChecks(
const std::size_t numSamplePoints);
291 std::unique_ptr<UnscaledSatfuncCheckPoint<Scalar>>
292 configureUnscaledCurveChecks(
const std::string& regionName,
293 const std::size_t numSamplePoints);
308 void configureScaledCurveChecks(
const UnscaledSatfuncCheckPoint<Scalar>& unscaledChecks,
309 const bool useImbibition,
310 const std::size_t numSamplePoints);
334 template <
typename Body>
335 void curveLoop(Body&& body);
347 template <
typename Body>
348 void curveLoop(Body&& body)
const;
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