Skip to content
Snippets Groups Projects
sfx_z0m_all_surface.f90 9.67 KiB
Newer Older
module sfx_z0m_all_surface
   !< @brief surface roughness parameterizations 

    use sfx_phys_const
    implicit none
    public :: get_dynamic_roughness_ch
    public :: get_dynamic_roughness_map
    public :: get_dynamic_roughness_ow
    public :: get_dynamic_roughness_fetch
    public :: get_dynamic_roughness_and
    ! --------------------------------------------------------------------------------


    

    ! --------------------------------------------------------------------------------
    real, parameter, private :: kappa = 0.40         !< von Karman constant [n/d]
    ! --------------------------------------------------------------------------------

    !< Charnock parameters
    !<     z0 = Re_visc_min * (nu / u_dyn) + gamma_c * (u_dyn^2 / g)
    ! --------------------------------------------------------------------------------
    
    real, parameter :: gamma_c = 0.0144
    real, parameter :: Re_visc_min = 0.111

    real, parameter :: h_charnock = 10.0
    real, parameter :: c1_charnock = log(h_charnock * (g / gamma_c))
    real, parameter :: c2_charnock = Re_visc_min * nu_air * c1_charnock
    real, parameter :: gamma_min = 0.01
    real, parameter :: gamma_max = 0.11
    real, parameter :: f_c = 100
    real, parameter :: eps = 1
    ! --------------------------------------------------------------------------------

    

   
 contains

    ! charnock roughness definition
    ! --------------------------------------------------------------------------------
    subroutine get_dynamic_roughness_ch(z0_m, u_dyn0, U, h, maxiters)
        ! ----------------------------------------------------------------------------
        real, intent(out) :: z0_m           !< aerodynamic roughness [m]
        real, intent(out) :: u_dyn0         !< dynamic velocity in neutral conditions [m/s]

        real, intent(in) :: h               !< constant flux layer height [m]
        real, intent(in) :: U               !< abs(wind speed) [m/s]
        integer, intent(in) :: maxiters     !< maximum number of iterations
        ! ----------------------------------------------------------------------------

        ! --- local variables
        real :: Uc                  ! wind speed at h_charnock [m/s]

        real :: a, b, c, c_min
        real :: f

        integer :: i, j
        ! ----------------------------------------------------------------------------

        Uc = U
        a = 0.0
        b = 25.0
        c_min = log(h_charnock) / kappa

        do i = 1, maxiters
            f = c1_charnock - 2.0 * log(Uc)
            do j = 1, maxiters
                c = (f + 2.0 * log(b)) / kappa
                ! looks like the check should use U10 instead of U
                !   but note that a1 should be set = 0 in cycle beforehand
                if (U <= 8.0e0) a = log(1.0 + c2_charnock * ((b / Uc)**3)) / kappa
                c = max(c - a, c_min)
                b = c
            end do
            z0_m = h_charnock * exp(-c * kappa)
            z0_m = max(z0_m, 0.000015e0)
            Uc = U * log(h_charnock / z0_m) / log(h / z0_m)
        end do

        ! --- define dynamic velocity in neutral conditions
        u_dyn0 = Uc / c

    end subroutine
    ! --------------------------------------------------------------------------------
    
      subroutine get_dynamic_roughness_ow(z0_m, u_dyn0, U, h, maxiters)
      !Owen 1964
        ! ----------------------------------------------------------------------------
        real, intent(out) :: z0_m           !< aerodynamic roughness [m]
        real, intent(out) :: u_dyn0         !< dynamic velocity in neutral conditions [m/s]

        real, intent(in) :: h               !< constant flux layer height [m]
        real, intent(in) :: U               !< abs(wind speed) [m/s]
        integer, intent(in) :: maxiters     !< maximum number of iterations
        ! ----------------------------------------------------------------------------

        ! --- local variables
        real :: Uc                  ! wind speed at h_charnock [m/s]

        real :: b1, b2, Cuz, betta_u, nu_m, C_z0,c
        real :: f

        integer :: i, j
        ! ----------------------------------------------------------------------------
        Uc=U
        C_z0=0.007
        betta_u=0.111
        nu_m=0.0000133
        b1=log(h*g/C_z0)
        b2=betta_u*nu_m*g/C_z0
        Cuz=25.0
        do i = 1, maxiters
        f = c1_charnock - 2.0 * log(Uc)
        c = (f + 2.0 * log(Cuz)) / kappa
        Cuz=(1.0/kappa)*(b1+log(U/Cuz)-log(b2+(U/Cuz)*(U/Cuz)))
        if(Cuz==0.0) exit
        z0_m=h*exp(-kappa*Cuz)
        end do
        

        u_dyn0 = Uc / c
        !  write(*,*) z0_m, 'oewn'
        end subroutine
     ! --------------------------------------------------------------------------------
   
        subroutine get_dynamic_roughness_fetch(z0_m, u_dyn0, U, depth, h, maxiters)
        ! ----------------------------------------------------------------------------
        real, intent(out) :: z0_m           !< aerodynamic roughness [m]
        real, intent(out) :: u_dyn0         !< dynamic velocity in neutral conditions [m/s]

        real, intent(in) :: U               !< abs(wind speed) [m/s]
        real, intent(in) :: depth           !< depth [m]
        real, intent(in) :: h               !< constant flux layer height [m]
        integer, intent(in) :: maxiters     !< maximum number of iterations
        ! ----------------------------------------------------------------------------

        ! --- local variables
        real :: Uc                  ! wind speed at h_charnock [m/s]

        real :: a, b, c, c_min
        real :: f

        real :: A_lake, B_lake, gamma_c, fetch, c1_charnock_lake, c2_charnock_lake

        integer :: i, j
        ! ----------------------------------------------------------------------------

        Uc = U
        a = 0.0
        b = 25.0
        c_min = log(h_charnock) / kappa
        
        fetch = 25.0 * depth !25.0 * depth

        !<     z0 = Re_visc_min * (nu / u_dyn) + gamma_c * (u_dyn^2 / g)
        !<     gamma_c = gamma_min + (gamma_max - gamma_min) * exp(-min(A_lake, B_lake))
        !<     А_lake = (fetch * g / U^2)^(1/3) / f_c
        !<     B_lake = eps (sqrt(depth * g)/U)

        do i = 1, maxiters
            A_lake = ((fetch * g / (U)**2)**(1/3)) / f_c
            B_lake = eps * (sqrt(depth * g)/U)
            gamma_c =  gamma_min + (gamma_max - gamma_min) * exp(-min(A_lake, B_lake))
            !write(*,*) A_lake
            !write(*,*) B_lake
            c1_charnock_lake = log(h_charnock * (g / gamma_c))
            c2_charnock_lake = Re_visc_min * nu_air * c1_charnock_lake
            f = c1_charnock_lake - 2.0 * log(Uc)
            do j = 1, maxiters
                c = (f + 2.0 * log(b)) / kappa
                if (U <= 8.0e0) a = log(1.0 + c2_charnock_lake * ((b / Uc)**3)) / kappa
                c = max(c - a, c_min)
                b = c
            end do
            z0_m = h_charnock * exp(-c * kappa)
            z0_m = max(z0_m, 0.000015e0)
            Uc = U * log(h_charnock / z0_m) / log(h / z0_m)
        end do

        ! --- define dynamic velocity in neutral conditions
        u_dyn0 = Uc / c

     end subroutine
   
   
