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 "NonlocalKernel.h" 11 : #include "Assembly.h" 12 : #include "MooseVariableFE.h" 13 : #include "Problem.h" 14 : #include "SubProblem.h" 15 : #include "SystemBase.h" 16 : #include "MooseMesh.h" 17 : 18 : #include "libmesh/threads.h" 19 : #include "libmesh/quadrature.h" 20 : 21 : InputParameters 22 42877 : NonlocalKernel::validParams() 23 : { 24 42877 : InputParameters params = Kernel::validParams(); 25 42877 : return params; 26 : } 27 : 28 44 : NonlocalKernel::NonlocalKernel(const InputParameters & parameters) : Kernel(parameters) 29 : { 30 44 : _mesh.errorIfDistributedMesh("NonlocalKernel"); 31 44 : mooseWarning("NonlocalKernel is a computationally expensive experimental capability used only " 32 : "for integral terms."); 33 44 : } 34 : 35 : void 36 238 : NonlocalKernel::computeJacobian() 37 : { 38 238 : prepareMatrixTag(_assembly, _var.number(), _var.number()); 39 238 : precalculateJacobian(); 40 1386 : for (_j = 0; _j < _phi.size(); 41 1148 : _j++) // looping order for _i & _j are reversed for performance improvement 42 : { 43 1148 : getUserObjectJacobian(_var.number(), _var.dofIndices()[_j]); 44 8820 : for (_i = 0; _i < _test.size(); _i++) 45 66024 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 46 58352 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpJacobian(); 47 : } 48 238 : accumulateTaggedLocalMatrix(); 49 : 50 238 : if (_has_diag_save_in) 51 : { 52 0 : unsigned int rows = _local_ke.m(); 53 0 : DenseVector<Number> diag(rows); 54 0 : for (unsigned int i = 0; i < rows; i++) 55 0 : diag(i) = _local_ke(i, i); 56 : 57 0 : for (const auto & var : _diag_save_in) 58 0 : var->sys().solution().add_vector(diag, var->dofIndices()); 59 0 : } 60 238 : } 61 : 62 : void 63 364 : NonlocalKernel::computeOffDiagJacobian(const unsigned int jvar_num) 64 : { 65 364 : if (jvar_num == _var.number()) 66 224 : computeJacobian(); 67 : else 68 : { 69 140 : const auto & jvar = getVariable(jvar_num); 70 : 71 140 : prepareMatrixTag(_assembly, _var.number(), jvar_num); 72 : 73 : // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting 74 : // on the displaced mesh 75 140 : const auto phi_size = jvar.dofIndices().size(); 76 : 77 140 : precalculateOffDiagJacobian(jvar_num); 78 420 : for (_j = 0; _j < phi_size; 79 280 : _j++) // looping order for _i & _j are reversed for performance improvement 80 : { 81 280 : getUserObjectJacobian(jvar_num, jvar.dofIndices()[_j]); 82 840 : for (_i = 0; _i < _test.size(); _i++) 83 1680 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 84 1120 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobian(jvar_num); 85 : } 86 140 : accumulateTaggedLocalMatrix(); 87 : } 88 364 : } 89 : 90 : void 91 196 : NonlocalKernel::computeNonlocalJacobian() 92 : { 93 196 : prepareMatrixTagNonlocal(_assembly, _var.number(), _var.number()); 94 : // compiling set of global IDs for the local DOFs on the element 95 196 : std::set<dof_id_type> local_dofindices(_var.dofIndices().begin(), _var.dofIndices().end()); 96 : // storing the global IDs for all the DOFs of the variable 97 196 : const std::vector<dof_id_type> & var_alldofindices = _var.allDofIndices(); 98 196 : unsigned int n_total_dofs = var_alldofindices.size(); 99 : 100 196 : precalculateJacobian(); 101 3472 : for (_k = 0; _k < n_total_dofs; 102 3276 : _k++) // looping order for _i & _k are reversed for performance improvement 103 : { 104 : // eliminating the local components 105 3276 : auto it = local_dofindices.find(var_alldofindices[_k]); 106 3276 : if (it == local_dofindices.end()) 107 : { 108 2212 : getUserObjectJacobian(_var.number(), var_alldofindices[_k]); 109 : // skip global DOFs that do not contribute to the jacobian 110 2212 : if (!globalDoFEnabled(_var, var_alldofindices[_k])) 111 0 : continue; 112 : 113 19404 : for (_i = 0; _i < _test.size(); _i++) 114 153720 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 115 136528 : _nonlocal_ke(_i, _k) += 116 136528 : _JxW[_qp] * _coord[_qp] * computeQpNonlocalJacobian(var_alldofindices[_k]); 117 : } 118 : } 119 196 : accumulateTaggedNonlocalMatrix(); 120 196 : } 121 : 122 : void 123 336 : NonlocalKernel::computeNonlocalOffDiagJacobian(unsigned int jvar) 124 : { 125 336 : if (jvar == _var.number()) 126 196 : computeNonlocalJacobian(); 127 : else 128 : { 129 140 : MooseVariableFEBase & jv = _sys.getVariable(_tid, jvar); 130 140 : prepareMatrixTagNonlocal(_assembly, _var.number(), jvar); 131 : // compiling set of global IDs for the local DOFs on the element 132 140 : std::set<dof_id_type> local_dofindices(jv.dofIndices().begin(), jv.dofIndices().end()); 133 : // storing the global IDs for all the DOFs of the variable 134 140 : const std::vector<dof_id_type> & jv_alldofindices = jv.allDofIndices(); 135 140 : unsigned int n_total_dofs = jv_alldofindices.size(); 136 : 137 140 : precalculateOffDiagJacobian(jvar); 138 560 : for (_k = 0; _k < n_total_dofs; 139 420 : _k++) // looping order for _i & _k are reversed for performance improvement 140 : { 141 : // eliminating the local components 142 420 : auto it = local_dofindices.find(jv_alldofindices[_k]); 143 420 : if (it == local_dofindices.end()) 144 : { 145 140 : getUserObjectJacobian(jvar, jv_alldofindices[_k]); 146 : // skip global DOFs that do not contribute to the jacobian 147 140 : if (!globalDoFEnabled(jv, jv_alldofindices[_k])) 148 0 : continue; 149 : 150 420 : for (_i = 0; _i < _test.size(); _i++) 151 840 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 152 1120 : _nonlocal_ke(_i, _k) += _JxW[_qp] * _coord[_qp] * 153 560 : computeQpNonlocalOffDiagJacobian(jvar, jv_alldofindices[_k]); 154 : } 155 : } 156 140 : accumulateTaggedNonlocalMatrix(); 157 140 : } 158 336 : }