24#ifndef DUMUX_PNM_NEWTON_CONSISTENCY_CHECKS
25#define DUMUX_PNM_NEWTON_CONSISTENCY_CHECKS
29#include <dune/common/exceptions.hh>
40template <
typename T,
typename U =
int>
41struct SaturationIndex {};
45struct SaturationIndex <T, decltype((void) T::saturationIdx, 0)>
47 static constexpr auto value = T::saturationIdx;
52struct SaturationIndex <T, decltype((void) T::switchIdx, 0)>
54 static constexpr auto value = T::switchIdx;
59struct SaturationIndex <T, decltype((void) T::s0Idx, 0)>
61 static constexpr auto value = T::s0Idx;
70template<
class Gr
idVariables,
class SolutionVector>
78 void performChecks(
const GridVariables& gridVariables,
const SolutionVector& uCurrentIter,
const SolutionVector& prevSol)
const
90 using Scalar =
typename SolutionVector::block_type::value_type;
91 const auto& problem = gridVariables.curGridVolVars().problem();
93 static const Scalar allowedSaturationChange = getParamFromGroup<Scalar>(problem.paramGroup(),
"Newton.AllowedSaturationChange", -1.0);
94 if (allowedSaturationChange < 0.0)
97 const auto& curGridVolVars = gridVariables.curGridVolVars();
98 const auto& prevGridVolVars = gridVariables.prevGridVolVars();
100 std::vector<bool> dofVisited(uCurrentIter.size(),
false);
102 auto fvGeometry =
localView(problem.gridGeometry());
103 auto curElemVolVars =
localView(curGridVolVars);
104 auto prevElemVolVars =
localView(prevGridVolVars);
106 for (
const auto& element : elements(problem.gridGeometry().gridView()))
108 fvGeometry.bindElement(element);
109 curElemVolVars.bindElement(element, fvGeometry, uCurrentIter);
110 prevElemVolVars.bindElement(element, fvGeometry, prevSol);
112 for (
const auto& scv : scvs(fvGeometry))
114 if (dofVisited[scv.dofIndex()])
116 dofVisited[scv.dofIndex()] =
true;
118 const Scalar satNew = curElemVolVars[scv].saturation(0);
119 const Scalar satOld = prevElemVolVars[scv].saturation(0);
121 const Scalar deltaS = abs(satNew - satOld);
122 static const bool considerRelativeShift = getParamFromGroup<bool>(problem.paramGroup(),
"Newton.SaturationChangeIsRelative",
false);
124 if (considerRelativeShift)
127 if (satOld > 1e-3 && deltaS / satOld > allowedSaturationChange)
128 DUNE_THROW(
NumericalProblem,
"Saturation change too high at dof " << scv.dofIndex() <<
", old sat. " << satOld <<
", new sat. " << satNew << std::endl);
130 else if (deltaS > allowedSaturationChange)
131 DUNE_THROW(
NumericalProblem,
"Saturation change too high at dof " << scv.dofIndex() <<
", old sat. " << satOld <<
", new sat. " << satNew << std::endl);
141 const auto& problem = gridVariables.curGridVolVars().problem();
143 static const bool doPlausibilityCheck = getParamFromGroup<bool>(problem.paramGroup(),
"Newton.PlausibilityCheck",
false);
144 if (!doPlausibilityCheck)
147 for (std::size_t i = 0; i < uCurrentIter.size(); ++i)
149 const auto& priVars = uCurrentIter[i];
150 using Indices =
typename GridVariables::VolumeVariables::Indices;
151 static constexpr auto saturationOrMoleFractionIndex = Detail::SaturationIndex<Indices>::value;
152 if (priVars[saturationOrMoleFractionIndex] < 0.0 || priVars[saturationOrMoleFractionIndex] > 1.0)
154 std::cout <<
"at dof " << i <<
": saturation " << priVars[saturationOrMoleFractionIndex] << std::endl;
166 const auto& invasionState = gridVariables.gridFluxVarsCache().invasionState();
167 invasionState.checkIfCapillaryPressureIsCloseToEntryPressure(uCurrentIter, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Free function to get the local view of a grid cache object.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: discretization/porenetwork/fvelementgeometry.hh:34
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Consistency checks for the PNM two-phase model.
Definition: newtonconsistencychecks.hh:72
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 treshol...
Definition: newtonconsistencychecks.hh:88
void checkIfValuesArePhysical(const GridVariables &gridVariables, const SolutionVector &uCurrentIter) const
Checks if the saturation is between zero and one.
Definition: newtonconsistencychecks.hh:139
void performChecks(const GridVariables &gridVariables, const SolutionVector &uCurrentIter, const SolutionVector &prevSol) const
Perform all checks.
Definition: newtonconsistencychecks.hh:78
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:163