https://mooseframework.inl.gov
HFEMDirichletBC.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 
10 #include "HFEMDirichletBC.h"
11 
13 
16 {
18  params.addParam<RealEigenVector>("value", "Value of the BC");
19  params.addCoupledVar("uhat", "The coupled variable");
20  params.addClassDescription("Imposes the Dirichlet BC with HFEM.");
21  return params;
22 }
23 
25  : LowerDIntegratedBC(parameters),
26  _value(isParamValid("value") ? getParam<Real>("value") : 0),
27  _uhat_var(isParamValid("uhat") ? getVar("uhat", 0) : nullptr),
28  _uhat(_uhat_var ? (_is_implicit ? &_uhat_var->slnLower() : &_uhat_var->slnLowerOld()) : nullptr)
29 {
30  if (_uhat_var)
31  {
32  for (const auto & id : _uhat_var->activeSubdomains())
33  if (_mesh.boundaryLowerDBlocks().count(id) == 0)
34  mooseDocumentedError("moose",
35  29151,
36  "'uhat' must be defined on lower-dimensional boundary subdomain '" +
38  "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe "
39  "check could be overly restrictive.");
40 
41  if (isParamValid("value"))
42  paramError("uhat", "'uhat' and 'value' can not be both provided");
43  }
44 }
45 
46 Real
48 {
49  return _lambda[_qp] * _test[_i][_qp];
50 }
51 
52 Real
54 {
55  if (_uhat)
56  return (_u[_qp] - (*_uhat)[_qp]) * _test_lambda[_i][_qp];
57  else
58  return (_u[_qp] - _value) * _test_lambda[_i][_qp];
59 }
60 
61 Real
63 {
64  return 0;
65 }
66 
67 Real
69 {
70  switch (type)
71  {
73  return _test_lambda[_i][_qp] * _phi[_j][_qp];
74 
76  return _phi_lambda[_j][_qp] * _test[_i][_qp];
77 
78  default:
79  break;
80  }
81 
82  return 0;
83 }
84 
85 Real
87  const MooseVariableFEBase & jvar)
88 {
89  if (_uhat_var && jvar.number() == _uhat_var->number() && type == Moose::LowerLower)
90  return -_test_lambda[_i][_qp] * _phi[_j][_qp];
91  else
92  return 0;
93 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
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
unsigned int number() const
Get variable number coming from libMesh.
virtual Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
HFEMDirichletBC(const InputParameters &parameters)
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 mooseDocumentedError(const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
Definition: MooseBase.h:273
This class provides an interface for common operations on field variables of both FE and FV types wit...
const VariableValue *const _uhat
Holds the coupled solution at the current quadrature point on the face.
const VariableTestValue & _test_lambda
test functions
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:90
const std::string & getSubdomainName(SubdomainID subdomain_id) const
Return the name of a block given an id.
Definition: MooseMesh.C:1763
unsigned int _qp
quadrature point index
const MooseVariable *const _uhat_var
Variable coupled in.
virtual Real computeLowerDQpResidual() override
Method for computing the Lower part of residual at quadrature points.
registerMooseObject("MooseApp", HFEMDirichletBC)
virtual Real computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &jvar) override
Method for computing an off-diagonal jacobian component at quadrature points.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
const VariableValue & _lambda
Holds the current solution at the current quadrature point on the face.
const std::set< SubdomainID > & boundaryLowerDBlocks() const
Definition: MooseMesh.h:1407
static InputParameters validParams()
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
const std::set< SubdomainID > & activeSubdomains() const
The subdomains the variable is active on.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
const VariablePhiValue & _phi_lambda
Shape functions.
ConstraintJacobianType
Definition: MooseTypes.h:797
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
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...
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
const Real _value
Boundary values.
virtual Real computeLowerDQpJacobian(Moose::ConstraintJacobianType type) override
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian at quadrature po...
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:104