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 "HFEMTestJump.h" 11 : 12 : registerMooseObject("MooseApp", HFEMTestJump); 13 : 14 : InputParameters 15 14336 : HFEMTestJump::validParams() 16 : { 17 14336 : InputParameters params = DGKernel::validParams(); 18 14336 : params.addRequiredCoupledVar("side_variable", "side variable to use as Lagrange multiplier"); 19 14336 : params.addClassDescription("Imposes constraints for HFEM with side-discontinuous variables."); 20 14336 : return params; 21 0 : } 22 : 23 37 : HFEMTestJump::HFEMTestJump(const InputParameters & parameters) 24 : : DGKernel(parameters), 25 37 : _lambda_var(*getVar("side_variable", 0)), 26 37 : _lambda(_is_implicit ? _lambda_var.sln() : _lambda_var.slnOld()), 27 37 : _phi_lambda(_lambda_var.phiFace()), 28 37 : _test_lambda(_lambda_var.phiFace()), 29 74 : _lambda_id(coupled("side_variable")) 30 : { 31 37 : } 32 : 33 : Real 34 2373480 : HFEMTestJump::computeQpResidual(Moose::DGResidualType type) 35 : { 36 : // Use normal vector at qp 0 to make solution depend on geometry 37 : // (well, geometry and quadrature rule, with curved boundaries...) 38 : // rather than element numbering 39 2373480 : const Real sign = (_normals[0] > Point()) ? 1 : -1; 40 : 41 2373480 : switch (type) 42 : { 43 1186740 : case Moose::Element: 44 1186740 : return -sign * _lambda[_qp] * _test[_i][_qp]; 45 : break; 46 : 47 1186740 : case Moose::Neighbor: 48 1186740 : return sign * _lambda[_qp] * _test_neighbor[_i][_qp]; 49 : break; 50 : } 51 0 : return 0; 52 : } 53 : 54 : Real 55 266880 : HFEMTestJump::computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar) 56 : { 57 266880 : if (jvar != _lambda_id) 58 0 : return 0; 59 : 60 266880 : const Real sign = (_normals[0] > Point()) ? 1 : -1; 61 : 62 266880 : switch (type) 63 : { 64 66720 : case Moose::ElementElement: 65 66720 : return -sign * _phi_lambda[_j][_qp] * _test[_i][_qp]; 66 : 67 66720 : case Moose::NeighborElement: 68 66720 : return sign * _phi_lambda[_j][_qp] * _test_neighbor[_i][_qp]; 69 : 70 66720 : case Moose::ElementNeighbor: 71 66720 : return 0; 72 : 73 66720 : case Moose::NeighborNeighbor: 74 66720 : return 0; 75 : 76 0 : default: 77 0 : break; 78 : } 79 : 80 0 : return 0; 81 : }