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

Kernel that adds contributions from a diffusion term of the turbulent variables, limited in the near-wall regions, discretized using the finite volume method to a linear system. More...

#include <LinearFVTurbulentDiffusion.h>

Inheritance diagram for LinearFVTurbulentDiffusion:
[legend]

Public Types

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

Public Member Functions

 LinearFVTurbulentDiffusion (const InputParameters &params)
 Class constructor. More...
 
virtual void initialSetup () override
 
virtual Real computeElemMatrixContribution () override
 
virtual Real computeNeighborMatrixContribution () override
 
virtual Real computeElemRightHandSideContribution () override
 
virtual Real computeNeighborRightHandSideContribution () override
 
virtual Real computeBoundaryMatrixContribution (const LinearFVBoundaryCondition &bc) override
 
virtual Real computeBoundaryRHSContribution (const LinearFVBoundaryCondition &bc) override
 
virtual void addMatrixContribution () override
 
virtual void addRightHandSideContribution () override
 
virtual bool hasFaceSide (const FaceInfo &fi, bool fi_elem_side) const override
 
virtual void setupFaceData (const FaceInfo *face_info)
 
void setCurrentFaceArea (const Real area)
 
virtual const MooseLinearVariableFV< Real > & variable () const override
 
void linkTaggedVectorsAndMatrices (const std::set< TagID > &vector_tags, const std::set< TagID > &matrix_tags)
 
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 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
 
void setRandomResetFrequency (ExecFlagType exec_flag)
 
unsigned long getRandomLong () const
 
Real getRandomReal () const
 
unsigned int getSeed (std::size_t id)
 
unsigned int getMasterSeed () const
 
bool isNodal () const
 
ExecFlagType getResetOnTime () const
 
void setRandomDataPointer (RandomData *random_data)
 
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
 
const std::vector< SubdomainName > & blocks () const
 
unsigned int numBlocks () const
 
virtual const std::set< SubdomainID > & blockIDs () const
 
unsigned int blocksMaxDimension () const
 
bool hasBlocks (const SubdomainName &name) const
 
bool hasBlocks (const std::vector< SubdomainName > &names) const
 
bool hasBlocks (const std::set< SubdomainName > &names) const
 
bool hasBlocks (SubdomainID id) const
 
bool hasBlocks (const std::vector< SubdomainID > &ids) const
 
bool hasBlocks (const std::set< SubdomainID > &ids) const
 
bool isBlockSubset (const std::set< SubdomainID > &ids) const
 
bool isBlockSubset (const std::vector< SubdomainID > &ids) const
 
bool hasBlockMaterialProperty (const std::string &prop_name)
 
const std::set< SubdomainID > & meshBlockIDs () const
 
virtual bool blockRestricted () const
 
virtual void checkVariable (const MooseVariableFieldBase &variable) 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
 
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 std::string deduceFunctorName (const std::string &name, const InputParameters &params)
 
static void setRMParamsAdvection (const InputParameters &obj_params, InputParameters &rm_params, const unsigned short conditional_extended_layers)
 
static void setRMParamsDiffusion (const InputParameters &obj_params, InputParameters &rm_params, const unsigned short conditional_extended_layers)
 

Public Attributes

const ConsoleStream _console
 

Protected Member Functions

std::string deduceFunctorName (const std::string &name) const
 
Real computeFluxMatrixContribution ()
 
Real computeFluxRHSContribution ()
 
