https://mooseframework.inl.gov
FunctorDirichletBC.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 "FunctorDirichletBC.h"
11 
13 
16 {
18 
19  params.addClassDescription("Imposes the Dirichlet boundary condition "
20  "$u(t,\\vec{x})=h(t,\\vec{x})$, "
21  "where $h$ is a functor and can have complex dependencies.");
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 
27  return params;
28 }
29 
31  : ADDirichletBCBaseTempl<Real>(parameters),
32  _functor(getFunctor<ADReal>("functor")),
33  _coef(getFunctor<ADReal>("coefficient"))
34 {
35 }
36 
37 ADReal
39 {
41  const Moose::StateArg time_arg = Moose::currentState();
42  return _coef(space_arg, time_arg) * _functor(space_arg, time_arg);
43 }
registerMooseObject("MooseApp", FunctorDirichletBC)
static const std::set< SubdomainID > undefined_subdomain_connection
A static member that can be used when the connection of a node to subdomains is unknown.
FunctorDirichletBC(const InputParameters &parameters)
Base class for automatic differentiation Dirichlet BCs.
const Moose::Functor< ADReal > & _functor
The functor value to impose.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
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...
Dirichlet boundary condition with functor inputs.
static InputParameters validParams()
static InputParameters validParams()
const Moose::Functor< ADReal > & _coef
Coefficient.
virtual ADReal computeQpValue() override
Compute the value of the Dirichlet BC at the current quadrature point.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Node *const & _current_node
current node being processed
Definition: ADNodalBC.h:42
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...
State argument for evaluating functors.
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...
StateArg currentState()