Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
LinearFVAdvection.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 
10 #include "LinearFVAdvection.h"
11 #include "Assembly.h"
12 #include "SubProblem.h"
14 
16 
19 {
21  params.addClassDescription("Represents the matrix and right hand side contributions of an "
22  "advection term in a partial differential equation.");
23  params.addRequiredParam<RealVectorValue>("velocity", "Constant advection velocity");
25  return params;
26 }
27 
29  : LinearFVFluxKernel(params), _velocity(getParam<RealVectorValue>("velocity"))
30 
31 {
32  Moose::FV::setInterpolationMethod(*this, _advected_interp_method, "advected_interp_method");
33 }
34 
35 void
37 {
38  for (const auto bc : _var.getBoundaryConditionMap())
39  if (!dynamic_cast<const LinearFVAdvectionDiffusionBC *>(bc.second))
40  mooseError(
41  bc.second->type(), " is not a compatible boundary condition with ", this->type(), "!");
42 }
43 
44 Real
46 {
47  const Real face_flux = _velocity * _current_face_info->normal();
48  const auto interp_coeffs =
50  return interp_coeffs.first * face_flux * _current_face_area;
51 }
52 
53 Real
55 {
56  const Real face_flux = _velocity * _current_face_info->normal();
57  const auto interp_coeffs =
59  return interp_coeffs.second * face_flux * _current_face_area;
60 }
61 
62 Real
64 {
65  return 0.0;
66 }
67 
68 Real
70 {
71  return 0.0;
72 }
73 
74 Real
76 {
77  const auto * const adv_bc = static_cast<const LinearFVAdvectionDiffusionBC *>(&bc);
78  mooseAssert(adv_bc, "This should be a valid BC!");
79 
80  const auto boundary_value_matrix_contrib = adv_bc->computeBoundaryValueMatrixContribution();
81 
82  // We support internal boundaries too so we have to make sure the normal points always outward
83  const auto factor = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM) ? 1.0 : -1.0;
84 
85  return boundary_value_matrix_contrib * factor * (_velocity * _current_face_info->normal()) *
87 }
88 
89 Real
91 {
92  const auto * const adv_bc = static_cast<const LinearFVAdvectionDiffusionBC *>(&bc);
93  mooseAssert(adv_bc, "This should be a valid BC!");
94 
95  // We support internal boundaries too so we have to make sure the normal points always outward
96  const auto factor = (_current_face_type == FaceInfo::VarFaceNeighbors::ELEM ? 1.0 : -1.0);
97 
98  const auto boundary_value_rhs_contrib = adv_bc->computeBoundaryValueRHSContribution();
99  return -boundary_value_rhs_contrib * factor * (_velocity * _current_face_info->normal()) *
101 }
Base class for boundary conditions for linear FV systems.
const RealVectorValue _velocity
Constant advecting velocity vector.
Kernel that adds contributions from an advection term discretized using the finite volume method to a...
std::pair< Real, Real > interpCoeffs(const InterpMethod m, const FaceInfo &fi, const bool one_is_elem, const T &face_flux=0.0)
Produce the interpolation coefficients in the equation:
Definition: MathFVUtils.h:114
Finite volume kernel that contributes approximations of discretized face flux terms to the matrix and...
MooseLinearVariableFV< Real > & _var
Reference to the linear finite volume variable.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
FaceInfo::VarFaceNeighbors _current_face_type
Face ownership information for the current face.
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...
registerMooseObject("MooseApp", LinearFVAdvection)
Moose::FV::InterpMethod _advected_interp_method
The interpolation method to use for the advected quantity.
LinearFVAdvection(const InputParameters &params)
Class constructor.
Base class for boundary conditions that are valid for advection diffusion problems.
static InputParameters validParams()
const FaceInfo * _current_face_info
Pointer to the face info we are operating on right now.
virtual Real computeNeighborMatrixContribution() override
Computes the system matrix contribution from the neighbor side on an internal face.
InputParameters advectedInterpolationParameter()
Definition: MathFVUtils.C:68
virtual Real computeNeighborRightHandSideContribution() override
Computes the right hand side contribution from the neighbor side on an internal face.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
virtual Real computeElemMatrixContribution() override
Computes the system matrix contribution from an element side on an internal 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
virtual Real computeBoundaryMatrixContribution(const LinearFVBoundaryCondition &bc) override
Computes the matrix contribution from a boundary face.
virtual Real computeBoundaryRHSContribution(const LinearFVBoundaryCondition &bc) override
Computes the right hand side contribution from a boundary face.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
const std::unordered_map< BoundaryID, LinearFVBoundaryCondition * > & getBoundaryConditionMap()
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
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...
Real _current_face_area
The current, coordinate system specific face area.
static InputParameters validParams()
bool setInterpolationMethod(const MooseObject &obj, Moose::FV::InterpMethod &interp_method, const std::string &param_name)
Sets one interpolation method.
Definition: MathFVUtils.C:110
virtual Real computeElemRightHandSideContribution() override
Computes the right hand side contribution from the element side on an internal face.