Moose::FaceArg singleSidedFaceArg (const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
 
void requestVariableCellGradient (const std::string &variable_name)
 
virtual void addUserObjectDependencyHelper (const UserObject &) const
 
virtual void addPostprocessorDependencyHelper (const PostprocessorName &) const
 
virtual void addVectorPostprocessorDependencyHelper (const VectorPostprocessorName &) const
 
T & declareRestartableData (const std::string &data_name, Args &&... args)
 
ManagedValue< T > declareManagedRestartableDataWithContext (const std::string &data_name, void *context, Args &&... args)
 
const T & getRestartableData (const std::string &data_name) const
 
T & declareRestartableDataWithContext (const std::string &data_name, void *context, Args &&... args)
 
T & declareRecoverableData (const std::string &data_name, Args &&... args)
 
T & declareRestartableDataWithObjectName (const std::string &data_name, const std::string &object_name, Args &&... args)
 
T & declareRestartableDataWithObjectNameWithContext (const std::string &data_name, const std::string &object_name, void *context, Args &&... args)
 
std::string restartableName (const std::string &data_name) 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 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 bool hasBlockMaterialPropertyHelper (const std::string &prop_name)
 
void initializeBlockRestrictable (const MooseObject *moose_object)
 
Moose::CoordinateSystemType getBlockCoordSystem ()
 
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)
 
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 ()
 

Protected Attributes

const Moose::Functor< Real > & _scaling_coeff
 The functor for the scaling coefficient for the diffusion term. More...
 
const std::vector< BoundaryName > & _wall_boundary_names
 Wall boundaries. More...
 
std::unordered_set< const Elem * > _wall_bounded
 List for wall bounded elements. More...
 
const Moose::Functor< Real > & _diffusion_coeff
 
const bool _use_nonorthogonal_correction
 
Real _flux_matrix_contribution
 
Real _flux_rhs_contribution
 
const FaceInfo_current_face_info
 
Real _current_face_area
 
FaceInfo::VarFaceNeighbors _current_face_type
 
bool _cached_matrix_contribution
 
bool _cached_rhs_contribution
 
const bool _force_boundary_execution
 
DenseVector< dof_id_type_dof_indices
 
DenseMatrix< Real_matrix_contribution
 
DenseVector< Real_rhs_contribution
 
MooseLinearVariableFV< Real > & _var
 
const unsigned int _var_num
 
const unsigned int _sys_num
 
FEProblemBase_fe_problem
 
SystemBase_sys
 
libMesh::LinearImplicitSystem_linear_system
 
const THREAD_ID _tid
 
MooseMesh_mesh
 
std::vector< NumericVector< Number > *> _vectors
 
std::vector< SparseMatrix< Number > *> _matrices
 
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
 
MooseApp_restartable_app
 
const std::string _restartable_system_name
 
const THREAD_ID _restartable_tid
 
const bool _restartable_read_only
 
FEProblemBase_mci_feproblem
 
SubProblem_subproblem
 
DenseVector< Number_local_re
 
DenseMatrix< Number_local_ke
 
DenseMatrix< Number_nonlocal_ke
 
const MaterialData_blk_material_data
 
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

Kernel that adds contributions from a diffusion term of the turbulent variables, limited in the near-wall regions, discretized using the finite volume method to a linear system.

Definition at line 19 of file LinearFVTurbulentDiffusion.h.

Constructor & Destructor Documentation

◆ LinearFVTurbulentDiffusion()

LinearFVTurbulentDiffusion::LinearFVTurbulentDiffusion ( const InputParameters params)

Class constructor.

Parameters
paramsThe InputParameters for the kernel.

Definition at line 35 of file LinearFVTurbulentDiffusion.C.

36  : LinearFVDiffusion(params),
37  _scaling_coeff(getFunctor<Real>("scaling_coeff")),
38  _wall_boundary_names(getParam<std::vector<BoundaryName>>("walls"))
39 {
42 }
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
MooseLinearVariableFV< Real > & _var
const Moose::Functor< Real > & _scaling_coeff
The functor for the scaling coefficient for the diffusion term.
LinearFVDiffusion(const InputParameters &params)
const T & getParam(const std::string &name) const
const bool _use_nonorthogonal_correction

Member Function Documentation

◆ addMatrixContribution()

void LinearFVTurbulentDiffusion::addMatrixContribution ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 129 of file LinearFVTurbulentDiffusion.C.

