https://mooseframework.inl.gov
PCNSFVHLLC.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 "PCNSFVHLLC.h"
11 #include "NS.h"
13 
15 
18 {
20  params.addRequiredParam<UserObjectName>(NS::fluid, "Fluid properties userobject");
21  params.set<unsigned short>("ghost_layers") = 2;
22  return params;
23 }
24 
26  : FVFluxKernel(params),
27  _fluid(UserObjectInterface::getUserObject<SinglePhaseFluidProperties>(NS::fluid)),
28  _specific_internal_energy_elem(getADMaterialProperty<Real>(NS::specific_internal_energy)),
29  _specific_internal_energy_neighbor(
30  getNeighborADMaterialProperty<Real>(NS::specific_internal_energy)),
31  _rho_et_elem(getADMaterialProperty<Real>(NS::total_energy_density)),
32  _rho_et_neighbor(getNeighborADMaterialProperty<Real>(NS::total_energy_density)),
33  _vel_elem(getADMaterialProperty<RealVectorValue>(NS::velocity)),
34  _vel_neighbor(getNeighborADMaterialProperty<RealVectorValue>(NS::velocity)),
35  _speed_elem(getADMaterialProperty<Real>(NS::speed)),
36  _speed_neighbor(getNeighborADMaterialProperty<Real>(NS::speed)),
37  _rho_elem(getADMaterialProperty<Real>(NS::density)),
38  _rho_neighbor(getNeighborADMaterialProperty<Real>(NS::density)),
39  _pressure_elem(getADMaterialProperty<Real>(NS::pressure)),
40  _pressure_neighbor(getNeighborADMaterialProperty<Real>(NS::pressure)),
41  _eps_elem(getMaterialProperty<Real>(NS::porosity)),
42  _eps_neighbor(getNeighborMaterialProperty<Real>(NS::porosity))
43 {
44 }
45 
46 std::array<ADReal, 3>
47 PCNSFVHLLC::waveSpeed(const ADReal & rho_elem,
48  const ADRealVectorValue & vel_elem,
49  const ADReal & e_elem,
50  const Real eps_elem,
51  const ADReal & rho_neighbor,
52  const ADRealVectorValue & vel_neighbor,
53  const ADReal & e_neighbor,
54  const Real eps_neighbor,
56  const ADRealVectorValue & normal)
57 {
58  const auto & rho1 = rho_elem;
59  const auto u1 = vel_elem.norm();
60  const auto q1 = normal * vel_elem;
61  const auto v1 = 1.0 / rho1;
62  const auto & e1 = e_elem;
63  const auto et1 = e1 + 0.5 * u1 * u1;
64  const auto p1 = fluid.p_from_v_e(v1, e1);
65  const auto ht1 = et1 + p1 / rho1;
66  const auto c1 = fluid.c_from_v_e(v1, e1);
67  const auto eps1 = eps_elem;
68 
69  const auto & rho2 = rho_neighbor;
70  const auto u2 = vel_neighbor.norm();
71  const auto q2 = normal * vel_neighbor;
72  const auto v2 = 1.0 / rho2;
73  const auto & e2 = e_neighbor;
74  const auto et2 = e2 + 0.5 * u2 * u2;
75  const auto p2 = fluid.p_from_v_e(v2, e2);
76  const auto ht2 = et2 + p2 / rho2;
77  const auto c2 = fluid.c_from_v_e(v2, e2);
78  const auto eps2 = eps_neighbor;
79 
80  // compute Roe-averaged variables
81  const auto sqrt_rho1 = std::sqrt(rho1);
82  const auto sqrt_rho2 = std::sqrt(rho2);
83  const auto u_roe = (sqrt_rho1 * u1 + sqrt_rho2 * u2) / (sqrt_rho1 + sqrt_rho2);
84  const auto q_roe = (sqrt_rho1 * q1 + sqrt_rho2 * q2) / (sqrt_rho1 + sqrt_rho2);
85  const auto e_roe = (sqrt_rho1 * e1 + sqrt_rho2 * e2) / (sqrt_rho1 + sqrt_rho2);
86  const auto rho_roe = std::sqrt(rho1 * rho2);
87  const auto v_roe = 1.0 / rho_roe;
88  const auto c_roe = fluid.c_from_v_e(v_roe, e_roe);
89 
90  // compute wave speeds
91  // I may want to change the estimate of these wave speeds!
92  auto SL = std::min(q1 - c1, q_roe - c_roe);
93  auto SR = std::max(q2 + c2, q_roe + c_roe);
94  auto SM = (eps2 * rho2 * q2 * (SR - q2) - eps1 * rho1 * q1 * (SL - q1) + eps1 * p1 - eps2 * p2) /
95  (eps2 * rho2 * (SR - q2) - eps1 * rho1 * (SL - q1));
96 
97  return {{std::move(SL), std::move(SM), std::move(SR)}};
98 }
99 
100 ADReal
102 {
105  auto wave_speeds = waveSpeed(_rho_elem[_qp],
106  _vel_elem[_qp],
108  _eps_elem[_qp],
113  _fluid,
114  _normal);
115  _SL = std::move(wave_speeds[0]);
116  _SM = std::move(wave_speeds[1]);
117  _SR = std::move(wave_speeds[2]);
118  if (_SL >= 0)
119  return fluxElem();
120  else if (_SR <= 0)
121  return fluxNeighbor();
122  else
123  {
124  if (_SM >= 0)
125  {
127  return fluxElem() + _SL * (f * hllcElem() - conservedVariableElem());
128  }
129  else
130  {
131  ADReal f =
134  }
135  }
136  mooseError("Should never get here");
137 }
const MaterialProperty< Real > & _eps_elem
porosities left == elem, right == neighbor
Definition: PCNSFVHLLC.h:105
static const std::string total_energy_density
Definition: NS.h:65
const ADMaterialProperty< RealVectorValue > & _vel_neighbor
Definition: PCNSFVHLLC.h:86
auto norm() const -> decltype(std::norm(T()))
static const std::string speed
Definition: NS.h:143
const ADMaterialProperty< RealVectorValue > & _vel_elem
velocities left == elem, right == neighbor
Definition: PCNSFVHLLC.h:85
const ADMaterialProperty< Real > & _specific_internal_energy_elem
internal energies left == elem, right == neighbor
Definition: PCNSFVHLLC.h:75
const RealVectorValue & normal() const
T & set(const std::string &name, bool quiet_mode=false)
static const std::string density
Definition: NS.h:33
auto raw_value(const Eigen::Map< T > &in)
ADReal _SL
the wave speeds
Definition: PCNSFVHLLC.h:61
ADReal _SM
Definition: PCNSFVHLLC.h:62
static const std::string fluid
Definition: NS.h:87
RealVectorValue _normal
virtual ADReal fluxElem()=0
flux functions on elem & neighbor, i.e. standard left/right values of F
static const std::string specific_internal_energy
Definition: NS.h:62
void addRequiredParam(const std::string &name, const std::string &doc_string)
ADReal _normal_speed_neighbor
Definition: PCNSFVHLLC.h:68
virtual ADReal hllcElem()=0
HLLC modifications to flux for elem & neighbor, see Toro.
static const std::string porosity
Definition: NS.h:104
const ADMaterialProperty< Real > & _specific_internal_energy_neighbor
Definition: PCNSFVHLLC.h:76
static InputParameters validParams()
const ADMaterialProperty< Real > & _rho_elem
densities left == elem, right == neighbor
Definition: PCNSFVHLLC.h:95
virtual ADReal fluxNeighbor()=0
Real f(Real x)
Test function for Brents method.
ADReal _SR
Definition: PCNSFVHLLC.h:63
const unsigned int _qp
virtual ADReal computeQpResidual() override
Definition: PCNSFVHLLC.C:101
Common class for single phase fluid properties.
static std::array< ADReal, 3 > waveSpeed(const ADReal &rho_elem, const ADRealVectorValue &vel_elem, const ADReal &e_elem, Real eps_elem, const ADReal &rho_neighbor, const ADRealVectorValue &vel_neighbor, const ADReal &e_neighbor, Real eps_neighbor, const SinglePhaseFluidProperties &fluid, const ADRealVectorValue &normal)
helper function for computing wave speed
Definition: PCNSFVHLLC.C:47
virtual ADReal conservedVariableNeighbor()=0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual ADReal hllcNeighbor()=0
static const std::string pressure
Definition: NS.h:56
void mooseError(Args &&... args) const
static InputParameters validParams()
Definition: PCNSFVHLLC.C:17
static const std::string velocity
Definition: NS.h:45
virtual ADReal conservedVariableElem()=0
const ADMaterialProperty< Real > & _rho_neighbor
Definition: PCNSFVHLLC.h:96
const MaterialProperty< Real > & _eps_neighbor
Definition: PCNSFVHLLC.h:106
PCNSFVHLLC(const InputParameters &params)
Definition: PCNSFVHLLC.C:25
ADReal _normal_speed_elem
speeds normal to the interface
Definition: PCNSFVHLLC.h:67
const SinglePhaseFluidProperties & _fluid
fluid properties
Definition: PCNSFVHLLC.h:72