www.mooseframework.org
FunctionPenaltyDirichletBC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 
14 
17 {
19  params.addClassDescription(
20  "Enforces a (possibly) time and space-dependent MOOSE Function Dirichlet boundary condition "
21  "in a weak sense by penalizing differences between the current "
22  "solution and the Dirichlet data.");
23  params.addRequiredParam<Real>("penalty", "Penalty scalar");
24  params.addRequiredParam<FunctionName>("function", "Forcing function");
25 
26  return params;
27 }
28 
30  : IntegratedBC(parameters), _func(getFunction("function")), _p(getParam<Real>("penalty"))
31 {
32 }
33 
34 Real
36 {
37  return _p * _test[_i][_qp] * (-_func.value(_t, _q_point[_qp]) + _u[_qp]);
38 }
39 
40 Real
42 {
43  return _p * _phi[_j][_qp] * _test[_i][_qp];
44 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
registerMooseObject("MooseApp", FunctionPenaltyDirichletBC)
static InputParameters validParams()
Definition: IntegratedBC.C:22
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 VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:90
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 Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
const MooseArray< Point > & _q_point
active quadrature points
FunctionPenaltyDirichletBC(const InputParameters &parameters)
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
A different approach to applying Dirichlet BCs.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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...
static InputParameters validParams()
Factory constructor, takes parameters so that all derived classes can be built using the same constru...
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:41
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:104