https://mooseframework.inl.gov
LinearFVAdvectionDiffusionScalarSymmetryBC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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("Adds a symmetry boundary condition for a scalar quantity.");
19  params.addParam<bool>(
20  "use_two_term_expansion",
21  false,
22  "If an approximate linear expansion should be used to compute the face value.");
23  return params;
24 }
25 
27  const InputParameters & parameters)
28  : LinearFVAdvectionDiffusionBC(parameters),
29  _two_term_expansion(getParam<bool>("use_two_term_expansion"))
30 {
33 }
34 
35 Real
37 {
38  // We allow internal boundaries too so we need to check which side we are on
42 
43  // By default we approximate the boundary value with the neighboring cell value
44  auto boundary_value = _var.getElemValue(*elem_info, determineState());
45 
46  // If we request linear extrapolation, we add the gradient term as well. We make sure
47  // that the zero normal gradient is respected (by subtracting the normal component).
49  {
50  const auto cell_gradient = _var.gradSln(*elem_info);
51  const auto & normal = _current_face_info->normal();
52  const auto d_cf = computeCellToFaceVector();
53 
54  const auto correction_vector = (d_cf - (d_cf * normal) * normal);
55  boundary_value += cell_gradient * correction_vector;
56  }
57 
58  return boundary_value;
59 }
60 
61 Real
63 {
64  // We don't have this on a symmetry plane
65  return 0.0;
66 }
67 
68 Real
70 {
71  // No matter if we have a one-term or two-term expansion we will always
72  // have a contribution to the matrix
73  return 1.0;
74 }
75 
76 Real
78 {
79  // If we request linear extrapolation, we add the gradient term to the right hand
80  // side.
82  {
83  // We allow internal boundaries too so we need to check which side we are on
84  const auto & elem_info = _current_face_type == FaceInfo::VarFaceNeighbors::ELEM
87 
88  const auto & d_cf = computeCellToFaceVector();
89  const auto & normal = _current_face_info->normal();
90  const auto correction_vector = (d_cf - (d_cf * normal) * normal);
91  const auto & cell_gradient = _var.gradSln(*elem_info);
92 
93  return cell_gradient * correction_vector;
94  }
95 
96  return 0.0;
97 }
98 
99 Real
101 {
102  // Nothing to add here, considering that we have a symmetry condition
103  return 0.0;
104 }
105 
106 Real
108 {
109  // Nothing to add here, considering that we have a symmetry condition
110  return 0.0;
111 }
RealVectorValue computeCellToFaceVector() const
Computes the vector connecting the cell and boundary face centers.
const bool _two_term_expansion
Switch for enabling linear extrapolation for the boundary face value.
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
LinearFVAdvectionDiffusionScalarSymmetryBC(const InputParameters &parameters)
Class constructor.
virtual Real computeBoundaryValueMatrixContribution() const override
Computes the boundary value&#39;s contribution to the linear system matrix.
Base class for boundary conditions that are valid for advection diffusion problems.
virtual Real computeBoundaryValue() const override
Computes the boundary value of this object.
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&#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.
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.
virtual Real computeBoundaryGradientRHSContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system right hand side.
virtual Real computeBoundaryNormalGradient() const override
Computes the normal gradient (often used in diffusion terms) on the boundary.
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 computeBoundaryGradientMatrixContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system matrix.
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...
void computeCellGradients()
Switch to request cell gradient computations.
registerMooseObject("MooseApp", LinearFVAdvectionDiffusionScalarSymmetryBC)
Class implementing a symmetry boundary condition for scalar quantities.