130 {
131  // Coumputing bounding map
132  const Elem * elem = _current_face_info->elemPtr();
133  const auto bounded_elem = _wall_bounded.find(elem) != _wall_bounded.end();
134  const Elem * neighbor = _current_face_info->neighborPtr();
135  const auto bounded_neigh = _wall_bounded.find(neighbor) != _wall_bounded.end();
136 
137  // If we are on an internal face, we populate the four entries in the system matrix
138  // which touch the face
140  {
141  // The dof ids of the variable corresponding to the element and neighbor
144 
145  // Compute the entries which will go to the neighbor (offdiagonal) and element
146  // (diagonal).
147  const auto elem_matrix_contribution = computeElemMatrixContribution();
148  const auto neighbor_matrix_contribution = computeNeighborMatrixContribution();
149 
150  // Populate matrix
151  if (hasBlocks(_current_face_info->elemInfo()->subdomain_id()) && !(bounded_elem))
152  {
153  _matrix_contribution(0, 0) = elem_matrix_contribution;
154  _matrix_contribution(0, 1) = neighbor_matrix_contribution;
155  }
156 
157  if (hasBlocks(_current_face_info->neighborInfo()->subdomain_id()) && !(bounded_neigh))
158  {
159  _matrix_contribution(1, 0) = -elem_matrix_contribution;
160  _matrix_contribution(1, 1) = -neighbor_matrix_contribution;
161  }
162  (*_linear_system.matrix).add_matrix(_matrix_contribution, _dof_indices.get_values());
163  }
164  // We are at a block boundary where the variable is not defined on one of the adjacent cells.
165  // We check if we have a boundary condition here
168  {
169  mooseAssert(_current_face_info->boundaryIDs().size() == 1,
170  "We should only have one boundary on every face.");
171 
172  LinearFVBoundaryCondition * bc_pointer =
174 
175  if (bc_pointer || _force_boundary_execution)
176  {
177  if (bc_pointer)
179  const auto matrix_contribution = computeBoundaryMatrixContribution(*bc_pointer);
180 
181  // We allow internal (for the mesh) boundaries too, so we have to check on which side we
182  // are on (assuming that this is a boundary for the variable)
183  if ((_current_face_type == FaceInfo::VarFaceNeighbors::ELEM) && !(bounded_elem))
184  {
185  const auto dof_id_elem = _current_face_info->elemInfo()->dofIndices()[_sys_num][_var_num];
186  (*_linear_system.matrix).add(dof_id_elem, dof_id_elem, matrix_contribution);
187  }
188  else if ((_current_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) && !(bounded_neigh))
189  {
190  const auto dof_id_neighbor =
192  (*_linear_system.matrix).add(dof_id_neighbor, dof_id_neighbor, matrix_contribution);
193  }
194  }
195  }
196 }
const unsigned int _var_num
virtual Real computeBoundaryMatrixContribution(const LinearFVBoundaryCondition &bc) override
const std::set< BoundaryID > & boundaryIDs() const
libMesh::LinearImplicitSystem & _linear_system
const ElemInfo * neighborInfo() const
void setupFaceData(const FaceInfo *face_info, const FaceInfo::VarFaceNeighbors face_type)
MooseLinearVariableFV< Real > & _var
const ElemInfo * elemInfo() const
LinearFVBoundaryCondition * getBoundaryCondition(const BoundaryID bd_id) const
virtual Real computeNeighborMatrixContribution() override
std::unordered_set< const Elem * > _wall_bounded
List for wall bounded elements.
FaceInfo::VarFaceNeighbors _current_face_type
const bool _force_boundary_execution
const FaceInfo * _current_face_info
const Elem * neighborPtr() const
virtual Real computeElemMatrixContribution() override
const Elem * elemPtr() const
DenseVector< dof_id_type > _dof_indices
const std::vector< std::vector< dof_id_type > > & dofIndices() const
SparseMatrix< Number > * matrix
bool hasBlocks(const SubdomainName &name) const
const unsigned int _sys_num
SubdomainID subdomain_id() const
DenseMatrix< Real > _matrix_contribution

