19 #include "MooseMesh.h"
21 #include "libmesh/quadrature.h"
27 InputParameters params = validParams<Material>();
29 params.addClassDescription(
"This is the base class all materials should use if you are trying to "
30 "use the Navier-Stokes Kernels.");
36 params.addRequiredCoupledVar(
NS::enthalpy,
"total enthalpy");
38 params.addRequiredCoupledVar(
NS::density,
"density");
43 params.addRequiredParam<UserObjectName>(
"fluid_properties",
44 "The name of the user object for fluid properties");
50 : Material(parameters),
51 _mesh_dimension(_mesh.dimension()),
53 _grad_v(_mesh_dimension >= 2 ? coupledGradient(
NS::
velocity_y) : _grad_zero),
54 _grad_w(_mesh_dimension == 3 ? coupledGradient(
NS::
velocity_z) : _grad_zero),
56 _viscous_stress_tensor(declareProperty<RealTensorValue>(
"viscous_stress_tensor")),
57 _thermal_conductivity(declareProperty<Real>(
"thermal_conductivity")),
61 _dynamic_viscosity(declareProperty<Real>(
"dynamic_viscosity")),
64 _calA(declareProperty<std::vector<RealTensorValue>>(
"calA")),
67 _calC(declareProperty<std::vector<RealTensorValue>>(
"calC")),
70 _calE(declareProperty<std::vector<std::vector<RealTensorValue>>>(
"calE")),
75 _v_vel(_mesh.dimension() >= 2 ? coupledValue(
NS::velocity_y) : _zero),
76 _w_vel(_mesh.dimension() == 3 ? coupledValue(
NS::
velocity_z) : _zero),
84 _rho_v(_mesh.dimension() >= 2 ? coupledValue(
NS::
momentum_y) : _zero),
85 _rho_w(_mesh.dimension() == 3 ? coupledValue(
NS::
momentum_z) : _zero),
91 _drhov_dt(_mesh.dimension() >= 2 ? coupledDot(
NS::
momentum_y) : _zero),
92 _drhow_dt(_mesh.dimension() == 3 ? coupledDot(
NS::
momentum_z) : _zero),
98 _grad_rho_v(_mesh.dimension() >= 2 ? coupledGradient(
NS::
momentum_y) : _grad_zero),
99 _grad_rho_w(_mesh.dimension() == 3 ? coupledGradient(
NS::
momentum_z) : _grad_zero),
103 _hsupg(declareProperty<Real>(
"hsupg")),
104 _tauc(declareProperty<Real>(
"tauc")),
105 _taum(declareProperty<Real>(
"taum")),
106 _taue(declareProperty<Real>(
"taue")),
107 _strong_residuals(declareProperty<std::vector<Real>>(
"strong_residuals")),
118 for (
unsigned int qp = 0; qp < _qrule->n_points(); ++qp)
125 grad_outer_u += grad_outer_u.transpose();
128 for (
unsigned int i = 0; i < 3; ++i)
132 for (
unsigned int i = 0; i < 3; ++i)
133 grad_outer_u(i, i) -= 2.0 / 3.0 * div_vel;
153 const Real Pr = 0.71;
259 _hsupg[qp] = _current_elem->hmin();
338 Real sqrt_term = 4. / _dt / _dt + velmag * velmag / h2;
341 _tauc[qp] = 1. / std::sqrt(sqrt_term);
342 _taum[qp] = 1. / std::sqrt(sqrt_term + visc_term * visc_term);
343 _taue[qp] = 1. / std::sqrt(sqrt_term + k_term * k_term);
370 RealVectorValue zero(0., 0., 0.);
373 Real velmag2 = vel.norm_sq();
391 for (
unsigned int i = 0; i < 3; ++i)
410 RealTensorValue calS =
_calC[qp][0] *
_calC[qp][0].transpose();
418 0.5 * (
_fp.
gamma() - 1.0) * velmag2;
419 _calA[qp][0] -= calS;
421 for (
unsigned int m = 1; m <= 3; ++m)
424 unsigned int m_local = m - 1;
431 _calA[qp][m] +=
_calC[qp][m_local].transpose();
442 for (
unsigned int k = 0; k < 3; ++k)
445 _calE[qp][k].resize(5);
449 RealTensorValue Ck_T =
_calC[qp][k].transpose();
452 _calE[qp][k][0].zero();
455 for (
unsigned int m = 1; m <= 3; ++m)
458 unsigned int m_local = m - 1;
461 _calE[qp][k][m].zero();
467 _calE[qp][k][4].zero();