https://mooseframework.inl.gov
ADNumericalFlux3EqnDGKernel.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 "THMIndicesVACE.h"
14 
15 registerMooseObject("ThermalHydraulicsApp", ADNumericalFlux3EqnDGKernel);
16 
19 {
21 
22  params.addClassDescription(
23  "Adds side fluxes for the 1-D, 1-phase, variable-area Euler equations");
24 
25  params.addRequiredCoupledVar("A_linear", "Cross-sectional area, linear");
26  params.addRequiredCoupledVar("rhoA", "Conserved variable: rho*A");
27  params.addRequiredCoupledVar("rhouA", "Conserved variable: rho*u*A");
28  params.addRequiredCoupledVar("rhoEA", "Conserved variable: rho*E*A");
29 
30  params.addRequiredParam<UserObjectName>("numerical_flux", "Name of numerical flux user object");
31 
32  return params;
33 }
34 
36  : ADDGKernel(parameters),
37 
38  _A_elem(adCoupledValue("A_linear")),
39  _A_neig(adCoupledNeighborValue("A_linear")),
40  _rhoA1(getADMaterialProperty<Real>("rhoA")),
41  _rhouA1(getADMaterialProperty<Real>("rhouA")),
42  _rhoEA1(getADMaterialProperty<Real>("rhoEA")),
43  _p1(getADMaterialProperty<Real>("p")),
44  _rhoA2(getNeighborADMaterialProperty<Real>("rhoA")),
45  _rhouA2(getNeighborADMaterialProperty<Real>("rhouA")),
46  _rhoEA2(getNeighborADMaterialProperty<Real>("rhoEA")),
47  _p2(getNeighborADMaterialProperty<Real>("p")),
48  _numerical_flux(getUserObject<ADNumericalFlux3EqnBase>("numerical_flux")),
49  _rhoA_var(coupled("rhoA")),
50  _rhouA_var(coupled("rhouA")),
51  _rhoEA_var(coupled("rhoEA")),
52  _jmap(getIndexMapping()),
53  _equation_index(_jmap.at(_var.number()))
54 {
55 }
56 
57 ADReal
59 {
60  // construct the left and right solution vectors from the reconstructed solution
61  std::vector<ADReal> U1 = {_rhoA1[_qp], _rhouA1[_qp], _rhoEA1[_qp], _A_elem[_qp]};
62  std::vector<ADReal> U2 = {_rhoA2[_qp], _rhouA2[_qp], _rhoEA2[_qp], _A_neig[_qp]};
63 
64  const Real nLR_dot_d = _current_side * 2 - 1.0;
65 
66  const std::vector<ADReal> & flux_elem =
67  _numerical_flux.getFlux(_current_side, _current_elem->id(), true, U1, U2, nLR_dot_d);
68  const std::vector<ADReal> & flux_neig =
69  _numerical_flux.getFlux(_current_side, _current_elem->id(), false, U1, U2, nLR_dot_d);
70 
71  ADReal re = 0.0;
72  switch (type)
73  {
74  case Moose::Element:
75  re = flux_elem[_equation_index] * _test[_i][_qp];
76  break;
77  case Moose::Neighbor:
78  re = -flux_neig[_equation_index] * _test_neighbor[_i][_qp];
79  break;
80  }
81  return re;
82 }
83 
84 std::map<unsigned int, unsigned int>
86 {
87  std::map<unsigned int, unsigned int> jmap;
88  jmap.insert(std::pair<unsigned int, unsigned int>(_rhoA_var, THMVACE1D::MASS));
89  jmap.insert(std::pair<unsigned int, unsigned int>(_rhouA_var, THMVACE1D::MOMENTUM));
90  jmap.insert(std::pair<unsigned int, unsigned int>(_rhoEA_var, THMVACE1D::ENERGY));
91 
92  return jmap;
93 }
static InputParameters validParams()
std::map< unsigned int, unsigned int > getIndexMapping() const
Creates the mapping of coupled variable index to index in Euler system.
const ADMaterialProperty< Real > & _rhouA1
const ADVariableValue & _A_elem
Area.
const ADNumericalFlux3EqnBase & _numerical_flux
Numerical flux user object.
virtual ADReal computeQpResidual(Moose::DGResidualType type)
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
void addRequiredParam(const std::string &name, const std::string &doc_string)
ADNumericalFlux3EqnDGKernel(const InputParameters &parameters)
registerMooseObject("ThermalHydraulicsApp", ADNumericalFlux3EqnDGKernel)
unsigned int _qp
const ADMaterialProperty< Real > & _rhouA2
const ADMaterialProperty< Real > & _rhoEA2
const VariableTestValue & _test_neighbor
const ADMaterialProperty< Real > & _rhoA1
const std::string & type() const
Base class for computing numerical fluxes for FlowModelSinglePhase.
const unsigned int _equation_index
index within the Euler system of the equation upon which this kernel acts
const ADMaterialProperty< Real > & _rhoEA1
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
const ADMaterialProperty< Real > & _rhoA2
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Elem *const & _current_elem
const unsigned int & _current_side
void addClassDescription(const std::string &doc_string)
const VariableTestValue & _test
Adds side fluxes for the 1-D, 1-phase, variable-area Euler equations.