https://mooseframework.inl.gov
HSCoupler2D2DRadiationUserObject.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 "MeshAlignment2D2D.h"
13 #include "THMUtils.h"
14 #include "HeatConductionNames.h"
15 
17 
20 {
22 
23  params.addRequiredParam<std::vector<Real>>("emissivities", "Emissivities of each surface");
24  params.addRequiredParam<std::vector<std::vector<Real>>>(
25  "view_factors", "View factors between each surface, as a matrix");
26  params.addRequiredParam<bool>(
27  "include_environment",
28  "Whether or not to include an environment surrounding all of the surfaces.");
29  params.addParam<Real>("T_environment", 0.0, "Environment temperature.");
30  params.addRequiredParam<UserObjectName>(
31  "temperature_uo",
32  "StoreVariableByElemIDSideUserObjects containing the temperature values for each element");
33  params.addRequiredParam<MeshAlignment2D2D *>("mesh_alignment", "Mesh alignment object");
34 
35  params.addClassDescription("Computes heat fluxes for HSCoupler2D2D.");
36 
37  return params;
38 }
39 
41  const InputParameters & parameters)
42  : SideUserObject(parameters),
43 
44  _emissivities(getParam<std::vector<Real>>("emissivities")),
45  _view_factors(getParam<std::vector<std::vector<Real>>>("view_factors")),
46  _include_environment(getParam<bool>("include_environment")),
47  _T_environment(getParam<Real>("T_environment")),
48  _temperature_uo(getUserObject<StoreVariableByElemIDSideUserObject>("temperature_uo")),
49  _mesh_alignment(*getParam<MeshAlignment2D2D *>("mesh_alignment")),
50  _n_hs(_emissivities.size()),
51  _n_surfaces(_include_environment ? _n_hs + 1 : _n_hs)
52 {
53 }
54 
55 void
57 {
58  _elem_id_to_heat_flux.clear();
59 }
60 
61 void
63 {
64  using std::pow;
65 
66  // store element IDs corresponding to this side
67  std::vector<dof_id_type> elem_ids;
68  elem_ids.push_back(_current_elem->id());
69  const auto coupled_elem_ids = _mesh_alignment.getCoupledSecondaryElemIDs(elem_ids[0]);
70  elem_ids.insert(elem_ids.end(), coupled_elem_ids.begin(), coupled_elem_ids.end());
71 
72  // store temperatures
73  std::vector<std::vector<ADReal>> temperatures;
74  for (const auto & elem_id : elem_ids)
75  temperatures.push_back(_temperature_uo.getVariableValues(elem_id));
76 
77  // loop over the quadrature points and compute the heat fluxes
78  const auto n_qp = _qrule->n_points();
79  std::vector<std::vector<ADReal>> heat_fluxes(_n_hs, std::vector<ADReal>(n_qp, 0.0));
80  for (unsigned int qp = 0; qp < n_qp; qp++)
81  {
82  // form matrix and RHS
83  DenseMatrix<ADReal> matrix(_n_surfaces, _n_surfaces);
84  DenseVector<ADReal> emittances(_n_surfaces);
85  for (unsigned int i = 0; i < _n_surfaces; ++i)
86  {
87  if (_include_environment && i == _n_surfaces - 1) // environment surface
88  {
90  matrix(i, i) = 1.0;
91  }
92  else // heat structure surface
93  {
94  emittances(i) =
95  _emissivities[i] * HeatConduction::Constants::sigma * pow(temperatures[i][qp], 4);
96  matrix(i, i) = 1.0;
97  for (unsigned int j = 0; j < _n_surfaces; ++j)
98  matrix(i, j) -= (1 - _emissivities[i]) * _view_factors[i][j];
99  }
100  }
101 
102  // solve for radiosities
103  DenseVector<ADReal> radiosities(_n_surfaces);
104  matrix.lu_solve(emittances, radiosities);
105 
106  // compute heat fluxes
107  for (unsigned int i = 0; i < _n_hs; ++i)
108  heat_fluxes[i][qp] =
109  1.0 / (1.0 - _emissivities[i]) * (emittances(i) - _emissivities[i] * radiosities(i));
110  }
111 
112  // store the heat fluxes by element ID
113  for (unsigned int i = 0; i < _n_hs; ++i)
114  _elem_id_to_heat_flux[elem_ids[i]] = heat_fluxes[i];
115 }
116 
117 void
119 {
120  const auto & other_uo = static_cast<const HSCoupler2D2DRadiationUserObject &>(uo);
121  for (auto & it : other_uo._elem_id_to_heat_flux)
122  if (_elem_id_to_heat_flux.find(it.first) == _elem_id_to_heat_flux.end())
123  _elem_id_to_heat_flux[it.first] = it.second;
124  else
125  {
126  auto & existing = _elem_id_to_heat_flux[it.first];
127  for (const auto qp : index_range(existing))
128  existing[qp] += it.second[qp];
129  }
130 }
131 
132 void
134 {
136 }
137 
138 const std::vector<ADReal> &
140 {
141  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
142 
143  auto it = _elem_id_to_heat_flux.find(elem_id);
144  if (it != _elem_id_to_heat_flux.end())
145  return it->second;
146  else
147  mooseError(name(), ": Requested heat flux for element ", elem_id, " was not computed.");
148 }
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()
registerMooseObject("ThermalHydraulicsApp", HSCoupler2D2DRadiationUserObject)
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const std::vector< ADReal > & getHeatFlux(dof_id_type elem_id) const
Gets the heat fluxes for each quadrature point for a given element ID.
const bool _include_environment
Whether or not to include an environment surrounding all of the surfaces.
std::map< dof_id_type, std::vector< ADReal > > _elem_id_to_heat_flux
Map of the element ID to the heat flux.
const std::vector< Real > & _emissivities
Emissivities of each boundary.
Stores variable values at each quadrature point on a side by element ID.
const Parallel::Communicator & comm() const
virtual void threadJoin(const UserObject &uo) override
void addRequiredParam(const std::string &name, const std::string &doc_string)
const std::vector< std::vector< Real > > & _view_factors
View factors between each surface.
const std::string & name() const
Builds mapping between multiple 2D boundaries.
const unsigned int _n_hs
Number of heat structures.
HSCoupler2D2DRadiationUserObject(const InputParameters &parameters)
const StoreVariableByElemIDSideUserObject & _temperature_uo
User object containing the temperature values on the boundary.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
const Real & _T_environment
Environment temperature.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const QBase *const & _qrule
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
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
const unsigned int _n_surfaces
Number of surfaces.
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const std::vector< dof_id_type > & getCoupledSecondaryElemIDs(const dof_id_type &primary_elem_id) const
Gets the coupled secondary element IDs for a given primary element ID.
const Elem *const & _current_elem
Computes heat fluxes for HSCoupler2D2DRadiation.
MooseUnits pow(const MooseUnits &, int)
for(PetscInt i=0;i< nvars;++i)
auto index_range(const T &sizable)
const MeshAlignment2D2D & _mesh_alignment
Mesh alignment object.
uint8_t dof_id_type