https://mooseframework.inl.gov
NumericalFluxGasMixDGKernel.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 
12 #include "MooseVariable.h"
13 #include "THMIndicesGasMix.h"
14 #include "THMNames.h"
15 
16 registerMooseObject("ThermalHydraulicsApp", NumericalFluxGasMixDGKernel);
17 
20 {
22 
23  params.addClassDescription("Adds side fluxes from NumericalFluxGasMixBase objects.");
24 
25  params.addRequiredCoupledVar("A_linear", "Cross-sectional area, linear");
26  params.addRequiredCoupledVar("xirhoA", "Conserved variable xi*rho*A");
27  params.addRequiredCoupledVar("rhoA", "Conserved variable rho*A");
28  params.addRequiredCoupledVar("rhouA", "Conserved variable rho*u*A");
29  params.addRequiredCoupledVar("rhoEA", "Conserved variable rho*E*A");
30 
31  params.addRequiredParam<UserObjectName>("numerical_flux",
32  "Name of NumericalFluxGasMixBase object");
33 
34  return params;
35 }
36 
38  : ADDGKernel(parameters),
39 
40  _A_elem(adCoupledValue("A_linear")),
41  _A_neig(adCoupledNeighborValue("A_linear")),
42  _xirhoA1(getADMaterialProperty<Real>(THM::XIRHOA)),
43  _rhoA1(getADMaterialProperty<Real>(THM::RHOA)),
44  _rhouA1(getADMaterialProperty<Real>(THM::RHOUA)),
45  _rhoEA1(getADMaterialProperty<Real>(THM::RHOEA)),
46  _xirhoA2(getNeighborADMaterialProperty<Real>(THM::XIRHOA)),
47  _rhoA2(getNeighborADMaterialProperty<Real>(THM::RHOA)),
48  _rhouA2(getNeighborADMaterialProperty<Real>(THM::RHOUA)),
49  _rhoEA2(getNeighborADMaterialProperty<Real>(THM::RHOEA)),
50  _numerical_flux(getUserObject<NumericalFluxGasMixBase>("numerical_flux")),
51  _xirhoA_var(coupled("xirhoA")),
52  _rhoA_var(coupled("rhoA")),
53  _rhouA_var(coupled("rhouA")),
54  _rhoEA_var(coupled("rhoEA")),
55  _jmap(getIndexMapping()),
56  _equation_index(_jmap.at(_var.number()))
57 {
58 }
59 
60 ADReal
62 {
63  // construct the left and right solution vectors from the reconstructed solution
64  std::vector<ADReal> U1(THMGasMix1D::N_FLUX_INPUTS, 0.0);
70 
71  std::vector<ADReal> U2(THMGasMix1D::N_FLUX_INPUTS, 0.0);
77 
78  const Real nLR_dot_d = _current_side * 2 - 1.0;
79 
80  const std::vector<ADReal> & flux_elem =
81  _numerical_flux.getFlux(_current_side, _current_elem->id(), true, U1, U2, nLR_dot_d);
82  const std::vector<ADReal> & flux_neig =
83  _numerical_flux.getFlux(_current_side, _current_elem->id(), false, U1, U2, nLR_dot_d);
84 
85  ADReal re = 0.0;
86  switch (type)
87  {
88  case Moose::Element:
89  re = flux_elem[_equation_index] * _test[_i][_qp];
90  break;
91  case Moose::Neighbor:
92  re = -flux_neig[_equation_index] * _test_neighbor[_i][_qp];
93  break;
94  }
95  return re;
96 }
97 
98 std::map<unsigned int, unsigned int>
100 {
101  std::map<unsigned int, unsigned int> jmap;
102  jmap.insert(std::pair<unsigned int, unsigned int>(_xirhoA_var, THMGasMix1D::SPECIES));
103  jmap.insert(std::pair<unsigned int, unsigned int>(_rhoA_var, THMGasMix1D::MASS));
104  jmap.insert(std::pair<unsigned int, unsigned int>(_rhouA_var, THMGasMix1D::MOMENTUM));
105  jmap.insert(std::pair<unsigned int, unsigned int>(_rhoEA_var, THMGasMix1D::ENERGY));
106 
107  return jmap;
108 }
const ADMaterialProperty< Real > & _rhoA2
const ADMaterialProperty< Real > & _rhoEA1
Base class for computing numerical fluxes for FlowModelGasMix.
const ADMaterialProperty< Real > & _xirhoA1
registerMooseObject("ThermalHydraulicsApp", NumericalFluxGasMixDGKernel)
NumericalFluxGasMixDGKernel(const InputParameters &parameters)
const ADMaterialProperty< Real > & _rhoA1
static InputParameters validParams()
unsigned int _i
DGResidualType
virtual const std::vector< ADReal > & getFlux(const unsigned int iside, const dof_id_type ielem, bool res_side_is_left, const std::vector< ADReal > &UL_1d, const std::vector< ADReal > &UR_1d, Real nLR_dot_d) const
Gets the 1D flux vector for an element/side combination.
DualNumber< Real, DNDerivativeType, true > ADReal
std::map< unsigned int, unsigned int > getIndexMapping() const
Creates the mapping of coupled variable index to index in equation system.
void addRequiredParam(const std::string &name, const std::string &doc_string)
unsigned int _qp
const VariableTestValue & _test_neighbor
const ADMaterialProperty< Real > & _rhouA1
const std::string & type() const
const ADMaterialProperty< Real > & _rhouA2
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Adds side fluxes from NumericalFluxGasMix objects.
const NumericalFluxGasMixBase & _numerical_flux
Numerical flux user object.
const Elem *const & _current_elem
const unsigned int & _current_side
const unsigned int _equation_index
index within the equation system of the equation upon which this kernel acts
void addClassDescription(const std::string &doc_string)
static const unsigned int N_FLUX_INPUTS
Number of numerical flux function inputs.
const ADMaterialProperty< Real > & _xirhoA2
virtual ADReal computeQpResidual(Moose::DGResidualType type) override
const ADVariableValue & _A_elem
Area.
const ADMaterialProperty< Real > & _rhoEA2
const VariableTestValue & _test
static InputParameters validParams()