https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
LinearFVConvectiveHeatTransferBC Class Reference

Class describing a convective heat transfer between two domains. More...

#include <LinearFVConvectiveHeatTransferBC.h>

Inheritance diagram for LinearFVConvectiveHeatTransferBC:
[legend]

Public Types

enum  TEST_TYPE
 
enum  ResidualTagType { ResidualTagType::NonReference, ResidualTagType::Reference }
 
typedef DataFileName DataFileParameterType
 

Public Member Functions

 LinearFVConvectiveHeatTransferBC (const InputParameters &parameters)
 Class constructor. More...
 
virtual Real computeBoundaryValue () const override
 
virtual Real computeBoundaryNormalGradient () const override
 
virtual Real computeBoundaryValueMatrixContribution () const override
 
virtual Real computeBoundaryValueRHSContribution () const override
 
virtual Real computeBoundaryGradientMatrixContribution () const override
 
virtual Real computeBoundaryGradientRHSContribution () const override
 
virtual bool includesMaterialPropertyMultiplier () const override
 
virtual bool useBoundaryGradientExtrapolation () const
 
virtual bool hasFaceSide (const FaceInfo &fi, bool fi_elem_side) const override
 
const SubProblemsubProblem () const
 
const MooseLinearVariableFV< Real > & variable () const
 
void setupFaceData (const FaceInfo *face_info, const FaceInfo::VarFaceNeighbors face_type)
 
const FaceInfocurrentFaceInfo () const
 
FaceInfo::VarFaceNeighbors currentFaceType () const
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
virtual const std::set< BoundaryID > & boundaryIDs () const
 
const std::vector< BoundaryName > & boundaryNames () const
 
unsigned int numBoundaryIDs () const
 
bool hasBoundary (const BoundaryName &name) const
 
bool hasBoundary (const std::vector< BoundaryName > &names) const
 
bool hasBoundary (const BoundaryID &id) const
 
bool hasBoundary (const std::vector< BoundaryID > &ids, TEST_TYPE type=ALL) const
 
bool hasBoundary (const std::set< BoundaryID > &ids, TEST_TYPE type=ALL) const
 
bool isBoundarySubset (const std::set< BoundaryID > &ids) const
 
bool isBoundarySubset (const std::vector< BoundaryID > &ids) const
 
bool hasBoundaryMaterialProperty (const std::string &prop_name) const
 
virtual bool boundaryRestricted () const
 
const std::set< BoundaryID > & meshBoundaryIDs () const
 
virtual bool checkVariableBoundaryIntegrity () const
 
virtual void initialSetup ()
 
virtual void timestepSetup ()
 
virtual void jacobianSetup ()
 
virtual void residualSetup ()
 
virtual void subdomainSetup ()
 
virtual void customSetup (const ExecFlagType &)
 
const ExecFlagEnumgetExecuteOnEnum () const
 
const FunctiongetFunction (const std::string &name) const
 
const FunctiongetFunctionByName (const FunctionName &name) const
 
bool hasFunction (const std::string &param_name) const
 
bool hasFunctionByName (const FunctionName &name) const
 
UserObjectName getUserObjectName (const std::string &param_name) const
 
const T & getUserObject (const std::string &param_name, bool is_dependency=true) const
 
