https://mooseframework.inl.gov
PorousFlowBroadbridgeWhiteTest.C
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 #include "gtest/gtest.h"
11 
13 
14 class PorousFlowBroadbridgeWhiteTest : public ::testing::Test
15 {
16 protected:
17  Real _ep = 1e-8;
18  Real _c = 1.5;
19  Real _sn = .1;
20  Real _ss = .95;
21  Real _las = 2.0;
22  Real _kn = 0.05;
23  Real _ks = .8;
24 };
25 
27 {
28  EXPECT_NEAR(1.0, PorousFlowBroadbridgeWhite::effectiveSaturation(50, _c, _sn, _ss, _las), 1.0E-5);
29  const Real eff = PorousFlowBroadbridgeWhite::effectiveSaturation(-1.0, _c, _sn, _ss, _las);
30  const Real th = (eff - _sn) / (_ss - _sn);
31  const Real t1 = (1.0 / _c) * std::log((_c - th) / (_c - 1.0) / th);
32  const Real t2 = (th - 1.0) / th;
33  EXPECT_NEAR(-1.0, _las * (t2 - t1), 1.0E-5);
34 }
35 
37 {
38  EXPECT_NEAR(
39  0.0, PorousFlowBroadbridgeWhite::dEffectiveSaturation(50, 0.7, 0.5, 0.6, 1.0), 1.0E-5);
40  const Real fd = (PorousFlowBroadbridgeWhite::effectiveSaturation(-1.0 + _ep, _c, _sn, _ss, _las) -
41  PorousFlowBroadbridgeWhite::effectiveSaturation(-1.0, _c, _sn, _ss, _las)) /
42  _ep;
43  EXPECT_NEAR(
44  fd, PorousFlowBroadbridgeWhite::dEffectiveSaturation(-1.0, _c, _sn, _ss, _las), 1.0E-5);
45 }
46 
48 {
49  EXPECT_NEAR(
50  0.0, PorousFlowBroadbridgeWhite::d2EffectiveSaturation(50, 0.7, 0.5, 0.6, 1.0), 1.0E-5);
51  const Real fd =
52  (PorousFlowBroadbridgeWhite::dEffectiveSaturation(-1.0 + _ep, _c, _sn, _ss, _las) -
53  PorousFlowBroadbridgeWhite::dEffectiveSaturation(-1.0, _c, _sn, _ss, _las)) /
54  _ep;
55  EXPECT_NEAR(
56  fd, PorousFlowBroadbridgeWhite::d2EffectiveSaturation(-1.0, _c, _sn, _ss, _las), 1.0E-5);
57 }
58 
60 {
61  EXPECT_NEAR(
62  _kn, PorousFlowBroadbridgeWhite::relativePermeability(0.01, _c, _sn, _ss, _kn, _ks), 1.0E-5);
63  const Real sat = 0.5;
64  const Real th = (sat - _sn) / (_ss - _sn);
65  const Real expect = _kn + (_ks - _kn) * th * th * (_c - 1.0) / (_c - th);
66  EXPECT_NEAR(expect,
67  PorousFlowBroadbridgeWhite::relativePermeability(sat, _c, _sn, _ss, _kn, _ks),
68  1.0E-5);
69  EXPECT_NEAR(
70  _ks, PorousFlowBroadbridgeWhite::relativePermeability(0.99, _c, _sn, _ss, _kn, _ks), 1.0E-5);
71 }
72 
74 {
75  EXPECT_NEAR(
76  0.0, PorousFlowBroadbridgeWhite::dRelativePermeability(0.01, _c, _sn, _ss, _kn, _ks), 1.0E-5);
77  const Real sat = 0.4;
78  const Real fd =
79  (PorousFlowBroadbridgeWhite::relativePermeability(sat + _ep, _c, _sn, _ss, _kn, _ks) -
80  PorousFlowBroadbridgeWhite::relativePermeability(sat, _c, _sn, _ss, _kn, _ks)) /
81  _ep;
82  EXPECT_NEAR(
83  fd, PorousFlowBroadbridgeWhite::dRelativePermeability(sat, _c, _sn, _ss, _kn, _ks), 1.0E-5);
84  EXPECT_NEAR(
85  0.0, PorousFlowBroadbridgeWhite::dRelativePermeability(0.99, _c, _sn, _ss, _kn, _ks), 1.0E-5);
86 }
87 
89 {
90  EXPECT_NEAR(0.0,
91  PorousFlowBroadbridgeWhite::d2RelativePermeability(0.01, _c, _sn, _ss, _kn, _ks),
92  1.0E-5);
93  const Real sat = 0.6;
94  const Real fd =
95  (PorousFlowBroadbridgeWhite::dRelativePermeability(sat + _ep, _c, _sn, _ss, _kn, _ks) -
96  PorousFlowBroadbridgeWhite::dRelativePermeability(sat, _c, _sn, _ss, _kn, _ks)) /
97  _ep;
98  EXPECT_NEAR(
99  fd, PorousFlowBroadbridgeWhite::d2RelativePermeability(sat, _c, _sn, _ss, _kn, _ks), 1.0E-5);
100  EXPECT_NEAR(0.0,
101  PorousFlowBroadbridgeWhite::d2RelativePermeability(0.99, _c, _sn, _ss, _kn, _ks),
102  1.0E-5);
103 }
104 
106 {
107  ADReal sat = 0.3;
108  sat.derivatives() = {};
109  Moose::derivInsert(sat.derivatives(), 0, 1.0);
110 
111  const auto adrelperm =
112  PorousFlowBroadbridgeWhite::relativePermeability(sat, _c, _sn, _ss, _kn, _ks);
113 
114  const auto relperm =
115  PorousFlowBroadbridgeWhite::relativePermeability(sat.value(), _c, _sn, _ss, _kn, _ks);
116  const auto drelperm =
117  PorousFlowBroadbridgeWhite::dRelativePermeability(sat.value(), _c, _sn, _ss, _kn, _ks);
118 
119  EXPECT_NEAR(adrelperm.value(), relperm, 1.0E-5);
120  EXPECT_NEAR(adrelperm.derivatives()[0], drelperm, 1.0E-5);
121 }
Real d2RelativePermeability(Real s, Real c, Real sn, Real ss, Real kn, Real ks)
Second derivative of relative permeability with respect to saturation.
Real dEffectiveSaturation(Real pc, Real c, Real sn, Real ss, Real las)
Derivative of effective saturation wrt capillary pressure.
TEST_F(PorousFlowBroadbridgeWhiteTest, sat)
T relativePermeability(const T &s, Real c, Real sn, Real ss, Real kn, Real ks)
Relative permeability as a function of saturation.
Real d2EffectiveSaturation(Real pc, Real c, Real sn, Real ss, Real las)
Second derivative of effective saturation wrt capillary pressure.
DualNumber< Real, DNDerivativeType, true > ADReal
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...
Real dRelativePermeability(Real s, Real c, Real sn, Real ss, Real kn, Real ks)
Derivative of relative permeability with respect to saturation.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)