20                              "$\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, "    21                              "where $h$ is a functor.");
    23   params.
addRequiredParam<MooseFunctorName>(
"functor", 
"The functor to impose");
    25       "coefficient", 1.0, 
"An optional functor coefficient to multiply the imposed functor");
    26   params.
addParam<
bool>(
"flux_is_inward",
    28                         "Set to true if a positive evaluation of the provided functor corresponds "    29                         "to a flux in the inward direction; else the outward direction");
    36     _functor(getFunctor<
ADReal>(
"functor")),
    37     _coef(getFunctor<
ADReal>(
"coefficient")),
    38     _sign(getParam<bool>(
"flux_is_inward") ? -1.0 : 1.0)
 
Neumann boundary condition with functor inputs. 
const Elem *const  & _current_elem
current element 
unsigned int _i
i-th, j-th index for enumerating test and shape functions 
static InputParameters validParams()
DualNumber< Real, DNDerivativeType, true > ADReal
static InputParameters validParams()
unsigned int _qp
quadrature point index 
virtual ADReal computeQpResidual() override
Compute this IntegratedBC's contribution to the residual at the current quadrature point...
const MooseArray< Point > & _q_point
active quadrature points 
const Moose::Functor< ADReal > & _functor
The functor to impose. 
FunctorNeumannBC(const InputParameters ¶meters)
const Moose::Functor< ADReal > & _coef
Coefficient. 
Base class for deriving any boundary condition of a integrated type. 
const QBase *const  & _qrule
active quadrature rule 
const unsigned int & _current_side
current side of the current element 
const ADTemplateVariableTestValue< T > & _test
test function values (in QPs) 
registerMooseObject("MooseApp", FunctorNeumannBC)
const Real _sign
Sign to apply to flux. 
Argument for requesting functor evaluation at quadrature point locations on an element side...