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

Kernel that adds contributions to the source and the sink of the turbulent kinetic energy dissipation discretized using the finite volume method to a linear system. More...

#include <LinearFVTKEDSourceSink.h>

Inheritance diagram for LinearFVTKEDSourceSink:
[legend]

Public Types

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

Public Member Functions

 LinearFVTKEDSourceSink (const InputParameters &params)
 Class constructor. More...
 
virtual void initialSetup () override
 
virtual Real computeMatrixContribution () override
 
virtual Real computeRightHandSideContribution () override
 
virtual void addMatrixContribution () override
 
virtual void addRightHandSideContribution () override
 
virtual void setCurrentElemInfo (const ElemInfo *elem_info)
 
void setCurrentElemVolume (const Real volume)
 
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)
 
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
 
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 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 unsigned int _dim
 The dimension of the simulation. More...
 
const Moose::Functor< Real > & _u_var
 x-velocity More...
 
const Moose::Functor< Real > * _v_var
 y-velocity More...
 
const Moose::Functor< Real > * _w_var
 z-velocity More...
 
const Moose::Functor< Real > & _k
 Turbulent kinetic energy. More...
 
const Moose::Functor< Real > & _rho
 Density. More...
 
const Moose::Functor< Real > & _mu
 Dynamic viscosity. More...
 
const Moose::Functor< Real > & _mu_t
 Turbulent dynamic viscosity. More...
 
const std::vector< BoundaryName > & _wall_boundary_names
 Wall boundaries. More...
 
const bool _linearized_model
 If the user wants to use the linearized model. More...
 
NS::WallTreatmentEnum _wall_treatment
 Method used for wall treatment. More...
 
const Real _C1_eps
 Value of the first epsilon closure coefficient. More...
 
const Real _C2_eps
 Value of the second epsilon closure coefficient. More...
 
const Real _C_mu
 C_mu constant. More...
 
const Real _C_pl
 Production Limiter Constant. More...
 
const ElemInfo_current_elem_info
 
Real _current_elem_volume
 
dof_id_type _dof_id
 
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
 
std::unordered_set< const Elem * > _wall_bounded
 Maps for wall treatment. More...
 
std::map< const Elem *, std::vector< Real > > _dist
 
std::map< const Elem *, std::vector< const FaceInfo * > > _face_infos
 

Detailed Description

Kernel that adds contributions to the source and the sink of the turbulent kinetic energy dissipation discretized using the finite volume method to a linear system.

Definition at line 20 of file LinearFVTKEDSourceSink.h.

Constructor & Destructor Documentation

◆ LinearFVTKEDSourceSink()

LinearFVTKEDSourceSink::LinearFVTKEDSourceSink ( const InputParameters params)

Class constructor.

Parameters
paramsThe InputParameters for the kernel.

Definition at line 51 of file LinearFVTKEDSourceSink.C.

