version 3.10-dev
thresholdcapillarypressures.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//
7
14#ifndef DUMUX_PNM_THRESHOLD_CAPILLARY_PRESSURES_HH
15#define DUMUX_PNM_THRESHOLD_CAPILLARY_PRESSURES_HH
16
17#include <cmath>
18
19#include <dune/common/exceptions.hh>
21
23
30template <class Scalar>
31Scalar pcSnapoffRegularShape(const Scalar surfaceTension,
32 const Scalar contactAngle,
33 const Scalar inscribedRadius,
34 const Throat::Shape shape) noexcept
35{
36 using std::cos;
37 using std::sin;
38
39 const Scalar cornerHalfAngle = Throat::cornerHalfAngles<Scalar>(shape)[0];
40
41 // check if snap-off is possible with such contact angle and corner half-angle
42 // if the criteriun is not fulfilled, return the lowest negative value to make sure that snap-off will not happen
43 if (cornerHalfAngle + contactAngle >= 0.5*M_PI)
44 return std::numeric_limits<Scalar>::lowest();
45
46 const Scalar cosContactAngle = cos(contactAngle);
47 const Scalar sinContactAngle = sin(contactAngle);
48
49 // Get tangent of corner half-angle
50 const Scalar tanCornerHalfAngle = [&]{ switch (shape){ // optimized tan(corner half-angle) for certain shapes
51 case Throat::Shape::equilateralTriangle: return 0.577;
52 case Throat::Shape::square: return 1.0;
53 default: using std::tan; return tan(cornerHalfAngle);
54 }}();
55
56 return surfaceTension / inscribedRadius * (cosContactAngle - sinContactAngle * tanCornerHalfAngle);
57}
58
59
64template <class Scalar>
65Scalar pcSnapoff(const Scalar surfaceTension,
66 const Scalar contactAngle,
67 const Scalar inscribedRadius,
68 const Throat::Shape shape)
69{
70 if (Throat::isRegularShape(shape)) // call the pc snap-off calculated for regular shapes(same corner angles and side length)
71 return pcSnapoffRegularShape(surfaceTension, contactAngle, inscribedRadius, shape);
72 else
73 DUNE_THROW(Dune::NotImplemented, "Pc snap-off is not implemented for this irregular shape");
74}
75
81template<class Scalar>
82Scalar pcEntry(const Scalar surfaceTension,
83 const Scalar contactAngle,
84 const Scalar inscribedRadius,
85 const Scalar shapeFactor) noexcept
86{
87 using std::sin;
88 using std::cos;
89 using std::sqrt;
90
91 const Scalar cosContactAngle = cos(contactAngle);
92 const Scalar sinContactAngle = sin(contactAngle);
93
94 const Scalar D = M_PI - 3*contactAngle + 3*sinContactAngle*cosContactAngle - (cosContactAngle*cosContactAngle) / (4*shapeFactor);
95 const Scalar F = (1 + sqrt(1 + 4*shapeFactor*D / (cosContactAngle*cosContactAngle))) / (1 + 2*sqrt(M_PI*shapeFactor));
96 return surfaceTension / inscribedRadius * cosContactAngle * (1 + 2*sqrt(M_PI*shapeFactor)) * F;
97}
98
99} // namespace Dumux::PoreNetwork::ThresholdCapillaryPressures
100
101#endif
Definition: thresholdcapillarypressures.hh:22
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
Scalar pcSnapoffRegularShape(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius, const Throat::Shape shape) noexcept
The snap-off capillary pressure of a pore throat with regular cross section shape (with same corner a...
Definition: thresholdcapillarypressures.hh:31
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:165
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:26
Scalar shapeFactor(Shape shape, const Scalar inscribedRadius)
Returns the value of the shape factor for a given shape.
Definition: throatproperties.hh:150
bool isRegularShape(Shape shape)
Returns if a shape is regular.
Definition: throatproperties.hh:182
This file contains functions related to calculate pore-throat properties.