version 3.10-dev
invasionstate.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_2P_INVASIONSTATE_HH
13#define DUMUX_PNM_2P_INVASIONSTATE_HH
14
15#include <vector>
16#include <type_traits>
17#include <dune/common/std/type_traits.hh>
20
21namespace Dumux::PoreNetwork {
22
27template<class P>
29{
30 using Problem = P;
31
32 template <class T>
33 using GlobalCapillaryPressureDetector = decltype(std::declval<T>().globalCapillaryPressure());
34
35 template<class T>
36 static constexpr bool hasGlobalCapillaryPressure()
37 { return Dune::Std::is_detected<GlobalCapillaryPressureDetector, T>::value; }
38
39 enum class EventType {invasion, snapOff, none};
40
41public:
42
43 TwoPInvasionState(const Problem& problem) : problem_(problem)
44 {
45 // initialize the invasion state
46 invadedCurrentIteration_.resize(problem.gridGeometry().gridView().size(0));
47 invadedPreviousTimeStep_.resize(problem.gridGeometry().gridView().size(0));
48
49 for (auto&& element : elements(problem.gridGeometry().gridView()))
50 {
51 const auto eIdx = problem.gridGeometry().elementMapper().index(element);
52 invadedCurrentIteration_[eIdx] = problem.initialInvasionState(element);
53 invadedPreviousTimeStep_[eIdx] = invadedCurrentIteration_[eIdx];
54 }
55
56 numThroatsInvaded_ = std::count(invadedCurrentIteration_.begin(), invadedCurrentIteration_.end(), true);
57 verbose_ = getParamFromGroup<bool>(problem.paramGroup(), "InvasionState.Verbosity", true);
58 restrictToGlobalCapillaryPressure_ = getParamFromGroup<bool>(problem.paramGroup(), "InvasionState.RestrictInvasionToGlobalCapillaryPressure", false);
59
60 if constexpr (hasGlobalCapillaryPressure<Problem>())
61 {
62 if (restrictToGlobalCapillaryPressure_)
63 std::cout << "\n *** Invasion behavior is restricted by a global capillary pressure defined in the problem! *** \n" << std::endl;
64 else
65 std::cout << "\n *** WARNING: global capillary pressure defined in the problem but InvasionState.RestrictInvasionToGlobalCapillaryPressure is set to false.\n"
66 << " Invasion behavior will NOT be restricted! ***\n" << std::endl;
67 }
68 }
69
71 template<class Element>
72 bool invaded(const Element& element) const
73 {
74 const auto eIdx = problem_.gridGeometry().elementMapper().index(element);
75 return invadedCurrentIteration_[eIdx];
76 }
77
79 std::size_t numThroatsInvaded() const
80 { return numThroatsInvaded_; }
81
83 template<class SolutionVector, class GridVolumeVariables, class GridFluxVariablesCache>
84 bool update(const SolutionVector& sol, const GridVolumeVariables& gridVolVars, GridFluxVariablesCache& gridFluxVarsCache)
85 {
86 hasChangedInCurrentIteration_ = false;
87 auto fvGeometry = localView(problem_.gridGeometry());
88 auto elemVolVars = localView(gridVolVars);
89 auto elemFluxVarsCache = localView(gridFluxVarsCache);
90 for (auto&& element : elements(problem_.gridGeometry().gridView()))
91 {
92 fvGeometry.bindElement(element);
93 elemVolVars.bind(element, fvGeometry, sol);
94 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
95
96 for (auto&& scvf : scvfs(fvGeometry))
97 {
98 // checks if invasion or snap-off occurred after Newton iteration step
99 if (const auto invasionResult = invasionSwitch_(element, elemVolVars, elemFluxVarsCache[scvf]); invasionResult)
100 {
101 hasChangedInCurrentIteration_ = true;
102 if constexpr (GridFluxVariablesCache::cachingEnabled)
103 {
104 const auto eIdx = problem_.gridGeometry().elementMapper().index(element);
105 gridFluxVarsCache.cache(eIdx, scvf.index()).update(problem_, element, fvGeometry, elemVolVars, scvf, invadedCurrentIteration_[eIdx]);
106 }
107 }
108 }
109 }
110 numThroatsInvaded_ = std::count(invadedCurrentIteration_.begin(), invadedCurrentIteration_.end(), true);
111 return hasChangedInCurrentIteration_;
112 }
113
115 void reset()
116 {
117 hasChangedInCurrentIteration_ = false;
118 invadedCurrentIteration_ = invadedPreviousTimeStep_;
119 }
120
122 bool hasChanged() const
123 { return hasChangedComparedToPreviousTimestep_; }
124
127 { return hasChangedInCurrentIteration_; }
128
130 void advance()
131 {
132 hasChangedComparedToPreviousTimestep_ = (invadedPreviousTimeStep_ != invadedCurrentIteration_);
133 invadedPreviousTimeStep_ = invadedCurrentIteration_;
134 }
135
136 template<class SolutionVector, class GridVolumeVariables, class GridFluxVariablesCache>
138 const GridVolumeVariables& gridVolVars,
139 const GridFluxVariablesCache& gridFluxVarsCache) const
140 {
141 using Scalar = typename SolutionVector::block_type::value_type;
142 static const Scalar accuracyCriterion = getParamFromGroup<Scalar>(problem_.paramGroup(), "InvasionState.AccuracyCriterion", -1.0);
143
144 if (accuracyCriterion < 0.0)
145 return;
146
147 auto fvGeometry = localView(problem_.gridGeometry());
148 auto elemVolVars = localView(gridVolVars);
149 auto elemFluxVarsCache = localView(gridFluxVarsCache);
150 for (auto&& element : elements(problem_.gridGeometry().gridView()))
151 {
152 // Only consider throats which have been invaded during the current time step
153 const auto eIdx = problem_.gridGeometry().elementMapper().index(element);
154 if (!invadedCurrentIteration_[eIdx] || invadedPreviousTimeStep_[eIdx] == invadedCurrentIteration_[eIdx])
155 continue;
156
157 fvGeometry.bindElement(element);
158 elemVolVars.bind(element, fvGeometry, sol);
159 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
160
161 for (auto&& scvf : scvfs(fvGeometry))
162 {
163 // checks if pc is close enough to the entry pressure value
164 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
165
166 using std::max;
167 const Scalar pc = max(elemVolVars[0].capillaryPressure(), elemVolVars[1].capillaryPressure());
168
169 if (pc < accuracyCriterion * fluxVarsCache.pcEntry())
170 DUNE_THROW(NumericalProblem, "At element " << eIdx << ": pc " << pc << " too far away form pcEntry " << fluxVarsCache.pcEntry());
171 }
172 }
173 }
174
175private:
176
178 template<class Element, class ElementVolumeVariables, class FluxVariablesCache>
179 auto invasionSwitch_(const Element& element,
180 const ElementVolumeVariables& elemVolVars,
181 const FluxVariablesCache& fluxVarsCache)
182
183 {
184 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
185 const auto& gridGeometry = problem_.gridGeometry();
186 const auto& spatialParams = problem_.spatialParams();
187 const auto eIdx = gridGeometry.elementMapper().index(element);
188 bool invadedBeforeSwitch = invadedCurrentIteration_[eIdx];
189 bool invadedAfterSwitch = invadedBeforeSwitch;
190
191 // Result type, containing the local scv index of the pore from which the invasion/snap-off occurred
192 // Evaluates to 'false' if no invasion/snap-off occurred
193 struct Result
194 {
195 std::uint8_t localScvIdxWithCriticalPc;
196 Scalar criticalPc;
197 EventType event = EventType::none;
198
199 operator bool() const
200 { return event != EventType::none; }
201 };
202
203 // Block non-wetting phase flux out of the outlet
204 static const auto blockNonwettingPhase = getParamFromGroup<std::vector<int>>(problem_.paramGroup(), "InvasionState.BlockNonwettingPhaseAtThroatLabel", std::vector<int>{});
205 if (!blockNonwettingPhase.empty() && std::find(blockNonwettingPhase.begin(), blockNonwettingPhase.end(), gridGeometry.throatLabel(eIdx)) != blockNonwettingPhase.end())
206 {
207 invadedCurrentIteration_[eIdx] = false;
208 return Result{}; // nothing happened
209 }
210
211 //Determine whether throat gets invaded or snap-off occurs
212 const std::array<Scalar, 2> pc = { elemVolVars[0].capillaryPressure(), elemVolVars[1].capillaryPressure() };
213 const auto pcMax = std::max_element(pc.begin(), pc.end());
214 const Scalar pcEntry = fluxVarsCache.pcEntry();
215 const Scalar pcSnapoff = fluxVarsCache.pcSnapoff();
216
217 // check if there is a user-specified global capillary pressure which needs to be obeyed
218 if (maybeRestrictToGlobalCapillaryPressure_(pcEntry))
219 {
220 if (*pcMax > pcEntry)
221 {
222 std::cout << "Throat " << eIdx << " would have been invaded by pc of " << *pcMax << "but a global capillary pressure restricion was set in the problem.";
223 std::cout << ". pcEntry: " << spatialParams.pcEntry(element, elemVolVars) << std::endl;
224 }
225
226 invadedCurrentIteration_[eIdx] = false;
227 return Result{}; //nothing happened
228 }
229
230 if (*pcMax > pcEntry)
231 invadedAfterSwitch = true;
232 else if (*pcMax <= pcSnapoff)
233 invadedAfterSwitch = false;
234
235 invadedCurrentIteration_[eIdx] = invadedAfterSwitch;
236
237 if (invadedBeforeSwitch == invadedAfterSwitch)
238 return Result{}; // nothing happened
239 else
240 {
241 Result result;
242 result.localScvIdxWithCriticalPc = std::distance(pc.begin(), pcMax);
243 result.criticalPc = *pcMax;
244 result.event = !invadedBeforeSwitch && invadedAfterSwitch ? EventType::invasion : EventType::snapOff;
245
246 if (verbose_)
247 {
248 const auto wPhaseIdx = spatialParams.template wettingPhase<typename ElementVolumeVariables::VolumeVariables::FluidSystem>(element, elemVolVars);
249 const std::array sw = { elemVolVars[0].saturation(wPhaseIdx), elemVolVars[1].saturation(wPhaseIdx) };
250 const auto vIdx = gridGeometry.gridView().indexSet().subIndex(element, result.localScvIdxWithCriticalPc, 1);
251 if (result.event == EventType::invasion)
252 {
253 std::cout << "Throat " << eIdx << " was invaded from pore " << vIdx << " :";
254 std::cout << " pc: " << *pcMax;
255 std::cout << ", pcEntry: " << spatialParams.pcEntry(element, elemVolVars);
256 std::cout << ", sw: " << sw[result.localScvIdxWithCriticalPc] << std::endl;
257 }
258 else
259 {
260 std::cout << "Snap-off occurred at throat " << eIdx << " from pore " << vIdx << " :";
261 std::cout << " pc: " << *pcMax;
262 std::cout << ", pcSnapoff: " << spatialParams.pcSnapoff(element, elemVolVars);
263 std::cout << ", sw: " << sw[result.localScvIdxWithCriticalPc] << std::endl;
264 }
265 }
266
267 return result;
268 }
269 }
270
273 template<class Scalar>
274 bool maybeRestrictToGlobalCapillaryPressure_(const Scalar pcEntry) const
275 {
276 if constexpr (hasGlobalCapillaryPressure<Problem>())
277 return restrictToGlobalCapillaryPressure_ && (pcEntry > problem_.globalCapillaryPressure());
278 else
279 return false;
280 }
281
282 std::vector<bool> invadedCurrentIteration_;
283 std::vector<bool> invadedPreviousTimeStep_;
284 bool hasChangedInCurrentIteration_ = false;
285 bool hasChangedComparedToPreviousTimestep_ = false;
286 std::size_t numThroatsInvaded_;
287 bool verbose_;
288 bool restrictToGlobalCapillaryPressure_;
289
290 const Problem& problem_;
291};
292
293} // end namespace Dumux::PoreNetwork
294
295#endif
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
This class updates the invasion state for the two-phase PNM.
Definition: invasionstate.hh:29
void checkIfCapillaryPressureIsCloseToEntryPressure(const SolutionVector &sol, const GridVolumeVariables &gridVolVars, const GridFluxVariablesCache &gridFluxVarsCache) const
Definition: invasionstate.hh:137
bool hasChangedInCurrentIteration() const
Return whether an invasion or snap-off occurred anywhere during the current Newton iteration.
Definition: invasionstate.hh:126
bool invaded(const Element &element) const
Return whether a given throat is invaded or not.
Definition: invasionstate.hh:72
bool hasChanged() const
Return whether an invasion or snap-off occurred anywhere. Can be used, e.g., for output file writing ...
Definition: invasionstate.hh:122
bool update(const SolutionVector &sol, const GridVolumeVariables &gridVolVars, GridFluxVariablesCache &gridFluxVarsCache)
Update the invasion state of all throats. This is done after each Newton step by a call from the Newt...
Definition: invasionstate.hh:84
TwoPInvasionState(const Problem &problem)
Definition: invasionstate.hh:43
void advance()
This is called after the Newton method has successfully finished one time step.
Definition: invasionstate.hh:130
std::size_t numThroatsInvaded() const
Return the number of currently invaded throats.
Definition: invasionstate.hh:79
void reset()
Restore the old invasion state after a Newton iteration has failed.
Definition: invasionstate.hh:115
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
Defines labels for pores and throats.
std::string capillaryPressure() noexcept
I/O name of capillary pressure.
Definition: name.hh:123
Scalar pcSnapoff(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius, const Throat::Shape shape)
The snap-off capillary pressure of a pore throat It checks if the cross section shape of the throat i...
Definition: thresholdcapillarypressures.hh:65
Scalar pcEntry(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius, const Scalar shapeFactor) noexcept
The entry capillary pressure of a pore throat.
Definition: thresholdcapillarypressures.hh:82
Definition: discretization/porenetwork/fvelementgeometry.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.