52  : LinearFVElementalKernel(params),
54  _u_var(getFunctor<Real>("u")),
55  _v_var(params.isParamValid("v") ? &(getFunctor<Real>("v")) : nullptr),
56  _w_var(params.isParamValid("w") ? &(getFunctor<Real>("w")) : nullptr),
57  _k(getFunctor<Real>(NS::TKE)),
58  _rho(getFunctor<Real>(NS::density)),
59  _mu(getFunctor<Real>(NS::mu)),
60  _mu_t(getFunctor<Real>(NS::mu_t)),
61  _wall_boundary_names(getParam<std::vector<BoundaryName>>("walls")),
62  _linearized_model(getParam<bool>("linearized_model")),
63  _wall_treatment(getParam<MooseEnum>("wall_treatment").getEnum<NS::WallTreatmentEnum>()),
64  _C1_eps(getParam<Real>("C1_eps")),
65  _C2_eps(getParam<Real>("C2_eps")),
66  _C_mu(getParam<Real>("C_mu")),
67  _C_pl(getParam<Real>("C_pl"))
68 {
69  if (_dim >= 2 && !_v_var)
70  paramError("v", "In two or more dimensions, the v velocity must be supplied!");
71 
72  if (_dim >= 3 && !_w_var)
73  paramError("w", "In three or more dimensions, the w velocity must be supplied!");
74 
75  // Strain tensor term requires velocity gradients;
76  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(&_u_var))
77  requestVariableCellGradient(getParam<MooseFunctorName>("u"));
78  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(_v_var))
79  requestVariableCellGradient(getParam<MooseFunctorName>("v"));
80  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(_w_var))
81  requestVariableCellGradient(getParam<MooseFunctorName>("w"));
82 }
virtual MooseMesh & mesh()=0
static const std::string mu_t
Definition: NS.h:125
SubProblem & _subproblem
const bool _linearized_model
If the user wants to use the linearized model.
const Moose::Functor< Real > * _v_var
y-velocity
void requestVariableCellGradient(const std::string &variable_name)
const Moose::Functor< Real > & _u_var
x-velocity
static const std::string density
Definition: NS.h:33
static const std::string TKE
Definition: NS.h:176
WallTreatmentEnum
Wall treatment options.
Definition: NS.h:182
NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
const Real _C2_eps
Value of the second epsilon closure coefficient.
const Moose::Functor< Real > & _rho
Density.
const Real _C_pl
Production Limiter Constant.
const unsigned int _dim
The dimension of the simulation.
virtual unsigned int dimension() const
static const std::string mu
Definition: NS.h:123
const T & getParam(const std::string &name) const
void paramError(const std::string &param, Args... args) const
const Real _C1_eps
Value of the first epsilon closure coefficient.
LinearFVElementalKernel(const InputParameters &params)
const Moose::Functor< Real > & _mu_t
Turbulent dynamic viscosity.
const Moose::Functor< Real > & _mu
Dynamic viscosity.
const Moose::Functor< Real > & _k
Turbulent kinetic energy.
const Real _C_mu
C_mu constant.
const Moose::Functor< Real > * _w_var
z-velocity
bool isParamValid(const std::string &name) const

Member Function Documentation

◆ computeMatrixContribution()

Real LinearFVTKEDSourceSink::computeMatrixContribution ( )
overridevirtual

Implements LinearFVElementalKernel.

Definition at line 95 of file LinearFVTKEDSourceSink.C.

96 {
97  if (_wall_bounded.find(_current_elem_info->elem()) != _wall_bounded.end())
98  // TKED value for near wall element will be directly assigned for this cell
99  return 1.0;
100  else
101  {
102  // Convenient definitions
103  const auto state = determineState();
104  const auto elem_arg = makeElemArg(_current_elem_info->elem());
105  const Real rho = _rho(elem_arg, state);
106  const Real TKE = _k(elem_arg, state);
107  const auto epsilon = _var.getElemValue(*_current_elem_info, state);
108 
109  // Compute destruction
110  const auto destruction = _C2_eps * rho * epsilon / TKE;
111 
112  // Assign to matrix (term gets multiplied by TKED)
113  return destruction * _current_elem_volume;
114  }
115 }
const ElemInfo * _current_elem_info
Moose::StateArg determineState() const
const Elem * elem() const
MooseLinearVariableFV< Real > & _var
static const std::string TKE
Definition: NS.h:176
const Real _C2_eps
Value of the second epsilon closure coefficient.
const Moose::Functor< Real > & _rho
Density.
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
const Moose::Functor< Real > & _k
Turbulent kinetic energy.
std::unordered_set< const Elem * > _wall_bounded
Maps for wall treatment.

◆ computeRightHandSideContribution()

Real LinearFVTKEDSourceSink::computeRightHandSideContribution ( )
overridevirtual

Implements LinearFVElementalKernel.

Definition at line 118 of file LinearFVTKEDSourceSink.C.

