https://mooseframework.inl.gov
ConservativeAdvectionBC.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 #include "Function.h"
12 
15 
16 template <bool is_ad>
19 {
21  params.addParam<MaterialPropertyName>(
22  "velocity_mat_prop",
23  "Velocity vector as a material property. Should be provided when we want the velocity value "
24  "to be determined implicitly (e.g. we don't have a Dirichlet condition)");
25  params.addParam<FunctionName>("velocity_function",
26  "Function describing the values of velocity on the boundary.");
27  params.addClassDescription(
28  "Boundary condition for advection when it is integrated by parts. Supports Dirichlet "
29  "(inlet-like) and implicit (outlet-like) conditions.");
30  params.addParam<MaterialPropertyName>("advected_quantity",
31  "An optional material property to be advected. If not "
32  "supplied, then the variable will be used.");
33  params.addParam<FunctionName>("primal_dirichlet_value",
34  "The value of the primal variable on the boundary.");
35  params.addParam<MaterialPropertyName>(
36  "primal_coefficient",
37  1.0,
38  "If a primal Dirichlet value is supplied, then a coefficient may be optionally multiplied "
39  "that multiples the Dirichlet value");
40  return params;
41 }
42 
45 {
47 }
48 
49 template <bool is_ad>
51  const InputParameters & parameters)
52  : GenericIntegratedBC<is_ad>(parameters),
53  _velocity_mat_prop(isParamValid("velocity_mat_prop")
54  ? &this->template getGenericMaterialProperty<RealVectorValue, is_ad>(
55  "velocity_mat_prop")
56  : nullptr),
57  _velocity_function(isParamValid("velocity_function") ? &getFunction("velocity_function")
58  : nullptr),
59  _user_supplied_adv_quant(isParamValid("advected_quantity")),
60  _adv_quant(
61  _user_supplied_adv_quant
62  ? this->template getGenericMaterialProperty<Real, is_ad>("advected_quantity").get()
63  : _u),
64  _primal_dirichlet(
65  isParamValid("primal_dirichlet_value") ? &getFunction("primal_dirichlet_value") : nullptr),
66  _primal_coeff(this->template getGenericMaterialProperty<Real, is_ad>("primal_coefficient"))
67 {
68  if (isParamSetByUser("primal_coefficient") && !_primal_dirichlet)
69  paramError("primal_coefficient",
70  "This parameter should only be provided when 'primal_dirichlet_value' is provided");
71  if (static_cast<bool>(_primal_dirichlet) + isParamValid("advected_quantity") > 1)
72  mooseError("Only one of 'primal_dirichlet_value' or 'advected_quantity' should be provided");
73  if (static_cast<bool>(_velocity_mat_prop) + static_cast<bool>(_velocity_function) != 1)
74  mooseError("Exactly one of 'velocity_mat_prop' or 'velocity_function' should be provided");
75 }
76 
78  : ConservativeAdvectionBCTempl<false>(parameters)
79 {
80 }
81 
82 template <bool is_ad>
85 {
86  const auto vdotn =
87  (_velocity_mat_prop
88  ? (*_velocity_mat_prop)[_qp]
89  : GenericRealVectorValue<is_ad>(_velocity_function->vectorValue(_t, _q_point[_qp]))) *
90  this->_normals[_qp];
91 
92  if (_primal_dirichlet)
93  return _test[_i][_qp] * vdotn * _primal_dirichlet->value(_t, _q_point[_qp]) *
94  _primal_coeff[_qp];
95  else
96  return _test[_i][_qp] * vdotn * _adv_quant[_qp];
97 }
98 
99 Real
101 {
103  return _test[_i][_qp] *
105  ? (*_velocity_mat_prop)[_qp]
107  this->_normals[_qp] * _phi[_j][_qp];
108  return 0.0;
109 }
110 
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
Moose::GenericType< Real, is_ad > GenericReal
Definition: MooseTypes.h:649
const bool _user_supplied_adv_quant
Flag to check if user has supplied an advective quantity or not.
const MooseArray< Point > & _normals
normals at quadrature points
Definition: IntegratedBC.h:85
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:439
ConservativeAdvectionBCTempl(const InputParameters &parameters)
static InputParameters validParams()
T * get(const std::unique_ptr< T > &u)
The MooseUtils::get() specializations are used to support making forwards-compatible code changes fro...
Definition: MooseUtils.h:1133
static InputParameters validParams()
static InputParameters validParams()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _i
i-th, j-th index for enumerating test and shape functions
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:90
registerMooseObject("MooseApp", ConservativeAdvectionBC)
unsigned int _qp
quadrature point index
const MooseArray< Point > & _q_point
active quadrature points
virtual GenericReal< is_ad > computeQpResidual() override
Method for computing the residual at quadrature points.
const Function *const _velocity_function
The velocity as a function.
const Function *const _primal_dirichlet
Dirichlet value for the primal variable.
Moose::GenericType< RealVectorValue, is_ad > GenericRealVectorValue
Definition: MooseTypes.h:653
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:271
virtual RealVectorValue vectorValue(Real t, const Point &p) const
Override this to evaluate the vector function at a point (t,x,y,z), by default this returns a zero ve...
Definition: Function.C:96
ConservativeAdvectionBC(const InputParameters &parameters)
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...
const GenericMaterialProperty< RealVectorValue, is_ad > *const _velocity_mat_prop
The velocity as a material property.
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:199
virtual Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
bool isParamSetByUser(const std::string &name) const
Test if the supplied parameter is set by a user, as opposed to not set or set to default.
Definition: MooseBase.h:205
A boundary condition for when the advection term is integrated by parts.