https://mooseframework.inl.gov
ADConservativeAdvectionBC.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 
14 
17 {
19  params.addParam<MaterialPropertyName>(
20  "velocity_mat_prop",
21  "Velocity vector as a material property. Should be provided when we want the velocity value "
22  "to be determined implicitly (e.g. we don't have a Dirichlet condition)");
23  params.addParam<FunctionName>("velocity_function",
24  "Function describing the values of velocity on the boundary.");
25  params.addClassDescription(
26  "Boundary condition for advection when it is integrated by parts. Supports Dirichlet "
27  "(inlet-like) and implicit (outlet-like) conditions.");
28  params.addParam<MaterialPropertyName>("advected_quantity",
29  "An optional material property to be advected. If not "
30  "supplied, then the variable will be used.");
31  params.addParam<FunctionName>("primal_dirichlet_value",
32  "The value of the primal variable on the boundary.");
33  params.addParam<MaterialPropertyName>(
34  "primal_coefficient",
35  1,
36  "If a primal Dirichlet value is supplied, then a coefficient may be optionally multiplied "
37  "that multiples the Dirichlet value");
38  return params;
39 }
40 
42  : ADIntegratedBC(parameters),
43  _velocity_mat_prop(isParamValid("velocity_mat_prop")
44  ? &getADMaterialProperty<RealVectorValue>("velocity_mat_prop")
45  : nullptr),
46  _velocity_function(isParamValid("velocity_function") ? &getFunction("velocity_function")
47  : nullptr),
48  _adv_quant(isParamValid("advected_quantity")
49  ? getADMaterialProperty<Real>("advected_quantity").get()
50  : _u),
51  _primal_dirichlet(
52  isParamValid("primal_dirichlet_value") ? &getFunction("primal_dirichlet_value") : nullptr),
53  _primal_coeff(getADMaterialProperty<Real>("primal_coefficient"))
54 {
55  if (isParamSetByUser("primal_coefficient") && !_primal_dirichlet)
56  paramError("primal_coefficient",
57  "This parameter should only be provided when 'primal_dirichlet_value' is provided");
58  if (static_cast<bool>(_primal_dirichlet) + isParamValid("advected_quantity") > 1)
59  mooseError("Only one of 'primal_dirichlet_value' or 'advected_quantity' should be provided");
60  if (static_cast<bool>(_velocity_mat_prop) + static_cast<bool>(_velocity_function) != 1)
61  mooseError("Exactly one of 'velocity_mat_prop' or 'velocity_function' should be provided");
62 }
63 
64 ADReal
66 {
67  const auto vdotn =
68  (_velocity_mat_prop ? (*_velocity_mat_prop)[_qp]
70  _normals[_qp];
71 
73  return _test[_i][_qp] * vdotn * _primal_dirichlet->value(_t, _q_point[_qp]) *
75  else
76  return _test[_i][_qp] * vdotn * _adv_quant[_qp];
77 }
const Function *const _velocity_function
The velocity as a function.
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:435
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:1155
const Function *const _primal_dirichlet
Dirichlet value for the primal variable.
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
virtual ADReal computeQpResidual() override
Compute this IntegratedBC&#39;s contribution to the residual at the current quadrature point...
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
static InputParameters validParams()
ADConservativeAdvectionBC(const InputParameters &parameters)
unsigned int _qp
quadrature point index
const MooseArray< Point > & _q_point
active quadrature points
registerMooseObject("MooseApp", ADConservativeAdvectionBC)
const ADMaterialProperty< Real > & _primal_coeff
Coefficient for multiplying the primal Dirichlet value.
Base class for deriving any boundary condition of a integrated type.
const ADMaterialProperty< RealVectorValue > *const _velocity_mat_prop
The velocity as a material property.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const ADTemplateVariableTestValue< T > & _test
test function values (in QPs)
const MooseArray< ADReal > & _adv_quant
The advected quantity.
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:267
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
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...
static InputParameters validParams()
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
virtual Real value(Real t, const Point &p) const
Override this to evaluate the scalar function at point (t,x,y,z), by default this returns zero...
Definition: Function.C:44
libMesh::VectorValue< ADReal > ADRealVectorValue
AD typedefs.
Definition: MooseTypes.h:369
A boundary condition for when the advection term is integrated by parts.
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:201
const MooseArray< ADPoint > & _normals
normals at quadrature points