◆ addRightHandSideContribution()

void LinearFVTurbulentDiffusion::addRightHandSideContribution ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 199 of file LinearFVTurbulentDiffusion.C.

200 {
201  // Coumputing bounding map
202  const Elem * elem = _current_face_info->elemPtr();
203  const auto bounded_elem = _wall_bounded.find(elem) != _wall_bounded.end();
204  const Elem * neighbor = _current_face_info->neighborPtr();
205  const auto bounded_neigh = _wall_bounded.find(neighbor) != _wall_bounded.end();
206 
207  // If we are on an internal face, we populate the two entries in the right hand side
208  // which touch the face
210  {
211  // The dof ids of the variable corresponding to the element and neighbor
214 
215  // Compute the entries which will go to the neighbor and element positions.
216  const auto elem_rhs_contribution = computeElemRightHandSideContribution();
217  const auto neighbor_rhs_contribution = computeNeighborRightHandSideContribution();
218 
219  // Populate right hand side
221  _rhs_contribution(0) = elem_rhs_contribution;
223  _rhs_contribution(1) = neighbor_rhs_contribution;
224 
226  .add_vector(_rhs_contribution.get_values().data(), _dof_indices.get_values());
227  }
228  // We are at a block boundary where the variable is not defined on one of the adjacent cells.
229  // We check if we have a boundary condition here
232  {
233  mooseAssert(_current_face_info->boundaryIDs().size() == 1,
234  "We should only have one boundary on every face.");
235  LinearFVBoundaryCondition * bc_pointer =
237 
238  if (bc_pointer || _force_boundary_execution)
239  {
240  if (bc_pointer)
242 
243  const auto rhs_contribution = computeBoundaryRHSContribution(*bc_pointer);
244 
245  // We allow internal (for the mesh) boundaries too, so we have to check on which side we
246  // are on (assuming that this is a boundary for the variable)
247  if ((_current_face_type == FaceInfo::VarFaceNeighbors::ELEM) && !(bounded_elem))
248  {
249  const auto dof_id_elem = _current_face_info->elemInfo()->dofIndices()[_sys_num][_var_num];
250  (*_linear_system.rhs).add(dof_id_elem, rhs_contribution);
251  }
252  else if ((_current_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) && !(bounded_neigh))
253  {
254  const auto dof_id_neighbor =
256  (*_linear_system.rhs).add(dof_id_neighbor, rhs_contribution);
257  }
258  }
259  }
260 }
const unsigned int _var_num
const std::set< BoundaryID > & boundaryIDs() const
libMesh::LinearImplicitSystem & _linear_system
const ElemInfo * neighborInfo() const
void setupFaceData(const FaceInfo *face_info, const FaceInfo::VarFaceNeighbors face_type)
MooseLinearVariableFV< Real > & _var
NumericVector< Number > * rhs
const ElemInfo * elemInfo() const
LinearFVBoundaryCondition * getBoundaryCondition(const BoundaryID bd_id) const
virtual Real computeBoundaryRHSContribution(const LinearFVBoundaryCondition &bc) override
std::unordered_set< const Elem * > _wall_bounded
List for wall bounded elements.
FaceInfo::VarFaceNeighbors _current_face_type
const bool _force_boundary_execution
DenseVector< Real > _rhs_contribution
const FaceInfo * _current_face_info
const Elem * neighborPtr() const
virtual Real computeElemRightHandSideContribution() override
const Elem * elemPtr() const
DenseVector< dof_id_type > _dof_indices
const std::vector< std::vector< dof_id_type > > & dofIndices() const
virtual Real computeNeighborRightHandSideContribution() override
bool hasBlocks(const SubdomainName &name) const
const unsigned int _sys_num
SubdomainID subdomain_id() const

