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 "ADIntegratedBC.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 "ADUtils.h" 19 : 20 : #include "libmesh/quadrature.h" 21 : 22 : template <typename T> 23 : InputParameters 24 58265 : ADIntegratedBCTempl<T>::validParams() 25 : { 26 58265 : InputParameters params = IntegratedBCBase::validParams(); 27 58265 : params += ADFunctorInterface::validParams(); 28 58265 : return params; 29 0 : } 30 : 31 : template <typename T> 32 1583 : ADIntegratedBCTempl<T>::ADIntegratedBCTempl(const InputParameters & parameters) 33 : : IntegratedBCBase(parameters), 34 : MooseVariableInterface<T>(this, 35 : false, 36 : "variable", 37 : Moose::VarKindType::VAR_SOLVER, 38 : std::is_same<T, Real>::value ? Moose::VarFieldType::VAR_FIELD_STANDARD 39 : : Moose::VarFieldType::VAR_FIELD_VECTOR), 40 : ADFunctorInterface(this), 41 3166 : _var(*this->mooseVariable()), 42 1583 : _normals(_assembly.adNormals()), 43 1583 : _ad_q_points(_assembly.adQPointsFace()), 44 1583 : _test(_var.phiFace()), 45 1583 : _grad_test(_var.adGradPhiFace()), 46 1583 : _u(_var.adSln()), 47 1583 : _grad_u(_var.adGradSln()), 48 1583 : _ad_JxW(_assembly.adJxWFace()), 49 1583 : _ad_coord(_assembly.adCoordTransformation()), 50 1583 : _phi(_assembly.phi(_var)), 51 9498 : _use_displaced_mesh(getParam<bool>("use_displaced_mesh")) 52 : { 53 1583 : _subproblem.haveADObjects(true); 54 : 55 1583 : addMooseVariableDependency(this->mooseVariable()); 56 : 57 1583 : _save_in.resize(_save_in_strings.size()); 58 1583 : _diag_save_in.resize(_diag_save_in_strings.size()); 59 : 60 1583 : for (unsigned int i = 0; i < _save_in_strings.size(); i++) 61 : { 62 0 : MooseVariable * var = &_subproblem.getStandardVariable(_tid, _save_in_strings[i]); 63 : 64 0 : if (var->feType() != _var.feType()) 65 0 : paramError( 66 : "save_in", 67 : "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ", 68 0 : moose::internal::incompatVarMsg(*var, _var)); 69 0 : _save_in[i] = var; 70 0 : var->sys().addVariableToZeroOnResidual(_save_in_strings[i]); 71 0 : addMooseVariableDependency(var); 72 : } 73 : 74 1583 : _has_save_in = _save_in.size() > 0; 75 : 76 1583 : for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++) 77 : { 78 0 : MooseVariable * var = &_subproblem.getStandardVariable(_tid, _diag_save_in_strings[i]); 79 : 80 0 : if (var->feType() != _var.feType()) 81 0 : paramError( 82 : "diag_save_in", 83 : "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ", 84 0 : moose::internal::incompatVarMsg(*var, _var)); 85 : 86 0 : _diag_save_in[i] = var; 87 0 : var->sys().addVariableToZeroOnJacobian(_diag_save_in_strings[i]); 88 0 : addMooseVariableDependency(var); 89 : } 90 : 91 1583 : _has_diag_save_in = _diag_save_in.size() > 0; 92 1583 : } 93 : 94 : template <typename T> 95 : void 96 303819 : ADIntegratedBCTempl<T>::computeResidual() 97 : { 98 303819 : _residuals.resize(_test.size(), 0); 99 2701543 : for (auto & r : _residuals) 100 2397724 : r = 0; 101 : 102 303819 : precalculateResidual(); 103 1134871 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 104 : { 105 831052 : const auto jxw_c = _JxW[_qp] * _coord[_qp]; 106 7651556 : for (_i = 0; _i < _test.size(); _i++) 107 6820504 : _residuals[_i] += jxw_c * raw_value(computeQpResidual()); 108 : } 109 : 110 303819 : addResiduals(_assembly, _residuals, _var.dofIndices(), _var.scalingFactor()); 111 : 112 303819 : if (_has_save_in) 113 0 : for (unsigned int i = 0; i < _save_in.size(); i++) 114 0 : _save_in[i]->sys().solution().add_vector(_residuals.data(), _save_in[i]->dofIndices()); 115 303819 : } 116 : 117 : template <typename T> 118 : void 119 26429 : ADIntegratedBCTempl<T>::computeResidualsForJacobian() 120 : { 121 26429 : if (_residuals_and_jacobians.size() != _test.size()) 122 776 : _residuals_and_jacobians.resize(_test.size(), 0); 123 160063 : for (auto & r : _residuals_and_jacobians) 124 133634 : r = 0; 125 : 126 26429 : precalculateResidual(); 127 26429 : if (_use_displaced_mesh) 128 126 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 129 : { 130 63 : _r = _ad_JxW[_qp]; 131 63 : _r *= _ad_coord[_qp]; 132 189 : for (_i = 0; _i < _test.size(); _i++) 133 126 : _residuals_and_jacobians[_i] += _r * computeQpResidual(); 134 : } 135 : else 136 85760 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 137 : { 138 59394 : const auto jxw_c = _JxW[_qp] * _coord[_qp]; 139 384814 : for (_i = 0; _i < _test.size(); _i++) 140 325420 : _residuals_and_jacobians[_i] += jxw_c * computeQpResidual(); 141 : } 142 26429 : } 143 : 144 : template <typename T> 145 : void 146 2272 : ADIntegratedBCTempl<T>::computeResidualAndJacobian() 147 : { 148 2272 : computeResidualsForJacobian(); 149 6816 : addResidualsAndJacobian( 150 2272 : _assembly, _residuals_and_jacobians, _var.dofIndices(), _var.scalingFactor()); 151 2272 : } 152 : 153 : template <typename T> 154 : void 155 0 : ADIntegratedBCTempl<T>::computeJacobian() 156 : { 157 0 : computeADJacobian(); 158 : 159 0 : if (_has_diag_save_in && !_sys.computingScalingJacobian()) 160 0 : mooseError("_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use " 161 : "the tagging system instead."); 162 0 : } 163 : 164 : template <typename T> 165 : void 166 24157 : ADIntegratedBCTempl<T>::computeADJacobian() 167 : { 168 24157 : computeResidualsForJacobian(); 169 24157 : addJacobian(_assembly, _residuals_and_jacobians, _var.dofIndices(), _var.scalingFactor()); 170 24157 : } 171 : 172 : template <typename T> 173 : void 174 26592 : ADIntegratedBCTempl<T>::computeOffDiagJacobian(const unsigned int jvar) 175 : { 176 : // Only need to do this once because AD does all the derivatives at once 177 26592 : if (jvar == _var.number()) 178 24157 : computeADJacobian(); 179 26592 : } 180 : 181 : template <typename T> 182 : void 183 0 : ADIntegratedBCTempl<T>::computeOffDiagJacobianScalar(unsigned int /*jvar*/) 184 : { 185 0 : } 186 : 187 : template class ADIntegratedBCTempl<Real>; 188 : template class ADIntegratedBCTempl<RealVectorValue>;