https://mooseframework.inl.gov
LinearFVAdvectionDiffusionFunctorNeumannBC.C
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 
11 
13 
16 {
18  params.addClassDescription(
19  "Adds a fixed diffusive flux BC which can be used for the assembly of linear "
20  "finite volume system and whose normal face gradient values are determined "
21  "using a functor. This kernel is only designed to work with advection-diffusion problems.");
22  params.addRequiredParam<MooseFunctorName>(
23  "functor", "The diffusive flux functor for this boundary condition.");
24  params.addParam<MooseFunctorName>("diffusion_coeff", 1.0, "The diffusion coefficient.");
25  return params;
26 }
27 
29  const InputParameters & parameters)
30  : LinearFVAdvectionDiffusionBC(parameters),
31  _functor(getFunctor<Real>("functor")),
32  _diffusion_coeff(getFunctor<Real>("diffusion_coeff"))
33 {
35 }
36 
37 Real
39 {
40  const auto face_arg = makeCDFace(*_current_face_info);
45  const auto d_cf = computeCellToFaceVector();
46  // For non-orthogonal meshes we compute an extra correction vector to increase order accuracy
47  // correction_vector is a vector orthogonal to the boundary normal
48  const auto correction_vector =
49  (d_cf - (d_cf * _current_face_info->normal()) * _current_face_info->normal());
50  return _var.getElemValue(*elem_info, determineState()) +
53  _var.gradSln(*elem_info) * correction_vector;
54 }
55 
56 Real
58 {
59  const auto face_arg = makeCDFace(*_current_face_info);
61  _diffusion_coeff(face_arg, determineState());
62 }
63 
64 Real
66 {
67  return 1.0;
68 }
69 
70 Real
72 {
73  const auto face_arg = makeCDFace(*_current_face_info);
77 
78  // Fetch the boundary value from the provided functor.
80  const auto d_cf = computeCellToFaceVector();
81  // For non-orthogonal meshes we compute an extra correction vector to increase order accuracy
82  // correction_vector is a vector orthogonal to the boundary normal
83  const auto correction_vector =
84  (d_cf - (d_cf * _current_face_info->normal()) * _current_face_info->normal());
87  _var.gradSln(*elem_info) * correction_vector;
88 }
89 
90 Real
92 {
93  return 0.0;
94 }
95 
96 Real
98 {
100 }
RealVectorValue computeCellToFaceVector() const
Computes the vector connecting the cell and boundary face centers.
Real computeCellToFaceDistance() const
Compute the distance between the cell center and the face.
LinearFVAdvectionDiffusionFunctorNeumannBC(const InputParameters &parameters)
Class constructor.
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
const ElemInfo * neighborInfo() const
Definition: FaceInfo.h:86
virtual Real computeBoundaryGradientMatrixContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system matrix.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const ElemInfo * elemInfo() const
Definition: FaceInfo.h:85
registerMooseObject("MooseApp", LinearFVAdvectionDiffusionFunctorNeumannBC)
Real distance(const Point &p)
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
Base class for boundary conditions that are valid for advection diffusion problems.
virtual Real computeBoundaryValueMatrixContribution() const override
Computes the boundary value&#39;s contribution to the linear system matrix.
virtual Real computeBoundaryNormalGradient() const override
Computes the normal gradient (often used in diffusion terms) on the boundary.
FaceInfo::VarFaceNeighbors _current_face_type
Face ownership information for the current face.
Class implementing a Neumann boundary condition for linear finite volume variables.
virtual Real computeBoundaryValue() const override
Computes the boundary value of this object.
const Point & normal() const
Returns the unit normal vector for the face oriented outward from the face&#39;s elem element...
Definition: FaceInfo.h:68
MooseLinearVariableFV< Real > & _var
Reference to the linear finite volume variable object.
virtual Real computeBoundaryValueRHSContribution() const override
Computes the boundary value&#39;s contribution to the linear system right hand side.
const Moose::Functor< Real > & _diffusion_coeff
The functor for the diffusion coefficient.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const FaceInfo * _current_face_info
Pointer to the face info we are operating on right now.
const Moose::Functor< Real > & _functor
The functor for this BC (can be variable, function, etc)
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
Get the solution value for the provided element and seed the derivative for the corresponding dof ind...
const VectorValue< Real > gradSln(const ElemInfo &elem_info) const
Get the variable gradient at a cell center.
virtual Real computeBoundaryGradientRHSContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system right hand side.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
Determine the single sided face argument when evaluating a functor on a face.
Moose::FaceArg makeCDFace(const FaceInfo &fi, const bool correct_skewness=false) const
Make a functor face argument with a central differencing limiter, e.g.
void computeCellGradients()
Switch to request cell gradient computations.