www.mooseframework.org
Functions
PorousFlowBroadbridgeWhite Namespace Reference

Broadbridge-White version of relative permeability, and effective saturation as a function of capillary pressure. More...

Functions

Real LambertW (Real z)
 Provides the Lambert W function, which satisfies W(z) * exp(W(z)) = z. More...
 
Real effectiveSaturation (Real pc, Real c, Real sn, Real ss, Real las)
 Effective saturation as a function of capillary pressure If pc>=0 this will yield 1, otherwise it will yield <1. More...
 
Real dEffectiveSaturation (Real pc, Real c, Real sn, Real ss, Real las)
 Derivative of effective saturation wrt capillary pressure. More...
 
Real d2EffectiveSaturation (Real pc, Real c, Real sn, Real ss, Real las)
 Second derivative of effective saturation wrt capillary pressure. More...
 
Real relativePermeability (Real s, Real c, Real sn, Real ss, Real kn, Real ks)
 Relative permeability as a function of saturation. More...
 
Real dRelativePermeability (Real s, Real c, Real sn, Real ss, Real kn, Real ks)
 Derivative of relative permeability with respect to saturation. More...
 
Real d2RelativePermeability (Real s, Real c, Real sn, Real ss, Real kn, Real ks)
 Second derivative of relative permeability with respect to saturation. More...
 

Detailed Description

Broadbridge-White version of relative permeability, and effective saturation as a function of capillary pressure.

