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