https://mooseframework.inl.gov
LinearFVAdvectionDiffusionFunctorDirichletBC.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 dirichlet BC which can be used for the assembly of linear "
20  "finite volume system and whose face values are determined using a functor. This kernel is "
21  "only designed to work with advection-diffusion problems.");
22  params.addRequiredParam<MooseFunctorName>("functor", "The functor for this boundary condition.");
23  return params;
24 }
25 
27  const InputParameters & parameters)
28  : LinearFVAdvectionDiffusionBC(parameters), _functor(getFunctor<Real>("functor"))
29 {
30 }
31 
32 Real
34 {
36 }
37 
38 Real
40 {
46  raw_value(_var(elem_arg, determineState()))) /
47  distance;
48 }
49 
50 Real
52 {
53  // Ths will not contribute to the matrix from the value considering that
54  // the value is independent of the solution.
55  return 0.0;
56 }
57 
58 Real
60 {
61  // Fetch the boundary value from the provided functor.
63 }
64 
65 Real
67 {
68  // The implicit term from the central difference approximation of the normal
69  // gradient.
70  return 1.0 / computeCellToFaceDistance();
71 }
72 
73 Real
75 {
76  // The boundary term from the central difference approximation of the
77  // normal gradient.
80 }
virtual Real computeBoundaryValue() const override
Computes the boundary value of this object.
Real computeCellToFaceDistance() const
Compute the distance between the cell center and the face.
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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...
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
Helper method to create an elemental argument for a functor that includes whether to perform skewness...
virtual Real computeBoundaryGradientRHSContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system right hand side.
registerMooseObject("MooseApp", LinearFVAdvectionDiffusionFunctorDirichletBC)
Base class for boundary conditions that are valid for advection diffusion problems.
const Elem * neighborPtr() const
Definition: FaceInfo.h:84
FaceInfo::VarFaceNeighbors _current_face_type
Face ownership information for the current face.
virtual Real computeBoundaryGradientMatrixContribution() const override
Computes the boundary gradient&#39;s contribution to the linear system matrix.
MooseLinearVariableFV< Real > & _var
Reference to the linear finite volume variable object.
const Elem * elemPtr() const
Definition: FaceInfo.h:82
const Moose::Functor< Real > & _functor
The functor for this BC (can be variable, function, etc)
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.
LinearFVAdvectionDiffusionFunctorDirichletBC(const InputParameters &parameters)
Class constructor.
virtual Real computeBoundaryNormalGradient() const override
Computes the normal gradient (often used in diffusion terms) on the boundary.
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...
virtual Real computeBoundaryValueMatrixContribution() const override
Computes the boundary value&#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.
Class implementing a Dirichlet boundary condition for linear finite volume variables.
virtual Real computeBoundaryValueRHSContribution() const override
Computes the boundary value&#39;s contribution to the linear system right hand side.