◆ computeBoundaryMatrixContribution()

Real LinearFVTurbulentDiffusion::computeBoundaryMatrixContribution ( const LinearFVBoundaryCondition bc)
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 81 of file LinearFVTurbulentDiffusion.C.

Referenced by addMatrixContribution().

82 {
83  const auto * const diff_bc = static_cast<const LinearFVAdvectionDiffusionBC *>(&bc);
84  mooseAssert(diff_bc, "This should be a valid BC!");
85 
86  auto grad_contrib = diff_bc->computeBoundaryGradientMatrixContribution() * _current_face_area;
87  // If the boundary condition does not include the diffusivity contribution then
88  // add it here.
89  if (!diff_bc->includesMaterialPropertyMultiplier())
90  {
91  const auto face_arg = singleSidedFaceArg(_current_face_info);
92  grad_contrib *= _diffusion_coeff(face_arg, determineState());
93  }
94 
95  return grad_contrib;
96 }
const Moose::Functor< Real > & _diffusion_coeff
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
Moose::StateArg determineState() const
const FaceInfo * _current_face_info

◆ computeBoundaryRHSContribution()

Real LinearFVTurbulentDiffusion::computeBoundaryRHSContribution ( const LinearFVBoundaryCondition bc)
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 99 of file LinearFVTurbulentDiffusion.C.

Referenced by addRightHandSideContribution().

100 {
101  const auto * const diff_bc = static_cast<const LinearFVAdvectionDiffusionBC *>(&bc);
102  mooseAssert(diff_bc, "This should be a valid BC!");
103 
104  const auto face_arg = singleSidedFaceArg(_current_face_info);
105  auto grad_contrib = diff_bc->computeBoundaryGradientRHSContribution() * _current_face_area;
106 
107  // If the boundary condition does not include the diffusivity contribution then
108  // add it here.
109  if (!diff_bc->includesMaterialPropertyMultiplier())
110  grad_contrib *= _diffusion_coeff(face_arg, determineState());
111 
112  // We add the nonorthogonal corrector for the face here. Potential idea: we could do
113  // this in the boundary condition too. For now, however, we keep it like this.
115  {
116  const auto correction_vector =
119 
120  grad_contrib += _diffusion_coeff(face_arg, determineState()) *
121  _var.gradSln(*_current_face_info->elemInfo()) * correction_vector *
123  }
124 
125  return grad_contrib;
126 }
const Moose::Functor< Real > & _diffusion_coeff
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
Moose::StateArg determineState() const
MooseLinearVariableFV< Real > & _var
const ElemInfo * elemInfo() const
const FaceInfo * _current_face_info
const Point & normal() const
const Point & eCN() const
const VectorValue< Real > gradSln(const ElemInfo &elem_info) const
const bool _use_nonorthogonal_correction

◆ computeElemMatrixContribution()

Real LinearFVTurbulentDiffusion::computeElemMatrixContribution ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 53 of file LinearFVTurbulentDiffusion.C.

Referenced by addMatrixContribution().

54 {
55  const auto face_arg = makeCDFace(*_current_face_info);
57 }
Moose::StateArg determineState() const
Real computeFluxMatrixContribution()
const Moose::Functor< Real > & _scaling_coeff
The functor for the scaling coefficient for the diffusion term.
const FaceInfo * _current_face_info
Moose::FaceArg makeCDFace(const FaceInfo &fi, const bool correct_skewness=false) const

◆ computeElemRightHandSideContribution()

Real LinearFVTurbulentDiffusion::computeElemRightHandSideContribution ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 67 of file LinearFVTurbulentDiffusion.C.

Referenced by addRightHandSideContribution().