Valid for Kn small. P Broadbridge, I White `‘Constant rate rainfall infiltration: A versatile nonlinear model, 1 Analytical solution’'. Water Resources Research 24 (1988) 145–154.

Function Documentation

◆ d2EffectiveSaturation()

Real PorousFlowBroadbridgeWhite::d2EffectiveSaturation ( Real  pc,
Real  c,
Real  sn,
Real  ss,
Real  las 
)

Second derivative of effective saturation wrt capillary pressure.

Parameters
pccapillary pressure
cBW's C parameter
snBW's S_n parameter
ssBW's S_s parameter
lasBW's lambda_s parameter
Returns
second derivative of effective saturation wrt capillary pressure

Definition at line 114 of file PorousFlowBroadbridgeWhite.C.

115 {
116  if (pressure >= 0.0)
117  return 0.0;
118  const Real x = (c - 1.0) * std::exp(c - 1.0 - c * pressure / las);
119  const Real lamw = LambertW(x);
120  return -(ss - sn) * Utility::pow<3>(c) / Utility::pow<2>(las) * lamw * (1.0 - 2.0 * lamw) /
121  Utility::pow<5>(1.0 + lamw);
122 }

Referenced by PorousFlowCapillaryPressureBW::d2EffectiveSaturation().

◆ d2RelativePermeability()

Real PorousFlowBroadbridgeWhite::d2RelativePermeability ( Real  s,
Real  c,
Real  sn,
Real  ss,
Real  kn,
Real  ks 
)

Second derivative of relative permeability with respect to saturation.

Parameters
ssaturation
cBW's C parameter
snBW's S_n parameter
ssBW's S_s parameter
knBW's K_n parameter
ksBW's K_s parameter
Returns
second derivative of relative permeability wrt saturation

Definition at line 155 of file PorousFlowBroadbridgeWhite.C.

156 {
157  if (s <= sn)
158  return 0.0;
159 
160  if (s >= ss)
161  return 0.0;
162 
163  const Real coef = (ks - kn) * (c - 1.0);
164  const Real th = (s - sn) / (ss - sn);
165  const Real krelpp = coef * (2.0 / (c - th) + 4.0 * th / Utility::pow<2>(c - th) +
166  2.0 * Utility::pow<2>(th) / Utility::pow<3>(c - th));
167  return krelpp / Utility::pow<2>(ss - sn);
168 }

◆ dEffectiveSaturation()

Real PorousFlowBroadbridgeWhite::dEffectiveSaturation ( Real  pc,
Real  c,
Real  sn,
Real  ss,
Real  las 
)

Derivative of effective saturation wrt capillary pressure.

Parameters
pccapillary pressure
cBW's C parameter
snBW's S_n parameter
ssBW's S_s parameter
lasBW's lambda_s parameter
Returns
derivative of effective saturation wrt capillary pressure

Definition at line 104 of file PorousFlowBroadbridgeWhite.C.

105 {
106  if (pressure >= 0.0)
107  return 0.0;
108  const Real x = (c - 1.0) * std::exp(c - 1.0 - c * pressure / las);
109  const Real lamw = LambertW(x);
110  return (ss - sn) * Utility::pow<2>(c) / las * lamw / Utility::pow<3>(1.0 + lamw);
111 }

Referenced by PorousFlowCapillaryPressureBW::dEffectiveSaturation().

◆ dRelativePermeability()

Real PorousFlowBroadbridgeWhite::dRelativePermeability ( Real  s,
Real  c,
Real  sn,
Real  ss,
Real  kn,
Real  ks 
)

Derivative of relative permeability with respect to saturation.

Parameters
ssaturation
cBW's C parameter
snBW's S_n parameter
ssBW's S_s parameter
knBW's K_n parameter
ksBW's K_s parameter
Returns
derivative of relative permeability wrt saturation

Definition at line 140 of file PorousFlowBroadbridgeWhite.C.

141 {
142  if (s <= sn)
143  return 0.0;
144 
145  if (s >= ss)
146  return 0.0;
147 
148  const Real coef = (ks - kn) * (c - 1.0);
149  const Real th = (s - sn) / (ss - sn);
150  const Real krelp = coef * (2.0 * th / (c - th) + Utility::pow<2>(th) / Utility::pow<2>(c - th));
151  return krelp / (ss - sn);
152 }

Referenced by PorousFlowRelativePermeabilityBW::dRelativePermeability().

◆ effectiveSaturation()

Real PorousFlowBroadbridgeWhite::effectiveSaturation ( Real  pc,
Real  c,
Real  sn,
Real  ss,
Real  las 
)

Effective saturation as a function of capillary pressure If pc>=0 this will yield 1, otherwise it will yield <1.

Parameters
pccapillary pressure
cBW's C parameter
snBW's S_n parameter
ssBW's S_s parameter
lasBW's lambda_s parameter
Returns
effective saturation

Definition at line 94 of file PorousFlowBroadbridgeWhite.C.

95 {
96  if (pressure >= 0.0)
97  return 1.0;
98  const Real x = (c - 1.0) * std::exp(c - 1.0 - c * pressure / las);
99  const Real th = c / (1.0 + LambertW(x)); // use branch 0 for positive x
100  return sn + (ss - sn) * th;
101 }

Referenced by PorousFlowCapillaryPressureBW::effectiveSaturation().

◆ LambertW()

Real PorousFlowBroadbridgeWhite::LambertW ( Real  z)

Provides the Lambert W function, which satisfies W(z) * exp(W(z)) = z.

Parameters
zz value, which must be positive
Returns
W

Definition at line 16 of file PorousFlowBroadbridgeWhite.C.

17 {
18  /* Lambert W function.
19  Was ~/C/LambertW.c written K M Briggs Keith dot Briggs at bt dot com 97 May 21.
20  Revised KMB 97 Nov 20; 98 Feb 11, Nov 24, Dec 28; 99 Jan 13; 00 Feb 23; 01 Apr 09
21 
22  Computes Lambert W function, principal branch.
23  See LambertW1.c for -1 branch.
24 
25  Returned value W(z) satisfies W(z)*exp(W(z))=z
26  test data...
27  W(1)= 0.5671432904097838730
28  W(2)= 0.8526055020137254914
29  W(20)=2.2050032780240599705
30  To solve (a+b*R)*exp(-c*R)-d=0 for R, use
31  R=-(b*W(-exp(-a*c/b)/b*d*c)+a*c)/b/c
32 
33  Test:
34  gcc -DTESTW LambertW.c -o LambertW -lm && LambertW
35  Library:
36  gcc -O3 -c LambertW.c
37 
38  Modified trivially by Andy to use MOOSE things
39  */
40  mooseAssert(z > 0, "LambertW function in RichardsSeff1BWsmall called with negative argument");
41 
42  const Real eps = 4.0e-16; //, em1=0.3678794411714423215955237701614608;
43  Real p, e, t, w;
44 
45  /* Uncomment this stuff is you ever need to call with a negative argument
46  if (z < -em1)
47  mooseError("LambertW: bad argument ", z, "\n");
48 
49  if (0.0 == z)
50  return 0.0;
51  if (z < -em1+1e-4)
52  {
53  // series near -em1 in sqrt(q)
54  Real q=z+em1,r=std::sqrt(q),q2=q*q,q3=q2*q;
55  return
56  -1.0
57  +2.331643981597124203363536062168*r
58  -1.812187885639363490240191647568*q
59  +1.936631114492359755363277457668*r*q
60  -2.353551201881614516821543561516*q2
61  +3.066858901050631912893148922704*r*q2
62  -4.175335600258177138854984177460*q3
63  +5.858023729874774148815053846119*r*q3
64  -8.401032217523977370984161688514*q3*q; // error approx 1e-16
65  }
66  */
67  /* initial approx for iteration... */
68  if (z < 1.0)
69  {
70  /* series near 0 */
71  p = std::sqrt(2.0 * (2.7182818284590452353602874713526625 * z + 1.0));
72  w = -1.0 + p * (1.0 + p * (-0.333333333333333333333 + p * 0.152777777777777777777777));
73  }
74  else
75  w = std::log(z); /* asymptotic */
76  if (z > 3.0)
77  w -= std::log(w); /* useful? */
78  for (unsigned int i = 0; i < 10; ++i)
79  {
80  /* Halley iteration */
81  e = std::exp(w);
82  t = w * e - z;
83  p = w + 1.0;
84  t /= e * p - 0.5 * (p + 1.0) * t / p;
85  w -= t;
86  if (std::abs(t) < eps * (1.0 + std::abs(w)))
87  return w; /* rel-abs error */
88  }
89  /* should never get here */
90  mooseError("LambertW: No convergence at z= ", z, "\n");
91 }

Referenced by d2EffectiveSaturation(), dEffectiveSaturation(), and effectiveSaturation().

◆ relativePermeability()

Real PorousFlowBroadbridgeWhite::relativePermeability ( Real  s,
Real  c,
Real  sn,
Real  ss,
Real  kn,
Real  ks 
)

Relative permeability as a function of saturation.

Parameters
ssaturation
cBW's C parameter
snBW's S_n parameter
ssBW's S_s parameter
knBW's K_n parameter
ksBW's K_s parameter
Returns
relative permeability

Definition at line 125 of file PorousFlowBroadbridgeWhite.C.

126 {
127  if (s <= sn)
128  return kn;
129 
130  if (s >= ss)
131  return ks;
132 
133  const Real coef = (ks - kn) * (c - 1.0);
134  const Real th = (s - sn) / (ss - sn);
135  const Real krel = kn + coef * Utility::pow<2>(th) / (c - th);
136  return krel;
137 }

Referenced by PorousFlowRelativePermeabilityBW::relativePermeability().

PorousFlowBroadbridgeWhite::LambertW
Real LambertW(Real z)
Provides the Lambert W function, which satisfies W(z) * exp(W(z)) = z.
Definition: PorousFlowBroadbridgeWhite.C:16
NS::pressure
const std::string pressure
Definition: NS.h:25