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

An advection kernel that implements the advection term for the turbulent variables limited for the first cells near the wall. More...

#include <LinearFVTurbulentAdvection.h>

Inheritance diagram for LinearFVTurbulentAdvection:
[legend]

Public Types

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

Public Member Functions

 LinearFVTurbulentAdvection (const InputParameters &params)
 
virtual void addMatrixContribution () override
 
virtual void addRightHandSideContribution () override
 
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 setupFaceData (const FaceInfo *face_info) override
 
virtual bool hasFaceSide (const FaceInfo &fi, bool fi_elem_side) const override
 
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
 
bool isKokkosObject (IsKokkosObjectKey &&) const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
const std::string & name () const
 
std::string typeAndName () const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
MooseObjectName uniqueName () const
 
const InputParametersparameters () const
 
const hit::Node * getHitNode () const
 
bool hasBase () const
 
const std::string & getBase () 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 &name) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) 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
 
std::string messagePrefix (const bool hit_prefix=true) const
 
std::string errorPrefix (const std::string &) const
 
void mooseError (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
void callMooseError (std::string msg, const bool with_prefix, const hit::Node *node=nullptr) 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 void callMooseError (MooseApp *const app, const InputParameters &params, std::string msg, const bool with_prefix, const hit::Node *node)
 
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
 

Static Public Attributes

static const std::string type_param
 
static const std::string name_param
 
static const std::string unique_name_param
 
static const std::string app_param
 
static const std::string moose_base_param
 
static const std::string kokkos_object_param
 

Protected Member Functions

std::string deduceFunctorName (const std::string &name) const
 
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 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)
 
void initializeKokkosBlockRestrictable (const Moose::Kokkos::Mesh *mesh)
 
Moose::CoordinateSystemType getBlockCoordSystem ()
 
KOKKOS_FUNCTION dof_id_type numKokkosBlockElements () const
 
KOKKOS_FUNCTION dof_id_type numKokkosBlockNodes () const
 
KOKKOS_FUNCTION dof_id_type numKokkosBlockSides () const
 
KOKKOS_FUNCTION ContiguousElementID kokkosBlockElementID (ThreadID tid) const
 
KOKKOS_FUNCTION ContiguousElementID kokkosBlockNodeID (ThreadID tid) const
 
KOKKOS_FUNCTION auto kokkosBlockElementSideID (ThreadID tid) const
 
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 RhieChowMassFlux_mass_flux_provider
 The Rhie-Chow user object that provides us with the face velocity. More...
 
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
 
Factory_factory
 
ActionFactory_action_factory
 
const std::string & _type
 
const std::string & _name
 
const InputParameters_pars
 
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
 

Private Attributes

std::pair< Real, Real_advected_interp_coeffs
 Container for the current advected interpolation coefficients on the face to make sure we don't compute it multiple times for different terms. More...
 
Real _mass_face_flux
 Container for the mass flux on the face which will be reused in the advection term's matrix and right hand side contribution. More...
 
Moose::FV::InterpMethod _advected_interp_method
 The interpolation method to use for the advected quantity. More...
 
const std::vector< BoundaryName > & _wall_boundary_names
 Wall boundaries. More...
 
std::unordered_set< const Elem * > _wall_bounded
 List for wall bounded elements. More...
 

Detailed Description

An advection kernel that implements the advection term for the turbulent variables limited for the first cells near the wall.

Definition at line 20 of file LinearFVTurbulentAdvection.h.

Constructor & Destructor Documentation

◆ LinearFVTurbulentAdvection()

LinearFVTurbulentAdvection::LinearFVTurbulentAdvection ( const InputParameters params)

Definition at line 37 of file LinearFVTurbulentAdvection.C.

