Line data Source code
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 : 12 : registerMooseObject("MooseApp", ArrayHFEMDirichletBC); 13 : 14 : InputParameters 15 14365 : ArrayHFEMDirichletBC::validParams() 16 : { 17 14365 : InputParameters params = ArrayLowerDIntegratedBC::validParams(); 18 14365 : params.addParam<RealEigenVector>("value", "Value of the BC"); 19 14365 : params.addCoupledVar("uhat", "The coupled variable"); 20 14365 : params.addClassDescription("Imposes the Dirichlet BC with HFEM."); 21 14365 : return params; 22 0 : } 23 : 24 52 : ArrayHFEMDirichletBC::ArrayHFEMDirichletBC(const InputParameters & parameters) 25 : : ArrayLowerDIntegratedBC(parameters), 26 156 : _value(isParamValid("value") ? getParam<RealEigenVector>("value") 27 104 : : RealEigenVector::Zero(_count)), 28 52 : _uhat_var(isParamValid("uhat") ? getArrayVar("uhat", 0) : nullptr), 29 104 : _uhat(_uhat_var ? (_is_implicit ? &_uhat_var->slnLower() : &_uhat_var->slnLowerOld()) : nullptr) 30 : { 31 52 : if (_value.size() != _count) 32 0 : paramError( 33 0 : "value", "Number of values must equal number of variable components (", _count, ")."); 34 : 35 52 : if (_uhat_var) 36 : { 37 52 : for (const auto & id : _uhat_var->activeSubdomains()) 38 26 : if (_mesh.boundaryLowerDBlocks().count(id) == 0) 39 0 : mooseDocumentedError("moose", 40 : 29151, 41 0 : "'uhat' must be defined on lower-dimensional boundary subdomain '" + 42 0 : _mesh.getSubdomainName(id) + 43 : "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe " 44 : "check could be overly restrictive."); 45 : 46 26 : if (_uhat_var->count() != _count) 47 0 : paramError("uhat", 48 : "The number of components must be equal to the number of " 49 : "components of 'variable'"); 50 : 51 26 : if (isParamValid("value")) 52 0 : paramError("uhat", "'uhat' and 'value' can not be both provided"); 53 : } 54 52 : } 55 : 56 : void 57 30720 : ArrayHFEMDirichletBC::computeQpResidual(RealEigenVector & residual) 58 : { 59 30720 : residual += _lambda[_qp] * _test[_i][_qp]; 60 30720 : } 61 : 62 : void 63 3072 : ArrayHFEMDirichletBC::computeLowerDQpResidual(RealEigenVector & r) 64 : { 65 3072 : if (_uhat) 66 1536 : r += (_u[_qp] - (*_uhat)[_qp]) * _test_lambda[_i][_qp]; 67 : else 68 1536 : r += (_u[_qp] - _value) * _test_lambda[_i][_qp]; 69 3072 : } 70 : 71 : RealEigenVector 72 32256 : ArrayHFEMDirichletBC::computeLowerDQpJacobian(Moose::ConstraintJacobianType type) 73 : { 74 32256 : RealEigenVector r = RealEigenVector::Zero(_count); 75 32256 : switch (type) 76 : { 77 15360 : case Moose::LowerPrimary: 78 30720 : return RealEigenVector::Constant(_count, _test_lambda[_i][_qp] * _phi[_j][_qp]); 79 : 80 15360 : case Moose::PrimaryLower: 81 30720 : return RealEigenVector::Constant(_count, _phi_lambda[_j][_qp] * _test[_i][_qp]); 82 : 83 1536 : default: 84 1536 : break; 85 : } 86 : 87 1536 : return r; 88 32256 : } 89 : 90 : RealEigenMatrix 91 36480 : ArrayHFEMDirichletBC::computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type, 92 : const MooseVariableFEBase & jvar) 93 : { 94 36480 : if (_uhat_var && jvar.number() == _uhat_var->number() && type == Moose::LowerLower) 95 : { 96 384 : RealEigenVector v = RealEigenVector::Constant(_count, -_test_lambda[_i][_qp] * _phi[_j][_qp]); 97 384 : RealEigenMatrix t = RealEigenMatrix::Zero(_var.count(), _var.count()); 98 384 : t.diagonal() = v; 99 384 : return t; 100 384 : } 101 : else 102 36096 : return ArrayLowerDIntegratedBC::computeLowerDQpOffDiagJacobian(type, jvar); 103 : }