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 "CoupledTiedValueConstraint.h" 11 : 12 : // MOOSE includes 13 : #include "MooseVariableFE.h" 14 : #include "SystemBase.h" 15 : 16 : #include "libmesh/sparse_matrix.h" 17 : 18 : registerMooseObject("MooseApp", CoupledTiedValueConstraint); 19 : 20 : InputParameters 21 14293 : CoupledTiedValueConstraint::validParams() 22 : { 23 14293 : InputParameters params = NodeFaceConstraint::validParams(); 24 14293 : params.addClassDescription( 25 : "Requires the value of two variables to be the consistent on both sides of an interface."); 26 14293 : params.addParam<Real>("scaling", 1, "scaling factor to be applied to constraint equations"); 27 14293 : params.set<bool>("use_displaced_mesh") = true; 28 14293 : return params; 29 0 : } 30 : 31 12 : CoupledTiedValueConstraint::CoupledTiedValueConstraint(const InputParameters & parameters) 32 : : NodeFaceConstraint(parameters), 33 12 : _scaling(getParam<Real>("scaling")), 34 24 : _residual_copy(_sys.residualGhosted()) 35 : { 36 12 : } 37 : 38 : Real 39 32 : CoupledTiedValueConstraint::computeQpSecondaryValue() 40 : { 41 32 : return _u_primary[_qp]; 42 : } 43 : 44 : Real 45 480 : CoupledTiedValueConstraint::computeQpResidual(Moose::ConstraintType type) 46 : { 47 480 : Real scaling_factor = _var.scalingFactor(); 48 480 : Real secondary_resid = 0; 49 480 : Real retVal = 0; 50 480 : switch (type) 51 : { 52 96 : case Moose::Secondary: 53 96 : retVal = (_u_secondary[_qp] - _u_primary[_qp]) * _test_secondary[_i][_qp] * _scaling; 54 96 : break; 55 384 : case Moose::Primary: 56 384 : secondary_resid = 57 384 : _residual_copy(_current_node->dof_number(0, _var.number(), 0)) / scaling_factor; 58 384 : retVal = secondary_resid * _test_primary[_i][_qp]; 59 384 : break; 60 0 : default: 61 0 : break; 62 : } 63 480 : return retVal; 64 : } 65 : 66 : Real 67 1600 : CoupledTiedValueConstraint::computeQpJacobian(Moose::ConstraintJacobianType type) 68 : { 69 1600 : Real scaling_factor = _var.scalingFactor(); 70 1600 : Real secondary_jac = 0; 71 1600 : Real retVal = 0; 72 1600 : switch (type) 73 : { 74 320 : case Moose::SecondarySecondary: 75 320 : retVal = _phi_secondary[_j][_qp] * _test_secondary[_i][_qp] * _scaling; 76 320 : break; 77 0 : case Moose::SecondaryPrimary: 78 0 : retVal = 0; 79 0 : break; 80 1280 : case Moose::PrimarySecondary: 81 : secondary_jac = 82 1280 : (*_jacobian)(_current_node->dof_number(0, _var.number(), 0), _connected_dof_indices[_j]); 83 1280 : retVal = secondary_jac * _test_primary[_i][_qp] / scaling_factor; 84 1280 : break; 85 0 : case Moose::PrimaryPrimary: 86 0 : retVal = 0; 87 0 : break; 88 0 : default: 89 0 : mooseError("Unsupported type"); 90 : break; 91 : } 92 1600 : return retVal; 93 : } 94 : 95 : Real 96 1280 : CoupledTiedValueConstraint::computeQpOffDiagJacobian(Moose::ConstraintJacobianType type, 97 : unsigned int jvar) 98 : { 99 1280 : Real retVal = 0; 100 : 101 1280 : if (jvar == _primary_var_num) 102 : { 103 1280 : switch (type) 104 : { 105 0 : case Moose::SecondarySecondary: 106 0 : retVal = 0; 107 0 : break; 108 256 : case Moose::SecondaryPrimary: 109 256 : retVal = -_phi_primary[_j][_qp] * _test_secondary[_i][_qp] * _scaling; 110 256 : break; 111 0 : case Moose::PrimarySecondary: 112 0 : retVal = 0; 113 0 : break; 114 1024 : case Moose::PrimaryPrimary: 115 1024 : retVal = 0; 116 1024 : break; 117 0 : default: 118 0 : mooseError("Unsupported type"); 119 : break; 120 : } 121 : } 122 : 123 1280 : return retVal; 124 : }