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 "HFEMDirichletBC.h" 11 : 12 : registerMooseObject("MooseApp", HFEMDirichletBC); 13 : 14 : InputParameters 15 14406 : HFEMDirichletBC::validParams() 16 : { 17 14406 : InputParameters params = LowerDIntegratedBC::validParams(); 18 14406 : params.addParam<RealEigenVector>("value", "Value of the BC"); 19 14406 : params.addCoupledVar("uhat", "The coupled variable"); 20 14406 : params.addClassDescription("Imposes the Dirichlet BC with HFEM."); 21 14406 : return params; 22 0 : } 23 : 24 75 : HFEMDirichletBC::HFEMDirichletBC(const InputParameters & parameters) 25 : : LowerDIntegratedBC(parameters), 26 71 : _value(isParamValid("value") ? getParam<Real>("value") : 0), 27 71 : _uhat_var(isParamValid("uhat") ? getVar("uhat", 0) : nullptr), 28 146 : _uhat(_uhat_var ? (_is_implicit ? &_uhat_var->slnLower() : &_uhat_var->slnLowerOld()) : nullptr) 29 : { 30 71 : if (_uhat_var) 31 : { 32 156 : for (const auto & id : _uhat_var->activeSubdomains()) 33 85 : if (_mesh.boundaryLowerDBlocks().count(id) == 0) 34 0 : mooseDocumentedError("moose", 35 : 29151, 36 0 : "'uhat' must be defined on lower-dimensional boundary subdomain '" + 37 0 : _mesh.getSubdomainName(id) + 38 : "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe " 39 : "check could be overly restrictive."); 40 : 41 71 : if (isParamValid("value")) 42 0 : paramError("uhat", "'uhat' and 'value' can not be both provided"); 43 : } 44 71 : } 45 : 46 : Real 47 489600 : HFEMDirichletBC::computeQpResidual() 48 : { 49 489600 : return _lambda[_qp] * _test[_i][_qp]; 50 : } 51 : 52 : Real 53 25344 : HFEMDirichletBC::computeLowerDQpResidual() 54 : { 55 25344 : if (_uhat) 56 25344 : return (_u[_qp] - (*_uhat)[_qp]) * _test_lambda[_i][_qp]; 57 : else 58 0 : return (_u[_qp] - _value) * _test_lambda[_i][_qp]; 59 : } 60 : 61 : Real 62 4809600 : HFEMDirichletBC::computeQpJacobian() 63 : { 64 4809600 : return 0; 65 : } 66 : 67 : Real 68 502272 : HFEMDirichletBC::computeLowerDQpJacobian(Moose::ConstraintJacobianType type) 69 : { 70 502272 : switch (type) 71 : { 72 244800 : case Moose::LowerPrimary: 73 244800 : return _test_lambda[_i][_qp] * _phi[_j][_qp]; 74 : 75 244800 : case Moose::PrimaryLower: 76 244800 : return _phi_lambda[_j][_qp] * _test[_i][_qp]; 77 : 78 12672 : default: 79 12672 : break; 80 : } 81 : 82 12672 : return 0; 83 : } 84 : 85 : Real 86 257472 : HFEMDirichletBC::computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type, 87 : const MooseVariableFEBase & jvar) 88 : { 89 257472 : if (_uhat_var && jvar.number() == _uhat_var->number() && type == Moose::LowerLower) 90 12672 : return -_test_lambda[_i][_qp] * _phi[_j][_qp]; 91 : else 92 244800 : return 0; 93 : }