https://mooseframework.inl.gov
FVFluxKernel.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "FVKernel.h"
13 #include "MathFVUtils.h"
14 #include "NeighborCoupleable.h"
18 #include "FVFaceResidualObject.h"
19 #include "FaceArgInterface.h"
20 
21 class FaceInfo;
22 
30 class FVFluxKernel : public FVKernel,
32  public NeighborMooseVariableInterface<Real>,
34  public FVFaceResidualObject,
36 {
37 public:
39  FVFluxKernel(const InputParameters & params);
40 
41  void computeResidual() override;
42  void computeJacobian() override;
43  void computeResidualAndJacobian() override;
44  void computeResidual(const FaceInfo & fi) override;
45  void computeJacobian(const FaceInfo & fi) override;
46  void computeResidualAndJacobian(const FaceInfo & fi) override;
47 
48  const MooseVariableFV<Real> & variable() const override { return _var; }
49 
50  bool hasFaceSide(const FaceInfo & fi, const bool fi_elem_side) const override;
51 
52 protected:
57  virtual ADReal computeQpResidual() = 0;
58 
62  virtual ADReal gradUDotNormal(const Moose::StateArg & time, const bool correct_skewness) const;
63 
69  virtual bool skipForBoundary(const FaceInfo & fi) const;
70 
71  const RealVectorValue & normal() const { return _normal; }
72 
74 
75  const unsigned int _qp = 0;
76 
81 
86 
89  const FaceInfo * _face_info = nullptr;
90 
93 
97  bool onBoundary(const FaceInfo & fi) const;
98 
102  Moose::ElemArg elemArg(bool correct_skewness = false) const;
103 
107  Moose::ElemArg neighborArg(bool correct_skewness = false) const;
108 
119  const FaceInfo * fi = nullptr,
121  bool correct_skewness = false,
122  const Moose::StateArg * state_limiter = nullptr) const;
123 
128  bool avoidBoundary(const FaceInfo & fi) const;
129 
131  std::unordered_set<BoundaryID> _boundaries_to_force;
132 
133 private:
136 
138  std::unordered_set<BoundaryID> _boundaries_to_avoid;
139 };
bool hasFaceSide(const FaceInfo &fi, const bool fi_elem_side) const override
Definition: FVFluxKernel.C:292
virtual ADReal gradUDotNormal(const Moose::StateArg &time, const bool correct_skewness) const
Calculates and returns "grad_u dot normal" on the face to be used for diffusive terms.
Definition: FVFluxKernel.C:232
std::unordered_set< BoundaryID > _boundaries_to_avoid
Which boundaries/sidesets to prevent the execution of flux kernels on.
Definition: FVFluxKernel.h:138
std::unordered_set< BoundaryID > _boundaries_to_force
Which boundaries/sidesets to force the execution of flux kernels on.
Definition: FVFluxKernel.h:131
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
Moose::ElemArg elemArg(bool correct_skewness=false) const
Definition: FVFluxKernel.C:239
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
Definition: FVFluxKernel.C:286
const FaceInfo * _face_info
This is holds meta-data for geometric information relevant to the current face including elem+neighbo...
Definition: FVFluxKernel.h:89
const RealVectorValue & normal() const
Definition: FVFluxKernel.h:71
Interface class for a finite volume residual object whose residuals are based on faces.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
bool avoidBoundary(const FaceInfo &fi) const
Returns whether to avoid execution on a boundary.
Definition: FVFluxKernel.C:265
const MooseVariableFV< Real > & variable() const override
Returns the variable that this object operates on.
Definition: FVFluxKernel.h:48
RealVectorValue _normal
This is the outward unit normal vector for the face the kernel is currently operating on...
Definition: FVFluxKernel.h:85
const bool _force_boundary_execution
Whether to force execution of flux kernels on all external boundaries.
Definition: FVFluxKernel.h:135
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
const ADVariableValue & _u_elem
The elem solution value of the kernel&#39;s _var for the current face.
Definition: FVFluxKernel.h:78
void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: FVFluxKernel.C:280
FVFluxKernel(const InputParameters &params)
Definition: FVFluxKernel.C:49
virtual bool skipForBoundary(const FaceInfo &fi) const
Kernels are called even on boundaries in case one is for a variable with a dirichlet BC - in which ca...
Definition: FVFluxKernel.C:91
bool onBoundary(const FaceInfo &fi) const
Return whether the supplied face is on a boundary of this object&#39;s execution.
Definition: FVFluxKernel.C:80
Enhances MooseVariableInterface interface provide values from neighbor elements.
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
static InputParameters validParams()
Definition: FVFluxKernel.C:22
A structure defining a "face" evaluation calling argument for Moose functors.
An interface for producers of functor face arguments, e.g.
const unsigned int _qp
Definition: FVFluxKernel.h:75
VarFaceNeighbors
This enum is used to indicate which side(s) of a face a particular variable is defined on...
Definition: FaceInfo.h:48
LimiterType
Definition: Limiter.h:26
A structure that is used to evaluate Moose functors logically at an element/cell center.
const ADVariableValue & _u_neighbor
The neighbor solution value of the kernel&#39;s _var for the current face.
Definition: FVFluxKernel.h:80
FaceInfo::VarFaceNeighbors _face_type
The face type.
Definition: FVFluxKernel.h:92
FVKernel is a base class for all finite volume method kernels.
Definition: FVKernel.h:32
virtual ADReal computeQpResidual()=0
This is the primary function that must be implemented for flux kernel terms.
forward declarations
Moose::ElemArg neighborArg(bool correct_skewness=false) const
Definition: FVFluxKernel.C:246
This interface is designed for DGKernel, InternalSideUserObject, InterfaceUserObject, where material properties on a side of both its primary side (face) and its secondary side (neighbor) all required.
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi=nullptr, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false, const Moose::StateArg *state_limiter=nullptr) const
Determine the single sided face argument when evaluating a functor on a face.
Definition: FVFluxKernel.C:253
State argument for evaluating functors.
MooseVariableFV< Real > & _var
Definition: FVFluxKernel.h:73
void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: FVFluxKernel.C:274
FVFluxKernel is used for calculating residual contributions from numerical fluxes from surface integr...
Definition: FVFluxKernel.h:30