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 "ADPenaltyPeriodicSegmentalConstraint.h" 11 : 12 : registerMooseObject("MooseApp", ADPenaltyPeriodicSegmentalConstraint); 13 : 14 : InputParameters 15 14387 : ADPenaltyPeriodicSegmentalConstraint::validParams() 16 : { 17 14387 : InputParameters params = ADMortarScalarBase::validParams(); 18 14387 : params.addClassDescription( 19 : "ADPenaltyPeriodicSegmentalConstraint enforces macro-micro periodic conditions between " 20 : "secondary and primary sides of a mortar interface using a penalty approach " 21 : "(no Lagrange multipliers needed). Must be used alongside PenaltyEqualValueConstraint."); 22 14387 : params.renameCoupledVar("scalar_variable", "epsilon", "Primary coupled scalar variable"); 23 14387 : params.addRequiredCoupledVar("sigma", "Controlled scalar averaging variable"); 24 43161 : params.addParam<Real>( 25 : "penalty_value", 26 28774 : 1.0, 27 : "Penalty value used to impose a generalized force capturing the mortar constraint equation"); 28 : 29 14387 : return params; 30 0 : } 31 : 32 61 : ADPenaltyPeriodicSegmentalConstraint::ADPenaltyPeriodicSegmentalConstraint( 33 61 : const InputParameters & parameters) 34 : : DerivativeMaterialInterface<ADMortarScalarBase>(parameters), 35 61 : _temp_jump_global(), 36 61 : _tau_s(), 37 61 : _kappa_aux_ptr(getScalarVar("sigma", 0)), 38 61 : _ka_order(_kappa_aux_ptr->order()), 39 61 : _kappa_aux(coupledScalarValue("sigma")), 40 122 : _pen_scale(getParam<Real>("penalty_value")) 41 : { 42 61 : if (_kappa_aux_ptr->kind() != Moose::VarKindType::VAR_AUXILIARY) 43 0 : paramError("sigma", 44 : "Must assign auxiliary scalar variable to sigma, rather than nonlinear variable"); 45 61 : } 46 : 47 : // Compute the stability parameters to use for all quadrature points 48 : void 49 2336 : ADPenaltyPeriodicSegmentalConstraint::precalculateResidual() 50 : { 51 2336 : precalculateStability(); 52 2336 : } 53 : 54 : // Compute the temperature jump for current quadrature point 55 : void 56 5968 : ADPenaltyPeriodicSegmentalConstraint::initScalarQpResidual() 57 : { 58 5968 : precalculateMaterial(); 59 5968 : } 60 : 61 : ADReal 62 68480 : ADPenaltyPeriodicSegmentalConstraint::computeQpResidual(const Moose::MortarType mortar_type) 63 : { 64 : // Compute penalty parameter times x-jump times average heat flux 65 68480 : RealVectorValue dx(_phys_points_primary[_qp] - _phys_points_secondary[_qp]); 66 68480 : ADRealVectorValue kappa_vec(_kappa[0], 0, 0); 67 68480 : Moose::derivInsert(kappa_vec(0).derivatives(), _kappa_var_ptr->dofIndices()[0], 1); 68 68480 : if (_k_order == 2) 69 : { 70 27008 : kappa_vec(1) = _kappa[1]; 71 27008 : Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1); 72 : } 73 41472 : else if (_k_order == 3) 74 : { 75 41472 : kappa_vec(1) = _kappa[1]; 76 41472 : kappa_vec(2) = _kappa[2]; 77 41472 : Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1); 78 41472 : Moose::derivInsert(kappa_vec(2).derivatives(), _kappa_var_ptr->dofIndices()[2], 1); 79 : } 80 68480 : ADReal r = _tau_s * (kappa_vec * dx); 81 : 82 68480 : switch (mortar_type) 83 : { 84 34240 : case Moose::MortarType::Secondary: 85 34240 : r *= _test_secondary[_i][_qp]; 86 34240 : break; 87 34240 : case Moose::MortarType::Primary: 88 34240 : r *= -_test_primary[_i][_qp]; 89 34240 : break; 90 0 : case Moose::MortarType::Lower: 91 0 : return 0; 92 0 : default: 93 0 : return 0; 94 : } 95 68480 : return r; 96 68480 : } 97 : 98 : ADReal 99 14528 : ADPenaltyPeriodicSegmentalConstraint::computeScalarQpResidual() 100 : { 101 : // Stability/penalty term for residual of scalar variable 102 14528 : ADReal r = _tau_s * _temp_jump_global; 103 14528 : RealVectorValue dx(_phys_points_primary[_qp] - _phys_points_secondary[_qp]); 104 : 105 14528 : r *= -dx(_h); 106 : 107 14528 : ADRealVectorValue kappa_vec(_kappa[0], 0, 0); 108 14528 : RealVectorValue kappa_aux_vec(_kappa_aux[0], 0, 0); 109 14528 : Moose::derivInsert(kappa_vec(0).derivatives(), _kappa_var_ptr->dofIndices()[0], 1); 110 14528 : if (_k_order == 2) 111 : { 112 6752 : kappa_vec(1) = _kappa[1]; 113 6752 : Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1); 114 6752 : kappa_aux_vec(1) = _kappa_aux[1]; 115 : } 116 7776 : else if (_k_order == 3) 117 : { 118 7776 : kappa_vec(1) = _kappa[1]; 119 7776 : kappa_vec(2) = _kappa[2]; 120 7776 : Moose::derivInsert(kappa_vec(1).derivatives(), _kappa_var_ptr->dofIndices()[1], 1); 121 7776 : Moose::derivInsert(kappa_vec(2).derivatives(), _kappa_var_ptr->dofIndices()[2], 1); 122 7776 : kappa_aux_vec(1) = _kappa_aux[1]; 123 7776 : kappa_aux_vec(2) = _kappa_aux[2]; 124 : } 125 : 126 14528 : r += dx(_h) * _tau_s * (kappa_vec * dx); 127 14528 : r -= dx(_h) * (kappa_aux_vec * _normals[_qp]); 128 : 129 29056 : return r; 130 14528 : } 131 : 132 : void 133 2336 : ADPenaltyPeriodicSegmentalConstraint::precalculateStability() 134 : { 135 : // Example showing how the penalty could be loaded from some function 136 2336 : _tau_s = _pen_scale; 137 2336 : } 138 : 139 : // Compute temperature jump and flux average/jump 140 : void 141 5968 : ADPenaltyPeriodicSegmentalConstraint::precalculateMaterial() 142 : { 143 5968 : _temp_jump_global = (_u_primary[_qp] - _u_secondary[_qp]); 144 5968 : }