version 3.10-dev
newtonconsistencychecks.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_PNM_NEWTON_CONSISTENCY_CHECKS
13#define DUMUX_PNM_NEWTON_CONSISTENCY_CHECKS
14
15#include <vector>
16#include <iostream>
17#include <dune/common/exceptions.hh>
21
22namespace Dumux::PoreNetwork {
23
24#ifndef DOXYGEN
25namespace Detail {
26
27// Primary template
28template <typename T, typename U = int>
29struct SaturationIndex {};
30
31// Specialization for 2p
32template <typename T>
33struct SaturationIndex <T, decltype((void) T::saturationIdx, 0)>
34{
35 static constexpr auto value = T::saturationIdx;
36};
37
38// Specialization for 2pnc
39template <typename T>
40struct SaturationIndex <T, decltype((void) T::switchIdx, 0)>
41{
42 static constexpr auto value = T::switchIdx;
43};
44
45// Specialization for MPNC
46template <typename T>
47struct SaturationIndex <T, decltype((void) T::s0Idx, 0)>
48{
49 static constexpr auto value = T::s0Idx;
50};
51} // end namespace Detail
52#endif
53
58template<class GridVariables, class SolutionVector>
60{
61public:
62
66 void performChecks(const GridVariables& gridVariables, const SolutionVector& uCurrentIter, const SolutionVector& prevSol) const
67 {
68 checkRelativeSaturationShift(gridVariables, uCurrentIter, prevSol);
69 checkIfValuesArePhysical(gridVariables, uCurrentIter);
70 checkIfCapillaryPressureIsCloseToEntryPressure(gridVariables, uCurrentIter);
71 }
72
76 void checkRelativeSaturationShift(const GridVariables& gridVariables, const SolutionVector& uCurrentIter, const SolutionVector& prevSol) const
77 {
78 using Scalar = typename SolutionVector::block_type::value_type;
79 const auto& problem = gridVariables.curGridVolVars().problem();
80
81 static const Scalar allowedSaturationChange = getParamFromGroup<Scalar>(problem.paramGroup(), "Newton.AllowedSaturationChange", -1.0);
82 if (allowedSaturationChange < 0.0)
83 return;
84
85 const auto& curGridVolVars = gridVariables.curGridVolVars();
86 const auto& prevGridVolVars = gridVariables.prevGridVolVars();
87
88 std::vector<bool> dofVisited(uCurrentIter.size(), false);
89
90 auto fvGeometry = localView(problem.gridGeometry());
91 auto curElemVolVars = localView(curGridVolVars);
92 auto prevElemVolVars = localView(prevGridVolVars);
93
94 for (const auto& element : elements(problem.gridGeometry().gridView()))
95 {
96 fvGeometry.bindElement(element);
97 curElemVolVars.bindElement(element, fvGeometry, uCurrentIter);
98 prevElemVolVars.bindElement(element, fvGeometry, prevSol);
99
100 for (const auto& scv : scvs(fvGeometry))
101 {
102 if (dofVisited[scv.dofIndex()])
103 continue;
104 dofVisited[scv.dofIndex()] = true;
105
106 const Scalar satNew = curElemVolVars[scv].saturation(0);
107 const Scalar satOld = prevElemVolVars[scv].saturation(0);
108 using std::abs;
109 const Scalar deltaS = abs(satNew - satOld);
110 static const bool considerRelativeShift = getParamFromGroup<bool>(problem.paramGroup(), "Newton.SaturationChangeIsRelative", false);
111
112 if (considerRelativeShift)
113 {
114 // satOld mgiht be (close to) zero, so have to take care of this
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);
117 }
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);
120 }
121 }
122 }
123
127 void checkIfValuesArePhysical(const GridVariables& gridVariables, const SolutionVector& uCurrentIter) const
128 {
129 const auto& problem = gridVariables.curGridVolVars().problem();
130
131 static const bool doPlausibilityCheck = getParamFromGroup<bool>(problem.paramGroup(), "Newton.PlausibilityCheck", false);
132 if (!doPlausibilityCheck)
133 return;
134
135 for (std::size_t i = 0; i < uCurrentIter.size(); ++i)
136 {
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)
141 {
142 std::cout << "at dof " << i << ": saturation " << priVars[saturationOrMoleFractionIndex] << std::endl;
143 DUNE_THROW(NumericalProblem, "Saturation is below 0 or above 1");
144 }
145 }
146 }
147
151 void checkIfCapillaryPressureIsCloseToEntryPressure(const GridVariables& gridVariables, const SolutionVector& uCurrentIter) const
152 {
153 // this check is implemented in the invasion state itself
154 const auto& invasionState = gridVariables.gridFluxVarsCache().invasionState();
155 invasionState.checkIfCapillaryPressureIsCloseToEntryPressure(uCurrentIter, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
156 }
157};
158} // end namespace Dumux
159#endif
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.