Newer
Older

Evgeny Mortikov
committed
module sfx_surface

Evgeny Mortikov
committed
! modules used
! --------------------------------------------------------------------------------
use sfx_phys_const

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------
! directives list
! --------------------------------------------------------------------------------
implicit none
! --------------------------------------------------------------------------------

Evgeny Mortikov
committed
! public interface

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------
#if defined(INCLUDE_CXX)
public :: set_c_struct_sfx_surface_param_values
#endif

Evgeny Mortikov
committed
public :: get_charnock_roughness
public :: get_thermal_roughness
public :: get_dynamic_roughness_all
public :: get_dynamic_roughness_definition

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------

Evgeny Mortikov
committed
! *: 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
integer, public, parameter :: surface_snow = 3 !< snow covered surface
integer, public, parameter :: surface_forest = 4 !< snow covered surface
integer, public, parameter :: surface_user = 5 !< snow covered surface
character(len = 16), parameter :: surface_ocean_tag = 'ocean'
character(len = 16), parameter :: surface_land_tag = 'land'
character(len = 16), parameter :: surface_lake_tag = 'lake'
character(len = 16), parameter :: surface_snow_tag = 'snow'
character(len = 16), parameter :: surface_forest_tag = 'forest'
character(len = 16), parameter :: surface_user_tag = 'user'

Evgeny Mortikov
committed
integer, public, parameter :: z0m_ch = 0
integer, public, parameter :: z0m_fe = 1
integer, public, parameter :: z0m_ow = 2
integer, public, parameter :: z0m_map = 3
integer, public, parameter :: z0m_user = 4
character(len = 16), parameter :: z0m_tag_ch = 'charnock'
character(len = 16), parameter :: z0m_tag_fe = 'fetch'
character(len = 16), parameter :: z0m_tag_ow = 'owen'
character(len = 16), parameter :: z0m_tag_map = 'map'
character(len = 16), parameter :: z0m_tag_user = 'user'
integer, public, parameter :: ocean_z0m_id = z0m_ch !< ocean surface
integer, public, parameter :: land_z0m_id = z0m_map !< land surface
integer, public, parameter :: lake_z0m_id = z0m_fe !< lake surface
integer, public, parameter :: snow_z0m_id = z0m_ow !< snow covered surface
integer, public, parameter :: forest_z0m_id = z0m_map !< forest csurface
integer, public, parameter :: usersf_z0m_id = z0m_map !< user surface
real, public, parameter :: depth = 1.0

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

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------
!< Charnock parameters
!< z0 = Re_visc_min * (nu / u_dyn) + gamma_c * (u_dyn^2 / g)

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------
!real, parameter :: gamma_c = 0.0144
!real, parameter :: Re_visc_min = 0.111

Evgeny Mortikov
committed
!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 * (g / gamma_c)

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------

Evgeny Mortikov
committed
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

Evgeny Mortikov
committed
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)

Evgeny Mortikov
committed
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

Evgeny Mortikov
committed
contains
! surface type definition
! --------------------------------------------------------------------------------
function get_surface_id(tag) result(id)
implicit none
character(len=*), intent(in) :: tag
integer :: id
id = - 1
if (trim(tag) == trim(surface_ocean_tag)) then
id = surface_ocean
else if (trim(tag) == trim(surface_land_tag)) then
id = surface_land
else if (trim(tag) == trim(surface_lake_tag)) then
id = surface_lake
else if (trim(tag) == trim(surface_snow_tag)) then
id = surface_snow
else if (trim(tag) == trim(surface_forest_tag)) then
id = surface_forest
else if (trim(tag) == trim(surface_user_tag)) then
id = surface_user
end if
end function
function get_surface_tag(id) result(tag)
implicit none
integer :: id
character(len=:), allocatable :: tag
tag = 'undefined'
if (id == surface_ocean) then
tag = surface_ocean_tag
else if (id == surface_land) then
tag = surface_land_tag
else if (id == surface_lake) then
tag = surface_lake_tag
else if (id == surface_snow) then
tag = surface_snow_tag
else if (id == surface_forest) then
tag = surface_forest_tag
else if (id == surface_user) then
tag = surface_user_tag
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
function get_surface_z0m_id(tag_z0m) result(z0m_id)
implicit none
character(len=*), intent(in) :: tag_z0m
integer :: z0m_id
z0m_id = - 1
if (trim(tag_z0m) == trim(z0m_tag_ch)) then
z0m_id = z0m_ch
else if (trim(tag_z0m) == trim(z0m_tag_fe)) then
z0m_id = z0m_fe
else if (trim(tag_z0m) == trim(z0m_tag_ow)) then
z0m_id = z0m_ow
else if (trim(tag_z0m) == trim(z0m_tag_map)) then
z0m_id = z0m_map
else if (trim(tag_z0m) == trim(z0m_tag_user)) then
z0m_id = z0m_user
end if
end function
function get_surface_z0m_tag(z0m_id) result(tag_z0m)
implicit none
integer :: z0m_id
character(len=:), allocatable :: tag_z0m
tag_z0m = 'undefined'
if (z0m_id == z0m_ch) then
tag_z0m = z0m_tag_ch
else if (z0m_id == z0m_fe) then
tag_z0m = z0m_tag_fe
else if (z0m_id == z0m_ow) then
tag_z0m = z0m_tag_ow
else if (z0m_id == z0m_map) then
tag_z0m = z0m_tag_map
else if (z0m_id == z0m_user) then
tag_z0m = z0m_tag_user
end if
end function
#if defined(INCLUDE_CXX)
subroutine set_c_struct_sfx_surface_param_values(surface_param)
use sfx_data
implicit none
type (sfx_surface_param), intent(inout) :: surface_param
surface_param%surface_ocean = surface_ocean
surface_param%surface_land = surface_land
surface_param%surface_lake = surface_lake
surface_param%kappa = kappa
surface_param%gamma_c = gamma_c
surface_param%Re_visc_min = Re_visc_min
surface_param%h_charnock = h_charnock
surface_param%c1_charnock = c1_charnock
surface_param%c2_charnock = c2_charnock
surface_param%Re_rough_min = Re_rough_min
surface_param%B1_rough = B1_rough
surface_param%B2_rough = B2_rough
surface_param%B3_rough = B3_rough
surface_param%B4_rough = B4_rough
surface_param%B_max_lake = B_max_lake
surface_param%B_max_ocean = B_max_ocean
surface_param%B_max_land = B_max_land
end subroutine set_c_struct_sfx_surface_param_values
#endif

