https://mooseframework.inl.gov
PorousFlowVanGenuchten.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "MooseTypes.h"
13 #include "libmesh/utility.h"
14 
31 {
41 
50 
59 
69 Real capillaryPressure(Real seff, Real alpha, Real m, Real pc_max);
70 
80 Real dCapillaryPressure(Real seff, Real alpha, Real m, Real pc_max);
81 
91 Real d2CapillaryPressure(Real seff, Real alpha, Real m, Real pc_max);
92 
99 template <typename T>
100 T
101 relativePermeability(const T & seff, Real m)
102 {
103  if (MetaPhysicL::raw_value(seff) <= 0.0)
104  return 0.0;
105  else if (MetaPhysicL::raw_value(seff) >= 1.0)
106  return 1.0;
107 
108  using std::pow, std::sqrt;
109 
110  const T a = 1.0 - pow(seff, 1.0 / m);
111  const T b = 1.0 - pow(a, m);
112 
113  return sqrt(seff) * Utility::pow<2>(b);
114 }
115 
122 Real dRelativePermeability(Real seff, Real m);
123 
130 Real d2RelativePermeability(Real seff, Real m);
131 
138 template <typename T>
139 T
140 relativePermeabilityNW(const T & seff, Real m)
141 {
142  if (MetaPhysicL::raw_value(seff) <= 0.0)
143  return 0.0;
144  else if (MetaPhysicL::raw_value(seff) >= 1.0)
145  return 1.0;
146 
147  using std::pow, std::sqrt;
148 
149  const T a = pow(1.0 - seff, 1.0 / m);
150  const T b = pow(1.0 - a, 2.0 * m);
151 
152  return sqrt(seff) * b;
153 }
154 
161 Real dRelativePermeabilityNW(Real seff, Real m);
162 
170 Real d2RelativePermeabilityNW(Real seff, Real m);
171 
181 {
183  {
187  };
192 
195  S(0.0),
196  Pc(std::numeric_limits<Real>::max()),
197  dPc(std::numeric_limits<Real>::lowest()){};
198 
200  : strategy(strategy), S(S), Pc(Pc), dPc(dPc){};
201 };
202 
212 {
214  {
217  };
222 
225  S(1.0),
226  Pc(0.0),
227  dPc(std::numeric_limits<Real>::lowest()){};
228 
230  : strategy(strategy), S(S), Pc(Pc), dPc(dPc){};
231 };
232 
250  Real sl,
251  Real slmin,
252  Real sgrdel,
253  Real alpha,
254  Real n,
255  const LowCapillaryPressureExtension & low_ext = LowCapillaryPressureExtension(),
256  const HighCapillaryPressureExtension & high_ext = HighCapillaryPressureExtension());
257 
274  Real sl,
275  Real slmin,
276  Real sgrdel,
277  Real alpha,
278  Real n,
279  const LowCapillaryPressureExtension & low_ext = LowCapillaryPressureExtension(),
280  const HighCapillaryPressureExtension & high_ext = HighCapillaryPressureExtension());
281 
297  Real sl,
298  Real slmin,
299  Real sgrdel,
300  Real alpha,
301  Real n,
302  const LowCapillaryPressureExtension & low_ext = LowCapillaryPressureExtension(),
303  const HighCapillaryPressureExtension & high_ext = HighCapillaryPressureExtension());
304 
319 Real
320 saturationHys(Real pc,
321  Real slmin,
322  Real sgrdel,
323  Real alpha,
324  Real n,
325  const LowCapillaryPressureExtension & low_ext = LowCapillaryPressureExtension(),
326  const HighCapillaryPressureExtension & high_ext = HighCapillaryPressureExtension());
327 
341 Real
342 dsaturationHys(Real pc,
343  Real slmin,
344  Real sgrdel,
345  Real alpha,
346  Real n,
347  const LowCapillaryPressureExtension & low_ext = LowCapillaryPressureExtension(),
348  const HighCapillaryPressureExtension & high_ext = HighCapillaryPressureExtension());
349 
363 Real
364 d2saturationHys(Real pc,
365  Real slmin,
366  Real sgrdel,
367  Real alpha,
368  Real n,
369  const LowCapillaryPressureExtension & low_ext = LowCapillaryPressureExtension(),
370  const HighCapillaryPressureExtension & high_ext = HighCapillaryPressureExtension());
371 
395  Real slr,
396  Real sgrdel,
397  Real sgrmax,
398  Real sldel,
399  Real m,
400  Real upper_liquid_param,
401  Real y0,
402  Real y0p,
403  Real y1,
404  Real y1p);
405 
429  Real slr,
430  Real sgrdel,
431  Real sgrmax,
432  Real sldel,
433  Real m,
434  Real upper_liquid_param,
435  Real y0,
436  Real y0p,
437  Real y1,
438  Real y1p);
439 
460  Real slr,
461  Real sgrdel,
462  Real sgrmax,
463  Real sldel,
464  Real m,
465  Real gamma,
466  Real k_rg_max,
467  Real y0p);
468 
489  Real slr,
490  Real sgrdel,
491  Real sgrmax,
492  Real sldel,
493  Real m,
494  Real gamma,
495  Real k_rg_max,
496  Real y0p);
497 }
Real dcapillaryPressureHys(Real sl, Real slmin, Real sgrdel, Real alpha, Real n, const LowCapillaryPressureExtension &low_ext=LowCapillaryPressureExtension(), const HighCapillaryPressureExtension &high_ext=HighCapillaryPressureExtension())
Derivative of capillaryPressureHys with respect to sl.
Real d2capillaryPressureHys(Real sl, Real slmin, Real sgrdel, Real alpha, Real n, const LowCapillaryPressureExtension &low_ext=LowCapillaryPressureExtension(), const HighCapillaryPressureExtension &high_ext=HighCapillaryPressureExtension())
Second derivative of capillaryPressureHys with respect to sl.
Real drelativePermeabilityHys(Real sl, Real slr, Real sgrdel, Real sgrmax, Real sldel, Real m, Real upper_liquid_param, Real y0, Real y0p, Real y1, Real y1p)
Derivative of Hysteretic relative permeability for liquid, with respect to liquid saturation...
Real relativePermeabilityNWHys(Real sl, Real slr, Real sgrdel, Real sgrmax, Real sldel, Real m, Real gamma, Real k_rg_max, Real y0p)
Hysteretic relative permeability for gas.
Parameters associated with the extension of the hysteretic wetting capillary pressure function to hig...
T relativePermeabilityNW(const T &seff, Real m)
Relative permeability for a non-wetting phase as a function of effective saturation.
auto raw_value(const Eigen::Map< T > &in)
van Genuchten effective saturation, capillary pressure and relative permeability functions.
auto max(const L &left, const R &right)
Real capillaryPressureHys(Real sl, Real slmin, Real sgrdel, Real alpha, Real n, const LowCapillaryPressureExtension &low_ext=LowCapillaryPressureExtension(), const HighCapillaryPressureExtension &high_ext=HighCapillaryPressureExtension())
Hysteretic capillary pressure function (Eqn(1) of Doughty2007) with extensions (page5 and Fig1 of Dou...
Real drelativePermeabilityNWHys(Real sl, Real slr, Real sgrdel, Real sgrmax, Real sldel, Real m, Real gamma, Real k_rg_max, Real y0p)
Derivative of hysteretic relative permeability for gas with respect to the liquid saturation...
Parameters associated with the extension of the hysteretic capillary pressure function to low saturat...
Real d2RelativePermeabilityNW(Real seff, Real m)
Second derivative of relative permeability for a non-wetting phase with respect to effective saturati...
Real relativePermeabilityHys(Real sl, Real slr, Real sgrdel, Real sgrmax, Real sldel, Real m, Real upper_liquid_param, Real y0, Real y0p, Real y1, Real y1p)
Hysteretic relative permeability for liquid.
Real d2CapillaryPressure(Real seff, Real alpha, Real m, Real pc_max)
Second derivative of capillary pressure wrt effective saturation.
Real dRelativePermeability(Real seff, Real m)
Derivative of relative permeability with respect to effective saturation.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
T relativePermeability(const T &seff, Real m)
Relative permeability as a function of effective saturation.
Real effectiveSaturation(Real p, Real alpha, Real m)
Effective saturation as a function of porepressure.
LowCapillaryPressureExtension(const ExtensionStrategy &strategy, Real S, Real Pc, Real dPc)
HighCapillaryPressureExtension(const ExtensionStrategy &strategy, Real S, Real Pc, Real dPc)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
static const std::string alpha
Definition: NS.h:134
Real dCapillaryPressure(Real seff, Real alpha, Real m, Real pc_max)
Derivative of capillary pressure wrt effective saturation.
Real d2EffectiveSaturation(Real p, Real alpha, Real m)
Second derivative of effective saturation wrt porepressure.
Real dsaturationHys(Real pc, Real slmin, Real sgrdel, Real alpha, Real n, const LowCapillaryPressureExtension &low_ext=LowCapillaryPressureExtension(), const HighCapillaryPressureExtension &high_ext=HighCapillaryPressureExtension())
Derivative of Hysteretic saturation function with respect to pc.
Real d2RelativePermeability(Real seff, Real m)
Second derivative of relative permeability with respect to effective saturation.
Real d2saturationHys(Real pc, Real slmin, Real sgrdel, Real alpha, Real n, const LowCapillaryPressureExtension &low_ext=LowCapillaryPressureExtension(), const HighCapillaryPressureExtension &high_ext=HighCapillaryPressureExtension())
Second derivative of Hysteretic saturation function with respect to pc.
Real dEffectiveSaturation(Real p, Real alpha, Real m)
Derivative of effective saturation wrt porepressure.
Real dRelativePermeabilityNW(Real seff, Real m)
Derivative of relative permeability for a non-wetting phase with respect to effective saturation...
MooseUnits pow(const MooseUnits &, int)
Real saturationHys(Real pc, Real slmin, Real sgrdel, Real alpha, Real n, const LowCapillaryPressureExtension &low_ext=LowCapillaryPressureExtension(), const HighCapillaryPressureExtension &high_ext=HighCapillaryPressureExtension())
Hysteretic saturation function (Eqn(1) of Doughty2007) with extensions (page5 and Fig1 of Doughty2008...
Real capillaryPressure(Real seff, Real alpha, Real m, Real pc_max)
Capillary pressure as a function of effective saturation.