Skip to content
Snippets Groups Projects
Commit 1de9b15a authored by Evgeny Mortikov's avatar Evgeny Mortikov
Browse files

separating 'surface' module for aerodynamics and thermal roughness calculation

parent 563f198b
No related branches found
No related tags found
No related merge requests found
...@@ -4,7 +4,7 @@ gfortran -c -cpp -Wuninitialized srcF/sfx_phys_const.f90 ...@@ -4,7 +4,7 @@ gfortran -c -cpp -Wuninitialized srcF/sfx_phys_const.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_common.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_common.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_data.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_data.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_roughness.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_surface.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_log_param.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_log_param.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_log.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_log.f90
...@@ -13,5 +13,5 @@ gfortran -c -cpp -Wuninitialized srcF/sfx_esm_param.f90 ...@@ -13,5 +13,5 @@ gfortran -c -cpp -Wuninitialized srcF/sfx_esm_param.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_esm.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_esm.f90
gfortran -c -cpp -Wuninitialized srcF/sfx_main.f90 gfortran -c -cpp -Wuninitialized srcF/sfx_main.f90
gfortran -o sfx.exe sfx_main.o sfx_log.o sfx_log_param.o sfx_esm.o sfx_esm_param.o sfx_roughness.o sfx_data.o sfx_common.o sfx_phys_const.o gfortran -o sfx.exe sfx_main.o sfx_log.o sfx_log_param.o sfx_esm.o sfx_esm_param.o sfx_surface.o sfx_data.o sfx_common.o sfx_phys_const.o
rm *.o *.mod rm *.o *.mod
...@@ -10,7 +10,7 @@ ifeq ($(COMPILER),gnu) ...@@ -10,7 +10,7 @@ ifeq ($(COMPILER),gnu)
FC = gfortran FC = gfortran
endif endif
OBJ_F90 = sfx_phys_const.o sfx_common.o sfx_data.o sfx_roughness.o sfx_log_param.o sfx_log.o sfx_esm_param.o sfx_esm.o sfx_main.o OBJ_F90 = sfx_phys_const.o sfx_common.o sfx_data.o sfx_surface.o sfx_log_param.o sfx_log.o sfx_esm_param.o sfx_esm.o sfx_main.o
OBJ_F = OBJ_F =
OBJ = $(OBJ_F90) $(OBJ_F) OBJ = $(OBJ_F90) $(OBJ_F)
......
module sfx_data module sfx_data
!> @brief surface flux module data !> @brief surface flux module data
! macro defs.
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
! modules used ! modules used
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
......
! sfx model macro definitions ! sfx model macro definitions
!#define SFX_FORCE_DEPRECATED_CODE !#define SFX_FORCE_DEPRECATED_ESM_CODE
!#define SFX_CHECK_NAN !#define SFX_CHECK_NAN
\ No newline at end of file
...@@ -9,7 +9,7 @@ module sfx_esm ...@@ -9,7 +9,7 @@ module sfx_esm
use sfx_common use sfx_common
#endif #endif
use sfx_data use sfx_data
use sfx_roughness use sfx_surface
use sfx_esm_param use sfx_esm_param
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
...@@ -53,7 +53,7 @@ contains ...@@ -53,7 +53,7 @@ contains
! ---------------------------------------------------------------------------- ! ----------------------------------------------------------------------------
do i = 1, n do i = 1, n
#ifdef SFX_FORCE_DEPRECATED_CODE #ifdef SFX_FORCE_DEPRECATED_ESM_CODE
#else #else
meteo_cell = meteoDataType(& meteo_cell = meteoDataType(&
h = meteo%h(i), & h = meteo%h(i), &
...@@ -171,25 +171,10 @@ contains ...@@ -171,25 +171,10 @@ contains
u_dyn0 = U * kappa / log(h0_m) u_dyn0 = U * kappa / log(h0_m)
end if end if
! --- define B = log(z0_m / z0_h) ! --- define thermal roughness & B = log(z0_m / z0_h)
Re = u_dyn0 * z0_m / nu_air Re = u_dyn0 * z0_m / nu_air
if(Re <= Re_rough_min) then call get_thermal_roughness(z0_t, B, z0_m, Re, surface_type)
B = B1_rough * alog(B3_rough * Re) + B2_rough
else
! *: B4 takes into account Re value at z' ~ O(10) z0
B = B4_rough * (Re**B2_rough)
end if
! --- apply max restriction based on surface type
if (surface_type == surface_ocean) then
B = min(B, B_max_ocean)
endif
if (surface_type == surface_land) then
B = min(B, B_max_land)
end if
! --- define roughness [thermal]
z0_t = z0_m / exp(B)
! --- define relative height [thermal] ! --- define relative height [thermal]
h0_t = h / z0_t h0_t = h / z0_t
......
...@@ -12,9 +12,6 @@ module sfx_esm_param ...@@ -12,9 +12,6 @@ module sfx_esm_param
implicit none implicit none
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
! *: add lake surface & add surface type as input value
integer, public, parameter :: surface_ocean = 0 !> ocean surface
integer, public, parameter :: surface_land = 1 !> land surface
!> von Karman constant [n/d] !> von Karman constant [n/d]
real, parameter :: kappa = 0.40 real, parameter :: kappa = 0.40
...@@ -35,20 +32,4 @@ module sfx_esm_param ...@@ -35,20 +32,4 @@ module sfx_esm_param
!> --- max Ri-bulk value in stable case ( < 1 / beta_m ) !> --- max Ri-bulk value in stable case ( < 1 / beta_m )
real, parameter :: Rib_max = 0.9 / beta_m real, parameter :: Rib_max = 0.9 / beta_m
!> Re fully roughness minimum value [n/d]
real, parameter :: Re_rough_min = 16.3
!> roughness model coeff. [n/d]
!> --- transitional mode
!> B = log(z0_m / z0_t) = B1 * log(B3 * Re) + B2
real, parameter :: B1_rough = 5.0 / 6.0
real, parameter :: B2_rough = 0.45
real, parameter :: B3_rough = kappa * Pr_m
!> --- fully rough mode (Re > Re_rough_min)
!> B = B4 * Re^(B2)
real, parameter :: B4_rough =(0.14 * (30.0**B2_rough)) * (Pr_m**0.8)
real, parameter :: B_max_land = 2.0
real, parameter :: B_max_ocean = 8.0
end module sfx_esm_param end module sfx_esm_param
...@@ -9,7 +9,7 @@ module sfx_log ...@@ -9,7 +9,7 @@ module sfx_log
use sfx_common use sfx_common
#endif #endif
use sfx_data use sfx_data
use sfx_roughness use sfx_surface
use sfx_log_param use sfx_log_param
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
...@@ -160,25 +160,10 @@ contains ...@@ -160,25 +160,10 @@ contains
u_dyn0 = U * kappa / log(h0_m) u_dyn0 = U * kappa / log(h0_m)
end if end if
! --- define B = log(z0_m / z0_h) ! --- define thermal roughness & B = log(z0_m / z0_h)
Re = u_dyn0 * z0_m / nu_air Re = u_dyn0 * z0_m / nu_air
if(Re <= Re_rough_min) then call get_thermal_roughness(z0_t, B, z0_m, Re, surface_type)
B = B1_rough * alog(B3_rough * Re) + B2_rough
else
! *: B4 takes into account Re value at z' ~ O(10) z0
B = B4_rough * (Re**B2_rough)
end if
! --- apply max restriction based on surface type
if (surface_type == surface_ocean) then
B = min(B, B_max_ocean)
endif
if (surface_type == surface_land) then
B = min(B, B_max_land)
end if
! --- define roughness [thermal]
z0_t = z0_m / exp(B)
! --- define relative height [thermal] ! --- define relative height [thermal]
h0_t = h / z0_t h0_t = h / z0_t
......
...@@ -12,29 +12,10 @@ module sfx_log_param ...@@ -12,29 +12,10 @@ module sfx_log_param
implicit none implicit none
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
! *: add lake surface & add surface type as input value
integer, public, parameter :: surface_ocean = 0 !> ocean surface
integer, public, parameter :: surface_land = 1 !> land surface
!> von Karman constant [n/d] !> von Karman constant [n/d]
real, parameter :: kappa = 0.40 real, parameter :: kappa = 0.40
!> inverse Prandtl (turbulent) number in neutral conditions [n/d] !> inverse Prandtl (turbulent) number in neutral conditions [n/d]
real, parameter :: Pr_t_0_inv = 1.15 real, parameter :: Pr_t_0_inv = 1.15
!> Re fully roughness minimum value [n/d]
real, parameter :: Re_rough_min = 16.3
!> roughness model coeff. [n/d]
!> --- transitional mode
!> B = log(z0_m / z0_t) = B1 * log(B3 * Re) + B2
real, parameter :: B1_rough = 5.0 / 6.0
real, parameter :: B2_rough = 0.45
real, parameter :: B3_rough = kappa * Pr_m
!> --- fully rough mode (Re > Re_rough_min)
!> B = B4 * Re^(B2)
real, parameter :: B4_rough =(0.14 * (30.0**B2_rough)) * (Pr_m**0.8)
real, parameter :: B_max_land = 2.0
real, parameter :: B_max_ocean = 8.0
end module sfx_log_param end module sfx_log_param
#include "sfx_def.fi" #include "sfx_def.fi"
module sfx_roughness module sfx_surface
!> @brief surface roughness parameterizations !> @brief surface roughness parameterizations
! modules used ! modules used
...@@ -14,13 +14,20 @@ module sfx_roughness ...@@ -14,13 +14,20 @@ module sfx_roughness
private private
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
! public interface
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
real, parameter :: kappa = 0.40 !> von Karman constant [n/d] public :: get_charnock_roughness
public :: get_thermal_roughness
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
! public interface
! *: add surface type as input value
integer, public, parameter :: surface_ocean = 0 !> ocean surface
integer, public, parameter :: surface_land = 1 !> land surface
integer, public, parameter :: surface_lake = 2 !> lake surface
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
public :: get_charnock_roughness real, parameter :: kappa = 0.40 !> von Karman constant [n/d]
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
!> Charnock parameters !> Charnock parameters
...@@ -34,6 +41,23 @@ module sfx_roughness ...@@ -34,6 +41,23 @@ module sfx_roughness
real, parameter :: c2_charnock = Re_visc_min * nu_air * c1_charnock real, parameter :: c2_charnock = Re_visc_min * nu_air * c1_charnock
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
!> Re fully roughness minimum value [n/d]
real, parameter :: Re_rough_min = 16.3
!> roughness model coeff. [n/d]
!> --- transitional mode
!> B = log(z0_m / z0_t) = B1 * log(B3 * Re) + B2
real, parameter :: B1_rough = 5.0 / 6.0
real, parameter :: B2_rough = 0.45
real, parameter :: B3_rough = kappa * Pr_m
!> --- fully rough mode (Re > Re_rough_min)
!> B = B4 * Re^(B2)
real, parameter :: B4_rough =(0.14 * (30.0**B2_rough)) * (Pr_m**0.8)
real, parameter :: B_max_land = 2.0
real, parameter :: B_max_ocean = 8.0
real, parameter :: B_max_lake = 8.0
contains contains
! charnock roughness definition ! charnock roughness definition
...@@ -83,4 +107,43 @@ contains ...@@ -83,4 +107,43 @@ contains
end subroutine end subroutine
! -------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------
end module sfx_roughness ! thermal roughness definition
! --------------------------------------------------------------------------------
subroutine get_thermal_roughness(z0_t, B, &
z0_m, Re, surface_type)
! ----------------------------------------------------------------------------
real, intent(out) :: z0_t !> thermal roughness [m]
real, intent(out) :: B !> = log(z0_m / z0_t) [n/d]
real, intent(in) :: z0_m !> aerodynamic roughness [m]
real, intent(in) :: Re !> roughness Reynolds number [n/d]
integer, intent(in) :: surface_type !> = [ocean] || [land] || [lake]
! ----------------------------------------------------------------------------
! --- local variables
! ----------------------------------------------------------------------------
! --- define B = log(z0_m / z0_t)
if (Re <= Re_rough_min) then
B = B1_rough * alog(B3_rough * Re) + B2_rough
else
! *: B4 takes into account Re value at z' ~ O(10) z0
B = B4_rough * (Re**B2_rough)
end if
! --- apply max restriction based on surface type
if (surface_type == surface_ocean) then
B = min(B, B_max_ocean)
else if (surface_type == surface_lake) then
B = min(B, B_max_lake)
else if (surface_type == surface_land) then
B = min(B, B_max_land)
end if
! --- define roughness [thermal]
z0_t = z0_m / exp(B)
end subroutine
! --------------------------------------------------------------------------------
end module sfx_surface
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment