www.mooseframework.org
PorousFlowRelativePermeabilityBC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 
11 #include "PorousFlowBrooksCorey.h"
12 
15 
16 template <bool is_ad>
19 {
21  params.addRequiredParam<Real>("lambda", "The Brooks-Corey exponent of the phase");
22  params.addParam<bool>("nw_phase", false, "Set true if this is the non-wetting phase");
23  params.addClassDescription("Brooks-Corey relative permeability");
24  return params;
25 }
26 
27 template <bool is_ad>
29  const InputParameters & parameters)
30  : PorousFlowRelativePermeabilityBaseTempl<is_ad>(parameters),
31  _lambda(this->template getParam<Real>("lambda")),
32  _is_nonwetting(this->template getParam<bool>("nw_phase"))
33 {
34 }
35 
36 template <bool is_ad>
39 {
40  if (_is_nonwetting)
42  else
44 }
45 
46 template <bool is_ad>
47 Real
49 {
50  if (_is_nonwetting)
52  else
54 }
55 
Real dRelativePermeabilityNW(Real seff, Real lambda)
Derivative of relative permeability of the non-wetting phase wrt to effective saturation.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
T relativePermeabilityW(const T &seff, Real lambda)
Relative permeability of the wetting phase as a function of effective saturation. ...
virtual GenericReal< is_ad > relativePermeability(GenericReal< is_ad > seff) const override
Relative permeability equation (must be overriden in derived class)
T relativePermeabilityNW(const T &seff, Real lambda)
Relative permeability of the non-wetting phase as a function of effective saturation.
virtual Real dRelativePermeability(Real seff) const override
Derivative of relative permeability with respect to effective saturation.
void addRequiredParam(const std::string &name, const std::string &doc_string)
registerMooseObject("PorousFlowApp", PorousFlowRelativePermeabilityBC)
Base class for PorousFlow relative permeability materials.
Real dRelativePermeabilityW(Real seff, Real lambda)
Derivative of relative permeability of the wetting phase wrt to effective saturation.
PorousFlowRelativePermeabilityBCTempl(const InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Material to calculate Brooks-Corey relative permeability of an arbitrary phase given the effective sa...
void addClassDescription(const std::string &doc_string)
typename Moose::GenericType< Real, is_ad > GenericReal