119 {
120  if (_wall_bounded.find(_current_elem_info->elem()) != _wall_bounded.end())
121  {
122  // Convenient definitions
123  const auto state = determineState();
124  const auto elem_arg = makeElemArg(_current_elem_info->elem());
125  const Real rho = _rho(elem_arg, state);
126  const Real mu = _mu(elem_arg, state);
127  const Real TKE = _k(elem_arg, state);
128 
129  // Convenient variables
130  Real destruction = 0.0;
131  std::vector<Real> y_plus_vec, velocity_grad_norm_vec;
132  Real tot_weight = 0.0;
133 
134  // Get velocity vector
135  RealVectorValue velocity(_u_var(elem_arg, state));
136  if (_v_var)
137  velocity(1) = (*_v_var)(elem_arg, state);
138  if (_w_var)
139  velocity(2) = (*_w_var)(elem_arg, state);
140 
141  // Get near wall faceInfo and distances from cell center to every wall
142  const auto & face_info_vec = libmesh_map_find(_face_infos, _current_elem_info->elem());
143  const auto & distance_vec = libmesh_map_find(_dist, _current_elem_info->elem());
144  mooseAssert(distance_vec.size(), "Should have found a distance vector");
145  mooseAssert(distance_vec.size() == face_info_vec.size(),
146  "Should be as many distance vectors as face info vectors");
147 
148  // Update y+ and wall face cell
149  for (unsigned int i = 0; i < distance_vec.size(); i++)
150  {
151  const auto distance = distance_vec[i];
152  mooseAssert(distance > 0, "Should be at a non-zero distance");
153 
154  Real y_plus;
155  if (_wall_treatment == NS::WallTreatmentEnum::NEQ) // Non-equilibrium / Non-iterative
156  y_plus = distance * std::sqrt(std::sqrt(_C_mu) * TKE) * rho / mu;
157  else // Equilibrium / Iterative
158  {
159  const auto parallel_speed = NS::computeSpeed<Real>(
160  velocity - velocity * face_info_vec[i]->normal() * face_info_vec[i]->normal());
161  y_plus = NS::findyPlus<Real>(mu, rho, std::max(parallel_speed, 1e-10), distance);
162  }
163 
164  y_plus_vec.push_back(y_plus);
165  tot_weight += 1.0;
166  }
167 
168  // Compute near wall epsilon value
169  for (const auto i : index_range(y_plus_vec))
170  {
171  const auto y_plus = y_plus_vec[i];
172 
173  if (y_plus < 11.25)
174  destruction += 2.0 * TKE * mu / rho / Utility::pow<2>(distance_vec[i]) / tot_weight;
175  else
176  destruction += std::pow(_C_mu, 0.75) * std::pow(TKE, 1.5) /
177  (NS::von_karman_constant * distance_vec[i]) / tot_weight;
178  }
179 
180  // Assign the computed value of TKED for element near the wall
181  return destruction;
182  }
183  else
184  {
185  // Convenient definitions
186  const auto state = determineState();
187  const auto elem_arg = makeElemArg(_current_elem_info->elem());
188  const Real rho = _rho(elem_arg, state);
189  const Real TKE = _k(elem_arg, state);
190  const Real TKED = _var.getElemValue(*_current_elem_info, state);
191 
192  // Compute production of TKE
193  const auto symmetric_strain_tensor_sq_norm =
195  Real production = _mu_t(elem_arg, state) * symmetric_strain_tensor_sq_norm;
196 
197  // Limit TKE production (needed for flows with stagnation zones)
198  const Real production_limit = _C_pl * rho * TKED;
199  production = std::min(production, production_limit);
200 
201  // Compute production - recasted with mu_t definition to avoid division by epsilon
202  const auto production_epsilon = _C1_eps * production * TKED / TKE;
203 
204  // Assign to matrix (term gets multiplied by TKED)
205  return production_epsilon * _current_elem_volume;
206  }
207 }
const ElemInfo * _current_elem_info
static constexpr Real von_karman_constant
Definition: NS.h:191
const Moose::Functor< Real > * _v_var
y-velocity
Moose::StateArg determineState() const
const Elem * elem() const
const Moose::Functor< Real > & _u_var
x-velocity
std::map< const Elem *, std::vector< Real > > _dist
MooseLinearVariableFV< Real > & _var
static const std::string TKE
Definition: NS.h:176
NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
const Moose::Functor< Real > & _rho
Density.
Real distance(const Point &p)
const Real _C_pl
Production Limiter Constant.
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
template Real findyPlus< Real >(const Real &mu, const Real &rho, const Real &u, Real dist)
static const std::string mu
Definition: NS.h:123
std::map< const Elem *, std::vector< const FaceInfo * > > _face_infos
const Real _C1_eps
Value of the first epsilon closure coefficient.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
const Moose::Functor< Real > & _mu_t
Turbulent dynamic viscosity.
static const std::string TKED
Definition: NS.h:177
static const std::string velocity
Definition: NS.h:45
const Moose::Functor< Real > & _mu
Dynamic viscosity.
template Real computeSpeed< Real >(const libMesh::VectorValue< Real > &velocity)
const Moose::Functor< Real > & _k
Turbulent kinetic energy.
const Real _C_mu
C_mu constant.
MooseUnits pow(const MooseUnits &, int)
std::unordered_set< const Elem * > _wall_bounded
Maps for wall treatment.
auto index_range(const T &sizable)
const Moose::Functor< Real > * _w_var
z-velocity
template Real computeShearStrainRateNormSquared< Real >(const Moose::Functor< Real > &u, const Moose::Functor< Real > *v, const Moose::Functor< Real > *w, const Moose::ElemArg &elem_arg, const Moose::StateArg &state)

