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