www.mooseframework.org
FunctionPenaltyDirichletBC.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 /****************************************************************/
15 #include "Function.h"
16 
17 template <>
20 {
22  params.addRequiredParam<Real>("penalty", "Penalty scalar");
23  params.addRequiredParam<FunctionName>("function", "Forcing function");
24 
25  return params;
26 }
27 
29  : IntegratedBC(parameters), _func(getFunction("function")), _p(getParam<Real>("penalty"))
30 {
31 }
32 
33 Real
35 {
36  return _p * _test[_i][_qp] * (-_func.value(_t, _q_point[_qp]) + _u[_qp]);
37 }
38 
39 Real
41 {
42  return _p * _phi[_j][_qp] * _test[_i][_qp];
43 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:105
virtual Real value(Real t, const Point &p)
Override this to evaluate the scalar function at point (t,x,y,z), by default this returns zero...
Definition: Function.C:43
InputParameters validParams< IntegratedBC >()
Definition: IntegratedBC.C:28
InputParameters validParams< FunctionPenaltyDirichletBC >()
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...
virtual Real computeQpJacobian() override
FunctionPenaltyDirichletBC(const InputParameters &parameters)
Factory constructor, takes parameters so that all derived classes can be built using the same constru...
virtual Real computeQpResidual() override
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:33
unsigned int _qp
quadrature point index
Definition: IntegratedBC.h:83
unsigned int _j
Definition: IntegratedBC.h:93
const MooseArray< Point > & _q_point
active quadrature points
Definition: IntegratedBC.h:87
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