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 "ConservativeAdvectionBC.h" 11 : #include "Function.h" 12 : 13 : registerMooseObject("MooseApp", ConservativeAdvectionBC); 14 : registerMooseObject("MooseApp", ADConservativeAdvectionBC); 15 : 16 : template <bool is_ad> 17 : InputParameters 18 29640 : ConservativeAdvectionBCTempl<is_ad>::validParams() 19 : { 20 29640 : InputParameters params = GenericIntegratedBC<is_ad>::validParams(); 21 118560 : params.addParam<MaterialPropertyName>( 22 : "velocity_mat_prop", 23 : "Velocity vector as a material property. Should be provided when we want the velocity value " 24 : "to be determined implicitly (e.g. we don't have a Dirichlet condition)"); 25 118560 : params.addParam<FunctionName>("velocity_function", 26 : "Function describing the values of velocity on the boundary."); 27 59280 : params.addClassDescription( 28 : "Boundary condition for advection when it is integrated by parts. Supports Dirichlet " 29 : "(inlet-like) and implicit (outlet-like) conditions."); 30 118560 : params.addParam<MaterialPropertyName>("advected_quantity", 31 : "An optional material property to be advected. If not " 32 : "supplied, then the variable will be used."); 33 118560 : params.addParam<FunctionName>("primal_dirichlet_value", 34 : "The value of the primal variable on the boundary."); 35 59280 : params.addParam<MaterialPropertyName>( 36 : "primal_coefficient", 37 59280 : 1.0, 38 : "If a primal Dirichlet value is supplied, then a coefficient may be optionally multiplied " 39 : "that multiples the Dirichlet value"); 40 29640 : return params; 41 0 : } 42 : 43 : InputParameters 44 14841 : ConservativeAdvectionBC::validParams() 45 : { 46 14841 : return ConservativeAdvectionBCTempl<false>::validParams(); 47 : } 48 : 49 : template <bool is_ad> 50 78 : ConservativeAdvectionBCTempl<is_ad>::ConservativeAdvectionBCTempl( 51 : const InputParameters & parameters) 52 : : GenericIntegratedBC<is_ad>(parameters), 53 78 : _velocity_mat_prop(isParamValid("velocity_mat_prop") 54 156 : ? &this->template getGenericMaterialProperty<RealVectorValue, is_ad>( 55 : "velocity_mat_prop") 56 : : nullptr), 57 234 : _velocity_function(isParamValid("velocity_function") ? &getFunction("velocity_function") 58 : : nullptr), 59 156 : _user_supplied_adv_quant(isParamValid("advected_quantity")), 60 78 : _adv_quant( 61 78 : _user_supplied_adv_quant 62 78 : ? this->template getGenericMaterialProperty<Real, is_ad>("advected_quantity").get() 63 : : _u), 64 78 : _primal_dirichlet( 65 234 : isParamValid("primal_dirichlet_value") ? &getFunction("primal_dirichlet_value") : nullptr), 66 234 : _primal_coeff(this->template getGenericMaterialProperty<Real, is_ad>("primal_coefficient")) 67 : { 68 234 : if (isParamSetByUser("primal_coefficient") && !_primal_dirichlet) 69 0 : paramError("primal_coefficient", 70 : "This parameter should only be provided when 'primal_dirichlet_value' is provided"); 71 234 : if (static_cast<bool>(_primal_dirichlet) + isParamValid("advected_quantity") > 1) 72 0 : mooseError("Only one of 'primal_dirichlet_value' or 'advected_quantity' should be provided"); 73 78 : if (static_cast<bool>(_velocity_mat_prop) + static_cast<bool>(_velocity_function) != 1) 74 0 : mooseError("Exactly one of 'velocity_mat_prop' or 'velocity_function' should be provided"); 75 78 : } 76 : 77 50 : ConservativeAdvectionBC::ConservativeAdvectionBC(const InputParameters & parameters) 78 50 : : ConservativeAdvectionBCTempl<false>(parameters) 79 : { 80 50 : } 81 : 82 : template <bool is_ad> 83 : GenericReal<is_ad> 84 86736 : ConservativeAdvectionBCTempl<is_ad>::computeQpResidual() 85 : { 86 3024 : const auto vdotn = 87 86736 : (_velocity_mat_prop 88 49416 : ? (*_velocity_mat_prop)[_qp] 89 127080 : : GenericRealVectorValue<is_ad>(_velocity_function->vectorValue(_t, _q_point[_qp]))) * 90 86736 : this->_normals[_qp]; 91 : 92 86736 : if (_primal_dirichlet) 93 43368 : return _test[_i][_qp] * vdotn * _primal_dirichlet->value(_t, _q_point[_qp]) * 94 44880 : _primal_coeff[_qp]; 95 : else 96 43368 : return _test[_i][_qp] * vdotn * _adv_quant[_qp]; 97 3024 : } 98 : 99 : Real 100 4896 : ConservativeAdvectionBC::computeQpJacobian() 101 : { 102 4896 : if (!_user_supplied_adv_quant && !_primal_dirichlet) 103 2448 : return _test[_i][_qp] * 104 2448 : (_velocity_mat_prop 105 4896 : ? (*_velocity_mat_prop)[_qp] 106 4896 : : RealVectorValue(_velocity_function->vectorValue(_t, _q_point[_qp]))) * 107 4896 : this->_normals[_qp] * _phi[_j][_qp]; 108 2448 : return 0.0; 109 : } 110 : 111 : template class ConservativeAdvectionBCTempl<false>; 112 : template class ConservativeAdvectionBCTempl<true>;