https://mooseframework.inl.gov
PenaltyDirichletBC.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 "PenaltyDirichletBC.h"
11 
13 
16 {
18  params.addRequiredParam<Real>("penalty", "Penalty scalar");
19  params.addParam<Real>("value", 0.0, "Boundary value of the variable");
20  params.declareControllable("value");
21  params.addClassDescription("Enforces a Dirichlet boundary condition "
22  "in a weak sense by penalizing differences between the current "
23  "solution and the Dirichlet data.");
24  return params;
25 }
26 
28  : IntegratedBC(parameters), _p(getParam<Real>("penalty")), _v(getParam<Real>("value"))
29 {
30 }
31 
32 Real
34 {
35  return _p * _test[_i][_qp] * (-_v + _u[_qp]);
36 }
37 
38 Real
40 {
41  return _p * _phi[_j][_qp] * _test[_i][_qp];
42 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
static InputParameters validParams()
static InputParameters validParams()
Definition: IntegratedBC.C:23
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
registerMooseObject("MooseApp", PenaltyDirichletBC)
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
Weakly enforce a Dirichlet BC using a penalty term.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
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...
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
void declareControllable(const std::string &name, std::set< ExecFlagType > execute_flags={})
Declare the given parameters as controllable.
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:104
PenaltyDirichletBC(const InputParameters &parameters)