Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
ADFunctionNeumannBC.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 "ADFunctionNeumannBC.h"
11 
12 #include "Function.h"
13 
15 
18 {
20  params.addRequiredParam<FunctionName>("function", "The function.");
21  params.addClassDescription("Imposes the integrated boundary condition "
22  "$\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, "
23  "where $h$ is a (possibly) time and space-dependent MOOSE Function.");
24  return params;
25 }
26 
28  : ADIntegratedBC(parameters), _func(getFunction("function"))
29 {
30 }
31 
32 ADReal
34 {
35  return -_test[_i][_qp] * _func.value(_t, _q_point[_qp]);
36 }
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 & _func
The function being used for setting the value.
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
static InputParameters validParams()
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...
unsigned int _qp
quadrature point index
const MooseArray< Point > & _q_point
active quadrature points
Base class for deriving any boundary condition of a integrated type.
virtual ADReal computeQpResidual() override
Compute this IntegratedBC&#39;s contribution to the residual at the current quadrature point...
registerMooseObject("MooseApp", ADFunctionNeumannBC)
static InputParameters validParams()
const ADTemplateVariableTestValue< T > & _test
test function values (in QPs)
ADFunctionNeumannBC(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...
Boundary condition of a Neumann style whose value is computed by a user-defined function.
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