const T & getUserObjectByName (const UserObjectName &object_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBase (const std::string &param_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBaseByName (const UserObjectName &object_name, bool is_dependency=true) const
 
bool isImplicit ()
 
Moose::StateArg determineState () const
 
bool isDefaultPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessor (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessorByName (const PostprocessorName &name) const
 
std::size_t coupledPostprocessors (const std::string &param_name) const
 
const PostprocessorName & getPostprocessorName (const std::string &param_name, const unsigned int index=0) const
 
const VectorPostprocessorValuegetVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name, bool needs_broadcast) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
bool hasVectorPostprocessor (const std::string &param_name, const std::string &vector_name) const
 
bool hasVectorPostprocessor (const std::string &param_name) const
 
bool hasVectorPostprocessorByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
bool hasVectorPostprocessorByName (const VectorPostprocessorName &name) const
 
const VectorPostprocessorName & getVectorPostprocessorName (const std::string &param_name) const
 
PenetrationLocatorgetPenetrationLocator (const BoundaryName &primary, const BoundaryName &secondary, Order order)
 
PenetrationLocatorgetQuadraturePenetrationLocator (const BoundaryName &primary, const BoundaryName &secondary, Order order)
 
NearestNodeLocatorgetNearestNodeLocator (const BoundaryName &primary, const BoundaryName &secondary)
 
NearestNodeLocatorgetQuadratureNearestNodeLocator (const BoundaryName &primary, const BoundaryName &secondary)
 
bool requiresGeometricSearch () const
 
virtual void meshChanged ()
 
void useVectorTag (const TagName &tag_name, VectorTagsKey)
 
void useVectorTag (TagID tag_id, VectorTagsKey)
 
void useMatrixTag (const TagName &tag_name, MatrixTagsKey)
 
void useMatrixTag (TagID tag_id, MatrixTagsKey)
 
bool isVectorTagged ()
 
bool isMatrixTagged ()
 
bool hasVectorTags () const
 
const std::set< TagID > & getVectorTags (VectorTagsKey) const
 
const std::set< TagID > & getMatrixTags (MatrixTagsKey) const
 
MooseVariableBasemooseVariableBase () const
 
MooseVariableField< Real > & mooseVariableField ()
 
MooseVariableFE< Real > * mooseVariable () const
 
MooseVariableFV< Real > * mooseVariableFV () const
 
MooseLinearVariableFV< Real > * mooseLinearVariableFV () const
 
const std::set< MooseVariableFieldBase *> & getMooseVariableDependencies () const
 
std::set< MooseVariableFieldBase *> checkAllVariables (const DofObjectType &dof_object, const std::set< MooseVariableFieldBase * > &vars_to_omit={})
 
std::set< MooseVariableFieldBase *> checkVariables (const DofObjectType &dof_object, const std::set< MooseVariableFieldBase * > &vars_to_check)
 
void addMooseVariableDependency (MooseVariableFieldBase *var)
 
void addMooseVariableDependency (const std::vector< MooseVariableFieldBase * > &vars)
 
Moose::FaceArg makeFace (const FaceInfo &fi, const Moose::FV::LimiterType limiter_type, const bool elem_is_upwind, const bool correct_skewness=false, const Moose::StateArg *state_limiter=nullptr) const
 
Moose::FaceArg makeCDFace (const FaceInfo &fi, const bool correct_skewness=false) const
 
const DistributiongetDistribution (const std::string &name) const
 
const T & getDistribution (const std::string &name) const
 
const DistributiongetDistribution (const std::string &name) const
 
const T & getDistribution (const std::string &name) const
 
const DistributiongetDistributionByName (const DistributionName &name) const
 
const T & getDistributionByName (const std::string &name) const
 
const DistributiongetDistributionByName (const DistributionName &name) const
 
const T & getDistributionByName (const std::string &name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
bool isVectorPostprocessorDistributed (const std::string &param_name) const
 
bool isVectorPostprocessorDistributed (const std::string &param_name) const
 
bool isVectorPostprocessorDistributedByName (const VectorPostprocessorName &name) const
 
bool isVectorPostprocessorDistributedByName (const VectorPostprocessorName &name) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 
static bool restricted (const std::set< BoundaryID > &ids)
 
static std::string deduceFunctorName (const std::string &name, const InputParameters &params)
 

Public Attributes

 ALL
 
 ANY
 
const ConsoleStream _console
 

Protected Member Functions

std::string deduceFunctorName (const std::string &name) const
 
Real computeCellToFaceDistance () const
 
RealVectorValue computeCellToFaceVector () const
 
Moose::FaceArg singleSidedFaceArg (const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
 
bool hasBoundaryMaterialPropertyHelper (const std::string &prop_name) const
 
virtual void addUserObjectDependencyHelper (const UserObject &) const
 
virtual void addPostprocessorDependencyHelper (const PostprocessorName &) const
 
virtual void addVectorPostprocessorDependencyHelper (const VectorPostprocessorName &) const
 
void prepareVectorTag (Assembly &assembly, unsigned int ivar)
 
void prepareVectorTag (Assembly &assembly, unsigned int ivar, ResidualTagType tag_type)
 
void prepareVectorTag (Assembly &assembly, unsigned int ivar, ResidualTagType tag_type)
 
void prepareVectorTag (Assembly &assembly, unsigned int ivar, ResidualTagType tag_type)
 
void prepareVectorTagNeighbor (Assembly &assembly, unsigned int ivar)
 
void prepareVectorTagLower (Assembly &assembly, unsigned int ivar)
 
void prepareMatrixTag (Assembly &assembly, unsigned int ivar, unsigned int jvar)
 
void prepareMatrixTag (Assembly &assembly, unsigned int ivar, unsigned int jvar, DenseMatrix< Number > &k) const
 
void prepareMatrixTagNonlocal (Assembly &assembly, unsigned int ivar, unsigned int jvar)
 
void prepareMatrixTagNeighbor (Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::DGJacobianType type)
 
void prepareMatrixTagNeighbor (Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::DGJacobianType type, DenseMatrix< Number > &k) const
 
void prepareMatrixTagLower (Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::ConstraintJacobianType type)
 
void accumulateTaggedLocalResidual ()
 
void assignTaggedLocalResidual ()
 
void accumulateTaggedLocalMatrix ()
 
void accumulateTaggedLocalMatrix (Assembly &assembly, unsigned int ivar, unsigned int jvar, const DenseMatrix< Number > &k)
 
void accumulateTaggedLocalMatrix (Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::DGJacobianType type, const DenseMatrix< Number > &k)
 
void accumulateTaggedNonlocalMatrix ()
 
void assignTaggedLocalMatrix ()
 
void addResiduals (Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addResiduals (Assembly &assembly, const DenseVector< T > &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addResiduals (Assembly &assembly, const ADResidualsPacket &packet)
 
void addResidualsAndJacobian (Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addResidualsAndJacobian (Assembly &assembly, const ADResidualsPacket &packet)
 
void addJacobian (Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addJacobian (Assembly &assembly, const ADResidualsPacket &packet)
 
void addJacobian (Assembly &assembly, DenseMatrix< Real > &local_k, const std::vector< dof_id_type > &row_indices, const std::vector< dof_id_type > &column_indices, Real scaling_factor)
 
void addResidualsWithoutConstraints (Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addResidualsAndJacobianWithoutConstraints (Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addJacobianWithoutConstraints (Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
 
void addJacobianElement (Assembly &assembly, Real value, dof_id_type row_index, dof_id_type column_index, Real scaling_factor)
 
void setResidual (SystemBase &sys, const T &residual, MooseVariableFE< T > &var)
 
void setResidual (SystemBase &sys, Real residual, dof_id_type dof_index)
 
void setResidual (SystemBase &sys, SetResidualFunctor set_residual_functor)
 
virtual const OutputTools< Real >::VariableValuevalue ()
 
virtual const OutputTools< Real >::VariableValuevalueOld ()
 
virtual const OutputTools< Real >::VariableValuevalueOlder ()
 
virtual const OutputTools< Real >::VariableValuedot ()
 
virtual const OutputTools< Real >::VariableValuedotDot ()
 
virtual const OutputTools< Real >::VariableValuedotOld ()
 
virtual const OutputTools< Real >::VariableValuedotDotOld ()
 
virtual const VariableValuedotDu ()
 
virtual const VariableValuedotDotDu ()
 
virtual const OutputTools< Real >::VariableGradientgradient ()
 
virtual const OutputTools< Real >::VariableGradientgradientOld ()
 
virtual const OutputTools< Real >::VariableGradientgradientOlder ()
 
virtual const OutputTools< Real >::VariableSecondsecond ()
 
virtual const OutputTools< Real >::VariableSecondsecondOld ()
 
virtual const OutputTools< Real >::VariableSecondsecondOlder ()
 
virtual const OutputTools< Real >::VariableTestSecondsecondTest ()
 
virtual const OutputTools< Real >::VariableTestSecondsecondTestFace ()
 
virtual const OutputTools< Real >::VariablePhiSecondsecondPhi ()
 
virtual const OutputTools< Real >::VariablePhiSecondsecondPhiFace ()
 
const Moose::Functor< T > & getFunctor (const std::string &name)
 
const Moose::Functor< T > & getFunctor (const std::string &name, THREAD_ID tid)
 
const Moose::Functor< T > & getFunctor (const std::string &name, SubProblem &subproblem)
 
const Moose::Functor< T > & getFunctor (const std::string &name, SubProblem &subproblem, THREAD_ID tid)
 
bool isFunctor (const std::string &name) const
 
bool isFunctor (const std::string &name, const SubProblem &subproblem) const
 
Moose::ElemArg makeElemArg (const Elem *elem, bool correct_skewnewss=false) const
 
void checkFunctorSupportsSideIntegration (const std::string &name, bool qp_integration)
 

Protected Attributes

const Moose::Functor< Real > & _temp_fluid
 The fluid temperature, we use the functor form to enable situations when the user wants to supply a solution-independent form for this. More...
 
const Moose::Functor< Real > & _temp_solid
 The solid/wall temperature, we use the functor form to enable situations when the user wants to supply a solution-independent form for this. More...
 
const Moose::Functor< Real > & _htc
 The convective heat transfer coefficient. More...
 
bool _var_is_fluid
 Helper boolean to see if the variable we have is the fluid variable. More...
 
const Moose::Functor< Real > * _rhs_temperature
 The temperature which will contribute to the right hand side. More...
 
const THREAD_ID _tid
 
SubProblem_subproblem
 
MooseMesh_mesh
 
FEProblemBase_fv_problem
 
MooseLinearVariableFV< Real > & _var
 
SystemBase_sys
 
const FaceInfo_current_face_info
 
FaceInfo::VarFaceNeighbors _current_face_type
 
const unsigned int _var_num
 
const unsigned int _sys_num
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
const ExecFlagEnum_execute_enum
 
const ExecFlagType_current_execute_flag
 
const InputParameters_ti_params
 
FEProblemBase_ti_feproblem
 
bool _is_implicit
 
Real_t
 
const Real_t_old
 
int_t_step
 
Real_dt
 
Real_dt_old
 
bool _is_transient
 
GeometricSearchData_geometric_search_data
 
bool _requires_geometric_search
 
FEProblemBase_mci_feproblem
 
DenseVector< Number_local_re
 
DenseMatrix< Number_local_ke
 
DenseMatrix< Number_nonlocal_ke
 
bool _nodal
 
MooseVariableFE< Real > * _variable
 
MooseVariableFV< Real > * _fv_variable
 
MooseLinearVariableFV< Real > * _linear_fv_variable
 
MooseVariableField< Real > * _field_variable
 
Assembly_mvi_assembly
 
const Parallel::Communicator & _communicator
 

Detailed Description

Class describing a convective heat transfer between two domains.

The heat flux is described by: h * (T_solid - T_liquid), where h is the heat transfer coefficient, while T_solid and T_liquid denote the solid and liquid temperatures, respectively.

Definition at line 21 of file LinearFVConvectiveHeatTransferBC.h.

Constructor & Destructor Documentation

◆ LinearFVConvectiveHeatTransferBC()

LinearFVConvectiveHeatTransferBC::LinearFVConvectiveHeatTransferBC ( const InputParameters parameters)

Class constructor.

Parameters
parametersThe InputParameters for the object

Definition at line 26 of file LinearFVConvectiveHeatTransferBC.C.

29  _temp_fluid(getFunctor<Real>(NS::T_fluid)),
30  _temp_solid(getFunctor<Real>(NS::T_solid)),
31  _htc(getFunctor<Real>("h")),
32  _var_is_fluid("wraps_" + _var.name() == _temp_fluid.functorName() ||
33  "wraps_" + _var.name() + "_raw_value" == _temp_fluid.functorName())
34 {
35  // We determine which one is the source variable
36  if (_var_is_fluid)
38  else
40 }
const Moose::Functor< Real > & _temp_solid
The solid/wall temperature, we use the functor form to enable situations when the user wants to suppl...
static const std::string T_solid
Definition: NS.h:107
const Moose::Functor< Real > * _rhs_temperature
The temperature which will contribute to the right hand side.
const std::string & name() const override
LinearFVAdvectionDiffusionBC(const InputParameters &parameters)
static const std::string T_fluid
Definition: NS.h:106
MooseLinearVariableFV< Real > & _var
const Moose::Functor< Real > & _htc
The convective heat transfer coefficient.
const InputParameters & parameters() const
bool _var_is_fluid
Helper boolean to see if the variable we have is the fluid variable.
const Moose::Functor< Real > & _temp_fluid
The fluid temperature, we use the functor form to enable situations when the user wants to supply a s...

Member Function Documentation

◆ computeBoundaryGradientMatrixContribution()

Real LinearFVConvectiveHeatTransferBC::computeBoundaryGradientMatrixContribution ( ) const
overridevirtual

Implements LinearFVAdvectionDiffusionBC.

Definition at line 97 of file LinearFVConvectiveHeatTransferBC.C.

98 {
99  const auto face = singleSidedFaceArg(_current_face_info);
100  const auto state = determineState();
101 
102  // We just put the heat transfer coefficient on the diagonal (multiplication with the
103  // surface area is taken care of in the kernel).
104  return _htc(face, state);
105 }
Moose::StateArg determineState() const
const FaceInfo * _current_face_info
const Moose::Functor< Real > & _htc
The convective heat transfer coefficient.
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const

◆ computeBoundaryGradientRHSContribution()

Real LinearFVConvectiveHeatTransferBC::computeBoundaryGradientRHSContribution ( ) const
overridevirtual

Implements LinearFVAdvectionDiffusionBC.

Definition at line 108 of file LinearFVConvectiveHeatTransferBC.C.

109 {
110  const auto state = determineState();
112 
113  // We check where the functor contributing to the right hand side lives. We do this
114  // because this functor lives on the domain where the variable of this kernel doesn't.
115  if (_rhs_temperature->hasFaceSide(*_current_face_info, true))
116  face.face_side = _current_face_info->elemPtr();
117  else
118  face.face_side = _current_face_info->neighborPtr();
119 
120  return _htc(face, state) * (*_rhs_temperature)(face, state);
121 }
Moose::StateArg determineState() const
const Moose::Functor< Real > * _rhs_temperature
The temperature which will contribute to the right hand side.
const Elem * neighborPtr() const
const Elem * elemPtr() const
const FaceInfo * _current_face_info
const Moose::Functor< Real > & _htc
The convective heat transfer coefficient.
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const

◆ computeBoundaryNormalGradient()

Real LinearFVConvectiveHeatTransferBC::computeBoundaryNormalGradient ( ) const
overridevirtual

Implements LinearFVAdvectionDiffusionBC.

Definition at line 53 of file LinearFVConvectiveHeatTransferBC.C.

54 {
55  const auto face = singleSidedFaceArg(_current_face_info);
56  const auto state = determineState();
57 
58  const auto elem_info = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM)
61 
62  const auto neighbor_info = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM)
65 
66  const auto fluid_side_elem_info = _var_is_fluid ? elem_info : neighbor_info;
67 
68  // All this fuss is just for cases when we have an internal boundary, then the flux will change
69  // signs depending on which side of the face we are at.
70  const auto multiplier = _current_face_info->normal() * (_current_face_info->faceCentroid() -
71  fluid_side_elem_info->centroid()) >
72  0
73  ? 1
74  : -1;
75 
76  return multiplier * _htc(face, state) * (_temp_fluid(face, state) - _temp_solid(face, state));
77 }
const Moose::Functor< Real > & _temp_solid
The solid/wall temperature, we use the functor form to enable situations when the user wants to suppl...
Moose::StateArg determineState() const
const ElemInfo * neighborInfo() const
const Point & faceCentroid() const
const ElemInfo * elemInfo() const
FaceInfo::VarFaceNeighbors _current_face_type
const Point & normal() const
const FaceInfo * _current_face_info
const Moose::Functor< Real > & _htc
The convective heat transfer coefficient.
bool _var_is_fluid
Helper boolean to see if the variable we have is the fluid variable.
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
const Moose::Functor< Real > & _temp_fluid
The fluid temperature, we use the functor form to enable situations when the user wants to supply a s...

◆ computeBoundaryValue()

Real LinearFVConvectiveHeatTransferBC::computeBoundaryValue ( ) const
overridevirtual

Implements LinearFVAdvectionDiffusionBC.

Definition at line 43 of file LinearFVConvectiveHeatTransferBC.C.

44 {
45  const auto elem_info = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM)
48 
49  return _var.getElemValue(*elem_info, determineState());
50 }
Moose::StateArg determineState() const
const ElemInfo * neighborInfo() const
const ElemInfo * elemInfo() const
FaceInfo::VarFaceNeighbors _current_face_type
MooseLinearVariableFV< Real > & _var
const FaceInfo * _current_face_info
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const

◆ computeBoundaryValueMatrixContribution()

Real LinearFVConvectiveHeatTransferBC::computeBoundaryValueMatrixContribution ( ) const
overridevirtual

Implements LinearFVAdvectionDiffusionBC.

Definition at line 80 of file LinearFVConvectiveHeatTransferBC.C.

81 {
82  // We approximate the face value with the cell value here.
83  // TODO: we can extend this to a 2-term expansion at some point when the need arises.
84  return 1.0;
85 }

◆ computeBoundaryValueRHSContribution()

Real LinearFVConvectiveHeatTransferBC::computeBoundaryValueRHSContribution ( ) const
overridevirtual

Implements LinearFVAdvectionDiffusionBC.

Definition at line 88 of file LinearFVConvectiveHeatTransferBC.C.

89 {
90  // We approximate the face value with the cell value, we
91  // don't need to add anything to the right hand side.
92  // TODO: we can extend this to a 2-term expansion at some point when the need arises.
93  return 0.0;
94 }

◆ includesMaterialPropertyMultiplier()

virtual bool LinearFVConvectiveHeatTransferBC::includesMaterialPropertyMultiplier ( ) const
inlineoverridevirtual

Reimplemented from LinearFVAdvectionDiffusionBC.

Definition at line 44 of file LinearFVConvectiveHeatTransferBC.h.

44 { return true; }

◆ validParams()

InputParameters LinearFVConvectiveHeatTransferBC::validParams ( )
static

Definition at line 16 of file LinearFVConvectiveHeatTransferBC.C.

17 {
19  params.addRequiredParam<MooseFunctorName>(NS::T_fluid, "The fluid temperature variable");
20  params.addRequiredParam<MooseFunctorName>(NS::T_solid, "The solid/wall temperature variable");
21  params.addRequiredParam<MooseFunctorName>("h", "The convective heat transfer coefficient");
22  params.addClassDescription("Class describing a convective heat transfer between two domains.");
23  return params;
24 }
static const std::string T_solid
Definition: NS.h:107
void addRequiredParam(const std::string &name, const std::string &doc_string)
static const std::string T_fluid
Definition: NS.h:106
void addClassDescription(const std::string &doc_string)
static InputParameters validParams()

Member Data Documentation

◆ _htc

const Moose::Functor<Real>& LinearFVConvectiveHeatTransferBC::_htc
protected

◆ _rhs_temperature

const Moose::Functor<Real>* LinearFVConvectiveHeatTransferBC::_rhs_temperature
protected

The temperature which will contribute to the right hand side.

When this is the fluid domain, the solid temperature will go to the right hand side. When it is the solid domain, the fluid temperature will contribute to the right hand side.

Definition at line 65 of file LinearFVConvectiveHeatTransferBC.h.

Referenced by computeBoundaryGradientRHSContribution(), and LinearFVConvectiveHeatTransferBC().

◆ _temp_fluid

const Moose::Functor<Real>& LinearFVConvectiveHeatTransferBC::_temp_fluid
protected

The fluid temperature, we use the functor form to enable situations when the user wants to supply a solution-independent form for this.

Definition at line 49 of file LinearFVConvectiveHeatTransferBC.h.

Referenced by computeBoundaryNormalGradient(), and LinearFVConvectiveHeatTransferBC().

◆ _temp_solid

const Moose::Functor<Real>& LinearFVConvectiveHeatTransferBC::_temp_solid
protected

The solid/wall temperature, we use the functor form to enable situations when the user wants to supply a solution-independent form for this.

Definition at line 53 of file LinearFVConvectiveHeatTransferBC.h.

Referenced by computeBoundaryNormalGradient(), and LinearFVConvectiveHeatTransferBC().

◆ _var_is_fluid

bool LinearFVConvectiveHeatTransferBC::_var_is_fluid
protected

Helper boolean to see if the variable we have is the fluid variable.

Definition at line 59 of file LinearFVConvectiveHeatTransferBC.h.

Referenced by computeBoundaryNormalGradient(), and LinearFVConvectiveHeatTransferBC().


The documentation for this class was generated from the following files: