https://mooseframework.inl.gov
HSCoupler2D3DUserObject.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 "MeshAlignment2D3D.h"
13 #include "THMUtils.h"
14 #include "HeatTransferModels.h"
15 #include "Function.h"
16 
17 registerMooseObject("ThermalHydraulicsApp", HSCoupler2D3DUserObject);
18 
21 {
23 
24  params.addRequiredCoupledVar("temperature", "Temperature");
25  params.addRequiredParam<Real>("radius_2d", "Radius of the 2D heat structure boundary [m]");
26  params.addParam<FunctionName>(
27  "emissivity_2d",
28  0,
29  "Emissivity of the 2D heat structure boundary as a function of temperature [K]");
30  params.addParam<FunctionName>(
31  "emissivity_3d",
32  0,
33  "Emissivity of the 3D heat structure boundary as a function of temperature [K]");
34  params.addRequiredParam<FunctionName>("gap_thickness",
35  "Gap thickness [m] as a function of temperature [K]");
36  params.addRequiredParam<FunctionName>(
37  "gap_thermal_conductivity",
38  "Gap thermal conductivity [W/(m-K)] as a function of temperature [K]");
39  params.addParam<FunctionName>(
40  "gap_htc", 0, "Gap heat transfer coefficient [W/(m^2-K)] as a function of temperature [K]");
41  params.addRequiredParam<UserObjectName>(
42  "temperature_2d_uo",
43  "StoreVariableByElemIDSideUserObject containing the temperature values on the 2D boundary");
44  params.addRequiredParam<MeshAlignment2D3D *>("mesh_alignment", "Mesh alignment object");
45 
46  params.addClassDescription("Computes heat fluxes for HSCoupler2D3D.");
47 
48  return params;
49 }
50 
52  : SideUserObject(parameters),
53 
54  _T_3d(adCoupledValue("temperature")),
55  _r_2d(getParam<Real>("radius_2d")),
56  _emissivity_2d_fn(getFunction("emissivity_2d")),
57  _emissivity_3d_fn(getFunction("emissivity_3d")),
58  _include_radiation(isParamSetByUser("emissivity_2d") && isParamSetByUser("emissivity_3d")),
59  _gap_thickness_fn(getFunction("gap_thickness")),
60  _k_gap_fn(getFunction("gap_thermal_conductivity")),
61  _htc_gap_fn(getFunction("gap_htc")),
62  _temperature_2d_uo(getUserObject<StoreVariableByElemIDSideUserObject>("temperature_2d_uo")),
63  _mesh_alignment(*getParam<MeshAlignment2D3D *>("mesh_alignment"))
64 {
66 }
67 
68 void
70 {
71  _elem_id_to_heat_flux.clear();
72 }
73 
74 void
76 {
77  const auto elem_id_3d = _current_elem->id();
78  const auto elem_id_2d = _mesh_alignment.getCoupledPrimaryElemID(elem_id_3d);
79 
81  const auto n_qp_3d = _qrule->n_points();
82 
83  const auto & T_2d_values = _temperature_2d_uo.getVariableValues(elem_id_2d);
84  const auto & area_2d = _mesh_alignment.getPrimaryArea(elem_id_2d);
85 
86  std::vector<ADReal> heat_flux_3d(n_qp_3d, 0.0);
87  std::vector<ADReal> heat_flux_2d(n_qp_2d, 0.0);
88  for (unsigned int qp_3d = 0; qp_3d < n_qp_3d; qp_3d++)
89  {
90  const auto qp_2d = _mesh_alignment.getCoupledPrimaryElemQpIndex(elem_id_3d, qp_3d);
91 
92  const auto & T_2d = T_2d_values[qp_2d];
93  const auto & T_3d = _T_3d[qp_3d];
94  const auto T_gap = 0.5 * (T_2d + T_3d);
95 
96  const auto gap_thickness = evaluateTemperatureFunction(_gap_thickness_fn, T_gap);
97  const auto k_gap = evaluateTemperatureFunction(_k_gap_fn, T_gap);
98 
99  if (!MooseUtils::absoluteFuzzyGreaterThan(gap_thickness, 0))
100  mooseError("Gap thickness must be > 0.");
101 
102  const auto r_3d = _r_2d + gap_thickness;
103 
104  const auto heat_flux_cond =
106  auto heat_flux = heat_flux_cond;
107 
108  const auto htc = evaluateTemperatureFunction(_htc_gap_fn, T_gap);
109  heat_flux += htc * (T_2d - T_3d);
110 
111  if (_include_radiation)
112  {
113  const auto emissivity_2d = evaluateTemperatureFunction(_emissivity_2d_fn, T_2d);
114  const auto emissivity_3d = evaluateTemperatureFunction(_emissivity_3d_fn, T_3d);
115 
117  _r_2d, r_3d, emissivity_2d, emissivity_3d, T_2d, T_3d);
118  }
119 
120  heat_flux_3d[qp_3d] = heat_flux;
121  heat_flux_2d[qp_2d] -= _JxW[qp_3d] * _coord[qp_3d] * heat_flux / area_2d[qp_2d];
122  }
123 
124  // Store values in maps
125  _elem_id_to_heat_flux[elem_id_3d] = heat_flux_3d;
126  auto it = _elem_id_to_heat_flux.find(elem_id_2d);
127  if (it == _elem_id_to_heat_flux.end())
128  _elem_id_to_heat_flux[elem_id_2d] = heat_flux_2d;
129  else
130  {
131  auto & heat_flux_2d_existing = _elem_id_to_heat_flux[elem_id_2d];
132  for (const auto qp_2d : index_range(heat_flux_2d_existing))
133  heat_flux_2d_existing[qp_2d] += heat_flux_2d[qp_2d];
134  }
135 }
136 
137 void
139 {
140  const auto & other_uo = static_cast<const HSCoupler2D3DUserObject &>(uo);
141  for (auto & it : other_uo._elem_id_to_heat_flux)
142  if (_elem_id_to_heat_flux.find(it.first) == _elem_id_to_heat_flux.end())
143  _elem_id_to_heat_flux[it.first] = it.second;
144  else
145  {
146  auto & existing = _elem_id_to_heat_flux[it.first];
147  for (const auto qp : index_range(existing))
148  existing[qp] += it.second[qp];
149  }
150 }
151 
152 void
154 {
156 }
157 
158 const std::vector<ADReal> &
160 {
161  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
162 
163  auto it = _elem_id_to_heat_flux.find(elem_id);
164  if (it != _elem_id_to_heat_flux.end())
165  return it->second;
166  else
167  mooseError(name(), ": Requested heat flux for element ", elem_id, " was not computed.");
168 }
169 
170 ADReal
172 {
173  ADReal f = fn.value(T);
174  f.derivatives() = T.derivatives() * fn.timeDerivative(T.value());
175 
176  return f;
177 }
HSCoupler2D3DUserObject(const InputParameters &parameters)
virtual void execute() override
const Function & _emissivity_3d_fn
Emissivity of the 3D heat structure boundary as a function of temperature.
const std::vector< ADReal > & getVariableValues(dof_id_type elem_id) const
Gets the variable values at each quadrature point on the provided element.
static InputParameters validParams()
Assembly & _assembly
const bool _include_radiation
Include radiation?
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
virtual Real timeDerivative(Real t, const Point &p) const
const Function & _gap_thickness_fn
Gap thickness as a function of temperature.
MeshAlignment2D3D & _mesh_alignment
Mesh alignment object.
const Function & _htc_gap_fn
Gap heat transfer coefficient as a function of temperature.
const ADVariableValue & _T_3d
Coupled heated perimeter variable.
Stores variable values at each quadrature point on a side by element ID.
const Parallel::Communicator & comm() const
DualNumber< Real, DNDerivativeType, true > ADReal
virtual const std::string & name() const
void addRequiredParam(const std::string &name, const std::string &doc_string)
const MooseArray< Real > & _JxW
unsigned int getCoupledPrimaryElemQpIndex(const dof_id_type &secondary_elem_id, const unsigned int &secondary_qp) const
Gets the quadrature point index on the primary element corresponding to the quadrature point index on...
const Function & _k_gap_fn
Gap thermal conductivity as a function of temperature.
Real f(Real x)
Test function for Brents method.
virtual void finalize() override
auto cylindricalGapRadiationHeatFlux(const T1 &r_inner, const T2 &r_outer, const T3 &emiss_inner, const T4 &emiss_outer, const T5 &T_inner, const T6 &T_outer, const Real &sigma=HeatConduction::Constants::sigma)
Computes the radiation heat flux across a cylindrical gap.
const MooseArray< Real > & _coord
const std::vector< ADReal > & getHeatFlux(dof_id_type elem_id) const
Gets the heat fluxes for each quadrature point for a given element ID.
virtual void initialize() override
void buildCoupledElemQpIndexMapSecondary(Assembly &assembly)
Builds the map used for getting the coupled quadrature point index.
Computes heat fluxes for HSCoupler2D3D.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
ADReal evaluateTemperatureFunction(const Function &fn, const ADReal &T) const
Evaluates a function of temperature.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< Real > & getPrimaryArea(const dof_id_type primary_elem_id) const
Gets the area for each quadrature point on a primary element.
auto cylindricalGapConductionHeatFlux(const T1 &k_gap, const T2 &r_inner, const T3 &r_outer, const T4 &T_inner, const T5 &T_outer)
Computes the conduction heat flux across a cylindrical gap.
const Real & _r_2d
Radius of the 2D heat structure boundary.
const QBase *const & _qrule
unsigned int getPrimaryNumberOfQuadraturePoints() const
Gets the number of quadrature points for faces on the primary boundary.
void allGatherADVectorMapSum(const Parallel::Communicator &comm, std::map< dof_id_type, std::vector< ADReal >> &this_map)
Parallel gather of a map of DoF ID to AD vector.
Definition: THMUtils.C:58
dof_id_type getCoupledPrimaryElemID(const dof_id_type &secondary_elem_id) const
Gets the coupled primary element ID for a given secondary element ID.
const StoreVariableByElemIDSideUserObject & _temperature_2d_uo
User object containing the temperature values on the 2D boundary.
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
std::map< dof_id_type, std::vector< ADReal > > _elem_id_to_heat_flux
Map of the element ID to the heat flux.
virtual void threadJoin(const UserObject &uo) override
const Elem *const & _current_elem
const Function & _emissivity_2d_fn
Emissivity of the 2D heat structure boundary as a function of temperature.
virtual Real value(Real t, const Point &p) const
registerMooseObject("ThermalHydraulicsApp", HSCoupler2D3DUserObject)
Builds mapping between a 2D boundary and a 3D boundary.
for(PetscInt i=0;i< nvars;++i)
auto index_range(const T &sizable)
bool absoluteFuzzyGreaterThan(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
uint8_t dof_id_type
static InputParameters validParams()