Generic gap heat transfer model, with h_gap = h_conduction + h_contact + h_radiation.
More...
#include <GapHeatTransfer.h>
Generic gap heat transfer model, with h_gap = h_conduction + h_contact + h_radiation.
Definition at line 24 of file GapHeatTransfer.h.
◆ GapHeatTransfer()
GapHeatTransfer::GapHeatTransfer |
( |
const InputParameters & |
parameters | ) |
|
Definition at line 90 of file GapHeatTransfer.C.
91 : IntegratedBC(parameters),
92 _gap_geometry_type(declareRestartableData<GapConductance::GAP_GEOMETRY>(
"gap_geometry_type",
97 getParam<std::string>(
"appended_property_name"))),
99 "gap_conductance" + getParam<std::string>(
"appended_property_name") +
"_dT")),
100 _min_gap(getParam<Real>(
"min_gap")),
102 _max_gap(getParam<Real>(
"max_gap")),
112 : &getQuadraturePenetrationLocator(
113 parameters.get<BoundaryName>(
"paired_boundary"),
114 getParam<std::vector<BoundaryName>>(
"boundary")[0],
115 Utility::string_to_enum<Order>(parameters.get<MooseEnum>(
"order")))),
117 _p1(declareRestartableData<Point>(
"cylinder_axis_point_1", Point(0, 1, 0))),
118 _p2(declareRestartableData<Point>(
"cylinder_axis_point_2", Point(0, 0, 0)))
120 if (isParamValid(
"displacements"))
123 for (
unsigned int i = 0; i < coupledComponents(
"displacements"); ++i)
129 if (isParamValid(
"disp_x"))
131 if (isParamValid(
"disp_y"))
133 if (isParamValid(
"disp_z"))
141 if (!parameters.isParamValid(
"paired_boundary"))
142 mooseError(std::string(
"No 'paired_boundary' provided for ") + _name);
146 if (!isCoupled(
"gap_distance"))
147 mooseError(std::string(
"No 'gap_distance' provided for ") + _name);
149 if (!isCoupled(
"gap_temp"))
150 mooseError(std::string(
"No 'gap_temp' provided for ") + _name);
◆ computeGapValues()
void GapHeatTransfer::computeGapValues |
( |
| ) |
|
|
protectedvirtual |
Definition at line 288 of file GapHeatTransfer.C.
298 Node * qnode = _mesh.getQuadratureNode(_current_elem, _current_side, _qp);
311 const Elem * slave_side = pinfo->_side;
312 std::vector<std::vector<Real>> & slave_side_phi = pinfo->_side_phi;
313 _gap_temp = _variable->getValue(slave_side, slave_side_phi);
316 if (tangential_tolerance != 0.0)
318 _edge_multiplier = 1.0 - pinfo->_tangential_distance / tangential_tolerance;
326 mooseWarning(
"No gap value information found for node ",
Referenced by computeQpJacobian(), computeQpOffDiagJacobian(), and computeQpResidual().
◆ computeQpJacobian()
Real GapHeatTransfer::computeQpJacobian |
( |
| ) |
|
|
overrideprotectedvirtual |
◆ computeQpOffDiagJacobian()
Real GapHeatTransfer::computeQpOffDiagJacobian |
( |
unsigned |
jvar | ) |
|
|
overrideprotectedvirtual |
Definition at line 208 of file GapHeatTransfer.C.
215 unsigned int coupled_component;
217 for (coupled_component = 0; coupled_component <
_disp_vars.size(); ++coupled_component)
257 const Point & normal(_normals[_qp]);
259 const Real dgap =
dgapLength(-normal(coupled_component));
263 return _test[_i][_qp] * dRdx * _phi[_j][_qp];
◆ computeQpResidual()
Real GapHeatTransfer::computeQpResidual |
( |
| ) |
|
|
overrideprotectedvirtual |
Definition at line 166 of file GapHeatTransfer.C.
179 Threads::spin_mutex::scoped_lock lock(Threads::spin_mutex);
181 _slave_flux->add(_var.dofIndices()[_i], slave_flux);
184 return _test[_i][_qp] * grad_t;
◆ computeSlaveFluxContribution()
Real GapHeatTransfer::computeSlaveFluxContribution |
( |
Real |
grad_t | ) |
|
|
protectedvirtual |
◆ dgapLength()
Real GapHeatTransfer::dgapLength |
( |
Real |
normalComponent | ) |
const |
|
protectedvirtual |
◆ gapLength()
Real GapHeatTransfer::gapLength |
( |
| ) |
const |
|
protectedvirtual |
◆ initialSetup()
void GapHeatTransfer::initialSetup |
( |
| ) |
|
|
overridevirtual |
Definition at line 155 of file GapHeatTransfer.C.
158 _assembly.coordSystem(),
159 _fe_problem.getAxisymmetricRadialCoord(),
◆ validParams()
InputParameters GapHeatTransfer::validParams |
( |
| ) |
|
|
static |
Definition at line 27 of file GapHeatTransfer.C.
30 params.addClassDescription(
"Transfers heat across a gap between two "
31 "surfaces dependent on the gap geometry specified.");
32 params.addParam<std::string>(
33 "appended_property_name",
"",
"Name appended to material properties to make them unique");
36 params.addParam<Real>(
"min_gap", 1.0e-6,
"A minimum gap size");
37 params.addParam<Real>(
"max_gap", 1.0e6,
"A maximum gap size");
38 params.addRangeCheckedParam<
unsigned int>(
39 "min_gap_order", 0,
"min_gap_order<=1",
"Order of the Taylor expansion below min_gap");
42 MooseEnum coord_types(
"default XYZ cyl",
"default");
43 params.addDeprecatedParam<MooseEnum>(
46 "Gap calculation type (default or XYZ).",
47 "The functionality of this parameter is replaced by 'gap_geometry_type'.");
49 MooseEnum gap_geom_types(
"PLATE CYLINDER SPHERE");
50 params.addParam<MooseEnum>(
"gap_geometry_type",
52 "Gap calculation type. Choices are: " + gap_geom_types.getRawNames());
54 params.addParam<RealVectorValue>(
"cylinder_axis_point_1",
55 "Start point for line defining cylindrical axis");
56 params.addParam<RealVectorValue>(
"cylinder_axis_point_2",
57 "End point for line defining cylindrical axis");
58 params.addParam<RealVectorValue>(
"sphere_origin",
"Origin for sphere geometry");
61 params.addParam<
bool>(
"quadrature",
63 "Whether or not to do Quadrature point based gap heat "
64 "transfer. If this is true then gap_distance and "
65 "gap_temp should NOT be provided (and will be "
66 "ignored) however paired_boundary IS then required.");
67 params.addParam<BoundaryName>(
"paired_boundary",
"The boundary to be penetrated");
69 MooseEnum orders(AddVariableAction::getNonlinearVariableOrders());
70 params.addParam<MooseEnum>(
"order", orders,
"The finite element order");
72 params.addParam<
bool>(
73 "warnings",
false,
"Whether to output warning messages concerning nodes not being found");
75 params.addCoupledVar(
"disp_x",
"The x displacement");
76 params.addCoupledVar(
"disp_y",
"The y displacement");
77 params.addCoupledVar(
"disp_z",
"The z displacement");
81 "The displacements appropriate for the simulation geometry and coordinate system");
84 params.addCoupledVar(
"gap_distance",
"Distance across the gap");
85 params.addCoupledVar(
"gap_temp",
"Temperature on the other side of the gap");
◆ _disp_vars
std::vector<unsigned int> GapHeatTransfer::_disp_vars |
|
protected |
◆ _edge_multiplier
Real GapHeatTransfer::_edge_multiplier |
|
protected |
◆ _gap_conductance
const MaterialProperty<Real>& GapHeatTransfer::_gap_conductance |
|
protected |
◆ _gap_conductance_dT
const MaterialProperty<Real>& GapHeatTransfer::_gap_conductance_dT |
|
protected |
◆ _gap_distance
Real GapHeatTransfer::_gap_distance |
|
protected |
◆ _gap_distance_value
const VariableValue& GapHeatTransfer::_gap_distance_value |
|
protected |
◆ _gap_geometry_type
◆ _gap_temp
Real GapHeatTransfer::_gap_temp |
|
protected |
◆ _gap_temp_value
const VariableValue& GapHeatTransfer::_gap_temp_value |
|
protected |
◆ _has_info
bool GapHeatTransfer::_has_info |
|
protected |
◆ _max_gap
const Real GapHeatTransfer::_max_gap |
|
protected |
◆ _min_gap
const Real GapHeatTransfer::_min_gap |
|
protected |
◆ _min_gap_order
const unsigned int GapHeatTransfer::_min_gap_order |
|
protected |
◆ _p1
Point& GapHeatTransfer::_p1 |
|
protected |
◆ _p2
Point& GapHeatTransfer::_p2 |
|
protected |
◆ _penetration_locator
PenetrationLocator* GapHeatTransfer::_penetration_locator |
|
protected |
◆ _quadrature
const bool GapHeatTransfer::_quadrature |
|
protected |
◆ _r1
Real GapHeatTransfer::_r1 |
|
protected |
◆ _r2
Real GapHeatTransfer::_r2 |
|
protected |
◆ _radius
Real GapHeatTransfer::_radius |
|
protected |
◆ _slave_flux
NumericVector<Number>* GapHeatTransfer::_slave_flux |
|
protected |
◆ _warnings
const bool GapHeatTransfer::_warnings |
|
protected |
The documentation for this class was generated from the following files:
static void setGapGeometryParameters(const InputParameters ¶ms, const Moose::CoordinateSystemType coord_sys, unsigned int axisymmetric_radial_coord, GAP_GEOMETRY &gap_geometry_type, Point &p1, Point &p2)
static void computeGapRadii(const GAP_GEOMETRY gap_geometry_type, const Point ¤t_point, const Point &p1, const Point &p2, const Real &gap_distance, const Point ¤t_normal, Real &r1, Real &r2, Real &radius)