38  : LinearFVFluxKernel(params),
39  _mass_flux_provider(getUserObject<RhieChowMassFlux>("rhie_chow_user_object")),
40  _advected_interp_coeffs(std::make_pair<Real, Real>(0, 0)),
41  _mass_face_flux(0.0),
42  _wall_boundary_names(getParam<std::vector<BoundaryName>>("walls"))
43 {
44  Moose::FV::setInterpolationMethod(*this, _advected_interp_method, "advected_interp_method");
45 }
const T & getParam(const std::string &name) const
Real _mass_face_flux
Container for the mass flux on the face which will be reused in the advection term&#39;s matrix and right...
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
std::pair< Real, Real > _advected_interp_coeffs
Container for the current advected interpolation coefficients on the face to make sure we don&#39;t compu...
Moose::FV::InterpMethod _advected_interp_method
The interpolation method to use for the advected quantity.
LinearFVFluxKernel(const InputParameters &params)
const RhieChowMassFlux & _mass_flux_provider
The Rhie-Chow user object that provides us with the face velocity.
bool setInterpolationMethod(const MooseObject &obj, Moose::FV::InterpMethod &interp_method, const std::string &param_name)

Member Function Documentation

◆ addMatrixContribution()

void LinearFVTurbulentAdvection::addMatrixContribution ( )
overridevirtual

Reimplemented from LinearFVFluxKernel.

Definition at line 56 of file LinearFVTurbulentAdvection.C.

57 {
58  // Coumputing bounding map
59  const Elem * elem = _current_face_info->elemPtr();
60  const auto bounded_elem = _wall_bounded.find(elem) != _wall_bounded.end();
61  const Elem * neighbor = _current_face_info->neighborPtr();
62  const auto bounded_neigh = _wall_bounded.find(neighbor) != _wall_bounded.end();
63 
64  // If we are on an internal face, we populate the four entries in the system matrix
65  // which touch the face
67  {
68  // The dof ids of the variable corresponding to the element and neighbor
71 
72  // Compute the entries which will go to the neighbor (offdiagonal) and element
73  // (diagonal).
74  const auto elem_matrix_contribution = computeElemMatrixContribution();
75  const auto neighbor_matrix_contribution = computeNeighborMatrixContribution();
76 
77  // Populate matrix
78  if (hasBlocks(_current_face_info->elemInfo()->subdomain_id()) && !(bounded_elem))
79  {
80  _matrix_contribution(0, 0) = elem_matrix_contribution;
81  _matrix_contribution(0, 1) = neighbor_matrix_contribution;
82  }
83 
84  if (hasBlocks(_current_face_info->neighborInfo()->subdomain_id()) && !(bounded_neigh))
85  {
86  _matrix_contribution(1, 0) = -elem_matrix_contribution;
87  _matrix_contribution(1, 1) = -neighbor_matrix_contribution;
88  }
89  (*_linear_system.matrix).add_matrix(_matrix_contribution, _dof_indices.get_values());
90  }
91  // We are at a block boundary where the variable is not defined on one of the adjacent cells.
92  // We check if we have a boundary condition here
95  {
96  mooseAssert(_current_face_info->boundaryIDs().size() == 1,
97  "We should only have one boundary on every face.");
98 
99  LinearFVBoundaryCondition * bc_pointer =
101 
102  if (bc_pointer || _force_boundary_execution)
103  {
104  if (bc_pointer)
106  const auto matrix_contribution = computeBoundaryMatrixContribution(*bc_pointer);
107 
108  // We allow internal (for the mesh) boundaries too, so we have to check on which side we
109  // are on (assuming that this is a boundary for the variable)
110  if ((_current_face_type == FaceInfo::VarFaceNeighbors::ELEM) && !(bounded_elem))
111  {
112  const auto dof_id_elem = _current_face_info->elemInfo()->dofIndices()[_sys_num][_var_num];
113  (*_linear_system.matrix).add(dof_id_elem, dof_id_elem, matrix_contribution);
114  }
115  else if ((_current_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) && !(bounded_neigh))
116  {
117  const auto dof_id_neighbor =
119  (*_linear_system.matrix).add(dof_id_neighbor, dof_id_neighbor, matrix_contribution);
120  }
121  }
122  }
123 }
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
FaceInfo::VarFaceNeighbors _current_face_type
const bool _force_boundary_execution
virtual Real computeNeighborMatrixContribution() override
const FaceInfo * _current_face_info
const Elem * neighborPtr() const
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
std::unordered_set< const Elem * > _wall_bounded
List for wall bounded elements.
const unsigned int _sys_num
virtual Real computeElemMatrixContribution() override
SubdomainID subdomain_id() const
DenseMatrix< Real > _matrix_contribution