◆ initialSetup()

void LinearFVTKEDSourceSink::initialSetup ( )
overridevirtual

Reimplemented from LinearFVElementalKernel.

Definition at line 85 of file LinearFVTKEDSourceSink.C.

86 {
92 }
SubProblem & _subproblem
std::map< const Elem *, std::vector< Real > > _dist
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...
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
virtual const std::set< SubdomainID > & blockIDs() const
std::map< const Elem *, std::vector< const FaceInfo * > > _face_infos
void getWallDistance(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::map< const Elem *, std::vector< Real >> &dist_map)
Map storing wall ditance for near-wall marked elements The map passed in dist_map gets cleared and re...
void getElementFaceArgs(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::map< const Elem *, std::vector< const FaceInfo *>> &face_info_map)
Map storing face arguments to wall bounded faces The map passed in face_info_map gets cleared and re-...
FEProblemBase & _fe_problem
virtual void initialSetup()
std::unordered_set< const Elem * > _wall_bounded
Maps for wall treatment.

◆ validParams()

InputParameters LinearFVTKEDSourceSink::validParams ( )
static

Definition at line 18 of file LinearFVTKEDSourceSink.C.

19 {
21  params.addClassDescription("Elemental kernel to compute the production and destruction "
22  " terms of turbulent kinetic energy dissipation (TKED).");
23  params.addRequiredParam<MooseFunctorName>("u", "The velocity in the x direction.");
24  params.addParam<MooseFunctorName>("v", "The velocity in the y direction.");
25  params.addParam<MooseFunctorName>("w", "The velocity in the z direction.");
26  params.addRequiredParam<MooseFunctorName>(NS::TKE, "Coupled turbulent kinetic energy.");
27  params.addRequiredParam<MooseFunctorName>(NS::density, "fluid density");
28  params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity.");
29  params.addRequiredParam<MooseFunctorName>(NS::mu_t, "Turbulent viscosity.");
30 
31  params.addParam<std::vector<BoundaryName>>(
32  "walls", {}, "Boundaries that correspond to solid walls.");
33  params.addParam<bool>(
34  "linearized_model",
35  true,
36  "Boolean to determine if the problem should be used in a linear or nonlinear solve");
37  MooseEnum wall_treatment("eq_newton eq_incremental eq_linearized neq", "neq");
38  params.addParam<MooseEnum>("wall_treatment",
39  wall_treatment,
40  "The method used for computing the wall functions "
41  "'eq_newton', 'eq_incremental', 'eq_linearized', 'neq'");
42 
43  params.addParam<Real>("C1_eps", 1.44, "First epsilon coefficient");
44  params.addParam<Real>("C2_eps", 1.92, "Second epsilon coefficient");
45  params.addParam<Real>("C_mu", 0.09, "Coupled turbulent kinetic energy closure.");
46  params.addParam<Real>("C_pl", 10.0, "Production limiter constant multiplier.");
47 
48  return params;
49 }
static const std::string mu_t
Definition: NS.h:125
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
static const std::string density
Definition: NS.h:33
static const std::string TKE
Definition: NS.h:176
void addRequiredParam(const std::string &name, const std::string &doc_string)
static const std::string mu
Definition: NS.h:123
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _C1_eps

