mus_turb_wallFunc_module Module

This module contains turbulent wall function type and routines to calculate friction velocity and stream-wise velocity component.

Constant parameters for Implicit equation solver


Uses

  • module~~mus_turb_wallfunc_module~~UsesGraph module~mus_turb_wallfunc_module mus_turb_wallFunc_module module~env_module env_module module~mus_turb_wallfunc_module->module~env_module module~aotus_module aotus_module module~mus_turb_wallfunc_module->module~aotus_module module~tem_aux_module tem_aux_module module~mus_turb_wallfunc_module->module~tem_aux_module module~tem_float_module tem_float_module module~mus_turb_wallfunc_module->module~tem_float_module module~tem_logging_module tem_logging_module module~mus_turb_wallfunc_module->module~tem_logging_module

Used by

  • module~~mus_turb_wallfunc_module~~UsedByGraph module~mus_turb_wallfunc_module mus_turb_wallFunc_module module~mus_bc_fluid_turbulent_module mus_bc_fluid_turbulent_module module~mus_bc_fluid_turbulent_module->module~mus_turb_wallfunc_module module~mus_bc_header_module mus_bc_header_module module~mus_bc_fluid_turbulent_module->module~mus_bc_header_module module~mus_bc_header_module->module~mus_turb_wallfunc_module module~mus_bndforce_module mus_bndForce_module module~mus_bndforce_module->module~mus_bc_header_module module~mus_bc_passivescalar_module mus_bc_passiveScalar_module module~mus_bc_passivescalar_module->module~mus_bc_header_module module~mus_bc_general_module mus_bc_general_module module~mus_bc_general_module->module~mus_bc_fluid_turbulent_module module~mus_bc_general_module->module~mus_bc_header_module module~mus_bc_fluid_noneqexpol_module mus_bc_fluid_nonEqExpol_module module~mus_bc_fluid_noneqexpol_module->module~mus_bc_header_module module~mus_bc_fluid_wall_module mus_bc_fluid_wall_module module~mus_bc_fluid_wall_module->module~mus_bc_header_module module~mus_weights_module mus_weights_module module~mus_weights_module->module~mus_bc_header_module module~mus_scheme_type_module mus_scheme_type_module module~mus_scheme_type_module->module~mus_bc_header_module module~mus_bc_nernstplanck_module mus_bc_nernstPlanck_module module~mus_bc_nernstplanck_module->module~mus_bc_header_module module~mus_construction_module mus_construction_module module~mus_construction_module->module~mus_bc_header_module module~mus_field_module mus_field_module module~mus_field_module->module~mus_bc_header_module module~mus_bc_species_module mus_bc_species_module module~mus_bc_species_module->module~mus_bc_header_module module~mus_bc_var_module mus_bc_var_module module~mus_bc_var_module->module~mus_bc_header_module module~mus_bc_fluid_module mus_bc_fluid_module module~mus_bc_fluid_module->module~mus_bc_header_module module~mus_bc_fluid_experimental_module mus_bc_fluid_experimental_module module~mus_bc_fluid_experimental_module->module~mus_bc_header_module module~mus_dynloadbal_module mus_dynLoadBal_module module~mus_dynloadbal_module->module~mus_bc_header_module module~mus_bc_poisson_module mus_bc_poisson_module module~mus_bc_poisson_module->module~mus_bc_header_module module~mus_connectivity_module mus_connectivity_module module~mus_connectivity_module->module~mus_bc_header_module

Contents


Variables

TypeVisibilityAttributesNameInitial
integer, private, parameter:: imEq_nIter =1000
real(kind=rk), private, parameter:: imEq_tol =1e-10

Constant parameters for Reichardt's law

real(kind=rk), private, parameter:: reivonKA =0.4_rk
real(kind=rk), private, parameter:: oneOvervonKA =1.0_rk/reivonKA
real(kind=rk), private, parameter:: reiC =7.8_rk
real(kind=rk), private, parameter:: reiB1 =11.0_rk
real(kind=rk), private, parameter:: reiB2 =3.0_rk

constants for Werner & Wengle model http://trace-portal.de/userguide/trace/page_wallFunctionWernerWengle.html

real(kind=rk), private, parameter:: ww_C_m =8.3_rk
real(kind=rk), private, parameter:: ww_m =1.0_rk/7.0_rk
real(kind=rk), private, parameter:: ww_uLmt_yPlus =ww_C_m**(1.0_rk/(1.0_rk-ww_m))
real(kind=rk), private, parameter:: ww_uLmt_vel =WW_uLmt_yPlus**2.0_rk
real(kind=rk), private, parameter:: ww_onePlusM =1.0_rk+ww_m
real(kind=rk), private, parameter:: ww_onePlusMInv =1.0_rk/ww_onePlusM
real(kind=rk), private, parameter:: ww_mDivOnePlusM =ww_m/ww_onePlusM

