https://mooseframework.inl.gov
PorousFlowFLACrelpermTest.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 
12 #include "PorousFlowFLACrelperm.h"
13 
14 const double eps = 1.0E-8;
15 
16 TEST(PorousFlowFLACrelpermTest, relperm)
17 {
18  EXPECT_NEAR(1.0, PorousFlowFLACrelperm::relativePermeability(1.0E30, 2.7), 1.0E-5);
19  EXPECT_NEAR(0.0, PorousFlowFLACrelperm::relativePermeability(-1.0, 2.7), 1.0E-5);
20  EXPECT_NEAR(0.111976072427008, PorousFlowFLACrelperm::relativePermeability(0.3, 2.7), 1.0E-5);
21  EXPECT_NEAR(0.208087549965399, PorousFlowFLACrelperm::relativePermeability(0.8, 12.7), 1.0E-5);
22 }
23 
24 TEST(PorousFlowFLACrelpermTest, drelperm)
25 {
26  Real fd;
27  EXPECT_NEAR(0.0, PorousFlowFLACrelperm::dRelativePermeability(1.0E30, 2.7), 1.0E-5);
28  EXPECT_NEAR(0.0, PorousFlowFLACrelperm::dRelativePermeability(-1.0, 2.7), 1.0E-5);
31  eps;
32  EXPECT_NEAR(fd, PorousFlowFLACrelperm::dRelativePermeability(0.3, 2.7), 1.0E-5);
35  eps;
36  EXPECT_NEAR(fd, PorousFlowFLACrelperm::dRelativePermeability(0.8, 0.65), 1.0E-5);
37 }
38 
39 TEST(PorousFlowFLACrelpermTest, d2relperm)
40 {
41  Real fd;
42  EXPECT_NEAR(0.0, PorousFlowFLACrelperm::d2RelativePermeability(1.0E30, 2.7), 1.0E-5);
43  EXPECT_NEAR(0.0, PorousFlowFLACrelperm::d2RelativePermeability(-1.0, 2.7), 1.0E-5);
46  eps;
47  EXPECT_NEAR(fd, PorousFlowFLACrelperm::d2RelativePermeability(0.3, 2.7), 1.0E-5);
50  eps;
51  EXPECT_NEAR(fd, PorousFlowFLACrelperm::d2RelativePermeability(0.8, 0.65), 1.0E-5);
52 }
53 
54 TEST(PorousFlowFLACrelpermTest, adrelperm)
55 {
56  ADReal sat = 0.3;
57  sat.derivatives() = {};
58  Moose::derivInsert(sat.derivatives(), 0, 1.0);
59 
60  const auto adrelperm = PorousFlowFLACrelperm::relativePermeability(sat, 2.7);
61 
62  const auto relperm = PorousFlowFLACrelperm::relativePermeability(sat.value(), 2.7);
63  const auto drelperm = PorousFlowFLACrelperm::dRelativePermeability(sat.value(), 2.7);
64 
65  EXPECT_NEAR(adrelperm.value(), relperm, 1.0E-5);
66  EXPECT_NEAR(adrelperm.derivatives()[0], drelperm, 1.0E-5);
67 }
Real d2RelativePermeability(Real seff, Real m)
Second derivative of relative permeability with respect to effective saturation.
T relativePermeability(const T &seff, Real m)
Relative permeability as a function of effective saturation.
DualNumber< Real, DNDerivativeType, true > ADReal
Real dRelativePermeability(Real seff, Real m)
Derivative of relative permeability with respect to effective saturation.
TEST(PorousFlowFLACrelpermTest, relperm)
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)
const double eps