www.mooseframework.org
PenaltyDirichletBC.C
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 #include "PenaltyDirichletBC.h"
15 #include "Function.h"
16 
17 template <>
20 {
22  params.addRequiredParam<Real>("penalty", "Penalty scalar");
23  params.addParam<Real>("value", 0.0, "Boundary value of the variable");
24  params.declareControllable("value");
25  params.addClassDescription("Enforces a Dirichlet boundary condition "
26  "in a weak sense by penalizing differences between the current "
27  "solution and the Dirichlet data.");
28  return params;
29 }
30 
32  : IntegratedBC(parameters), _p(getParam<Real>("penalty")), _v(getParam<Real>("value"))
33 {
34 }
35 
36 Real
38 {
39  return _p * _test[_i][_qp] * (-_v + _u[_qp]);
40 }
41 
42 Real
44 {
45  return _p * _phi[_j][_qp] * _test[_i][_qp];
46 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:105
InputParameters validParams< IntegratedBC >()
Definition: IntegratedBC.C:28
InputParameters validParams< PenaltyDirichletBC >()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:98
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...
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:33
void declareControllable(const std::string &name)
Declare the given parameters as controllable.
unsigned int _qp
quadrature point index
Definition: IntegratedBC.h:83
unsigned int _j
Definition: IntegratedBC.h:93
virtual Real computeQpJacobian() override
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 option parameter and a documentation string to the InputParameters object...
virtual Real computeQpResidual() override
unsigned int _i
i-th, j-th index for enumerating test and shape functions
Definition: IntegratedBC.h:93
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:112
PenaltyDirichletBC(const InputParameters &parameters)