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 "NodeElemConstraint.h" 11 : 12 : // MOOSE includes 13 : #include "Assembly.h" 14 : #include "MooseMesh.h" 15 : 16 : InputParameters 17 28792 : NodeElemConstraint::validParams() 18 : { 19 28792 : InputParameters params = NodeElemConstraintBase::validParams(); 20 28792 : return params; 21 : } 22 : 23 132 : NodeElemConstraint::NodeElemConstraint(const InputParameters & parameters) 24 : : NodeElemConstraintBase(parameters), 25 264 : _u_primary(_primary_var.slnNeighbor()), 26 132 : _u_secondary(_var.dofValues()), 27 132 : _grad_phi_primary(_assembly.gradPhiNeighbor(_primary_var)), 28 132 : _grad_test_primary(_var.gradPhiNeighbor()), 29 264 : _grad_u_primary(_primary_var.gradSlnNeighbor()) 30 : { 31 132 : } 32 : 33 : void 34 6560 : NodeElemConstraint::computeResidual() 35 : { 36 6560 : _qp = 0; 37 : 38 6560 : prepareVectorTagNeighbor(_assembly, _primary_var.number()); 39 49440 : for (_i = 0; _i < _test_primary.size(); _i++) 40 42880 : _local_re(_i) += computeQpResidual(Moose::Primary); 41 6560 : accumulateTaggedLocalResidual(); 42 : 43 6560 : prepareVectorTag(_assembly, _var.number()); 44 13120 : for (_i = 0; _i < _test_secondary.size(); _i++) 45 6560 : _local_re(_i) += computeQpResidual(Moose::Secondary); 46 6560 : accumulateTaggedLocalResidual(); 47 6560 : } 48 : 49 : void 50 3280 : NodeElemConstraint::computeJacobian() 51 : { 52 3280 : getConnectedDofIndices(_var.number()); 53 : 54 3280 : _Kee.resize(_test_secondary.size(), _connected_dof_indices.size()); 55 3280 : _Kne.resize(_test_primary.size(), _connected_dof_indices.size()); 56 : 57 6560 : for (_i = 0; _i < _test_secondary.size(); _i++) 58 : // Loop over the connected dof indices so we can get all the jacobian contributions 59 38560 : for (_j = 0; _j < _connected_dof_indices.size(); _j++) 60 35280 : _Kee(_i, _j) += computeQpJacobian(Moose::SecondarySecondary); 61 : 62 3280 : prepareMatrixTagNeighbor(_assembly, _var.number(), _primary_var.number(), Moose::ElementNeighbor); 63 3280 : if (_local_ke.m() && _local_ke.n()) 64 6560 : for (_i = 0; _i < _test_secondary.size(); _i++) 65 24720 : for (_j = 0; _j < _phi_primary.size(); _j++) 66 21440 : _local_ke(_i, _j) += computeQpJacobian(Moose::SecondaryPrimary); 67 3280 : accumulateTaggedLocalMatrix(); 68 : 69 24720 : for (_i = 0; _i < _test_primary.size(); _i++) 70 : // Loop over the connected dof indices so we can get all the jacobian contributions 71 272000 : for (_j = 0; _j < _connected_dof_indices.size(); _j++) 72 250560 : _Kne(_i, _j) += computeQpJacobian(Moose::PrimarySecondary); 73 : 74 3280 : prepareMatrixTagNeighbor( 75 3280 : _assembly, _primary_var.number(), _primary_var.number(), Moose::NeighborNeighbor); 76 3280 : if (_local_ke.m() && _local_ke.n()) 77 24720 : for (_i = 0; _i < _test_primary.size(); _i++) 78 173760 : for (_j = 0; _j < _phi_primary.size(); _j++) 79 152320 : _local_ke(_i, _j) += computeQpJacobian(Moose::PrimaryPrimary); 80 3280 : accumulateTaggedLocalMatrix(); 81 3280 : } 82 : 83 : void 84 320 : NodeElemConstraint::computeOffDiagJacobian(const unsigned int jvar_num) 85 : { 86 320 : getConnectedDofIndices(jvar_num); 87 : 88 320 : _Kee.resize(_test_secondary.size(), _connected_dof_indices.size()); 89 320 : _Kne.resize(_test_primary.size(), _connected_dof_indices.size()); 90 : 91 640 : for (_i = 0; _i < _test_secondary.size(); _i++) 92 : // Loop over the connected dof indices so we can get all the jacobian contributions 93 320 : for (_j = 0; _j < _connected_dof_indices.size(); _j++) 94 0 : _Kee(_i, _j) += computeQpOffDiagJacobian(Moose::SecondarySecondary, jvar_num); 95 : 96 320 : prepareMatrixTagNeighbor(_assembly, _var.number(), jvar_num, Moose::ElementNeighbor); 97 640 : for (_i = 0; _i < _test_secondary.size(); _i++) 98 1600 : for (_j = 0; _j < _phi_primary.size(); _j++) 99 1280 : _local_ke(_i, _j) += computeQpOffDiagJacobian(Moose::SecondaryPrimary, jvar_num); 100 320 : accumulateTaggedLocalMatrix(); 101 : 102 320 : if (_Kne.m() && _Kne.n()) 103 0 : for (_i = 0; _i < _test_primary.size(); _i++) 104 : // Loop over the connected dof indices so we can get all the jacobian contributions 105 0 : for (_j = 0; _j < _connected_dof_indices.size(); _j++) 106 0 : _Kne(_i, _j) += computeQpOffDiagJacobian(Moose::PrimarySecondary, jvar_num); 107 : 108 320 : prepareMatrixTagNeighbor(_assembly, _primary_var.number(), jvar_num, Moose::NeighborNeighbor); 109 1600 : for (_i = 0; _i < _test_primary.size(); _i++) 110 6400 : for (_j = 0; _j < _phi_primary.size(); _j++) 111 5120 : _local_ke(_i, _j) += computeQpOffDiagJacobian(Moose::PrimaryPrimary, jvar_num); 112 320 : accumulateTaggedLocalMatrix(); 113 320 : } 114 : 115 : void 116 3600 : NodeElemConstraint::getConnectedDofIndices(unsigned int var_num) 117 : { 118 3600 : MooseVariableFEBase & var = _sys.getVariable(0, var_num); 119 : 120 3600 : _connected_dof_indices.clear(); 121 3600 : std::set<dof_id_type> unique_dof_indices; 122 : 123 3600 : auto node_to_elem_pair = _node_to_elem_map.find(_current_node->id()); 124 : mooseAssert(node_to_elem_pair != _node_to_elem_map.end(), "Missing entry in node to elem map"); 125 3600 : const std::vector<dof_id_type> & elems = node_to_elem_pair->second; 126 : 127 : // Get the dof indices from each elem connected to the node 128 14640 : for (const auto & cur_elem : elems) 129 : { 130 11040 : std::vector<dof_id_type> dof_indices; 131 : 132 11040 : var.getDofIndices(_mesh.elemPtr(cur_elem), dof_indices); 133 : 134 70240 : for (const auto & dof : dof_indices) 135 59200 : unique_dof_indices.insert(dof); 136 11040 : } 137 : 138 38880 : for (const auto & dof : unique_dof_indices) 139 35280 : _connected_dof_indices.push_back(dof); 140 : 141 3600 : _phi_secondary.resize(_connected_dof_indices.size()); 142 : 143 3600 : const dof_id_type current_node_var_dof_index = _sys.getVariable(0, var_num).nodalDofIndex(); 144 : 145 : // Fill up _phi_secondary so that it is 1 when j corresponds to the dof associated with this node 146 : // and 0 for every other dof 147 : // This corresponds to evaluating all of the connected shape functions at _this_ node 148 3600 : _qp = 0; 149 38880 : for (unsigned int j = 0; j < _connected_dof_indices.size(); j++) 150 : { 151 35280 : _phi_secondary[j].resize(1); 152 : 153 35280 : if (_connected_dof_indices[j] == current_node_var_dof_index) 154 3280 : _phi_secondary[j][_qp] = 1.0; 155 : else 156 32000 : _phi_secondary[j][_qp] = 0.0; 157 : } 158 3600 : } 159 : 160 : Real 161 0 : NodeElemConstraint::computeQpJacobian(Moose::ConstraintJacobianType /*type*/) 162 : { 163 0 : mooseError("Derived classes must implement computeQpJacobian."); 164 : return 0; 165 : }