https://mooseframework.inl.gov
CNSFVHLLCBase.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 "CNSFVHLLCBase.h"
11 #include "NS.h"
12 #include "HLLCUserObject.h"
14 
15 namespace nms = NS;
17 
20 {
22  params.addRequiredParam<UserObjectName>(nms::fluid, "Fluid properties userobject");
23  return params;
24 }
25 
27  : FVFluxKernel(params),
28  _fluid(dynamic_cast<FEProblemBase *>(&_subproblem)
29  ->getUserObject<SinglePhaseFluidProperties>(nms::fluid)),
30  _specific_internal_energy_elem(getADMaterialProperty<Real>(nms::specific_internal_energy)),
31  _specific_internal_energy_neighbor(
32  getNeighborADMaterialProperty<Real>(nms::specific_internal_energy)),
33  _rho_et_elem(getADMaterialProperty<Real>(nms::total_energy_density)),
34  _rho_et_neighbor(getNeighborADMaterialProperty<Real>(nms::total_energy_density)),
35  _vel_elem(getADMaterialProperty<RealVectorValue>(nms::velocity)),
36  _vel_neighbor(getNeighborADMaterialProperty<RealVectorValue>(nms::velocity)),
37  _speed_elem(getADMaterialProperty<Real>(nms::speed)),
38  _speed_neighbor(getNeighborADMaterialProperty<Real>(nms::speed)),
39  _rho_elem(getADMaterialProperty<Real>(nms::density)),
40  _rho_neighbor(getNeighborADMaterialProperty<Real>(nms::density)),
41  _pressure_elem(getADMaterialProperty<Real>(nms::pressure)),
42  _pressure_neighbor(getNeighborADMaterialProperty<Real>(nms::pressure))
43 {
44 }
45 
48 {
49  return {_fluid,
50  _rho_elem[_qp],
52  _vel_elem[_qp],
56 }
57 
58 std::array<ADReal, 3>
59 CNSFVHLLCBase::waveSpeed(const HLLCData & hllc_data, const ADRealVectorValue & normal)
60 {
61  const ADReal & rho1 = hllc_data.rho_elem;
62  const ADReal u1 = hllc_data.vel_elem.norm();
63  const ADReal q1 = normal * hllc_data.vel_elem;
64  const ADReal v1 = 1.0 / rho1;
65  const ADReal & e1 = hllc_data.e_elem;
66  const ADReal E1 = e1 + 0.5 * u1 * u1;
67  const ADReal p1 = hllc_data.fluid.p_from_v_e(v1, e1);
68  const ADReal H1 = E1 + p1 / rho1;
69  const ADReal c1 = hllc_data.fluid.c_from_v_e(v1, e1);
70 
71  const ADReal & rho2 = hllc_data.rho_neighbor;
72  const ADReal u2 = hllc_data.vel_neighbor.norm();
73  const ADReal q2 = normal * hllc_data.vel_neighbor;
74  const ADReal v2 = 1.0 / rho2;
75  const ADReal & e2 = hllc_data.e_neighbor;
76  const ADReal E2 = e2 + 0.5 * u2 * u2;
77  const ADReal p2 = hllc_data.fluid.p_from_v_e(v2, e2);
78  const ADReal H2 = E2 + p2 / rho2;
79  const ADReal c2 = hllc_data.fluid.c_from_v_e(v2, e2);
80 
81  // compute Roe-averaged variables
82  const ADReal sqrt_rho1 = std::sqrt(rho1);
83  const ADReal sqrt_rho2 = std::sqrt(rho2);
84  const ADReal u_roe = (sqrt_rho1 * u1 + sqrt_rho2 * u2) / (sqrt_rho1 + sqrt_rho2);
85  const ADReal q_roe = (sqrt_rho1 * q1 + sqrt_rho2 * q2) / (sqrt_rho1 + sqrt_rho2);
86  const ADReal H_roe = (sqrt_rho1 * H1 + sqrt_rho2 * H2) / (sqrt_rho1 + sqrt_rho2);
87  const ADReal h_roe = H_roe - 0.5 * u_roe * u_roe;
88  const ADReal rho_roe = std::sqrt(rho1 * rho2);
89  const ADReal v_roe = 1.0 / rho_roe;
90  const ADReal e_roe = hllc_data.fluid.e_from_v_h(v_roe, h_roe);
91  const ADReal c_roe = hllc_data.fluid.c_from_v_e(v_roe, e_roe);
92 
93  // compute wave speeds
94  ADReal SL = std::min(q1 - c1, q_roe - c_roe);
95  ADReal SR = std::max(q2 + c2, q_roe + c_roe);
96  ADReal SM = (rho2 * q2 * (SR - q2) - rho1 * q1 * (SL - q1) + p1 - p2) /
97  (rho2 * (SR - q2) - rho1 * (SL - q1));
98 
99  return {{std::move(SL), std::move(SM), std::move(SR)}};
100 }
static InputParameters validParams()
Definition: CNSFVHLLCBase.C:19
static const std::string total_energy_density
Definition: NS.h:65
H1
static std::array< ADReal, 3 > waveSpeed(const HLLCData &hllc_data, const ADRealVectorValue &normal)
helper function for computing wave speed
Definition: CNSFVHLLCBase.C:59
const ADMaterialProperty< Real > & _specific_internal_energy_elem
internal energies left == elem, right == neighbor
Definition: CNSFVHLLCBase.h:78
const ADReal & e_elem
internal energies left == elem, right == neighbor
Definition: CNSFVHLLCBase.h:41
auto norm() const -> decltype(std::norm(T()))
const ADRealVectorValue & vel_elem
velocities left == elem, right == neighbor
Definition: CNSFVHLLCBase.h:36
const ADReal & e_neighbor
Definition: CNSFVHLLCBase.h:42
static const std::string speed
Definition: NS.h:143
const ADReal & rho_elem
densities left == elem, right == neighbor
Definition: CNSFVHLLCBase.h:31
const RealVectorValue & normal() const
CNSFVHLLCBase(const InputParameters &params)
Definition: CNSFVHLLCBase.C:26
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
HLLCData hllcData() const
Definition: CNSFVHLLCBase.C:47
DualNumber< Real, DNDerivativeType, true > ADReal
const ADMaterialProperty< Real > & _specific_internal_energy_neighbor
Definition: CNSFVHLLCBase.h:79
static const std::string specific_internal_energy
Definition: NS.h:62
void addRequiredParam(const std::string &name, const std::string &doc_string)
const ADMaterialProperty< Real > & _rho_neighbor
Definition: CNSFVHLLCBase.h:97
static InputParameters validParams()
const SinglePhaseFluidProperties & fluid
fluid properties
Definition: CNSFVHLLCBase.h:28
const unsigned int _qp
H2
Common class for single phase fluid properties.
const ADReal & rho_neighbor
Definition: CNSFVHLLCBase.h:32
const ADMaterialProperty< RealVectorValue > & _vel_neighbor
Definition: CNSFVHLLCBase.h:87
const SinglePhaseFluidProperties & _fluid
fluid properties
Definition: CNSFVHLLCBase.h:75
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string pressure
Definition: NS.h:56
const ADRealVectorValue & vel_neighbor
Definition: CNSFVHLLCBase.h:37
static const std::string velocity
Definition: NS.h:45
Helper structure for holding data necessary for computing HLLC fluxes.
Definition: CNSFVHLLCBase.h:25
const ADMaterialProperty< Real > & _rho_elem
densities left == elem, right == neighbor
Definition: CNSFVHLLCBase.h:96
const ADMaterialProperty< RealVectorValue > & _vel_elem
velocities left == elem, right == neighbor
Definition: CNSFVHLLCBase.h:86