https://mooseframework.inl.gov
ArrayPenaltyDirichletBC.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 
13 
16 {
18  params.addParam<Real>("penalty", 4, "Penalty scalar");
19  params.addRequiredParam<RealEigenVector>("value", "Boundary value of the array variable");
20  params.addClassDescription(
21  "Enforces a Dirichlet boundary condition "
22  "in a weak sense with $p(\\vec{u}^\\ast, \\vec{u} - \\vec{u}_0)$, where $p$ is the constant "
23  "scalar penalty; $\\vec{u}^\\ast$ is the test functions and $\\vec{u} - \\vec{u}_0$ is the "
24  "differences between the current solution and the Dirichlet data.");
25  return params;
26 }
27 
29  : ArrayIntegratedBC(parameters),
30  _p(getParam<Real>("penalty")),
31  _v(getParam<RealEigenVector>("value"))
32 {
33  if (_v.size() != _count)
34  paramError(
35  "value", "Number of 'values' must equal number of variable components (", _count, ").");
36 }
37 
38 void
40 {
41  residual = _p * _test[_i][_qp] * (_u[_qp] - _v);
42 }
43 
46 {
47  return RealEigenVector::Constant(_count, _p * _phi[_j][_qp] * _test[_i][_qp]);
48 }
const ArrayVariablePhiValue & _phi
shape function values (in QPs)
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
ArrayPenaltyDirichletBC(const InputParameters &parameters)
const RealEigenVector & _v
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 ArrayVariableValue & _u
the values of the unknown variable this BC is acting on
const unsigned int _count
Number of components of the array variable.
static InputParameters validParams()
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", ArrayPenaltyDirichletBC)
virtual void computeQpResidual(RealEigenVector &residual) override
Method for computing the residual at quadrature points, to be filled in residual. ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
const ArrayVariableTestValue & _test
test function values (in QPs)
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...
Base class for deriving any boundary condition of a integrated type.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
virtual RealEigenVector computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.