12#ifndef DUMUX_PNM_NEWTON_CONSISTENCY_CHECKS
13#define DUMUX_PNM_NEWTON_CONSISTENCY_CHECKS
17#include <dune/common/exceptions.hh>
28template <
typename T,
typename U =
int>
29struct SaturationIndex {};
33struct SaturationIndex <T, decltype((void) T::saturationIdx, 0)>
35 static constexpr auto value = T::saturationIdx;
40struct SaturationIndex <T, decltype((void) T::switchIdx, 0)>
42 static constexpr auto value = T::switchIdx;
47struct SaturationIndex <T, decltype((void) T::s0Idx, 0)>
49 static constexpr auto value = T::s0Idx;
58template<
class Gr
idVariables,
class SolutionVector>
66 void performChecks(
const GridVariables& gridVariables,
const SolutionVector& uCurrentIter,
const SolutionVector& prevSol)
const
78 using Scalar =
typename SolutionVector::block_type::value_type;
79 const auto& problem = gridVariables.curGridVolVars().problem();
81 static const Scalar allowedSaturationChange = getParamFromGroup<Scalar>(problem.paramGroup(),
"Newton.AllowedSaturationChange", -1.0);
82 if (allowedSaturationChange < 0.0)
85 const auto& curGridVolVars = gridVariables.curGridVolVars();
86 const auto& prevGridVolVars = gridVariables.prevGridVolVars();
88 std::vector<bool> dofVisited(uCurrentIter.size(),
false);
90 auto fvGeometry =
localView(problem.gridGeometry());
91 auto curElemVolVars =
localView(curGridVolVars);
92 auto prevElemVolVars =
localView(prevGridVolVars);
94 for (
const auto& element : elements(problem.gridGeometry().gridView()))
96 fvGeometry.bindElement(element);
97 curElemVolVars.bindElement(element, fvGeometry, uCurrentIter);
98 prevElemVolVars.bindElement(element, fvGeometry, prevSol);
100 for (
const auto& scv : scvs(fvGeometry))
102 if (dofVisited[scv.dofIndex()])
104 dofVisited[scv.dofIndex()] =
true;
106 const Scalar satNew = curElemVolVars[scv].saturation(0);
107 const Scalar satOld = prevElemVolVars[scv].saturation(0);
109 const Scalar deltaS = abs(satNew - satOld);
110 static const bool considerRelativeShift = getParamFromGroup<bool>(problem.paramGroup(),
"Newton.SaturationChangeIsRelative",
false);
112 if (considerRelativeShift)
115 if (satOld > 1e-3 && deltaS / satOld > allowedSaturationChange)
116 DUNE_THROW(
NumericalProblem,
"Saturation change too high at dof " << scv.dofIndex() <<
", old sat. " << satOld <<
", new sat. " << satNew << std::endl);
118 else if (deltaS > allowedSaturationChange)
119 DUNE_THROW(
NumericalProblem,
"Saturation change too high at dof " << scv.dofIndex() <<
", old sat. " << satOld <<
", new sat. " << satNew << std::endl);
129 const auto& problem = gridVariables.curGridVolVars().problem();
131 static const bool doPlausibilityCheck = getParamFromGroup<bool>(problem.paramGroup(),
"Newton.PlausibilityCheck",
false);
132 if (!doPlausibilityCheck)
135 for (std::size_t i = 0; i < uCurrentIter.size(); ++i)
137 const auto& priVars = uCurrentIter[i];
138 using Indices =
typename GridVariables::VolumeVariables::Indices;
139 static constexpr auto saturationOrMoleFractionIndex = Detail::SaturationIndex<Indices>::value;
140 if (priVars[saturationOrMoleFractionIndex] < 0.0 || priVars[saturationOrMoleFractionIndex] > 1.0)
142 std::cout <<
"at dof " << i <<
": saturation " << priVars[saturationOrMoleFractionIndex] << std::endl;
154 const auto& invasionState = gridVariables.gridFluxVarsCache().invasionState();
155 invasionState.checkIfCapillaryPressureIsCloseToEntryPressure(uCurrentIter, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
Consistency checks for the PNM two-phase model.
Definition: newtonconsistencychecks.hh:60
void checkRelativeSaturationShift(const GridVariables &gridVariables, const SolutionVector &uCurrentIter, const SolutionVector &prevSol) const
Checks if the relative shift of saturation between to consecutive time steps is below a given thresho...
Definition: newtonconsistencychecks.hh:76
void checkIfValuesArePhysical(const GridVariables &gridVariables, const SolutionVector &uCurrentIter) const
Checks if the saturation is between zero and one.
Definition: newtonconsistencychecks.hh:127
void performChecks(const GridVariables &gridVariables, const SolutionVector &uCurrentIter, const SolutionVector &prevSol) const
Perform all checks.
Definition: newtonconsistencychecks.hh:66
void checkIfCapillaryPressureIsCloseToEntryPressure(const GridVariables &gridVariables, const SolutionVector &uCurrentIter) const
Checks if the capillary pressure at pore from which a throat was invaded is sufficiently close to the...
Definition: newtonconsistencychecks.hh:151
Some exceptions thrown in DuMux
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
Free function to get the local view of a grid cache object.
Definition: discretization/porenetwork/fvelementgeometry.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.