3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
localrulesforplatonicbody.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 *****************************************************************************/
25#ifndef DUMUX_PNM_2P_LOCAL_RULES_FOR_PLATONIC_BODY_HH
26#define DUMUX_PNM_2P_LOCAL_RULES_FOR_PLATONIC_BODY_HH
27
28#include <algorithm>
29#include <cassert>
30#include <cmath>
31#include <dune/common/math.hh>
37
39
44template<class Scalar>
46{
47 PlatonicBodyParams() = default;
48
49 PlatonicBodyParams& setPoreInscribedRadius(Scalar r) { radius_ = r; return *this;}
50 PlatonicBodyParams& setPoreShape(Pore::Shape s) { shape_ = s; return *this;}
51 PlatonicBodyParams& setSurfaceTension(Scalar st) { surfaceTension_ = st; return *this;}
52
53 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
54 PlatonicBodyParams(const SpatialParams& spatialParams,
55 const Element& element,
56 const SubControlVolume& scv,
57 const ElemSol& elemSol)
58 : shape_(spatialParams.gridGeometry().poreGeometry()[scv.dofIndex()])
59 , radius_(spatialParams.poreInscribedRadius(element, scv, elemSol))
60 , surfaceTension_(spatialParams.surfaceTension(element, scv, elemSol))
61 {}
62
63 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
64 void update(const SpatialParams& spatialParams,
65 const Element& element,
66 const SubControlVolume& scv,
67 const ElemSol& elemSol)
68 {
69 const auto& gridGeometry = spatialParams.gridGeometry();
70 shape_ = gridGeometry.poreGeometry()[scv.dofIndex()];
71 radius_ = spatialParams.poreInscribedRadius(element, scv, elemSol);
72 surfaceTension_ = spatialParams.surfaceTension(element, scv, elemSol);
73 }
74
75 Pore::Shape poreShape() const { return shape_; }
76
77 Scalar poreInscribedRadius() const { return radius_; }
78
79 Scalar surfaceTension() const { return surfaceTension_; }
80
81 bool operator== (const PlatonicBodyParams& p) const
82 {
83 return Dune::FloatCmp::eq(radius_, p.radius_, 1e-6)
84 && Dune::FloatCmp::eq(surfaceTension_, p.surfaceTension_, 1e-6)
85 && shape_ == p.shape_;
86 }
87
88private:
89 Pore::Shape shape_;
90 Scalar radius_;
91 Scalar surfaceTension_;
92};
93
101template<Pore::Shape shape>
103{
104 template<class Scalar>
106
107 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
108 static auto makeParams(const SpatialParams& spatialParams,
109 const Element& element,
110 const SubControlVolume& scv,
111 const ElemSol& elemSol)
112 {
113 using Scalar = std::decay_t<decltype(spatialParams.poreInscribedRadius(element, scv, elemSol))>;
114 return Params<Scalar>(spatialParams, element, scv, elemSol);
115 }
116
123 template<class Scalar>
124 static Scalar pc(Scalar sw, const Params<Scalar>& params)
125 {
126 assert(isPlatonicBody_(params.poreShape()));
127
128 using std::clamp;
129 sw = clamp(sw, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
130
131 const Scalar poreRadius = params.poreInscribedRadius();
132 const Scalar sigma = params.surfaceTension();
133 // TODO incorporate contact angle!!!
134
135 using std::exp;
136 return 2.0*sigma / (poreRadius*(1.0 - exp(-expFactor_<Scalar>()*sw)));
137 }
138
145 template<class Scalar>
146 static Scalar sw(Scalar pc, const Params<Scalar>& params)
147 {
148 assert(isPlatonicBody_(params.poreShape()));
149
150 using std::max;
151 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
152
153 const Scalar poreRadius = params.poreInscribedRadius();
154 const Scalar sigma = params.surfaceTension();
155
156 using std::log;
157 return -1.0/expFactor_<Scalar>()* log(1.0 - 2.0*sigma/(poreRadius*pc));
158 }
159
167 template<class Scalar>
168 static Scalar dpc_dsw(Scalar sw, const Params<Scalar>& params)
169 {
170 assert(isPlatonicBody_(params.poreShape()));
171
172 using std::clamp;
173 sw = clamp(sw, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
174
175 const Scalar sigma = params.surfaceTension();
176 const Scalar poreRadius = params.poreInscribedRadius();
177 using std::exp;
178 const Scalar e = exp(expFactor_<Scalar>()*sw);
179 return -(2.0*expFactor_<Scalar>()*sigma*e) / (poreRadius*(1.0-e)*(1.0-e));
180 }
181
189 template<class Scalar>
190 static Scalar dsw_dpc(Scalar pc, const Params<Scalar>& params)
191 {
192 assert(isPlatonicBody_(params.poreShape()));
193
194 using std::max;
195 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
196
197 const Scalar sigma = params.surfaceTension();
198 const Scalar poreRadius = params.poreInscribedRadius();
199 return sigma / (expFactor_<Scalar>()*sigma*pc - 0.5*expFactor_<Scalar>()*poreRadius*pc*pc);
200 }
201
202private:
203
204 template<class Scalar>
205 static constexpr Scalar expFactor_()
206 {
207 if constexpr (shape == Pore::Shape::tetrahedron)
208 return 3.87;
209 else if constexpr (shape == Pore::Shape::cube)
210 return 6.83;
211 else if constexpr (shape == Pore::Shape::octahedron)
212 return 8.71;
213 else if constexpr (shape == Pore::Shape::dodecahedron)
214 return 22.87;
215 else if constexpr (shape == Pore::Shape::icosahedron)
216 return 24.11;
217 else
218 {
219 static_assert(AlwaysFalse<Scalar>::value, "Shape not supported");
220 return 0;
221 }
222 }
223
224 static constexpr bool isPlatonicBody_(Pore::Shape s)
225 {
226 return s == Pore::Shape::tetrahedron ||
227 s == Pore::Shape::cube ||
231 }
232};
233
234template<class Scalar, class BaseLaw>
236{
238 using BaseLawParams = typename BaseLaw::template Params<Scalar>;
239public:
240
245 {
246 linear, spline, powerLaw
247 };
248
249 template<class S>
250 struct Params
251 {
252 // export method type
261 { pcLowSw_ = pcLowSw; }
262
266 S pcLowSw() const
267 { return pcLowSw_; }
268
273 { pcHighSw_ = pcHighSw; }
274
281 S pcHighSw() const
282 { return pcHighSw_; }
283
288 { highSwRegularizationMethod_ = method; }
289
294 { return highSwRegularizationMethod_; }
295
296private:
297 S pcLowSw_ = 0.01;
298 S pcHighSw_ = 0.99;
300 };
301
302
303 template<class MaterialLaw>
304 void init(const MaterialLaw* m, const Params<Scalar>& p, const std::string& paramGroup = "")
305 {
306 initPcParameters_(m, p, paramGroup);
307 }
308
309 OptionalScalar<Scalar> pc(const Scalar sw) const
310 {
311 // make sure that the capilary pressure observes a
312 // derivative != 0 for 'illegal' saturations. This is
313 // required for example by newton solvers (if the
314 // derivative is calculated numerically) in order to get the
315 // saturation moving to the right direction if it
316 // temporarily is in an 'illegal' range.
317 if (sw < pcLowSw_)
318 return pcLowSwPcValue_() + pcDerivativeLowSw_() * (sw - pcLowSw_);
319
320 if (sw <= pcHighSw_)
321 return {}; // standard
322 else if (sw < 1.0) // regularized part below sw = 1.0
323 {
324 using std::pow;
325 if (highSwRegularizationMethod_ == HighSwRegularizationMethod::powerLaw)
326 return pcHighSwPcValue_() * pow(((1.0-sw)/(1.0-pcHighSw_)), 1.0/3.0);
327
328 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::linear)
329 {
330 const Scalar slope = -pcHighSwPcValue_() / (1.0 - pcHighSw_);
331 return pcHighSwPcValue_() + (sw - pcHighSw_) * slope;
332 }
333
334 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::spline)
335 return pcSpline_().eval(sw);
336
337 else
338 DUNE_THROW(Dune::NotImplemented, "Regularization not method not implemented");
339 }
340 else // regularized part above sw = 1.0
341 return pcDerivativeHighSwEnd_()*(sw - 1.0);
342 }
343
347 OptionalScalar<Scalar> sw(const Scalar pc) const
348 {
349 if (pc <= 0.0)
350 {
351 if (pcHighSw_ >= 1.0)
352 return 1.0; // no regularization for high sw
353 else
354 return pc/pcDerivativeHighSwEnd_() + 1.0;
355 }
356
357 // low saturation
358 if (pc > pcLowSwPcValue_())
359 return (pc - pcLowSwPcValue_())/pcDerivativeLowSw_() + pcLowSw_;
360
361 // high saturation
362 else if (pc <= pcHighSwPcValue_())
363 {
364 if (highSwRegularizationMethod_ == HighSwRegularizationMethod::powerLaw)
365 {
366 // invert power law
367 using Dune::power;
368 return power(pc/pcHighSwPcValue_(), 3) * (pcHighSw_ - 1.0) + 1.0;
369 }
370
371 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::linear)
372 return pc/pcDerivativeHighSwEnd_() + 1.0;
373
374 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::spline)
375 // invert spline
376 return pcSpline_().intersectInterval(pcHighSw_, 1.0, 0.0, 0.0, 0.0, pc);
377 }
378
379 // no regularization
380 return {};
381 }
382
387 {
388 if (sw <= pcLowSw_)
389 return pcDerivativeLowSw_();
390
391 else if (sw >= 1.0)
392 return pcDerivativeHighSwEnd_();
393
394 else if (sw > pcHighSw_)
395 {
396 if (highSwRegularizationMethod_ == HighSwRegularizationMethod::powerLaw)
397 {
398 using std::pow;
399 return pcHighSwPcValue_()/3.0 * 1.0 /((pcHighSw_-1.0) * pow((sw-1.0)/(pcHighSw_-1.0), 2.0/3.0));
400 }
401
402 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::linear)
403 return pcDerivativeHighSwEnd_();
404
405 else
406 return pcSpline_().evalDerivative(sw);
407 }
408
409 else
410 return {}; // no regularization
411 }
412
417 {
418 if (pc <= 0.0)
419 {
420 if (pcHighSw_ >= 1.0)
421 return 0.0;
422 else
423 return 1.0/pcDerivativeHighSwEnd_();
424 }
425
426 // derivative of the inverse of the function is one over derivative of the function
427 else if (pc <= pcHighSwPcValue_())
428 {
429 if (highSwRegularizationMethod_ == HighSwRegularizationMethod::powerLaw)
430 {
431 using Dune::power;
432 return (3.0*pcHighSw_ - 3.0) * power(pc, 2) / power(pcHighSwPcValue_(), 3);
433 }
434
435 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::linear)
436 return 1.0/pcDerivativeHighSwEnd_();
437
438 else
439 return 1.0/pcSpline_().evalDerivative(pcSpline_().intersectInterval(pcHighSw_, 1.0, 0.0, 0.0, 0.0, pc));
440 }
441
442 else if (pc >= pcLowSwPcValue_())
443 return 1.0/pcDerivativeLowSw_();
444
445 else
446 return {}; // no regularization
447 }
448
449 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
450 void updateParams(const SpatialParams& spatialParams,
451 const Element& element,
452 const SubControlVolume& scv,
453 const ElemSol& elemSol)
454 {}
455
456private:
457
458 template<class MaterialLaw>
459 void initPcParameters_(const MaterialLaw* m, const Params<Scalar>& params, const std::string& paramGroup)
460 {
461 highSwRegularizationMethod_ = params.highSwRegularizationMethod();
462
463 // maybe overwrite method using input
464 static const bool hasMethod = hasParamInGroup(paramGroup, "HighSwRegularizationMethod");
465 if (hasMethod)
466 {
467 // set method only once
468 static const HighSwRegularizationMethod inputmethod = [&]
469 {
470 static const auto input = getParamFromGroup<std::string>(paramGroup, "HighSwRegularizationMethod");
471 if (input == "Linear")
473 else if (input == "Spline")
475 else if (input == "PowerLaw")
477 else
478 DUNE_THROW(Dune::InvalidStateException, input << " is not a valid regularization method");
479 }();
480
481 highSwRegularizationMethod_ = inputmethod;
482 }
483
484 static const bool splineZeroSlope = getParamFromGroup<bool>(paramGroup, "HighSwSplineZeroSlope", true);
485 highSwSplineZeroSlope_ = splineZeroSlope;
486
487 // print name only once
488 [[maybe_unused]] static const bool printName = [&]
489 {
490 std::string name;
491 if (highSwRegularizationMethod_ == HighSwRegularizationMethod::linear)
492 name = "Linear";
493 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::powerLaw)
494 name = "PowerLaw";
495 else if (highSwRegularizationMethod_ == HighSwRegularizationMethod::spline)
496 name = "Spline";
497 else
498 DUNE_THROW(Dune::NotImplemented, "Regularization not method not implemented");
499
500 std::cout << "\n*****\nUsing " << name << " as regularization method for high Sw\n*****" << std::endl;
501 return true;
502 }();
503
504 using std::isnan;
505 static const auto pcLowSwInput = getParamFromGroup<Scalar>(paramGroup, "RegularizationLowSw", params.pcLowSw());
506 pcLowSw_ = !isnan(pcLowSwInput) ? pcLowSwInput : params.pcLowSw();
507
508 static const auto pcHighSwInput = getParamFromGroup<Scalar>(paramGroup, "RegularizationHighSw", std::numeric_limits<Scalar>::quiet_NaN());
509 pcHighSw_ = !isnan(pcHighSwInput) ? pcHighSwInput : params.pcHighSw();
510
511 baseLawParamsPtr_ = &m->basicParams();
512
513 static const auto highSwFixedSlopeInput = getParamFromGroup<Scalar>(paramGroup, "RegularizationHighSwFixedSlope", std::numeric_limits<Scalar>::quiet_NaN());
514 highSwFixedSlope_ = highSwFixedSlopeInput;
515
516 // Note: we do not pre-calculate all end-point values (e.g., pcLowSwPcValue_ and pcDerivativeLowSw_)
517 // as done in, e.g., VanGenuchten. This is because this law will generally be instantiated only as
518 // a temporary object since all pore bodies generally have different parameters.
519 // We still cache above-mentioned values, but only if they are actually needed (see below).
520 }
521
522 // the capillary pressure for the lower saturation threshold
523 Scalar pcLowSwPcValue_() const
524 {
525 // calculated value within first function call, used cached value otherwise
526 if (!optionalPcLowSwPcValue_)
527 optionalPcLowSwPcValue_ = BaseLaw::pc(pcLowSw_, *baseLawParamsPtr_);
528 return optionalPcLowSwPcValue_.value();
529 }
530
531 // dpc_dsw for the lower saturation threshold
532 Scalar pcDerivativeLowSw_() const
533 {
534 if (!optionalPcDerivativeLowSw_)
535 optionalPcDerivativeLowSw_ = BaseLaw::dpc_dsw(pcLowSw_, *baseLawParamsPtr_);
536 return optionalPcDerivativeLowSw_.value();
537 }
538
539 // the capillary pressure for the upper saturation threshold
540 Scalar pcHighSwPcValue_() const
541 {
542 if (!optionalPcHighSwPcValue_)
543 optionalPcHighSwPcValue_ = BaseLaw::pc(pcHighSw_, *baseLawParamsPtr_);
544 return optionalPcHighSwPcValue_.value();
545 }
546
547 // dpc_dsw at Sw = 1.0
548 Scalar pcDerivativeHighSwEnd_() const
549 {
550 using std::isnan;
551 if (!isnan(highSwFixedSlope_))
552 return highSwFixedSlope_;
553 else
554 return (0.0 - pcHighSwPcValue_()) / (1.0 - pcHighSw_);
555 }
556
557 const Spline<Scalar>& pcSpline_() const
558 {
559 if (!optionalPcSpline_)
560 {
561 const auto slopes = highSwSplineZeroSlope_ ? std::array{0.0, 0.0}
562 : std::array{BaseLaw::dpc_dsw(pcHighSw_, *baseLawParamsPtr_), pcDerivativeHighSwEnd_()};
563
564 optionalPcSpline_ = Spline<Scalar>(pcHighSw_, 1.0, // x0, x1
565 pcHighSwPcValue_(), 0, // y0, y1
566 slopes[0], slopes[1]); // m0, m1
567 }
568 return optionalPcSpline_.value();
569 }
570
571 Scalar pcLowSw_, pcHighSw_;
572 mutable OptionalScalar<Scalar> optionalPcLowSwPcValue_, optionalPcDerivativeLowSw_;
573 mutable OptionalScalar<Scalar> optionalPcHighSwPcValue_;
574 HighSwRegularizationMethod highSwRegularizationMethod_;
575 const BaseLawParams* baseLawParamsPtr_;
576 mutable std::optional<Spline<Scalar>> optionalPcSpline_;
577 bool highSwSplineZeroSlope_;
578 Scalar highSwFixedSlope_;
579};
580
585template<Pore::Shape shape, typename Scalar = double>
590
595template<Pore::Shape shape, typename Scalar = double>
597
598} // end namespace Dumux::FluidMatrix
599
600#endif
Provides 3rd order polynomial splines.
A wrapper that can either contain a valid Scalar or NaN.
Base classes for standard pore-local pc-Sw curves.
This file contains functions related to calculate pore-body properties.
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:374
Definition: localrulesforplatonicbody.hh:38
Shape
Collection of different pore-body shapes.
Definition: poreproperties.hh:35
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:176
T value() const
Definition: optionalscalar.hh:48
Template which always yields a false value.
Definition: typetraits.hh:36
A tag to turn off regularization and it's overhead.
Definition: noregularization.hh:34
The parameter type.
Definition: localrulesforplatonicbody.hh:46
PlatonicBodyParams & setPoreShape(Pore::Shape s)
Definition: localrulesforplatonicbody.hh:50
PlatonicBodyParams & setSurfaceTension(Scalar st)
Definition: localrulesforplatonicbody.hh:51
PlatonicBodyParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: localrulesforplatonicbody.hh:54
Pore::Shape poreShape() const
Definition: localrulesforplatonicbody.hh:75
bool operator==(const PlatonicBodyParams &p) const
Definition: localrulesforplatonicbody.hh:81
Scalar surfaceTension() const
Definition: localrulesforplatonicbody.hh:79
PlatonicBodyParams & setPoreInscribedRadius(Scalar r)
Definition: localrulesforplatonicbody.hh:49
Scalar poreInscribedRadius() const
Definition: localrulesforplatonicbody.hh:77
void update(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: localrulesforplatonicbody.hh:64
Implementation of the simplified pore-local capillary pressure-saturation curve for platonic bodies (...
Definition: localrulesforplatonicbody.hh:103
static Scalar dpc_dsw(Scalar sw, const Params< Scalar > &params)
The partial derivative of the capillary pressure w.r.t. the wetting phase saturation.
Definition: localrulesforplatonicbody.hh:168
static Scalar sw(Scalar pc, const Params< Scalar > &params)
The wetting-phase saturation of a pore body.
Definition: localrulesforplatonicbody.hh:146
static Scalar pc(Scalar sw, const Params< Scalar > &params)
The capillary pressure-saturation curve.
Definition: localrulesforplatonicbody.hh:124
static Scalar dsw_dpc(Scalar pc, const Params< Scalar > &params)
The partial derivative of the wetting phase saturation w.r.t. the capillary pressure.
Definition: localrulesforplatonicbody.hh:190
static auto makeParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: localrulesforplatonicbody.hh:108
OptionalScalar< Scalar > dsw_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: localrulesforplatonicbody.hh:416
OptionalScalar< Scalar > sw(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: localrulesforplatonicbody.hh:347
void init(const MaterialLaw *m, const Params< Scalar > &p, const std::string &paramGroup="")
Definition: localrulesforplatonicbody.hh:304
HighSwRegularizationMethod
The available options for regularizing the pc-SW curve at high wetting-phase saturations.
Definition: localrulesforplatonicbody.hh:245
OptionalScalar< Scalar > pc(const Scalar sw) const
Definition: localrulesforplatonicbody.hh:309
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: localrulesforplatonicbody.hh:450
OptionalScalar< Scalar > dpc_dsw(const Scalar sw) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: localrulesforplatonicbody.hh:386
void setpcHighSw(S pcHighSw)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: localrulesforplatonicbody.hh:272
void setpcLowSw(S pcLowSw)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: localrulesforplatonicbody.hh:260
S pcLowSw() const
Threshold saturation below which the capillary pressure is regularized.
Definition: localrulesforplatonicbody.hh:266
S pcHighSw() const
Threshold saturation above which the capillary pressure is regularized.
Definition: localrulesforplatonicbody.hh:281
void setHighSwRegularizationMethod(HighSwRegularizationMethod method)
Set the regularization method for high saturations.
Definition: localrulesforplatonicbody.hh:287
HighSwRegularizationMethod highSwRegularizationMethod() const
Return the regularization method for high saturations.
Definition: localrulesforplatonicbody.hh:293
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:42