68 {
69  const auto face_arg = makeCDFace(*_current_face_info);
71 }
Moose::StateArg determineState() const
const Moose::Functor< Real > & _scaling_coeff
The functor for the scaling coefficient for the diffusion term.
const FaceInfo * _current_face_info
Real computeFluxRHSContribution()
Moose::FaceArg makeCDFace(const FaceInfo &fi, const bool correct_skewness=false) const

◆ computeNeighborMatrixContribution()

Real LinearFVTurbulentDiffusion::computeNeighborMatrixContribution ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 60 of file LinearFVTurbulentDiffusion.C.

Referenced by addMatrixContribution().

61 {
62  const auto face_arg = makeCDFace(*_current_face_info);
64 }
Moose::StateArg determineState() const
Real computeFluxMatrixContribution()
const Moose::Functor< Real > & _scaling_coeff
The functor for the scaling coefficient for the diffusion term.
const FaceInfo * _current_face_info
Moose::FaceArg makeCDFace(const FaceInfo &fi, const bool correct_skewness=false) const

◆ computeNeighborRightHandSideContribution()

Real LinearFVTurbulentDiffusion::computeNeighborRightHandSideContribution ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 74 of file LinearFVTurbulentDiffusion.C.

Referenced by addRightHandSideContribution().

75 {
76  const auto face_arg = makeCDFace(*_current_face_info);
78 }
Moose::StateArg determineState() const
const Moose::Functor< Real > & _scaling_coeff
The functor for the scaling coefficient for the diffusion term.
const FaceInfo * _current_face_info
Real computeFluxRHSContribution()
Moose::FaceArg makeCDFace(const FaceInfo &fi, const bool correct_skewness=false) const

◆ initialSetup()

void LinearFVTurbulentDiffusion::initialSetup ( )
overridevirtual

Reimplemented from LinearFVDiffusion.

Definition at line 45 of file LinearFVTurbulentDiffusion.C.

46 {
50 }
SubProblem & _subproblem
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
void getWallBoundedElements(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::unordered_set< const Elem *> &wall_bounded)
Map marking wall bounded elements The map passed in wall_bounded_map gets cleared and re-populated...
virtual const std::set< SubdomainID > & blockIDs() const
std::unordered_set< const Elem * > _wall_bounded
List for wall bounded elements.
virtual void initialSetup() override
FEProblemBase & _fe_problem

◆ validParams()

InputParameters LinearFVTurbulentDiffusion::validParams ( )
static

Definition at line 20 of file LinearFVTurbulentDiffusion.C.

21 {
23  params.addClassDescription(
24  "Represents the matrix and right hand side contributions of a "
25  "diffusion term for a turbulent variable in a partial differential equation.");
26 
27  params.addParam<MooseFunctorName>(
28  "scaling_coeff", 1.0, "The scaling coefficient for the diffusion term.");
29 
30  params.addParam<std::vector<BoundaryName>>(
31  "walls", {}, "Boundaries that correspond to solid walls.");
32  return params;
33 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _scaling_coeff

const Moose::Functor<Real>& LinearFVTurbulentDiffusion::_scaling_coeff
protected

The functor for the scaling coefficient for the diffusion term.

Definition at line 50 of file LinearFVTurbulentDiffusion.h.

Referenced by computeElemMatrixContribution(), computeElemRightHandSideContribution(), computeNeighborMatrixContribution(), and computeNeighborRightHandSideContribution().

◆ _wall_boundary_names

const std::vector<BoundaryName>& LinearFVTurbulentDiffusion::_wall_boundary_names
protected

Wall boundaries.

Definition at line 53 of file LinearFVTurbulentDiffusion.h.

Referenced by initialSetup().

◆ _wall_bounded

std::unordered_set<const Elem *> LinearFVTurbulentDiffusion::_wall_bounded
protected

List for wall bounded elements.

Definition at line 56 of file LinearFVTurbulentDiffusion.h.

Referenced by addMatrixContribution(), addRightHandSideContribution(), and initialSetup().


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