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 "IntegratedBC.h" 11 : 12 : // MOOSE includes 13 : #include "Assembly.h" 14 : #include "SubProblem.h" 15 : #include "SystemBase.h" 16 : #include "MooseVariableFE.h" 17 : #include "MooseVariableScalar.h" 18 : #include "FEProblemBase.h" 19 : 20 : #include "libmesh/quadrature.h" 21 : 22 : InputParameters 23 541427 : IntegratedBC::validParams() 24 : { 25 541427 : InputParameters params = IntegratedBCBase::validParams(); 26 541427 : return params; 27 : } 28 : 29 7054 : IntegratedBC::IntegratedBC(const InputParameters & parameters) 30 : : IntegratedBCBase(parameters), 31 : MooseVariableInterface<Real>(this, 32 : false, 33 : "variable", 34 : Moose::VarKindType::VAR_SOLVER, 35 : Moose::VarFieldType::VAR_FIELD_STANDARD), 36 14098 : _var(*mooseVariable()), 37 7049 : _normals(_assembly.normals()), 38 7049 : _phi(_assembly.phiFace(_var)), 39 7049 : _grad_phi(_assembly.gradPhiFace(_var)), 40 7049 : _test(_var.phiFace()), 41 7049 : _grad_test(_var.gradPhiFace()), 42 7049 : _u(_is_implicit ? _var.sln() : _var.slnOld()), 43 14103 : _grad_u(_is_implicit ? _var.gradSln() : _var.gradSlnOld()) 44 : { 45 7049 : addMooseVariableDependency(mooseVariable()); 46 : 47 7049 : _save_in.resize(_save_in_strings.size()); 48 7049 : _diag_save_in.resize(_diag_save_in_strings.size()); 49 : 50 7153 : for (unsigned int i = 0; i < _save_in_strings.size(); i++) 51 : { 52 104 : MooseVariable * var = &_subproblem.getStandardVariable(_tid, _save_in_strings[i]); 53 : 54 104 : if (var->feType() != _var.feType()) 55 0 : paramError( 56 : "save_in", 57 : "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ", 58 0 : moose::internal::incompatVarMsg(*var, _var)); 59 104 : _save_in[i] = var; 60 104 : var->sys().addVariableToZeroOnResidual(_save_in_strings[i]); 61 104 : addMooseVariableDependency(var); 62 : } 63 : 64 7049 : _has_save_in = _save_in.size() > 0; 65 : 66 7101 : for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++) 67 : { 68 52 : MooseVariable * var = &_subproblem.getStandardVariable(_tid, _diag_save_in_strings[i]); 69 : 70 52 : if (var->feType() != _var.feType()) 71 0 : paramError( 72 : "diag_save_in", 73 : "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ", 74 0 : moose::internal::incompatVarMsg(*var, _var)); 75 : 76 52 : _diag_save_in[i] = var; 77 52 : var->sys().addVariableToZeroOnJacobian(_diag_save_in_strings[i]); 78 52 : addMooseVariableDependency(var); 79 : } 80 : 81 7049 : _has_diag_save_in = _diag_save_in.size() > 0; 82 7049 : } 83 : 84 : void 85 2054748 : IntegratedBC::computeResidual() 86 : { 87 2054748 : prepareVectorTag(_assembly, _var.number()); 88 : 89 2054748 : precalculateResidual(); 90 : 91 7759664 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 92 : { 93 5704916 : precalculateQpResidual(); 94 54272214 : for (_i = 0; _i < _test.size(); _i++) 95 48567298 : _local_re(_i) += _JxW[_qp] * _coord[_qp] * computeQpResidual(); 96 : } 97 : 98 2054748 : accumulateTaggedLocalResidual(); 99 : 100 2054748 : if (_has_save_in) 101 : { 102 1410 : libMesh::Threads::spin_mutex::scoped_lock lock(libMesh::Threads::spin_mtx); 103 3000 : for (unsigned int i = 0; i < _save_in.size(); i++) 104 1590 : _save_in[i]->sys().solution().add_vector(_local_re, _save_in[i]->dofIndices()); 105 1410 : } 106 2054748 : } 107 : 108 : void 109 238339 : IntegratedBC::computeJacobian() 110 : { 111 238339 : prepareMatrixTag(_assembly, _var.number(), _var.number()); 112 : 113 238339 : precalculateJacobian(); 114 : 115 981619 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 116 : { 117 743280 : precalculateQpJacobian(); 118 7075513 : for (_i = 0; _i < _test.size(); _i++) 119 94319652 : for (_j = 0; _j < _phi.size(); _j++) 120 87987419 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpJacobian(); 121 : } 122 : 123 238339 : accumulateTaggedLocalMatrix(); 124 : 125 238339 : if (_has_diag_save_in) 126 : { 127 178 : unsigned int rows = _local_ke.m(); 128 178 : DenseVector<Number> diag(rows); 129 890 : for (unsigned int i = 0; i < rows; i++) 130 712 : diag(i) = _local_ke(i, i); 131 : 132 178 : libMesh::Threads::spin_mutex::scoped_lock lock(libMesh::Threads::spin_mtx); 133 356 : for (unsigned int i = 0; i < _diag_save_in.size(); i++) 134 178 : _diag_save_in[i]->sys().solution().add_vector(diag, _diag_save_in[i]->dofIndices()); 135 178 : } 136 238339 : } 137 : 138 : void 139 121783 : IntegratedBC::computeOffDiagJacobian(const unsigned int jvar_num) 140 : { 141 121783 : const auto & jvar = getVariable(jvar_num); 142 : 143 121783 : if (jvar_num == _var.number()) 144 : { 145 108467 : computeJacobian(); 146 108467 : return; 147 : } 148 : 149 13316 : prepareMatrixTag(_assembly, _var.number(), jvar_num); 150 : 151 13316 : precalculateOffDiagJacobian(jvar_num); 152 : 153 : // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting 154 : // on the displaced mesh, so we obtain the variable on the proper system 155 13316 : auto phi_size = jvar.dofIndices().size(); 156 : mooseAssert(phi_size * jvar.count() == _local_ke.n(), 157 : "The size of the phi container does not match the number of local Jacobian columns"); 158 : 159 51834 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 160 : { 161 38518 : precalculateQpOffDiagJacobian(jvar); 162 403342 : for (_i = 0; _i < _test.size(); _i++) 163 969700 : for (_j = 0; _j < phi_size; _j++) 164 604876 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobian(jvar_num); 165 : } 166 : 167 13316 : accumulateTaggedLocalMatrix(); 168 : } 169 : 170 : void 171 72 : IntegratedBC::computeOffDiagJacobianScalar(const unsigned int jvar) 172 : { 173 72 : prepareMatrixTag(_assembly, _var.number(), jvar); 174 : 175 72 : MooseVariableScalar & jv = _sys.getScalarVariable(_tid, jvar); 176 198 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 177 594 : for (_i = 0; _i < _test.size(); _i++) 178 936 : for (_j = 0; _j < jv.order(); _j++) 179 468 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobianScalar(jvar); 180 : 181 72 : accumulateTaggedLocalMatrix(); 182 72 : } 183 : 184 : void 185 204 : IntegratedBC::computeResidualAndJacobian() 186 : { 187 204 : computeResidual(); 188 : 189 912 : for (const auto & [ivariable, jvariable] : _fe_problem.couplingEntries(_tid, _sys.number())) 190 : { 191 708 : const unsigned int ivar = ivariable->number(); 192 708 : const unsigned int jvar = jvariable->number(); 193 : 194 708 : if (ivar != _var.number()) 195 336 : continue; 196 : 197 372 : if (_is_implicit) 198 : { 199 372 : prepareShapes(jvar); 200 372 : computeOffDiagJacobian(jvar); 201 : } 202 : } 203 : 204 : /// TODO: add nonlocal Jacobians and scalar Jacobians 205 204 : }