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 "ScalarLagrangeMultiplier.h" 11 : 12 : // MOOSE includes 13 : #include "Assembly.h" 14 : #include "MooseVariableFE.h" 15 : #include "MooseVariableScalar.h" 16 : #include "SystemBase.h" 17 : 18 : #include "libmesh/quadrature.h" 19 : 20 : registerMooseObject("MooseApp", ScalarLagrangeMultiplier); 21 : 22 : InputParameters 23 14668 : ScalarLagrangeMultiplier::validParams() 24 : { 25 14668 : InputParameters params = Kernel::validParams(); 26 14668 : params.addClassDescription("This class is used to enforce integral of phi = V_0 with a " 27 : "Lagrange multiplier approach."); 28 14668 : params.addRequiredCoupledVar("lambda", "Lagrange multiplier variable"); 29 : 30 14668 : return params; 31 0 : } 32 : 33 214 : ScalarLagrangeMultiplier::ScalarLagrangeMultiplier(const InputParameters & parameters) 34 214 : : Kernel(parameters), _lambda_var(coupledScalar("lambda")), _lambda(coupledScalarValue("lambda")) 35 : { 36 206 : } 37 : 38 412 : ScalarLagrangeMultiplier::~ScalarLagrangeMultiplier() {} 39 : 40 : Real 41 513552 : ScalarLagrangeMultiplier::computeQpResidual() 42 : { 43 513552 : return _lambda[0] * _test[_i][_qp]; 44 : } 45 : 46 : void 47 89994 : ScalarLagrangeMultiplier::computeOffDiagJacobianScalar(unsigned int jvar) 48 : { 49 : // Note: Here we are assembling the contributions for _both_ Eq. (9) and (10) 50 : // in the detailed writeup on this problem [0]. This is important to highlight 51 : // because it is slightly different from the way things usually work in MOOSE 52 : // because we are computing Jacobian contributions _for a different equation_ 53 : // than the equation which this Kernel is assigned in the input file. We do 54 : // this for both simplicity (results in slightly less code) and efficiency: 55 : // the contribution is symmetric, so it can be computed once and used twice. 56 : // 57 : // [0]: https://github.com/idaholab/large_media/blob/master/framework/scalar_constraint_kernel.pdf 58 : 59 89994 : MooseVariableScalar & jv = _sys.getScalarVariable(_tid, jvar); 60 : 61 89994 : prepareMatrixTag(_assembly, _var.number(), jvar); 62 180798 : for (_i = 0; _i < _test.size(); _i++) 63 181608 : for (_j = 0; _j < jv.order(); _j++) 64 554940 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 65 464136 : _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobianScalar(jvar); 66 89994 : accumulateTaggedLocalMatrix(); 67 : 68 89994 : const auto ke_copy = _local_ke; 69 89994 : prepareMatrixTag(_assembly, jvar, _var.number()); 70 89994 : ke_copy.get_transpose(_local_ke); 71 89994 : accumulateTaggedLocalMatrix(); 72 89994 : } 73 : 74 : Real 75 464136 : ScalarLagrangeMultiplier::computeQpOffDiagJacobianScalar(unsigned int jvar) 76 : { 77 464136 : if (jvar == _lambda_var) 78 464136 : return _test[_i][_qp]; 79 : else 80 0 : return 0.; 81 : }