https://mooseframework.inl.gov
ADVolumeJunction1PhaseBC.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 "THMIndicesVACE.h"
13 
14 registerMooseObject("ThermalHydraulicsApp", ADVolumeJunction1PhaseBC);
15 
18 {
20 
21  params.addRequiredParam<unsigned int>("connection_index", "Index of the connected flow channel");
22  params.addRequiredParam<UserObjectName>("volume_junction_uo", "Volume junction user object name");
23 
24  params.addRequiredCoupledVar("A_elem", "Cross-sectional area, elemental");
25  params.addRequiredCoupledVar("A_linear", "Cross-sectional area, linear");
26 
27  params.addRequiredCoupledVar("rhoA", "Flow channel variable: rho*A");
28  params.addRequiredCoupledVar("rhouA", "Flow channel variable: rho*u*A");
29  params.addRequiredCoupledVar("rhoEA", "Flow channel variable: rho*E*A");
30 
31  params.addClassDescription(
32  "Adds boundary fluxes for flow channels connected to a 1-phase volume junction");
33 
34  return params;
35 }
36 
38  : ADOneDIntegratedBC(params),
39 
40  _connection_index(getParam<unsigned int>("connection_index")),
41  _volume_junction_uo(getUserObject<ADVolumeJunction1PhaseUserObject>("volume_junction_uo")),
42 
43  _A_elem(adCoupledValue("A_elem")),
44  _A_linear(adCoupledValue("A_linear")),
45 
46  _rhoA_jvar(coupled("rhoA")),
47  _rhouA_jvar(coupled("rhouA")),
48  _rhoEA_jvar(coupled("rhoEA")),
49 
50  _flow_channel_jvar_map(getFlowChannelIndexMapping()),
51  _equation_index(_flow_channel_jvar_map.at(_var.number()))
52 {
53 }
54 
55 ADReal
57 {
58  const auto & flux = _volume_junction_uo.getFlux(_connection_index);
59 
60  // Note that the ratio A_linear / A_elem is necessary because A_elem is passed
61  // to the flux function, but A_linear is to be used on the boundary.
62  return flux[_equation_index] * _A_linear[_qp] / _A_elem[_qp] * _normal * _test[_i][_qp];
63 }
64 
65 std::map<unsigned int, unsigned int>
67 {
68  std::map<unsigned int, unsigned int> jvar_map;
69  jvar_map.insert(std::pair<unsigned int, unsigned int>(_rhoA_jvar, THMVACE1D::MASS));
70  jvar_map.insert(std::pair<unsigned int, unsigned int>(_rhouA_jvar, THMVACE1D::MOMENTUM));
71  jvar_map.insert(std::pair<unsigned int, unsigned int>(_rhoEA_jvar, THMVACE1D::ENERGY));
72 
73  return jvar_map;
74 }
const ADVariableValue & _A_linear
Cross-sectional area, linear.
const ADVariableValue & _A_elem
Cross-sectional area, elemental.
const unsigned int _connection_index
Index of the connected flow channel.
const ADVolumeJunction1PhaseUserObject & _volume_junction_uo
1-phase volume junction user object
virtual ADReal computeQpResidual() override
std::map< unsigned int, unsigned int > getFlowChannelIndexMapping() const
Creates the mapping of coupled variable index to local equation system index for flow channel variabl...
DualNumber< Real, DNDerivativeType, true > ADReal
const unsigned int _rhoA_jvar
Flow channel rho*A coupled variable index.
void addRequiredParam(const std::string &name, const std::string &doc_string)
unsigned int _qp
registerMooseObject("ThermalHydraulicsApp", ADVolumeJunction1PhaseBC)
const std::vector< ADReal > & getFlux(const unsigned int &connection_index) const override
Gets the flux vector for a connection.
static InputParameters validParams()
Computes and caches flux and residual vectors for a 1-phase volume junction.
const unsigned int _rhoEA_jvar
Flow channel rho*E*A coupled variable index.
const unsigned int _equation_index
Index within local system of the equation upon which this object acts.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
const ADTemplateVariableTestValue< T > & _test
Adds boundary fluxes for flow channels connected to a 1-phase volume junction.
void addClassDescription(const std::string &doc_string)
const unsigned int _rhouA_jvar
Flow channel rho*u*A coupled variable index.
const Real _normal
Component of outward normals along 1-D direction.
ADVolumeJunction1PhaseBC(const InputParameters &params)
Base class for integrated boundary conditions for 1D problems in 3D space.
void ErrorVector unsigned int
static InputParameters validParams()