Line data Source code
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 "MooseError.h" 14 : #include "libmesh/utility.h" 15 : 16 : /** 17 : * Broadbridge-White version of relative permeability, 18 : * and effective saturation as a function of capillary pressure. 19 : * Valid for Kn small. 20 : * P Broadbridge, I White ``Constant rate rainfall 21 : * infiltration: A versatile nonlinear model, 1 Analytical solution''. 22 : * Water Resources Research 24 (1988) 145--154. 23 : */ 24 : namespace PorousFlowBroadbridgeWhite 25 : { 26 : /** 27 : * Provides the Lambert W function, which satisfies 28 : * W(z) * exp(W(z)) = z 29 : * @param z z value, which must be positive 30 : * @return W 31 : */ 32 : Real LambertW(Real z); 33 : 34 : /** 35 : * Effective saturation as a function of capillary pressure 36 : * If pc>=0 this will yield 1, otherwise it will yield <1 37 : * @param pc capillary pressure 38 : * @param c BW's C parameter 39 : * @param sn BW's S_n parameter 40 : * @param ss BW's S_s parameter 41 : * @param las BW's lambda_s parameter 42 : * @return effective saturation 43 : */ 44 : Real effectiveSaturation(Real pc, Real c, Real sn, Real ss, Real las); 45 : 46 : /** 47 : * Derivative of effective saturation wrt capillary pressure 48 : * @param pc capillary pressure 49 : * @param c BW's C parameter 50 : * @param sn BW's S_n parameter 51 : * @param ss BW's S_s parameter 52 : * @param las BW's lambda_s parameter 53 : * @return derivative of effective saturation wrt capillary pressure 54 : */ 55 : Real dEffectiveSaturation(Real pc, Real c, Real sn, Real ss, Real las); 56 : 57 : /** 58 : * Second derivative of effective saturation wrt capillary pressure 59 : * @param pc capillary pressure 60 : * @param c BW's C parameter 61 : * @param sn BW's S_n parameter 62 : * @param ss BW's S_s parameter 63 : * @param las BW's lambda_s parameter 64 : * @return second derivative of effective saturation wrt capillary pressure 65 : */ 66 : Real d2EffectiveSaturation(Real pc, Real c, Real sn, Real ss, Real las); 67 : 68 : /** 69 : * Relative permeability as a function of saturation 70 : * @param s saturation 71 : * @param c BW's C parameter 72 : * @param sn BW's S_n parameter 73 : * @param ss BW's S_s parameter 74 : * @param kn BW's K_n parameter 75 : * @param ks BW's K_s parameter 76 : * @return relative permeability 77 : */ 78 : template <typename T> 79 : T 80 2743111 : relativePermeability(const T & s, Real c, Real sn, Real ss, Real kn, Real ks) 81 : { 82 2743111 : if (MetaPhysicL::raw_value(s) <= sn) 83 0 : return kn; 84 : 85 2743110 : if (MetaPhysicL::raw_value(s) >= ss) 86 0 : return ks; 87 : 88 2718138 : const T coef = (ks - kn) * (c - 1.0); 89 2718139 : const T th = (s - sn) / (ss - sn); 90 2718138 : const T krel = kn + coef * Utility::pow<2>(th) / (c - th); 91 2718138 : return krel; 92 : } 93 : 94 : /** 95 : * Derivative of relative permeability with respect to saturation 96 : * @param s saturation 97 : * @param c BW's C parameter 98 : * @param sn BW's S_n parameter 99 : * @param ss BW's S_s parameter 100 : * @param kn BW's K_n parameter 101 : * @param ks BW's K_s parameter 102 : * @return derivative of relative permeability wrt saturation 103 : */ 104 : Real dRelativePermeability(Real s, Real c, Real sn, Real ss, Real kn, Real ks); 105 : 106 : /** 107 : * Second derivative of relative permeability with respect to saturation 108 : * @param s saturation 109 : * @param c BW's C parameter 110 : * @param sn BW's S_n parameter 111 : * @param ss BW's S_s parameter 112 : * @param kn BW's K_n parameter 113 : * @param ks BW's K_s parameter 114 : * @return second derivative of relative permeability wrt saturation 115 : */ 116 : Real d2RelativePermeability(Real s, Real c, Real sn, Real ss, Real kn, Real ks); 117 : }