subroutine get_dynamic_roughness_map(z0_m, u_dyn0, U, h, z0m_map)
      
        ! ----------------------------------------------------------------------------
        real, intent(out) :: z0_m           !< aerodynamic roughness [m]
        real, intent(out) :: u_dyn0         !< dynamic velocity in neutral conditions [m/s]

        real, intent(in) :: h               !< constant flux layer height [m]
        real, intent(in) :: z0m_map           !< aerodynamic roughness from map[m]
        real, intent(in) :: U               !< abs(wind speed) [m/s]
        ! ----------------------------------------------------------------------------
        real :: h0_m
        
        z0_m=z0m_map  
        h0_m = h / z0_m
        u_dyn0 = U * kappa / log(h0_m)
        write(*,*) z0_m, 'map'
        end subroutine
     ! --------------------------------------------------------------------------------

        subroutine get_dynamic_roughness_and(z0_m, u_dyn0, U, h, z0m_map)
      
            ! ----------------------------------------------------------------------------
            real, intent(out) :: z0_m           !< aerodynamic roughness [m]
            real, intent(out) :: u_dyn0         !< dynamic velocity in neutral conditions [m/s]
    
            real, intent(in) :: h               !< constant flux layer height [m]
            real, intent(in) :: z0m_map           !< aerodynamic roughness from map[m]
            real, intent(in) :: U               !< abs(wind speed) [m/s]
            ! ----------------------------------------------------------------------------
            real :: h0_m, u_star_prev, nu, g
            
            
            real :: tolerance
            integer :: max_iterations, iter
        
            
            nu = 1.7e-5  
            g = 9.81     
            u_dyn0 = 0.2 
            
            
            tolerance = 1.0e-5  
            max_iterations = 10
        
            
            do iter = 1, max_iterations
                
                u_star_prev = u_dyn0
        
                
                z0_m = (0.135 * nu / u_dyn0) + (0.035 * u_dyn0**2 / g) * &
                      (1.0 + exp(-((u_dyn0 - 0.18) / 0.1)**2))
        
                      h0_m = h / z0_m
                      u_dyn0 =  U * kappa / log(h0_m) 
        
                
                if (abs(u_dyn0 - u_star_prev) < tolerance) exit
            end do

            u_dyn0 = U * kappa / log(h0_m)
           
             
            end subroutine
         ! --------------------------------------------------------------------------------


   
    end module  sfx_z0m_all_surface