https://mooseframework.inl.gov
ArrayHFEMDirichletBC.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 "ArrayHFEMDirichletBC.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  : ArrayLowerDIntegratedBC(parameters),
26  _value(isParamValid("value") ? getParam<RealEigenVector>("value")
27  : RealEigenVector::Zero(_count)),
28  _uhat_var(isParamValid("uhat") ? getArrayVar("uhat", 0) : nullptr),
29  _uhat(_uhat_var ? (_is_implicit ? &_uhat_var->slnLower() : &_uhat_var->slnLowerOld()) : nullptr)
30 {
31  if (_value.size() != _count)
32  paramError(
33  "value", "Number of values must equal number of variable components (", _count, ").");
34 
35  if (_uhat_var)
36  {
37  for (const auto & id : _uhat_var->activeSubdomains())
38  if (_mesh.boundaryLowerDBlocks().count(id) == 0)
39  mooseDocumentedError("moose",
40  29151,
41  "'uhat' must be defined on lower-dimensional boundary subdomain '" +
43  "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe "
44  "check could be overly restrictive.");
45 
46  if (_uhat_var->count() != _count)
47  paramError("uhat",
48  "The number of components must be equal to the number of "
49  "components of 'variable'");
50 
51  if (isParamValid("value"))
52  paramError("uhat", "'uhat' and 'value' can not be both provided");
53  }
54 }
55 
56 void
58 {
59  residual += _lambda[_qp] * _test[_i][_qp];
60 }
61 
62 void
64 {
65  if (_uhat)
66  r += (_u[_qp] - (*_uhat)[_qp]) * _test_lambda[_i][_qp];
67  else
68  r += (_u[_qp] - _value) * _test_lambda[_i][_qp];
69 }
70 
73 {
74  RealEigenVector r = RealEigenVector::Zero(_count);
75  switch (type)
76  {
78  return RealEigenVector::Constant(_count, _test_lambda[_i][_qp] * _phi[_j][_qp]);
79 
81  return RealEigenVector::Constant(_count, _phi_lambda[_j][_qp] * _test[_i][_qp]);
82 
83  default:
84  break;
85  }
86 
87  return r;
88 }
89 
92  const MooseVariableFEBase & jvar)
93 {
94  if (_uhat_var && jvar.number() == _uhat_var->number() && type == Moose::LowerLower)
95  {
96  RealEigenVector v = RealEigenVector::Constant(_count, -_test_lambda[_i][_qp] * _phi[_j][_qp]);
97  RealEigenMatrix t = RealEigenMatrix::Zero(_var.count(), _var.count());
98  t.diagonal() = v;
99  return t;
100  }
101  else
103 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
virtual void computeQpResidual(RealEigenVector &residual) override
Method for computing the residual at quadrature points, to be filled in residual. ...
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
unsigned int number() const
Get variable number coming from libMesh.
const ArrayMooseVariable *const _uhat_var
Variable coupled in.
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
const ArrayVariableValue & _lambda
Holds the current solution at the current quadrature point on the face.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
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
const ArrayVariableValue & _u
the values of the unknown variable this BC is acting on
This class provides an interface for common operations on field variables of both FE and FV types wit...
const ArrayVariableValue *const _uhat
Holds the coupled solution at the current quadrature point on the face.
const unsigned int _count
Number of components of the array variable.
const ArrayVariableTestValue & _test_lambda
test functions
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 RealEigenVector _value
Boundary values.
virtual void computeLowerDQpResidual(RealEigenVector &residual) override
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
ArrayMooseVariable & _var
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
virtual RealEigenVector computeLowerDQpJacobian(Moose::ConstraintJacobianType type) override
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian at quadrature po...
virtual RealEigenMatrix computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &jvar) override
Method for computing an off-diagonal jacobian component at quadrature points.
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
const std::set< SubdomainID > & boundaryLowerDBlocks() const
Definition: MooseMesh.h:1407
ArrayHFEMDirichletBC(const InputParameters &parameters)
const ArrayVariablePhiValue & _phi_lambda
Shape functions.
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.
Base class for deriving any boundary condition of a integrated type.
const ArrayVariableTestValue & _test
test function values (in QPs)
ConstraintJacobianType
Definition: MooseTypes.h:797
virtual RealEigenMatrix computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Method for computing an off-diagonal jacobian component at quadrature points.
static InputParameters validParams()
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
registerMooseObject("MooseApp", ArrayHFEMDirichletBC)