https://mooseframework.inl.gov
LinearFVPressureFluxBC.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 
10 #include "LinearFVPressureFluxBC.h"
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 the H/A flux. This kernel is only designed to work with advection-diffusion "
22  "problems.");
23  params.addRequiredParam<MooseFunctorName>("HbyA_flux", "The total HbyA face flux value.");
24  params.addRequiredParam<MooseFunctorName>(
25  "Ainv", "The 1/A where A is the momentum system diagonal vector.");
26  return params;
27 }
28 
30  : LinearFVAdvectionDiffusionBC(parameters),
31  _HbyA_flux(getFunctor<Real>("HbyA_flux")),
32  _Ainv(getFunctor<RealVectorValue>("Ainv"))
33 {
34 }
35 
36 Real
38 {
39  const auto face_arg = makeCDFace(*_current_face_info);
44 
45  return _var.getElemValue(*elem_info, determineState()) -
47  std::max(_Ainv(face_arg, determineState())(0), 1e-8) *
48  distance; // We use the 0th component of Ainv. Components of Ainv are
49  // equal for most applications, and Ainv(0) has a value for 1D,2D,3D.
50 }
51 
52 Real
54 {
55  const auto face_arg = singleSidedFaceArg(_current_face_info);
56  return -_HbyA_flux(face_arg, determineState()) /
57  std::max(_Ainv(face_arg, determineState())(0), 1e-8);
58 }
59 
60 Real
62 {
63  return 1.0;
64 }
65 
66 Real
68 {
69  const auto face_arg = singleSidedFaceArg(_current_face_info);
71 
72  return -_HbyA_flux(face_arg, determineState()) /
73  std::max(_Ainv(face_arg, determineState())(0), 1e-8) * distance;
74 }
75 
76 Real
78 {
79  return 0.0;
80 }
81 
82 Real
84 {
86 }
static InputParameters validParams()
Real computeCellToFaceDistance() const
Moose::StateArg determineState() const
const ElemInfo * neighborInfo() const
virtual Real computeBoundaryNormalGradient() const override
const ElemInfo * elemInfo() const
Real distance(const Point &p)
void addRequiredParam(const std::string &name, const std::string &doc_string)
LinearFVPressureFluxBC(const InputParameters &parameters)
Class constructor.
const Moose::Functor< RealVectorValue > & _Ainv
The functor for the 1/A tensor serving as a diffusion coefficient.
FaceInfo::VarFaceNeighbors _current_face_type
Class implementing a flux boundary condition for linear finite volume pressure variables used in the ...
MooseLinearVariableFV< Real > & _var
registerMooseObject("NavierStokesApp", LinearFVPressureFluxBC)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeBoundaryValueRHSContribution() const override
const FaceInfo * _current_face_info
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
virtual Real computeBoundaryGradientMatrixContribution() const override
virtual Real computeBoundaryValueMatrixContribution() const override
void addClassDescription(const std::string &doc_string)
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false) const
virtual Real computeBoundaryGradientRHSContribution() const override
virtual Real computeBoundaryValue() const override
Moose::FaceArg makeCDFace(const FaceInfo &fi, const bool correct_skewness=false) const
static InputParameters validParams()
const Moose::Functor< Real > & _HbyA_flux
The H/A flux functor for this BC (can be variable, function, etc)