const Real LinearFVTKEDSourceSink::_C1_eps
protected

Value of the first epsilon closure coefficient.

Definition at line 70 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution().

◆ _C2_eps

const Real LinearFVTKEDSourceSink::_C2_eps
protected

Value of the second epsilon closure coefficient.

Definition at line 73 of file LinearFVTKEDSourceSink.h.

Referenced by computeMatrixContribution().

◆ _C_mu

const Real LinearFVTKEDSourceSink::_C_mu
protected

C_mu constant.

Definition at line 76 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution().

◆ _C_pl

const Real LinearFVTKEDSourceSink::_C_pl
protected

Production Limiter Constant.

Definition at line 79 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution().

◆ _dim

const unsigned int LinearFVTKEDSourceSink::_dim
protected

The dimension of the simulation.

Definition at line 39 of file LinearFVTKEDSourceSink.h.

Referenced by LinearFVTKEDSourceSink().

◆ _dist

std::map<const Elem *, std::vector<Real> > LinearFVTKEDSourceSink::_dist
protected

Definition at line 84 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution(), and initialSetup().

◆ _face_infos

std::map<const Elem *, std::vector<const FaceInfo *> > LinearFVTKEDSourceSink::_face_infos
protected

Definition at line 85 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution(), and initialSetup().

◆ _k

const Moose::Functor<Real>& LinearFVTKEDSourceSink::_k
protected

Turbulent kinetic energy.

Definition at line 49 of file LinearFVTKEDSourceSink.h.

Referenced by computeMatrixContribution(), and computeRightHandSideContribution().

◆ _linearized_model

const bool LinearFVTKEDSourceSink::_linearized_model
protected

If the user wants to use the linearized model.

Definition at line 64 of file LinearFVTKEDSourceSink.h.

◆ _mu

const Moose::Functor<Real>& LinearFVTKEDSourceSink::_mu
protected

Dynamic viscosity.

Definition at line 55 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution().

◆ _mu_t

const Moose::Functor<Real>& LinearFVTKEDSourceSink::_mu_t
protected

Turbulent dynamic viscosity.

Definition at line 58 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution().

◆ _rho

const Moose::Functor<Real>& LinearFVTKEDSourceSink::_rho
protected

Density.

Definition at line 52 of file LinearFVTKEDSourceSink.h.

Referenced by computeMatrixContribution(), and computeRightHandSideContribution().

◆ _u_var

const Moose::Functor<Real>& LinearFVTKEDSourceSink::_u_var
protected

x-velocity

Definition at line 42 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution(), and LinearFVTKEDSourceSink().

◆ _v_var

const Moose::Functor<Real>* LinearFVTKEDSourceSink::_v_var
protected

y-velocity

Definition at line 44 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution(), and LinearFVTKEDSourceSink().

◆ _w_var

const Moose::Functor<Real>* LinearFVTKEDSourceSink::_w_var
protected

z-velocity

Definition at line 46 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution(), and LinearFVTKEDSourceSink().

◆ _wall_boundary_names

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

Wall boundaries.

Definition at line 61 of file LinearFVTKEDSourceSink.h.

Referenced by initialSetup().

◆ _wall_bounded

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

Maps for wall treatment.

Definition at line 83 of file LinearFVTKEDSourceSink.h.

Referenced by computeMatrixContribution(), computeRightHandSideContribution(), and initialSetup().

◆ _wall_treatment

NS::WallTreatmentEnum LinearFVTKEDSourceSink::_wall_treatment
protected

Method used for wall treatment.

Definition at line 67 of file LinearFVTKEDSourceSink.h.

Referenced by computeRightHandSideContribution().


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