https://mooseframework.inl.gov
LinearFVAdvectionDiffusionFunctorRobinBCBase.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 
14 {
16  return params;
17 }
18 
20  const InputParameters & parameters)
21  : LinearFVAdvectionDiffusionBC(parameters)
22 {
23 }
24 
25 Real
27 {
28  const auto face = singleSidedFaceArg(_current_face_info);
30  "This should not be assigned on an internal face!");
31  const auto & elem_info = _current_face_type == FaceInfo::VarFaceNeighbors::ELEM
34  const auto state = determineState();
35 
36  const auto alpha = getAlpha(face, state);
37  const auto beta = getBeta(face, state);
38  const auto gamma = getGamma(face, state);
39 
40  const auto phi = _var.getElemValue(*elem_info, state);
41  const auto grad_phi = _var.gradSln(*elem_info);
42 
43  const auto & nhat = _current_face_info->normal();
44 
45  const auto d_cf = computeCellToFaceVector(); // vector from boundary cell centre to boundary face
46  const auto projection = d_cf * nhat;
47  const auto vc = d_cf - (projection * nhat);
48  return ((alpha * phi) + (alpha * grad_phi * vc) + (gamma * projection)) /
49  (alpha + (beta * projection));
50 }
51 
52 Real
54 {
55  const auto face = singleSidedFaceArg(_current_face_info);
56  const auto state = determineState();
57  const auto alpha = getAlpha(face, state);
58  mooseAssert(!MooseUtils::isZero(alpha), "Alpha should not be 0!");
59  const auto beta = getBeta(face, state);
60  const auto gamma = getGamma(face, state);
61  return (gamma - beta * computeBoundaryValue()) / alpha;
62 }
63 
64 // implicit terms for advection kernel
65 Real
67 {
68  const auto face = singleSidedFaceArg(_current_face_info);
69  const auto state = determineState();
70  const auto alpha = getAlpha(face, state);
71  const auto beta = getBeta(face, state);
72  const auto & nhat = _current_face_info->normal();
73 
74  return alpha / (alpha + (beta * computeCellToFaceVector() * nhat));
75 }
76 
77 // explicit terms for advection kernel
78 Real
80 {
81  const auto face = singleSidedFaceArg(_current_face_info);
82  const auto state = determineState();
84  "This should not be assigned on an internal face!");
85  const auto & elem_info = _current_face_type == FaceInfo::VarFaceNeighbors::ELEM
88 
89  const auto alpha = getAlpha(face, state);
90  const auto beta = getBeta(face, state);
91  const auto gamma = getGamma(face, state);
92 
93  const auto & grad_phi = _var.gradSln(*elem_info);
94 
95  const auto & nhat = _current_face_info->normal();
96 
97  const auto d_cf = computeCellToFaceVector(); // vector from boundary cell centre to boundary face
98  const auto projection = d_cf * nhat;
99  const auto vc = d_cf - (projection * nhat); // correction vector for non-orthogonal cells
100 
101  return (gamma * projection / (alpha + (beta * projection))) +
102  (alpha * grad_phi * vc / (alpha + (beta * projection)));
103 }
104 
105 // implicit terms for diffusion kernel
106 Real
108 {
109  const auto face = singleSidedFaceArg(_current_face_info);
110  const auto state = determineState();
111 
112  const auto alpha = getAlpha(face, state);
113  const auto beta = getBeta(face, state);
114 
115  const auto & nhat = _current_face_info->normal();
116 
117  return beta / (alpha + (beta * computeCellToFaceVector() * nhat));
118 }
119 
120 // explicit terms for diffusion kernel
121 Real
123 {
125  "This should not be assigned on an internal face!");
126  const auto & elem_info = _current_face_type == FaceInfo::VarFaceNeighbors::ELEM
129  const auto face = singleSidedFaceArg(_current_face_info);
130  const auto state = determineState();
131  const auto & grad_phi = _var.gradSln(*elem_info);
132 
133  const auto alpha = getAlpha(face, state);
134  const auto beta = getBeta(face, state);
135  const auto gamma = getGamma(face, state);
136 
137  const auto & nhat = _current_face_info->normal();
138 
139  const auto d_cf = computeCellToFaceVector(); // vector from boundary cell centre to boundary face
140  const auto projection = d_cf * nhat;
141  const auto vc = d_cf - (projection * nhat); // correction vector for non-orthogonal cells
142 
143  return (gamma / alpha) + (-beta * gamma * projection / alpha / (alpha + (beta * projection))) +
144  (-beta * grad_phi * vc / (alpha + (beta * projection)));
145 }
RealVectorValue computeCellToFaceVector() const
Computes the vector connecting the cell and boundary face centers.
virtual Real getGamma(Moose::FaceArg face, Moose::StateArg state) const =0
virtual Real computeBoundaryValueMatrixContribution() const override
Computes the boundary value's contribution to the linear system matrix.
virtual Real computeBoundaryNormalGradient() const override
Computes the normal gradient (often used in diffusion terms) on the boundary.
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
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const ElemInfo * elemInfo() const
Definition: FaceInfo.h:85
virtual Real getBeta(Moose::FaceArg face, Moose::StateArg state) const =0
Base class for boundary conditions that are valid for advection diffusion problems.
virtual Real computeBoundaryValue() const override
Computes the boundary value of this object.
virtual Real computeBoundaryGradientRHSContribution() const override
Computes the boundary gradient's contribution to the linear system right hand side.
virtual Real computeBoundaryValueRHSContribution() const override
Computes the boundary value's contribution to the linear system right hand side.
FaceInfo::VarFaceNeighbors _current_face_type
Face ownership information for the current face.
const Point & normal() const
Returns the unit normal vector for the face oriented outward from the face's elem element...
Definition: FaceInfo.h:68
virtual Real getAlpha(Moose::FaceArg face, Moose::StateArg state) const =0
Getter functions (consistent entry point for all derived classes)
MooseLinearVariableFV< Real > & _var
Reference to the linear finite volume variable object.
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.
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.
LinearFVAdvectionDiffusionFunctorRobinBCBase(const InputParameters &parameters)
Class constructor.
virtual Real computeBoundaryGradientMatrixContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system matrix.
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.