Loading [MathJax]/extensions/tex2jax.js
3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PNM_2P_INVASIONSTATE_HH
25#define DUMUX_PNM_2P_INVASIONSTATE_HH
26
27#include <vector>
28#include <type_traits>
29#include <dune/common/std/type_traits.hh>
32
33namespace Dumux::PoreNetwork {
34
39template<class P>
41{
42 using Problem = P;
43
44 template <class T>
45 using GlobalCapillaryPressureDetector = decltype(std::declval<T>().globalCapillaryPressure());
46
47 template<class T>
48 static constexpr bool hasGlobalCapillaryPressure()
49 { return Dune::Std::is_detected<GlobalCapillaryPressureDetector, T>::value; }
50
51 enum class EventType {invasion, snapOff, none};
52
53public:
54
55 TwoPInvasionState(const Problem& problem) : problem_(problem)
56 {
57 // initialize the invasion state
58 invadedCurrentIteration_.resize(problem.gridGeometry().gridView().size(0));
59 invadedPreviousTimeStep_.resize(problem.gridGeometry().gridView().size(0));
60
61 for (auto&& element : elements(problem.gridGeometry().gridView()))
62 {
63 const auto eIdx = problem.gridGeometry().elementMapper().index(element);
64 invadedCurrentIteration_[eIdx] = problem.initialInvasionState(element);
65 invadedPreviousTimeStep_[eIdx] = invadedCurrentIteration_[eIdx];
66 }
67
68 numThroatsInvaded_ = std::count(invadedCurrentIteration_.begin(), invadedCurrentIteration_.end(), true);
69 verbose_ = getParamFromGroup<bool>(problem.paramGroup(), "InvasionState.Verbosity", true);
70 restrictToGlobalCapillaryPressure_ = getParamFromGroup<bool>(problem.paramGroup(), "InvasionState.RestrictInvasionToGlobalCapillaryPressure", false);
71
72 if constexpr (hasGlobalCapillaryPressure<Problem>())
73 {
74 if (restrictToGlobalCapillaryPressure_)
75 std::cout << "\n *** Invasion behavior is restricted by a global capillary pressure defined in the problem! *** \n" << std::endl;
76 else
77 std::cout << "\n *** WARNING: global capillary pressure defined in the problem but InvasionState.RestrictInvasionToGlobalCapillaryPressure is set to false.\n"
78 << " Invasion behavior will NOT be restricted! ***\n" << std::endl;
79 }
80 }
81
83 template<class Element>
84 bool invaded(const Element& element) const
85 {
86 const auto eIdx = problem_.gridGeometry().elementMapper().index(element);
87 return invadedCurrentIteration_[eIdx];
88 }
89
91 std::size_t numThroatsInvaded() const
92 { return numThroatsInvaded_; }
93
95 template<class SolutionVector, class GridVolumeVariables, class GridFluxVariablesCache>
96 bool update(const SolutionVector& sol, const GridVolumeVariables& gridVolVars, GridFluxVariablesCache& gridFluxVarsCache)
97 {
98 hasChangedInCurrentIteration_ = false;
99 for (auto&& element : elements(problem_.gridGeometry().gridView()))
100 {
101 auto fvGeometry = localView(problem_.gridGeometry());
102 fvGeometry.bindElement(element);
103
104 auto elemVolVars = localView(gridVolVars);
105 elemVolVars.bind(element, fvGeometry, sol);
106
107 auto elemFluxVarsCache = localView(gridFluxVarsCache);
108 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
109
110 for (auto&& scvf : scvfs(fvGeometry))
111 {
112 // checks if invasion or snap-off occured after Newton iteration step
113 if (const auto invasionResult = invasionSwitch_(element, elemVolVars, elemFluxVarsCache[scvf]); invasionResult)
114 {
115 hasChangedInCurrentIteration_ = true;
116 if constexpr (GridFluxVariablesCache::cachingEnabled)
117 {
118 const auto eIdx = problem_.gridGeometry().elementMapper().index(element);
119 gridFluxVarsCache.cache(eIdx, scvf.index()).update(problem_, element, fvGeometry, elemVolVars, scvf, invadedCurrentIteration_[eIdx]);
120 }
121 }
122 }
123 }
124 numThroatsInvaded_ = std::count(invadedCurrentIteration_.begin(), invadedCurrentIteration_.end(), true);
125 return hasChangedInCurrentIteration_;
126 }
127
129 void reset()
130 {
131 hasChangedInCurrentIteration_ = false;
132 invadedCurrentIteration_ = invadedPreviousTimeStep_;
133 }
134
136 bool hasChanged() const
137 { return hasChangedComparedToPreviousTimestep_; }
138
141 { return hasChangedInCurrentIteration_; }
142
144 void advance()
145 {
146 hasChangedComparedToPreviousTimestep_ = (invadedPreviousTimeStep_ != invadedCurrentIteration_);
147 invadedPreviousTimeStep_ = invadedCurrentIteration_;
148 }
149
150 template<class SolutionVector, class GridVolumeVariables, class GridFluxVariablesCache>
152 const GridVolumeVariables& gridVolVars,
153 const GridFluxVariablesCache& gridFluxVarsCache) const
154 {
155 using Scalar = typename SolutionVector::block_type::value_type;
156 static const Scalar accuracyCriterion = getParamFromGroup<Scalar>(problem_.paramGroup(), "InvasionState.AccuracyCriterion", -1.0);
157
158 if (accuracyCriterion < 0.0)
159 return;
160
161 for (auto&& element : elements(problem_.gridGeometry().gridView()))
162 {
163 // Only consider throats which have been invaded during the current time step
164 const auto eIdx = problem_.gridGeometry().elementMapper().index(element);
165 if (!invadedCurrentIteration_[eIdx] || invadedPreviousTimeStep_[eIdx] == invadedCurrentIteration_[eIdx])
166 continue;
167
168 auto fvGeometry = localView(problem_.gridGeometry());
169 fvGeometry.bindElement(element);
170
171 auto elemVolVars = localView(gridVolVars);
172 elemVolVars.bind(element, fvGeometry, sol);
173
174 auto elemFluxVarsCache = localView(gridFluxVarsCache);
175 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
176
177 for (auto&& scvf : scvfs(fvGeometry))
178 {
179 // checks if pc is close enough to the entry pressure value
180 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
181
182 using std::max;
183 const Scalar pc = max(elemVolVars[0].capillaryPressure(), elemVolVars[1].capillaryPressure());
184
185 if (pc < accuracyCriterion * fluxVarsCache.pcEntry())
186 DUNE_THROW(NumericalProblem, "At element " << eIdx << ": pc " << pc << " too far away form pcEntry " << fluxVarsCache.pcEntry());
187 }
188 }
189 }
190
191private:
192
194 template<class Element, class ElementVolumeVariables, class FluxVariablesCache>
195 auto invasionSwitch_(const Element& element,
196 const ElementVolumeVariables& elemVolVars,
197 const FluxVariablesCache& fluxVarsCache)
198
199 {
200 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
201 const auto& gridGeometry = problem_.gridGeometry();
202 const auto& spatialParams = problem_.spatialParams();
203 const auto eIdx = gridGeometry.elementMapper().index(element);
204 bool invadedBeforeSwitch = invadedCurrentIteration_[eIdx];
205 bool invadedAfterSwitch = invadedBeforeSwitch;
206
207 // Result type, containing the local scv index of the pore from which the invasion/snap-off occurred
208 // Evaluates to 'false' if no invasion/snap-off occurred
209 struct Result
210 {
211 std::uint8_t localScvIdxWithCriticalPc;
212 Scalar criticalPc;
213 EventType event = EventType::none;
214
215 operator bool() const
216 { return event != EventType::none; }
217 };
218
219 // Block non-wetting phase flux out of the outlet
220 static const auto blockNonwettingPhase = getParamFromGroup<std::vector<int>>(problem_.paramGroup(), "InvasionState.BlockNonwettingPhaseAtThroatLabel", std::vector<int>{Labels::outlet});
221 if (!blockNonwettingPhase.empty() && std::find(blockNonwettingPhase.begin(), blockNonwettingPhase.end(), gridGeometry.throatLabel(eIdx)) != blockNonwettingPhase.end())
222 {
223 invadedCurrentIteration_[eIdx] = false;
224 return Result{}; // nothing happened
225 }
226
227 //Determine whether throat gets invaded or snap-off occurs
228 const std::array<Scalar, 2> pc = { elemVolVars[0].capillaryPressure(), elemVolVars[1].capillaryPressure() };
229 const auto pcMax = std::max_element(pc.begin(), pc.end());
230 const Scalar pcEntry = fluxVarsCache.pcEntry();
231 const Scalar pcSnapoff = fluxVarsCache.pcSnapoff();
232
233 // check if there is a user-specified global capillary pressure which needs to be obeyed
234 if (maybeRestrictToGlobalCapillaryPressure_(pcEntry))
235 {
236 if (*pcMax > pcEntry)
237 {
238 std::cout << "Throat " << eIdx << " would have been invaded by pc of " << *pcMax << "but a global capillary pressure restricion was set in the problem.";
239 std::cout << ". pcEntry: " << spatialParams.pcEntry(element, elemVolVars) << std::endl;
240 }
241
242 invadedCurrentIteration_[eIdx] = false;
243 return Result{}; //nothing happened
244 }
245
246 if (*pcMax > pcEntry)
247 invadedAfterSwitch = true;
248 else if (*pcMax <= pcSnapoff)
249 invadedAfterSwitch = false;
250
251 invadedCurrentIteration_[eIdx] = invadedAfterSwitch;
252
253 if (invadedBeforeSwitch == invadedAfterSwitch)
254 return Result{}; // nothing happened
255 else
256 {
257 Result result;
258 result.localScvIdxWithCriticalPc = std::distance(pc.begin(), pcMax);
259 result.criticalPc = *pcMax;
260 result.event = !invadedBeforeSwitch && invadedAfterSwitch ? EventType::invasion : EventType::snapOff;
261
262 if (verbose_)
263 {
264 const auto wPhaseIdx = spatialParams.template wettingPhase<typename ElementVolumeVariables::VolumeVariables::FluidSystem>(element, elemVolVars);
265 const std::array sw = { elemVolVars[0].saturation(wPhaseIdx), elemVolVars[1].saturation(wPhaseIdx) };
266 const auto vIdx = gridGeometry.gridView().indexSet().subIndex(element, result.localScvIdxWithCriticalPc, 1);
267 if (result.event == EventType::invasion)
268 {
269 std::cout << "Throat " << eIdx << " was invaded from pore " << vIdx << " :";
270 std::cout << " pc: " << *pcMax;
271 std::cout << ", pcEntry: " << spatialParams.pcEntry(element, elemVolVars);
272 std::cout << ", sw: " << sw[result.localScvIdxWithCriticalPc] << std::endl;
273 }
274 else
275 {
276 std::cout << "Snap-off occured at throat " << eIdx << " from pore " << vIdx << " :";
277 std::cout << " pc: " << *pcMax;
278 std::cout << ", pcSnapoff: " << spatialParams.pcSnapoff(element, elemVolVars);
279 std::cout << ", sw: " << sw[result.localScvIdxWithCriticalPc] << std::endl;
280 }
281 }
282
283 return result;
284 }
285 }
286
289 template<class Scalar>
290 bool maybeRestrictToGlobalCapillaryPressure_(const Scalar pcEntry) const
291 {
292 if constexpr (hasGlobalCapillaryPressure<Problem>())
293 return restrictToGlobalCapillaryPressure_ && (pcEntry > problem_.globalCapillaryPressure());
294 else
295 return false;
296 }
297
298 std::vector<bool> invadedCurrentIteration_;
299 std::vector<bool> invadedPreviousTimeStep_;
300 bool hasChangedInCurrentIteration_ = false;
301 bool hasChangedComparedToPreviousTimestep_ = false;
302 std::size_t numThroatsInvaded_;
303 bool verbose_;
304 bool restrictToGlobalCapillaryPressure_;
305
306 const Problem& problem_;
307};
308
309} // end namespace Dumux::PoreNetwork
310
311#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Defines labels for pores and throats.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
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:33
std::string capillaryPressure() noexcept
I/O name of capillary pressure.
Definition: name.hh:135
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
This class updates the invasion state for the two-phase PNM.
Definition: invasionstate.hh:41
void checkIfCapillaryPressureIsCloseToEntryPressure(const SolutionVector &sol, const GridVolumeVariables &gridVolVars, const GridFluxVariablesCache &gridFluxVarsCache) const
Definition: invasionstate.hh:151
bool hasChangedInCurrentIteration() const
Return whether an invasion or snap-off occurred anywhere during the current Newton iteration.
Definition: invasionstate.hh:140
bool invaded(const Element &element) const
Return whether a given throat is invaded or not.
Definition: invasionstate.hh:84
bool hasChanged() const
Return whether an invasion or snap-off occurred anywhere. Can be used, e.g., for output file writing ...
Definition: invasionstate.hh:136
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:96
TwoPInvasionState(const Problem &problem)
Definition: invasionstate.hh:55
void advance()
This is called after the Newton method has successfully finished one time step.
Definition: invasionstate.hh:144
std::size_t numThroatsInvaded() const
Return the number of currently invaded throats.
Definition: invasionstate.hh:91
void reset()
Restore the old invasion state after a Newton iteration has failed.
Definition: invasionstate.hh:129
static constexpr int outlet
Label for pores/throats on an outlet.
Definition: labels.hh:41