Evgeny Mortikov
committed
! --------------------------------------------------------------------------------
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
subroutine get_dynamic_roughness_definition(surface_type, ocean_z0m_id, land_z0m_id, lake_z0m_id, snow_z0m_id, &
forest_z0m_id, usersf_z0m_id, z0m_id)
! ----------------------------------------------------------------------------
integer, intent(out) :: z0m_id
integer, intent(in) :: surface_type
integer, intent(in) :: ocean_z0m_id
integer, intent(in) :: land_z0m_id
integer, intent(in) :: lake_z0m_id
integer, intent(in) :: snow_z0m_id
integer, intent(in) :: forest_z0m_id
integer, intent(in) :: usersf_z0m_id
! ---------------------------------------------------------------------------
if (surface_type == surface_ocean) then
z0m_id = ocean_z0m_id
else if (surface_type == surface_land) then
z0m_id = land_z0m_id
else if (surface_type == surface_snow) then
z0m_id = snow_z0m_id
else if (surface_type == surface_lake) then
z0m_id = lake_z0m_id
else if (surface_type == surface_forest) then
z0m_id = forest_z0m_id
else if (surface_type == surface_user) then
z0m_id = usersf_z0m_id
end if
end subroutine
! ----------------------------------------------------------------------------
subroutine get_dynamic_roughness_all(z0_m, u_dyn0, U, depth, h,&
maxiters, z0m_map, z0m_id)
! ----------------------------------------------------------------------------
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]
real, intent(in) :: z0m_map !< aerodynamic roughness from map [m]
integer, intent(in) :: maxiters !< maximum number of iterations
integer, intent(in) :: z0m_id
! ---------------------------------------------------------------------------
if (z0m_id == z0m_ch) then
call get_dynamic_roughness_ch(z0_m, u_dyn0, U, h, maxiters)
else if (z0m_id == z0m_fe) then
call get_dynamic_roughness_fetch(z0_m, u_dyn0, U, depth, h, maxiters)
else if (z0m_id == z0m_ow) then
call get_dynamic_roughness_ow(z0_m, u_dyn0, U, h, maxiters)
else if (z0m_id == z0m_map) then
call get_dynamic_roughness_map(z0_m, u_dyn0, U, h, z0m_map)
else if (z0m_id == z0m_user) then
write(*, *) 'z0m_user'
end if
end subroutine
! --------------------------------------------------------------------------------
subroutine get_charnock_roughness(z0_m, u_dyn0, U, h, maxiters)

Evgeny Mortikov
committed
! ----------------------------------------------------------------------------
real, intent(out) :: z0_m !< aerodynamic roughness [m]
real, intent(out) :: u_dyn0 !< dynamic velocity in neutral conditions [m/s]

Evgeny Mortikov
committed
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

Evgeny Mortikov
committed
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
! ----------------------------------------------------------------------------
! --- 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
! --------------------------------------------------------------------------------

Evgeny Mortikov
committed
! 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]

Evgeny Mortikov
committed
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]

Evgeny Mortikov
committed
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
! ----------------------------------------------------------------------------
! --- 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