https://mooseframework.inl.gov
ADVectorFunctionNeumannBC.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 #include "MooseTypes.h"
13 
15 
18 {
20  params.addParam<FunctionName>("function_x", 0, "The function for the x component");
21  params.addParam<FunctionName>("function_y", 0, "The function for the y component");
22  params.addParam<FunctionName>("function_z", 0, "The function for the z component");
23  params.addClassDescription(
24  "Imposes the integrated boundary condition "
25  "$\\frac{\\partial \\vec{u}}{\\partial n} = \\vec{h}$, "
26  "where $\\vec{h}$ is a (possibly) time and space-dependent MOOSE Function.");
27  return params;
28 }
29 
31  : ADVectorIntegratedBC(parameters),
32  _function_x(getFunction("function_x")),
33  _function_y(getFunction("function_y")),
34  _function_z(getFunction("function_z"))
35 {
36 }
37 
38 ADReal
40 {
44 
45  return -_test[_i][_qp] * func_u;
46 }
ADVectorFunctionNeumannBC(const InputParameters &parameters)
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 Function & _function_z
z component function
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
static InputParameters validParams()
unsigned int _qp
quadrature point index
virtual ADReal computeQpResidual() override
Compute this IntegratedBC&#39;s contribution to the residual at the current quadrature point...
const MooseArray< Point > & _q_point
active quadrature points
static InputParameters validParams()
Base class for deriving any boundary condition of a integrated type.
const Function & _function_y
y component function
Boundary condition of a Neumann style whose value is computed by a user-defined function for vector v...
registerMooseObject("MooseApp", ADVectorFunctionNeumannBC)
const ADTemplateVariableTestValue< T > & _test
test function values (in QPs)
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...
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
const Function & _function_x
x component function