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 "NonlocalIntegratedBC.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 14331 : NonlocalIntegratedBC::validParams() 23 : { 24 14331 : InputParameters params = IntegratedBC::validParams(); 25 14331 : return params; 26 : } 27 : 28 34 : NonlocalIntegratedBC::NonlocalIntegratedBC(const InputParameters & parameters) 29 34 : : IntegratedBC(parameters) 30 : { 31 34 : _mesh.errorIfDistributedMesh("NonlocalIntegratedBC"); 32 34 : mooseWarning("NonlocalIntegratedBC is a computationally expensive experimental capability used " 33 : "only for integral terms."); 34 34 : } 35 : 36 : void 37 438 : NonlocalIntegratedBC::computeJacobian() 38 : { 39 438 : prepareMatrixTag(_assembly, _var.number(), _var.number()); 40 2190 : for (_j = 0; _j < _phi.size(); 41 1752 : _j++) // looping order for _i & _j are reversed for performance improvement 42 : { 43 1752 : getUserObjectJacobian(_var.number(), _var.dofIndices()[_j]); 44 8760 : for (_i = 0; _i < _test.size(); _i++) 45 21024 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 46 14016 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpJacobian(); 47 : } 48 438 : accumulateTaggedLocalMatrix(); 49 : 50 438 : 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 438 : } 61 : 62 : void 63 868 : NonlocalIntegratedBC::computeOffDiagJacobian(const unsigned int jvar_num) 64 : { 65 868 : if (jvar_num == _var.number()) 66 434 : computeJacobian(); 67 : else 68 : { 69 434 : const auto & jvar = getVariable(jvar_num); 70 434 : prepareMatrixTag(_assembly, _var.number(), jvar_num); 71 : 72 : // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting 73 : // on the displaced mesh 74 434 : auto phi_size = jvar.dofIndices().size(); 75 : 76 2170 : for (_j = 0; _j < phi_size; 77 1736 : _j++) // looping order for _i & _j are reversed for performance improvement 78 : { 79 1736 : getUserObjectJacobian(jvar_num, jvar.dofIndices()[_j]); 80 8680 : for (_i = 0; _i < _test.size(); _i++) 81 20832 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 82 13888 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobian(jvar.number()); 83 : } 84 434 : accumulateTaggedLocalMatrix(); 85 : } 86 868 : } 87 : 88 : void 89 0 : NonlocalIntegratedBC::computeNonlocalJacobian() 90 : { 91 0 : prepareMatrixTagNonlocal(_assembly, _var.number(), _var.number()); 92 : // compiling set of global IDs for the local DOFs on the element 93 0 : std::set<dof_id_type> local_dofindices(_var.dofIndices().begin(), _var.dofIndices().end()); 94 : // storing the global IDs for all the DOFs of the variable 95 0 : const std::vector<dof_id_type> & var_alldofindices = _var.allDofIndices(); 96 0 : unsigned int n_total_dofs = var_alldofindices.size(); 97 : 98 0 : for (_k = 0; _k < n_total_dofs; 99 0 : _k++) // looping order for _i & _k are reversed for performance improvement 100 : { 101 : // eliminating the local components 102 0 : auto it = local_dofindices.find(var_alldofindices[_k]); 103 0 : if (it == local_dofindices.end()) 104 : { 105 0 : getUserObjectJacobian(_var.number(), var_alldofindices[_k]); 106 : // skip global DOFs that do not contribute to the jacobian 107 0 : if (!globalDoFEnabled(_var, var_alldofindices[_k])) 108 0 : continue; 109 : 110 0 : for (_i = 0; _i < _test.size(); _i++) 111 0 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 112 0 : _nonlocal_ke(_i, _k) += 113 0 : _JxW[_qp] * _coord[_qp] * computeQpNonlocalJacobian(var_alldofindices[_k]); 114 : } 115 : } 116 0 : accumulateTaggedNonlocalMatrix(); 117 0 : } 118 : 119 : void 120 434 : NonlocalIntegratedBC::computeNonlocalOffDiagJacobian(unsigned int jvar) 121 : { 122 434 : if (jvar == _var.number()) 123 0 : computeNonlocalJacobian(); 124 : else 125 : { 126 434 : MooseVariableFEBase & jv = _sys.getVariable(_tid, jvar); 127 434 : prepareMatrixTagNonlocal(_assembly, _var.number(), jvar); 128 : // compiling set of global IDs for the local DOFs on the element 129 434 : std::set<dof_id_type> local_dofindices(jv.dofIndices().begin(), jv.dofIndices().end()); 130 : // storing the global IDs for all the DOFs of the variable 131 434 : const std::vector<dof_id_type> & jv_alldofindices = jv.allDofIndices(); 132 434 : unsigned int n_total_dofs = jv_alldofindices.size(); 133 : 134 51380 : for (_k = 0; _k < n_total_dofs; 135 50946 : _k++) // looping order for _i & _k are reversed for performance improvement 136 : { 137 : // eliminating the local components 138 50946 : auto it = local_dofindices.find(jv_alldofindices[_k]); 139 50946 : if (it == local_dofindices.end()) 140 : { 141 49210 : getUserObjectJacobian(jvar, jv_alldofindices[_k]); 142 : // skip global DOFs that do not contribute to the jacobian 143 49210 : if (!globalDoFEnabled(jv, jv_alldofindices[_k])) 144 0 : continue; 145 : 146 246050 : for (_i = 0; _i < _test.size(); _i++) 147 590520 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 148 787360 : _nonlocal_ke(_i, _k) += _JxW[_qp] * _coord[_qp] * 149 393680 : computeQpNonlocalOffDiagJacobian(jvar, jv_alldofindices[_k]); 150 : } 151 : } 152 434 : accumulateTaggedNonlocalMatrix(); 153 434 : } 154 434 : }