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 "DGFunctionDiffusionDirichletBC.h" 11 : 12 : // MOOSE includes 13 : #include "Function.h" 14 : #include "MooseVariableFE.h" 15 : 16 : #include "libmesh/numeric_vector.h" 17 : #include "libmesh/utility.h" 18 : 19 : // C++ includes 20 : #include <cmath> 21 : 22 : registerMooseObject("MooseApp", DGFunctionDiffusionDirichletBC); 23 : 24 : InputParameters 25 15321 : DGFunctionDiffusionDirichletBC::validParams() 26 : { 27 15321 : InputParameters params = IntegratedBC::validParams(); 28 15321 : params.addClassDescription( 29 : "Diffusion Dirichlet boundary condition for discontinuous Galerkin method."); 30 15321 : params.addParam<Real>("value", 0.0, "The value the variable should have on the boundary"); 31 15321 : params.addRequiredParam<FunctionName>("function", "The forcing function."); 32 15321 : params.addRequiredParam<Real>("epsilon", "Epsilon"); 33 15321 : params.addRequiredParam<Real>("sigma", "Sigma"); 34 45963 : params.addParam<MaterialPropertyName>( 35 30642 : "diff", 1, "The diffusion (or thermal conductivity or viscosity) coefficient."); 36 : 37 15321 : return params; 38 0 : } 39 : 40 551 : DGFunctionDiffusionDirichletBC::DGFunctionDiffusionDirichletBC(const InputParameters & parameters) 41 : : IntegratedBC(parameters), 42 551 : _func(getFunction("function")), 43 551 : _epsilon(getParam<Real>("epsilon")), 44 551 : _sigma(getParam<Real>("sigma")), 45 1102 : _diff(getMaterialProperty<Real>("diff")) 46 : { 47 551 : } 48 : 49 : Real 50 1659858 : DGFunctionDiffusionDirichletBC::computeQpResidual() 51 : { 52 1659858 : const int elem_b_order = std::max(libMesh::Order(1), _var.order()); 53 : const Real h_elem = 54 1659858 : _current_elem_volume / _current_side_volume * 1.0 / Utility::pow<2>(elem_b_order); 55 : 56 1659858 : const Real fn = _func.value(_t, _q_point[_qp]); 57 1659858 : Real r = 0.0; 58 1659858 : r -= (_diff[_qp] * _grad_u[_qp] * _normals[_qp] * _test[_i][_qp]); 59 1659858 : r += _epsilon * (_u[_qp] - fn) * _diff[_qp] * _grad_test[_i][_qp] * _normals[_qp]; 60 1659858 : r += _sigma / h_elem * (_u[_qp] - fn) * _test[_i][_qp]; 61 : 62 1659858 : return r; 63 : } 64 : 65 : Real 66 864742 : DGFunctionDiffusionDirichletBC::computeQpJacobian() 67 : { 68 864742 : const int elem_b_order = std::max(libMesh::Order(1), _var.order()); 69 : const Real h_elem = 70 864742 : _current_elem_volume / _current_side_volume * 1.0 / Utility::pow<2>(elem_b_order); 71 : 72 864742 : Real r = 0.0; 73 864742 : r -= (_diff[_qp] * _grad_phi[_j][_qp] * _normals[_qp] * _test[_i][_qp]); 74 864742 : r += _epsilon * _phi[_j][_qp] * _diff[_qp] * _grad_test[_i][_qp] * _normals[_qp]; 75 864742 : r += _sigma / h_elem * _phi[_j][_qp] * _test[_i][_qp]; 76 : 77 864742 : return r; 78 : }