https://mooseframework.inl.gov
VectorPenaltyDirichletBC.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 
11 #include "Function.h"
12 
14 
17 {
19  params.addRequiredParam<Real>("penalty", "The penalty coefficient");
20  params.addParam<FunctionName>("x_exact_sln", 0, "The exact solution for the x component");
21  params.addParam<FunctionName>("y_exact_sln", 0, "The exact solution for the y component");
22  params.addParam<FunctionName>("z_exact_sln", 0, "The exact solution for the z component");
23  params.addClassDescription("Enforces a Dirichlet boundary condition for "
24  "vector nonlinear variables in a weak sense by "
25  "applying a penalty to the difference in the "
26  "current solution and the Dirichlet data.");
27  return params;
28 }
29 
31  : VectorIntegratedBC(parameters),
32  _penalty(getParam<Real>("penalty")),
33  _exact_x(getFunction("x_exact_sln")),
34  _exact_y(getFunction("y_exact_sln")),
35  _exact_z(getFunction("z_exact_sln"))
36 {
37 }
38 
39 Real
41 {
45 
46  return _penalty * _test[_i][_qp] * (_u[_qp] - u_exact);
47 }
48 
49 Real
51 {
52  return _penalty * _test[_i][_qp] * _phi[_j][_qp];
53 }
const Function & _exact_z
The exact solution for the z component.
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
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
Real _penalty
The penalty coefficient.
virtual Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
const MooseArray< Point > & _q_point
active quadrature points
Base class for deriving any boundary condition of a integrated type.
const Function & _exact_y
The exact solution for the y component.
VectorPenaltyDirichletBC(const InputParameters &parameters)
registerMooseObject("MooseApp", VectorPenaltyDirichletBC)
const VectorVariableTestValue & _test
test function values (in QPs)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Function & _exact_x
The exact solution for the x component.
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...
static InputParameters validParams()
const VectorVariableValue & _u
the values of the unknown variable this BC is acting on
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:44
Enforces a Dirichlet boundary condition for vector nonlinear variables in a weak sense by applying a ...
static InputParameters validParams()
const VectorVariablePhiValue & _phi
shape function values (in QPs)