https://mooseframework.inl.gov
VectorCurlPenaltyDirichletBC.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>("function",
21  "The boundary condition vector function, "
22  "use as an alternative to a component-wise specification");
23  params.addParam<FunctionName>("function_x", 0, "The function for the x component");
24  params.addParam<FunctionName>("function_y", 0, "The function for the y component");
25  params.addParam<FunctionName>("function_z", 0, "The function for the z component");
26  params.addClassDescription("Enforces a Dirichlet boundary condition for the curl of vector "
27  "nonlinear variables in a weak sense by applying a penalty to the "
28  "difference in the current solution and the Dirichlet data.");
29  return params;
30 }
31 
33  : VectorIntegratedBC(parameters),
34  _penalty(getParam<Real>("penalty")),
35  _function(isParamValid("function") ? &getFunction("function") : nullptr),
36  _function_x(getFunction("function_x")),
37  _function_y(getFunction("function_y")),
38  _function_z(getFunction("function_z"))
39 {
40 }
41 
42 Real
44 {
45  RealVectorValue u_exact;
46  if (_function)
47  u_exact = _function->vectorValue(_t, _q_point[_qp]);
48  else
49  u_exact = {_function_x.value(_t, _q_point[_qp]),
52  RealVectorValue Ncu = (_u[_qp] - u_exact).cross(_normals[_qp]);
53  return _penalty * Ncu * _test[_i][_qp].cross(_normals[_qp]);
54 }
55 
56 Real
58 {
59  return _penalty * _phi[_j][_qp].cross(_normals[_qp]) * _test[_i][_qp].cross(_normals[_qp]);
60 }
const Function & _function_z
The function 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
const Function *const _function
The vector function for the exact solution.
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
Enforces a Dirichlet boundary condition for the curl of vector nonlinear variables in a weak sense by...
const MooseArray< Point > & _q_point
active quadrature points
virtual Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
Base class for deriving any boundary condition of a integrated type.
VectorCurlPenaltyDirichletBC(const InputParameters &parameters)
const MooseArray< Point > & _normals
normals at quadrature points
const VectorVariableTestValue & _test
test function values (in QPs)
const Function & _function_x
The function for the x component.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Function & _function_y
The function for the y component.
virtual RealVectorValue vectorValue(Real t, const Point &p) const
Override this to evaluate the vector function at a point (t,x,y,z), by default this returns a zero ve...
Definition: Function.C:96
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()
Real _penalty
The penalty coefficient.
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
registerMooseObject("MooseApp", VectorCurlPenaltyDirichletBC)
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
const VectorVariablePhiValue & _phi
shape function values (in QPs)