https://mooseframework.inl.gov
HLLCUserObject.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 "HLLCUserObject.h"
11 #include "NS.h"
13 
15 
16 registerMooseObject("NavierStokesApp", HLLCUserObject);
17 
20 {
22  params.addRequiredParam<UserObjectName>(NS::fluid, "Fluid properties userobject");
23  params.addParam<MaterialPropertyName>(
25  "An optional parameter for specifying a porosity material property. If not specified "
26  "free-flow conditions will be assumed.");
27  params.addClassDescription(
28  "Computes free-flow wave speeds on internal sides, useful in HLLC contexts");
29  return params;
30 }
31 
33  : InternalSideUserObject(parameters),
34  _fluid(UserObjectInterface::getUserObject<SinglePhaseFluidProperties>(NS::fluid)),
35  _face_info(_mesh.faceInfo()),
36  _vel_elem(getADMaterialProperty<RealVectorValue>(NS::velocity)),
37  _vel_neighbor(getNeighborADMaterialProperty<RealVectorValue>(NS::velocity)),
38  _speed_elem(getADMaterialProperty<Real>(NS::speed)),
39  _speed_neighbor(getNeighborADMaterialProperty<Real>(NS::speed)),
40  _pressure_elem(getADMaterialProperty<Real>(NS::pressure)),
41  _pressure_neighbor(getNeighborADMaterialProperty<Real>(NS::pressure)),
42  _rho_elem(getADMaterialProperty<Real>(NS::density)),
43  _rho_neighbor(getNeighborADMaterialProperty<Real>(NS::density)),
44  _specific_internal_energy_elem(getADMaterialProperty<Real>(NS::specific_internal_energy)),
45  _specific_internal_energy_neighbor(
46  getNeighborADMaterialProperty<Real>(NS::specific_internal_energy)),
47  _eps_elem(isParamValid(NS::porosity) ? &getMaterialProperty<Real>(NS::porosity) : nullptr),
48  _eps_neighbor(isParamValid(NS::porosity) ? &getNeighborMaterialProperty<Real>(NS::porosity)
49  : nullptr)
50 {
51 }
52 
53 void
55 {
56  // This means that we have no finite volume variables, so we need to build the
57  // finite volume information on the mesh so that we can still use this object
58  if (!_face_info.size())
60 
61  for (unsigned int j = 0; j < _face_info.size(); ++j)
62  {
63  const Elem * elem = &_face_info[j]->elem();
64  unsigned int side = _face_info[j]->elemSideID();
65  side_type elem_side(elem, side);
66  _side_to_face_info[elem_side] = j;
67  }
68 }
69 
70 void
72 {
73  using std::sqrt, std::min, std::max;
74 
76 
77  const ADReal & rho1 = _rho_elem[_qp];
78  const ADReal & u1 = _speed_elem[_qp];
79  const ADReal q1 = fi.normal() * _vel_elem[_qp];
80  const ADReal v1 = 1.0 / rho1;
82  const ADReal E1 = e1 + 0.5 * u1 * u1;
83  const ADReal & p1 = _pressure_elem[_qp];
84  const ADReal H1 = E1 + p1 / rho1;
85  const ADReal c1 = _fluid.c_from_v_e(v1, e1);
86  const Real eps1 = _eps_elem ? (*_eps_elem)[_qp] : 1.;
87 
88  const ADReal & rho2 = _rho_neighbor[_qp];
89  const ADReal & u2 = _speed_neighbor[_qp];
90  const ADReal q2 = fi.normal() * _vel_neighbor[_qp];
91  const ADReal v2 = 1.0 / rho2;
93  const ADReal E2 = e2 + 0.5 * u2 * u2;
94  const ADReal & p2 = _pressure_neighbor[_qp];
95  const ADReal H2 = E2 + p2 / rho2;
96  const ADReal c2 = _fluid.c_from_v_e(v2, e2);
97  const Real eps2 = _eps_neighbor ? (*_eps_neighbor)[_qp] : 1.;
98 
99  // compute Roe-averaged variables
100  const ADReal sqrt_rho1 = sqrt(rho1);
101  const ADReal sqrt_rho2 = sqrt(rho2);
102  const ADReal u_roe = (sqrt_rho1 * u1 + sqrt_rho2 * u2) / (sqrt_rho1 + sqrt_rho2);
103  const ADReal q_roe = (sqrt_rho1 * q1 + sqrt_rho2 * q2) / (sqrt_rho1 + sqrt_rho2);
104  const ADReal H_roe = (sqrt_rho1 * H1 + sqrt_rho2 * H2) / (sqrt_rho1 + sqrt_rho2);
105  const ADReal h_roe = H_roe - 0.5 * u_roe * u_roe;
106  const ADReal rho_roe = sqrt(rho1 * rho2);
107  const ADReal v_roe = 1.0 / rho_roe;
108  const ADReal e_roe = _fluid.e_from_v_h(v_roe, h_roe);
109  const ADReal c_roe = _fluid.c_from_v_e(v_roe, e_roe);
110 
111  // compute wave speeds
112  const ADReal SL = min(q1 - c1, q_roe - c_roe);
113  const ADReal SR = max(q2 + c2, q_roe + c_roe);
114  const ADReal SM =
115  (eps2 * rho2 * q2 * (SR - q2) - eps1 * rho1 * q1 * (SL - q1) + eps1 * p1 - eps2 * p2) /
116  (eps2 * rho2 * (SR - q2) - eps1 * rho1 * (SL - q1));
117 
118  // store these results in _wave_speed
120  _wave_speed[elem_side] = {SL, SM, SR};
121 }
122 
123 bool
124 HLLCUserObject::hasData(const Elem * const elem, const unsigned int side) const
125 {
126  side_type elem_side(elem, side);
127  return (_wave_speed.find(elem_side) != _wave_speed.end());
128 }
129 
130 void
132 {
133  const auto & pps = static_cast<const HLLCUserObject &>(y);
134  for (auto & ws : pps._wave_speed)
135  {
136  const auto & it = _wave_speed.find(ws.first);
137  if (it == _wave_speed.end())
138  _wave_speed[ws.first] = ws.second;
139  else
140  for (unsigned int j = 0; j < 3; ++j)
141  _wave_speed[ws.first][j] = ws.second[j];
142  }
143 }
144 
145 std::vector<ADReal>
146 HLLCUserObject::waveSpeed(const Elem * elem, unsigned int side) const
147 {
148  side_type elem_side(elem, side);
149  auto it = _wave_speed.find(elem_side);
150  if (it == _wave_speed.end())
151  mooseError("what the heck (version 2)?");
152  return it->second;
153 }
154 
155 const FaceInfo &
156 HLLCUserObject::faceInfoHelper(const Elem * elem, unsigned int side) const
157 {
158  side_type elem_side(elem, side);
159  auto it = _side_to_face_info.find(elem_side);
160  if (it == _side_to_face_info.end())
161  mooseError("what the heck?");
162  return *_face_info[it->second];
163 }
const unsigned int & _current_side
const SinglePhaseFluidProperties & _fluid
fluid properties
H1
const ADMaterialProperty< Real > & _speed_elem
const ADMaterialProperty< Real > & _speed_neighbor
static const std::string speed
Definition: NS.h:143
const ADMaterialProperty< Real > & _pressure_neighbor
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const ADMaterialProperty< RealVectorValue > & _vel_elem
material properties computed by VarMat that Riemann solver needs
const ADMaterialProperty< Real > & _specific_internal_energy_elem
static InputParameters validParams()
static const std::string density
Definition: NS.h:33
auto raw_value(const Eigen::Map< T > &in)
static const std::string fluid
Definition: NS.h:87
const std::vector< double > y
virtual void initialSetup() override
registerMooseObject("NavierStokesApp", HLLCUserObject)
bool hasData(const Elem *elem, unsigned int side) const
Query whether this processor has data for the provided element and side.
std::map< side_type, std::vector< ADReal > > _wave_speed
data structure storing the wave speeds SL, SM, SR
virtual void threadJoin(const UserObject &y) override
static const std::string specific_internal_energy
Definition: NS.h:62
void addRequiredParam(const std::string &name, const std::string &doc_string)
auto max(const L &left, const R &right)
static const std::string porosity
Definition: NS.h:104
virtual void execute() override
const ADMaterialProperty< Real > & _rho_elem
const MaterialProperty< Real > *const _eps_elem
const MaterialProperty< Real > *const _eps_neighbor
static InputParameters validParams()
unsigned int _qp
quadrature point dummy
H2
Common class for single phase fluid properties.
const ADMaterialProperty< Real > & _specific_internal_energy_neighbor
const Point & normal() const
std::map< side_type, unsigned int > _side_to_face_info
face info lookup allows searching for face info entry from elem/side pair
const std::vector< const FaceInfo * > & _face_info
FV face info from MooseMesh.
HLLCUserObject(const InputParameters &parameters)
const Elem *const & _current_elem
std::vector< ADReal > waveSpeed(const Elem *elem, unsigned int side) const
accessor for the wave speed
std::pair< const Elem *, unsigned int > side_type
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const ADMaterialProperty< Real > & _pressure_elem
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
static const std::string pressure
Definition: NS.h:56
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
static const std::string velocity
Definition: NS.h:45
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const FaceInfo & faceInfoHelper(const Elem *elem, unsigned int side) const
helper function for returning the FaceInfo object for an elem/side pair
const ADMaterialProperty< RealVectorValue > & _vel_neighbor
auto min(const L &left, const R &right)
for(PetscInt i=0;i< nvars;++i)
const ADMaterialProperty< Real > & _rho_neighbor
void setupFiniteVolumeMeshData() const