◆ addRightHandSideContribution()

void LinearFVTurbulentAdvection::addRightHandSideContribution ( )
overridevirtual

Reimplemented from LinearFVFluxKernel.

Definition at line 126 of file LinearFVTurbulentAdvection.C.

127 {
128  // Coumputing bounding map
129  const Elem * elem = _current_face_info->elemPtr();
130  const auto bounded_elem = _wall_bounded.find(elem) != _wall_bounded.end();
131  const Elem * neighbor = _current_face_info->neighborPtr();
132  const auto bounded_neigh = _wall_bounded.find(neighbor) != _wall_bounded.end();
133 
134  // If we are on an internal face, we populate the two entries in the right hand side
135  // which touch the face
137  {
138  // The dof ids of the variable corresponding to the element and neighbor
141 
142  // Compute the entries which will go to the neighbor and element positions.
143  const auto elem_rhs_contribution = computeElemRightHandSideContribution();
144  const auto neighbor_rhs_contribution = computeNeighborRightHandSideContribution();
145 
146  // Populate right hand side
147  if (hasBlocks(_current_face_info->elemInfo()->subdomain_id()) && !(bounded_elem))
148  _rhs_contribution(0) = elem_rhs_contribution;
149  if (hasBlocks(_current_face_info->neighborInfo()->subdomain_id()) && !(bounded_neigh))
150  _rhs_contribution(1) = neighbor_rhs_contribution;
151 
153  .add_vector(_rhs_contribution.get_values().data(), _dof_indices.get_values());
154  }
155  // We are at a block boundary where the variable is not defined on one of the adjacent cells.
156  // We check if we have a boundary condition here
159  {
160  mooseAssert(_current_face_info->boundaryIDs().size() == 1,
161  "We should only have one boundary on every face.");
162  LinearFVBoundaryCondition * bc_pointer =
164 
165  if (bc_pointer || _force_boundary_execution)
166  {
167  if (bc_pointer)
169 
170  const auto rhs_contribution = computeBoundaryRHSContribution(*bc_pointer);
171 
172  // We allow internal (for the mesh) boundaries too, so we have to check on which side we
173  // are on (assuming that this is a boundary for the variable)
174  if ((_current_face_type == FaceInfo::VarFaceNeighbors::ELEM) && !(bounded_elem))
175  {
176  const auto dof_id_elem = _current_face_info->elemInfo()->dofIndices()[_sys_num][_var_num];
177  (*_linear_system.rhs).add(dof_id_elem, rhs_contribution);
178  }
179  else if ((_current_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) && !(bounded_neigh))
180  {
181  const auto dof_id_neighbor =
183  (*_linear_system.rhs).add(dof_id_neighbor, rhs_contribution);
184  }
185  }
186  }
187 }
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
FaceInfo::VarFaceNeighbors _current_face_type
const bool _force_boundary_execution
DenseVector< Real > _rhs_contribution
virtual Real computeBoundaryRHSContribution(const LinearFVBoundaryCondition &bc) override
const FaceInfo * _current_face_info
const Elem * neighborPtr() const
const Elem * elemPtr() const
DenseVector< dof_id_type > _dof_indices
const std::vector< std::vector< dof_id_type > > & dofIndices() const
virtual Real computeNeighborRightHandSideContribution() override
virtual Real computeElemRightHandSideContribution() override
bool hasBlocks(const SubdomainName &name) const
std::unordered_set< const Elem * > _wall_bounded
List for wall bounded elements.
const unsigned int _sys_num
SubdomainID subdomain_id() const

◆ computeBoundaryMatrixContribution()

Real LinearFVTurbulentAdvection::computeBoundaryMatrixContribution ( const LinearFVBoundaryCondition bc)
overridevirtual

Implements LinearFVFluxKernel.

Definition at line 214 of file LinearFVTurbulentAdvection.C.

