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 "MortarConstraint.h" 11 : 12 : // MOOSE includes 13 : #include "Assembly.h" 14 : #include "MooseVariable.h" 15 : #include "SystemBase.h" 16 : 17 : InputParameters 18 72093 : MortarConstraint::validParams() 19 : { 20 72093 : InputParameters params = MortarConstraintBase::validParams(); 21 72093 : return params; 22 : } 23 : 24 383 : MortarConstraint::MortarConstraint(const InputParameters & parameters) 25 : : MortarConstraintBase(parameters), 26 383 : _lambda_dummy(), 27 383 : _lambda(_var ? _var->slnLower() : _lambda_dummy), 28 383 : _u_secondary(_secondary_var.sln()), 29 383 : _u_primary(_primary_var.slnNeighbor()), 30 383 : _grad_u_secondary(_secondary_var.gradSln()), 31 383 : _grad_u_primary(_primary_var.gradSlnNeighbor()), 32 383 : _phi(nullptr), 33 383 : _grad_phi(nullptr) 34 : { 35 383 : } 36 : 37 : void 38 114284 : MortarConstraint::computeResidual(Moose::MortarType mortar_type) 39 : { 40 114284 : unsigned int test_space_size = 0; 41 114284 : switch (mortar_type) 42 : { 43 55556 : case Moose::MortarType::Secondary: 44 55556 : prepareVectorTag(_assembly, _secondary_var.number()); 45 55556 : test_space_size = _test_secondary.size(); 46 55556 : break; 47 : 48 55556 : case Moose::MortarType::Primary: 49 55556 : prepareVectorTagNeighbor(_assembly, _primary_var.number()); 50 55556 : test_space_size = _test_primary.size(); 51 55556 : break; 52 : 53 3172 : case Moose::MortarType::Lower: 54 : mooseAssert(_var, "LM variable is null"); 55 3172 : prepareVectorTagLower(_assembly, _var->number()); 56 3172 : test_space_size = _test.size(); 57 3172 : break; 58 : } 59 : 60 717584 : for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++) 61 9195988 : for (_i = 0; _i < test_space_size; _i++) 62 8592688 : _local_re(_i) += _JxW_msm[_qp] * _coord[_qp] * computeQpResidual(mortar_type); 63 : 64 114284 : accumulateTaggedLocalResidual(); 65 114284 : } 66 : 67 : void 68 66452 : MortarConstraint::computeJacobian(Moose::MortarType mortar_type) 69 : { 70 66452 : std::size_t test_space_size = 0; 71 : typedef Moose::ConstraintJacobianType JType; 72 : typedef Moose::MortarType MType; 73 : std::array<JType, 3> jacobian_types; 74 : 75 66452 : switch (mortar_type) 76 : { 77 32212 : case MType::Secondary: 78 32212 : test_space_size = _secondary_var.dofIndices().size(); 79 32212 : jacobian_types = { 80 : {JType::SecondarySecondary, JType::SecondaryPrimary, JType::SecondaryLower}}; 81 32212 : break; 82 : 83 32212 : case MType::Primary: 84 32212 : test_space_size = _primary_var.dofIndicesNeighbor().size(); 85 32212 : jacobian_types = {{JType::PrimarySecondary, JType::PrimaryPrimary, JType::PrimaryLower}}; 86 32212 : break; 87 : 88 2028 : case MType::Lower: 89 2028 : test_space_size = _var ? _var->dofIndicesLower().size() : 0; 90 2028 : jacobian_types = {{JType::LowerSecondary, JType::LowerPrimary, JType::LowerLower}}; 91 2028 : break; 92 : } 93 : 94 66452 : auto & ce = _assembly.couplingEntries(); 95 217864 : for (const auto & it : ce) 96 : { 97 151412 : MooseVariableFEBase & ivariable = *(it.first); 98 151412 : MooseVariableFEBase & jvariable = *(it.second); 99 : 100 151412 : unsigned int ivar = ivariable.number(); 101 151412 : unsigned int jvar = jvariable.number(); 102 : 103 151412 : switch (mortar_type) 104 : { 105 60532 : case MType::Secondary: 106 60532 : if (ivar != _secondary_var.number()) 107 67608 : continue; 108 37996 : break; 109 : 110 60532 : case MType::Primary: 111 60532 : if (ivar != _primary_var.number()) 112 22536 : continue; 113 37996 : break; 114 : 115 30348 : case MType::Lower: 116 30348 : if (!_var || _var->number() != ivar) 117 22536 : continue; 118 7812 : break; 119 : } 120 : 121 83804 : std::array<size_t, 3> shape_space_sizes{{jvariable.dofIndices().size(), 122 83804 : jvariable.dofIndicesNeighbor().size(), 123 167608 : jvariable.dofIndicesLower().size()}}; 124 : std::array<const VariablePhiValue *, 3> phis; 125 : std::array<const VariablePhiGradient *, 3> grad_phis; 126 : std::array<const VectorVariablePhiValue *, 3> vector_phis; 127 : std::array<const VectorVariablePhiGradient *, 3> vector_grad_phis; 128 83804 : if (jvariable.isVector()) 129 : { 130 0 : const auto & temp_var = static_cast<MooseVariableFE<RealVectorValue> &>(jvariable); 131 0 : vector_phis = {{&temp_var.phiFace(), &temp_var.phiFaceNeighbor(), &temp_var.phiLower()}}; 132 0 : vector_grad_phis = { 133 0 : {&temp_var.gradPhiFace(), &temp_var.gradPhiFaceNeighbor(), &temp_var.gradPhiLower()}}; 134 : } 135 : else 136 : { 137 83804 : const auto & temp_var = static_cast<MooseVariableFE<Real> &>(jvariable); 138 83804 : phis = {{&temp_var.phiFace(), &temp_var.phiFaceNeighbor(), &temp_var.phiLower()}}; 139 83804 : grad_phis = { 140 83804 : {&temp_var.gradPhiFace(), &temp_var.gradPhiFaceNeighbor(), &temp_var.gradPhiLower()}}; 141 : } 142 : 143 335216 : for (MooseIndex(3) type_index = 0; type_index < 3; ++type_index) 144 : { 145 251412 : const auto jacobian_type = jacobian_types[type_index]; 146 : 147 251412 : prepareMatrixTagLower(_assembly, ivar, jvar, jacobian_type); 148 : 149 : /// Set the proper phis 150 251412 : if (jvariable.isVector()) 151 : { 152 0 : _vector_phi = vector_phis[type_index]; 153 0 : _vector_grad_phi = vector_grad_phis[type_index]; 154 : } 155 : else 156 : { 157 251412 : _phi = phis[type_index]; 158 251412 : _grad_phi = grad_phis[type_index]; 159 : } 160 : 161 3614532 : for (_i = 0; _i < test_space_size; _i++) 162 30372544 : for (_j = 0; _j < shape_space_sizes[type_index]; _j++) 163 204922432 : for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++) 164 177913008 : _local_ke(_i, _j) += 165 177913008 : _JxW_msm[_qp] * _coord[_qp] * computeQpJacobian(jacobian_type, jvar); 166 251412 : accumulateTaggedLocalMatrix(); 167 : } 168 : } 169 66452 : }