Werner & Wengle two layer model does not cover buffer layer so three different layer model from Schmidt is used to compute stream-wise velocity at first element. Viscous sublayer < 5 buffer layer > 5 and < 30 log layer > 30

real(kind=rk), private, parameter:: sc_lLmt =5.0_rk
real(kind=rk), private, parameter:: sc_uLmt =30.0_rk
real(kind=rk), private, parameter:: sc_vonKarman =0.4_rk
real(kind=rk), private, parameter:: sc_a =((1.0_rk/sc_vonKarman)*log(30.0_rk)+0.2_rk)/log(6.0_rk)
real(kind=rk), private, parameter:: sc_b =5.0_rk-sc_a*log(5.0_rk)

Abstract Interfaces

abstract interface

Interface definition for the turbulent wall bc routines

  • private subroutine mus_proc_calcFricVel(velTau, velSW, distToBnd, viscKine, nElems)

    This abstract interface defines the interface to calculate turbulent wall friction velocity from given velocity and distance to boundary. All inputs and output are in lattice units.

    Arguments

    TypeIntentOptionalAttributesName
    real(kind=rk), intent(inout) :: velTau(:)

    Friction velocity computed from wall model. it is inout to provide velTau from previous timestep as initial velTau for fixed-point or Newton iteration solver

    real(kind=rk), intent(in) :: velSW(:)

    Stream-wise velocity component from which friction velocity is computed

    real(kind=rk), intent(in) :: distToBnd(:)

    Distance to the boundary in the discrete normal direction

    real(kind=rk), intent(in) :: viscKine(:)

    Kinematic viscosity

    integer, intent(in) :: nElems

    Number of elements in input and output arrays

abstract interface

Interface definition for the turbulent wall bc routines

  • private subroutine mus_proc_calcStreamWiseVel(velSW, velTau, distToBnd, viscKine, nElems)

    This abstract interface defines the interface to calculate stream-wise velocity component from friction velocity and distance to boundary. All inputs and output are in lattice units.

    Arguments

    TypeIntentOptionalAttributesName
    real(kind=rk), intent(out) :: velSW(:)

    Stream-wise velocity component from wall model

    real(kind=rk), intent(in) :: velTau(:)

    Friction velocity computd from wall model

    real(kind=rk), intent(in) :: distToBnd(:)

    Distance to the boundary in the discrete normai direction

    real(kind=rk), intent(in) :: viscKine(:)

    Kinematic viscosity

    integer, intent(in) :: nElems

    Number of elements in input and output arrays


Derived Types

type, public :: mus_turb_wallFunc_type

Contains function pointers to compute fuction velocity and stream-wise velocity component

Components

TypeVisibilityAttributesNameInitial
logical, private :: isActive =.false.

is true if wall function is active

character(len=labelLen), private :: wall_func

Wall model function

character(len=labelLen), private :: nonlinear_solver

Nonlinear solver type

real(kind=rk), private :: vonKarman =0.4_rk

Von-Karman constant. Default = 0.4_rk

logical, private :: useVanDriest =.true.

Use vanDriest damping function to damp turbulent viscosity

type(mus_turb_wallFunc_data_type), private, allocatable:: dataOnLvl(:)

Contains data computed in turbulent wall bc routine on each level

procedure(mus_proc_calcFricVel), private, pointer, nopass:: calcFricVel=> null()

Function pointer to compute friction velocity

procedure(mus_proc_calcStreamWiseVel), private, pointer, nopass:: calcStreamWiseVel=> null()

Function pointer to compute strean-wise velocity component

type, private :: mus_turb_wallFunc_data_type

Contains friction velocity and turbulent viscosity on boundary elements on each level

Components

TypeVisibilityAttributesNameInitial
real(kind=rk), private, allocatable:: tVisc(:)

Turbulent viscosity on boundary element computed using mixing length formulation in lattice unit nu_t = (vonKarmandistToBnd)*2 * |du/dy|

real(kind=rk), private, allocatable:: velTau(:)

Friction velocity on first neighbor in normal direction in lattice unit computed from wall model

real(kind=rk), private, allocatable:: distToBnd(:)

Distance to boundary from first fluid in normal direction in lattice unit.

real(kind=rk), private, allocatable:: neighDistToBnd(:)

Distance to boundary from first fluid neighbor in normal direction in lattice unit.


Subroutines

public subroutine mus_load_turb_wallFunc(me, conf, parent)

This routine loads wall model and nonlinear solver type for nonlinear equation

Arguments

TypeIntentOptionalAttributesName
type(mus_turb_wallFunc_type), intent(inout) :: me

Turbulent wall model type to fill assign wall model

type(flu_state) :: conf

