https://mooseframework.inl.gov
FunctorNeumannBC.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 "FunctorNeumannBC.h"
11 
13 
16 {
18 
19  params.addClassDescription("Imposes the integrated boundary condition "
20  "$\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, "
21  "where $h$ is a functor.");
22 
23  params.addRequiredParam<MooseFunctorName>("functor", "The functor to impose");
24  params.addParam<MooseFunctorName>(
25  "coefficient", 1.0, "An optional functor coefficient to multiply the imposed functor");
26  params.addParam<bool>("flux_is_inward",
27  true,
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");
30 
31  return params;
32 }
33 
35  : ADIntegratedBC(parameters),
36  _functor(getFunctor<ADReal>("functor")),
37  _coef(getFunctor<ADReal>("coefficient")),
38  _sign(getParam<bool>("flux_is_inward") ? -1.0 : 1.0)
39 {
40 }
41 
42 ADReal
44 {
46  return _sign * _coef(space_arg, Moose::currentState()) *
47  _functor(space_arg, Moose::currentState()) * _test[_i][_qp];
48 }
Neumann boundary condition with functor inputs.
const Elem *const & _current_elem
current element
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
static InputParameters validParams()
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
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
const Moose::Functor< ADReal > & _functor
The functor to impose.
FunctorNeumannBC(const InputParameters &parameters)
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)
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...
registerMooseObject("MooseApp", FunctorNeumannBC)
const Real _sign
Sign to apply to flux.
StateArg currentState()
Argument for requesting functor evaluation at quadrature point locations on an element side...