Referenced by addMatrixContribution().

215 {
216  const auto * const adv_bc = static_cast<const LinearFVAdvectionDiffusionBC *>(&bc);
217  mooseAssert(adv_bc, "This should be a valid BC!");
218 
219  const auto boundary_value_matrix_contrib = adv_bc->computeBoundaryValueMatrixContribution();
220 
221  // We support internal boundaries too so we have to make sure the normal points always outward
222  const auto factor = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM) ? 1.0 : -1.0;
223 
224  return boundary_value_matrix_contrib * factor * _mass_face_flux * _current_face_area;
225 }
Real _mass_face_flux
Container for the mass flux on the face which will be reused in the advection term&#39;s matrix and right...
FaceInfo::VarFaceNeighbors _current_face_type

◆ computeBoundaryRHSContribution()

Real LinearFVTurbulentAdvection::computeBoundaryRHSContribution ( const LinearFVBoundaryCondition bc)
overridevirtual

Implements LinearFVFluxKernel.

Definition at line 228 of file LinearFVTurbulentAdvection.C.

Referenced by addRightHandSideContribution().

229 {
230  const auto * const adv_bc = static_cast<const LinearFVAdvectionDiffusionBC *>(&bc);
231  mooseAssert(adv_bc, "This should be a valid BC!");
232 
233  // We support internal boundaries too so we have to make sure the normal points always outward
234  const auto factor = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM ? 1.0 : -1.0);
235 
236  const auto boundary_value_rhs_contrib = adv_bc->computeBoundaryValueRHSContribution();
237  return -boundary_value_rhs_contrib * factor * _mass_face_flux * _current_face_area;
238 }
Real _mass_face_flux
Container for the mass flux on the face which will be reused in the advection term&#39;s matrix and right...
FaceInfo::VarFaceNeighbors _current_face_type

◆ computeElemMatrixContribution()

Real LinearFVTurbulentAdvection::computeElemMatrixContribution ( )
overridevirtual

Implements LinearFVFluxKernel.

Definition at line 190 of file LinearFVTurbulentAdvection.C.

Referenced by addMatrixContribution().

191 {
193 }
Real _mass_face_flux
Container for the mass flux on the face which will be reused in the advection term&#39;s matrix and right...
std::pair< Real, Real > _advected_interp_coeffs
Container for the current advected interpolation coefficients on the face to make sure we don&#39;t compu...

◆ computeElemRightHandSideContribution()

Real LinearFVTurbulentAdvection::computeElemRightHandSideContribution ( )
overridevirtual

Implements LinearFVFluxKernel.

Definition at line 202 of file LinearFVTurbulentAdvection.C.

Referenced by addRightHandSideContribution().

203 {
204  return 0.0;
205 }

◆ computeNeighborMatrixContribution()

Real LinearFVTurbulentAdvection::computeNeighborMatrixContribution ( )
overridevirtual

Implements LinearFVFluxKernel.

Definition at line 196 of file LinearFVTurbulentAdvection.C.

Referenced by addMatrixContribution().

197 {
199 }
Real _mass_face_flux
Container for the mass flux on the face which will be reused in the advection term&#39;s matrix and right...
std::pair< Real, Real > _advected_interp_coeffs
Container for the current advected interpolation coefficients on the face to make sure we don&#39;t compu...

◆ computeNeighborRightHandSideContribution()

Real LinearFVTurbulentAdvection::computeNeighborRightHandSideContribution ( )
overridevirtual

Implements LinearFVFluxKernel.

Definition at line 208 of file LinearFVTurbulentAdvection.C.

Referenced by addRightHandSideContribution().

209 {
210  return 0.0;
211 }

◆ initialSetup()

void LinearFVTurbulentAdvection::initialSetup ( )
overridevirtual

Reimplemented from LinearFVFluxKernel.

Definition at line 48 of file LinearFVTurbulentAdvection.C.

49 {
53 }
SubProblem & _subproblem
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
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
FEProblemBase & _fe_problem
virtual void initialSetup()
std::unordered_set< const Elem * > _wall_bounded
List for wall bounded elements.

