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 "INSFVTurbulentViscosityWallFunction.h" 11 : #include "Function.h" 12 : #include "NavierStokesMethods.h" 13 : 14 : registerMooseObject("NavierStokesApp", INSFVTurbulentViscosityWallFunction); 15 : 16 : InputParameters 17 1916 : INSFVTurbulentViscosityWallFunction::validParams() 18 : { 19 1916 : InputParameters params = FVDirichletBCBase::validParams(); 20 1916 : params.addClassDescription("Adds Dirichlet BC for wall values of the turbulent viscosity."); 21 3832 : params.addRequiredParam<MooseFunctorName>("u", "The velocity in the x direction."); 22 3832 : params.addParam<MooseFunctorName>("v", "The velocity in the y direction."); 23 3832 : params.addParam<MooseFunctorName>("w", "The velocity in the z direction."); 24 1916 : params.addRequiredParam<MooseFunctorName>(NS::density, "Density"); 25 1916 : params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity."); 26 1916 : params.addRequiredParam<MooseFunctorName>(NS::mu_t, "The turbulent viscosity."); 27 1916 : params.addParam<MooseFunctorName>(NS::TKE, "The turbulent kinetic energy."); 28 3832 : params.addParam<Real>("C_mu", 0.09, "Coupled turbulent kinetic energy closure."); 29 : 30 3832 : MooseEnum wall_treatment("eq_newton eq_incremental eq_linearized neq", "neq"); 31 3832 : params.addParam<MooseEnum>( 32 : "wall_treatment", wall_treatment, "The method used for computing the wall functions"); 33 1916 : return params; 34 1916 : } 35 : 36 441 : INSFVTurbulentViscosityWallFunction::INSFVTurbulentViscosityWallFunction( 37 441 : const InputParameters & params) 38 : : FVDirichletBCBase(params), 39 441 : _dim(_subproblem.mesh().dimension()), 40 882 : _u_var(getFunctor<ADReal>("u")), 41 1323 : _v_var(params.isParamValid("v") ? &(getFunctor<ADReal>("v")) : nullptr), 42 441 : _w_var(params.isParamValid("w") ? &(getFunctor<ADReal>("w")) : nullptr), 43 441 : _rho(getFunctor<ADReal>(NS::density)), 44 441 : _mu(getFunctor<ADReal>(NS::mu)), 45 441 : _mu_t(getFunctor<ADReal>(NS::mu_t)), 46 441 : _k(getFunctor<ADReal>(NS::TKE)), 47 882 : _C_mu(getParam<Real>("C_mu")), 48 882 : _wall_treatment(getParam<MooseEnum>("wall_treatment").getEnum<NS::WallTreatmentEnum>()), 49 882 : _preserve_sparsity_pattern(_fv_problem.preserveMatrixSparsityPattern()) 50 : { 51 441 : } 52 : 53 : ADReal 54 389040 : INSFVTurbulentViscosityWallFunction::boundaryValue(const FaceInfo & fi, 55 : const Moose::StateArg & /* state */) const 56 : { 57 : // if on an internal face (internal to the mesh, but an external boundary of the flow area), 58 : // we have to select the element on which the flow variables / viscosity is defined 59 : const bool use_elem = (fi.faceType(_var_sys_numbers_pair) == FaceInfo::VarFaceNeighbors::ELEM); 60 : const Real wall_dist = std::abs( 61 389040 : ((use_elem ? fi.elemCentroid() : fi.neighborCentroid()) - fi.faceCentroid()) * fi.normal()); 62 389040 : const Elem * const elem_ptr = use_elem ? fi.elemPtr() : fi.neighborPtr(); 63 389040 : const auto elem_arg = makeElemArg(elem_ptr); 64 : 65 : // Get the previous non linear iteration values 66 : const auto old_state = Moose::StateArg(1, Moose::SolutionIterationType::Nonlinear); 67 389040 : const auto mu = _mu(elem_arg, old_state); 68 389040 : const auto rho = _rho(elem_arg, old_state); 69 : 70 : // Get the velocity vector 71 778080 : ADRealVectorValue velocity(_u_var(elem_arg, old_state)); 72 389040 : if (_v_var) 73 389040 : velocity(1) = (*_v_var)(elem_arg, old_state); 74 389040 : if (_w_var) 75 0 : velocity(2) = (*_w_var)(elem_arg, old_state); 76 : 77 : // Compute the velocity and direction of the velocity component that is parallel to the wall 78 : const auto parallel_speed = 79 1167120 : NS::computeSpeed<ADReal>(velocity - velocity * (fi.normal()) * (fi.normal())); 80 : 81 : // Switch for determining the near wall quantities 82 : // wall_treatment can be: "eq_newton eq_incremental eq_linearized neq" 83 : ADReal y_plus; 84 : ADReal mut_log; // turbulent log-layer viscosity 85 : ADReal mu_wall; // total wall viscosity to obtain the shear stress at the wall 86 : 87 389040 : if (_wall_treatment == NS::WallTreatmentEnum::EQ_NEWTON) 88 : { 89 : // Full Newton-Raphson solve to find the wall quantities from the law of the wall 90 341160 : const auto u_tau = NS::findUStar<ADReal>(mu, rho, parallel_speed, wall_dist); 91 341160 : y_plus = wall_dist * u_tau * rho / mu; 92 341160 : mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed; 93 341160 : mut_log = mu_wall - mu; 94 : } 95 47880 : else if (_wall_treatment == NS::WallTreatmentEnum::EQ_INCREMENTAL) 96 : { 97 : // Incremental solve on y_plus to get the near-wall quantities 98 21960 : y_plus = NS::findyPlus<ADReal>(mu, rho, std::max(parallel_speed, 1e-10), wall_dist); 99 0 : mu_wall = mu * (NS::von_karman_constant * y_plus / 100 65880 : std::log(std::max(NS::E_turb_constant * y_plus, 1 + 1e-4))); 101 21960 : mut_log = mu_wall - mu; 102 : } 103 25920 : else if (_wall_treatment == NS::WallTreatmentEnum::EQ_LINEARIZED) 104 : { 105 : // Linearized approximation to the wall function to find the near-wall quantities faster 106 12960 : const ADReal a_c = 1 / NS::von_karman_constant; 107 25920 : const ADReal b_c = 1 / NS::von_karman_constant * 108 38880 : (std::log(NS::E_turb_constant * std::max(wall_dist, 1.0) / mu) + 1.0); 109 12960 : const ADReal c_c = parallel_speed; 110 : 111 64800 : const auto u_tau = (-b_c + std::sqrt(std::pow(b_c, 2) + 4.0 * a_c * c_c)) / (2.0 * a_c); 112 12960 : y_plus = wall_dist * u_tau * rho / mu; 113 12960 : mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed; 114 12960 : mut_log = mu_wall - mu; 115 : } 116 12960 : else if (_wall_treatment == NS::WallTreatmentEnum::NEQ) 117 : { 118 : // Assign non-equilibrium wall function value 119 25920 : y_plus = std::pow(_C_mu, 0.25) * wall_dist * std::sqrt(_k(elem_arg, old_state)) * rho / mu; 120 0 : mu_wall = mu * (NS::von_karman_constant * y_plus / 121 38880 : std::log(std::max(NS::E_turb_constant * y_plus, 1.0 + 1e-4))); 122 12960 : mut_log = mu_wall - mu; 123 : } 124 : else 125 : mooseAssert(false, 126 : "For `INSFVTurbulentViscosityWallFunction` , wall treatment should not reach here"); 127 : 128 389040 : ADReal residual = 0; 129 : // To keep the same sparsity pattern for all y_plus 130 389040 : if (_preserve_sparsity_pattern) 131 778080 : residual = 0 * mut_log * y_plus; 132 : 133 389040 : if (y_plus <= 5.0) 134 : // sub-laminar layer 135 : residual += 0.0; 136 135408 : else if (y_plus >= 30.0) 137 : // log-layer 138 119466 : residual += std::max(mut_log, NS::mu_t_low_limit); 139 : else 140 : { 141 : // buffer layer 142 15942 : const auto blending_function = (y_plus - 5.0) / 25.0; 143 : // the blending depends on the mut_log at y+=30 144 15942 : const auto mut_log = mu * _mut_30; 145 15942 : residual += std::max(blending_function * mut_log, NS::mu_t_low_limit); 146 : } 147 389040 : return residual; 148 : }