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 "LinearFVAdvectionDiffusionFunctorDirichletBC.h" 11 : 12 : registerMooseObject("MooseApp", LinearFVAdvectionDiffusionFunctorDirichletBC); 13 : 14 : InputParameters 15 16513 : LinearFVAdvectionDiffusionFunctorDirichletBC::validParams() 16 : { 17 16513 : InputParameters params = LinearFVAdvectionDiffusionBC::validParams(); 18 16513 : params.addClassDescription( 19 : "Adds a dirichlet BC which can be used for the assembly of linear " 20 : "finite volume system and whose face values are determined using a functor. This kernel is " 21 : "only designed to work with advection-diffusion problems."); 22 16513 : params.addRequiredParam<MooseFunctorName>("functor", "The functor for this boundary condition."); 23 16513 : return params; 24 0 : } 25 : 26 1124 : LinearFVAdvectionDiffusionFunctorDirichletBC::LinearFVAdvectionDiffusionFunctorDirichletBC( 27 1124 : const InputParameters & parameters) 28 1124 : : LinearFVAdvectionDiffusionBC(parameters), _functor(getFunctor<Real>("functor")) 29 : { 30 1124 : } 31 : 32 : Real 33 229006 : LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryValue() const 34 : { 35 229006 : return _functor(singleSidedFaceArg(_current_face_info), determineState()); 36 : } 37 : 38 : Real 39 0 : LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryNormalGradient() const 40 : { 41 0 : const auto elem_arg = makeElemArg(_current_face_type == FaceInfo::VarFaceNeighbors::ELEM 42 0 : ? _current_face_info->elemPtr() 43 0 : : _current_face_info->neighborPtr()); 44 0 : const Real distance = computeCellToFaceDistance(); 45 0 : return (_functor(singleSidedFaceArg(_current_face_info), determineState()) - 46 0 : raw_value(_var(elem_arg, determineState()))) / 47 0 : distance; 48 : } 49 : 50 : Real 51 93660 : LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryValueMatrixContribution() const 52 : { 53 : // Ths will not contribute to the matrix from the value considering that 54 : // the value is independent of the solution. 55 93660 : return 0.0; 56 : } 57 : 58 : Real 59 93660 : LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryValueRHSContribution() const 60 : { 61 : // Fetch the boundary value from the provided functor. 62 93660 : return _functor(singleSidedFaceArg(_current_face_info), determineState()); 63 : } 64 : 65 : Real 66 63444 : LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryGradientMatrixContribution() const 67 : { 68 : // The implicit term from the central difference approximation of the normal 69 : // gradient. 70 63444 : return 1.0 / computeCellToFaceDistance(); 71 : } 72 : 73 : Real 74 63444 : LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryGradientRHSContribution() const 75 : { 76 : // The boundary term from the central difference approximation of the 77 : // normal gradient. 78 63444 : return _functor(singleSidedFaceArg(_current_face_info), determineState()) / 79 63444 : computeCellToFaceDistance(); 80 : }