◆ setupFaceData()

void LinearFVTurbulentAdvection::setupFaceData ( const FaceInfo face_info)
overridevirtual

Reimplemented from LinearFVFluxKernel.

Definition at line 241 of file LinearFVTurbulentAdvection.C.

242 {
244 
245  // Caching the mass flux on the face which will be reused in the advection term's matrix and right
246  // hand side contributions
248 
249  // Caching the interpolation coefficients so they will be reused for the matrix and right hand
250  // side terms
253 }
std::pair< Real, Real > interpCoeffs(const InterpMethod m, const FaceInfo &fi, const bool one_is_elem, const T &face_flux=0.0)
Real getMassFlux(const FaceInfo &fi) const
Get the face velocity times density (used in advection terms)
virtual void setupFaceData(const FaceInfo *face_info)
Real _mass_face_flux
Container for the mass flux on the face which will be reused in the advection term&#39;s matrix and right...
const FaceInfo * _current_face_info
std::pair< Real, Real > _advected_interp_coeffs
Container for the current advected interpolation coefficients on the face to make sure we don&#39;t compu...
Moose::FV::InterpMethod _advected_interp_method
The interpolation method to use for the advected quantity.
const RhieChowMassFlux & _mass_flux_provider
The Rhie-Chow user object that provides us with the face velocity.

◆ validParams()

InputParameters LinearFVTurbulentAdvection::validParams ( )
static

Definition at line 19 of file LinearFVTurbulentAdvection.C.

20 {
22  params.addClassDescription("Represents the matrix and right hand side contributions of an "
23  "advection term for a turbulence variable.");
24 
25  params.addRequiredParam<UserObjectName>(
26  "rhie_chow_user_object",
27  "The rhie-chow user-object which is used to determine the face velocity.");
28 
30 
31  params.addParam<std::vector<BoundaryName>>(
32  "walls", {}, "Boundaries that correspond to solid walls.");
33 
34  return params;
35 }
void addRequiredParam(const std::string &name, const std::string &doc_string)
static InputParameters validParams()
InputParameters advectedInterpolationParameter()
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _advected_interp_coeffs

std::pair<Real, Real> LinearFVTurbulentAdvection::_advected_interp_coeffs
private

Container for the current advected interpolation coefficients on the face to make sure we don't compute it multiple times for different terms.

Definition at line 53 of file LinearFVTurbulentAdvection.h.

Referenced by computeElemMatrixContribution(), computeNeighborMatrixContribution(), and setupFaceData().

◆ _advected_interp_method

Moose::FV::InterpMethod LinearFVTurbulentAdvection::_advected_interp_method
private

The interpolation method to use for the advected quantity.

Definition at line 60 of file LinearFVTurbulentAdvection.h.

Referenced by LinearFVTurbulentAdvection(), and setupFaceData().

◆ _mass_face_flux

Real LinearFVTurbulentAdvection::_mass_face_flux
private

Container for the mass flux on the face which will be reused in the advection term's matrix and right hand side contribution.

Definition at line 57 of file LinearFVTurbulentAdvection.h.

Referenced by computeBoundaryMatrixContribution(), computeBoundaryRHSContribution(), computeElemMatrixContribution(), computeNeighborMatrixContribution(), and setupFaceData().

◆ _mass_flux_provider

const RhieChowMassFlux& LinearFVTurbulentAdvection::_mass_flux_provider
protected

The Rhie-Chow user object that provides us with the face velocity.

Definition at line 48 of file LinearFVTurbulentAdvection.h.

Referenced by setupFaceData().

◆ _wall_boundary_names

const std::vector<BoundaryName>& LinearFVTurbulentAdvection::_wall_boundary_names
private

Wall boundaries.

Definition at line 63 of file LinearFVTurbulentAdvection.h.

Referenced by initialSetup().

◆ _wall_bounded

std::unordered_set<const Elem *> LinearFVTurbulentAdvection::_wall_bounded
private

List for wall bounded elements.

Definition at line 66 of file LinearFVTurbulentAdvection.h.

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


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