lua flu state

integer, intent(in) :: parent

bc parent handle

private subroutine assign_turb_wallFunc(me)

Assign wall models pointers to calculate friction velocity and stream-wise velocity component

Arguments

TypeIntentOptionalAttributesName
type(mus_turb_wallFunc_type), intent(inout) :: me

Turbulent wall model type to fill assign wall model

private subroutine fricVel_musker_fixedPoint(velTau, velSW, distToBnd, viscKine, nElems)

This routine computes friction velocity from Musker wall model profile using fixed-point iteration method

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(inout) :: velTau(:)

Friction velocity computed from wall model. it is inout to provide velTau from previous timestep as initial velTau for fixed-point or Newton iteration solver

real(kind=rk), intent(in) :: velSW(:)

Stream-wise velocity component from which friction velocity is computed

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine fricVel_musker_newton(velTau, velSW, distToBnd, viscKine, nElems)

This routine computes friction velocity from Musker wall model profile using Newton iteration method

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(inout) :: velTau(:)

Friction velocity computed from wall model. it is inout to provide velTau from previous timestep as initial velTau for fixed-point or Newton iteration solver

real(kind=rk), intent(in) :: velSW(:)

Stream-wise velocity component from which friction velocity is computed

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine streamWiseVel_musker(velSW, velTau, distToBnd, viscKine, nElems)

This routines computes streamWise velocity component from friction velocity and distance to boundary using Musker profile.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(out) :: velSW(:)

Stream-wise velocity component from wall model

real(kind=rk), intent(in) :: velTau(:)

Friction velocity computd from wall model

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine fricVel_reichardt_fixedPoint(velTau, velSW, distToBnd, viscKine, nElems)

This routine computes friction velocity from reichardt wall model profile using fixed-point iteration method

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(inout) :: velTau(:)

Friction velocity computed from wall model. it is inout to provide velTau from previous timestep as initial velTau for fixed-point or Newton iteration solver

real(kind=rk), intent(in) :: velSW(:)

Stream-wise velocity component from which friction velocity is computed

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine fricVel_reichardt_newton(velTau, velSW, distToBnd, viscKine, nElems)

This routine computes friction velocity from reichardt wall model profile using Newton iteration method

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(inout) :: velTau(:)

Friction velocity computed from wall model. it is inout to provide velTau from previous timestep as initial velTau for fixed-point or Newton iteration solver

real(kind=rk), intent(in) :: velSW(:)

Stream-wise velocity component from which friction velocity is computed

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine streamWiseVel_reichardt(velSW, velTau, distToBnd, viscKine, nElems)

This routines computes streamWise velocity component from friction velocity and distance to boundary using reichardt profile.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(out) :: velSW(:)

Stream-wise velocity component from wall model

real(kind=rk), intent(in) :: velTau(:)

Friction velocity computd from wall model

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine fricVel_powerLaw(velTau, velSW, distToBnd, viscKine, nElems)

This routine computes friction velocity from explicit two equation Werner and Wengle wall model. Haussmann, Marc; BARRETO, Alejandro CLARO; KOUYI, Gislain LIPEME; Rivière, Nicolas; Nirschl, Hermann; Krause, Mathias J. (2019): Large-eddy simulation coupled with wall models for turbulent channel flows at high Reynolds numbers with a lattice Boltzmann method — Application to Coriolis mass flowmeter. In Computers & Mathematics with Applications 78 (10), pp. 3285–3302. DOI: 10.1016/j.camwa.2019.04.033.

Read more…

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(inout) :: velTau(:)

Friction velocity computed from wall model. it is inout to provide velTau from previous timestep as initial velTau for fixed-point or Newton iteration solver

real(kind=rk), intent(in) :: velSW(:)

Stream-wise velocity component from which friction velocity is computed

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays

private subroutine streamWiseVel_powerLaw(velSW, velTau, distToBnd, viscKine, nElems)

This routines computes streamWise velocity component from friction velocity and distance to boundary using two equation Werner and Wengle profile. Ref to following paper for Schmitt three layer equations. Haussmann, M. et al. (2019) ‘Large-eddy simulation coupled with wall models for turbulent channel flows at high Reynolds numbers with a lattice Boltzmann method — Application to Coriolis mass flowmeter’, Computers & Mathematics with Applications. Elsevier Ltd, 78(10), pp. 3285–3302.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(out) :: velSW(:)

Stream-wise velocity component from wall model

real(kind=rk), intent(in) :: velTau(:)

Friction velocity computd from wall model

real(kind=rk), intent(in) :: distToBnd(:)

Distance to the boundary in the discrete normai direction

real(kind=rk), intent(in) :: viscKine(:)

Kinematic viscosity

integer, intent(in) :: nElems

Number of elements in input and output arrays