Newer
Older

Evgeny Mortikov
committed
module sfx_data
! modules used
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
! directives list
! --------------------------------------------------------------------------------
implicit none
private
! --------------------------------------------------------------------------------
! public interface
! --------------------------------------------------------------------------------
public :: allocate_meteo_vec, deallocate_meteo_vec
public :: allocate_sfx_vec, deallocate_sfx_vec

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

Evgeny Mortikov
committed
type, public :: meteoDataType
real :: h !< constant flux layer height [m]
real :: U !< abs(wind speed) at 'h' [m/s]
real :: dT !< difference between potential temperature at 'h' and at surface [K]
real :: Tsemi !< semi-sum of potential temperature at 'h' and at surface [K]
real :: dQ !< difference between humidity at 'h' and at surface [g/g]
real :: z0_m !< surface aerodynamic roughness (should be < 0 for water bodies surface)

Evgeny Mortikov
committed
end type
!> @brief meteorological input for surface flux calculation
!> &details using arrays as input

Evgeny Mortikov
committed
type, public :: meteoDataVecType
real, allocatable :: h(:) !< constant flux layer height [m]
real, allocatable :: U(:) !< abs(wind speed) at 'h' [m/s]
real, allocatable :: dT(:) !< difference between potential temperature at 'h' and at surface [K]
real, allocatable :: Tsemi(:) !< semi-sum of potential temperature at 'h' and at surface [K]
real, allocatable :: dQ(:) !< difference between humidity at 'h' and at surface [g/g]
real, allocatable :: z0_m(:) !< surface aerodynamic roughness (should be < 0 for water bodies surface)

Evgeny Mortikov
committed
end type
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------

Evgeny Mortikov
committed
type, public :: sfxDataType
real :: zeta !< = z/L [n/d]
real :: Rib !< bulk Richardson number [n/d]
real :: Re !< Reynolds number [n/d]
real :: B !< = log(z0_m / z0_h) [n/d]
real :: z0_m !< aerodynamic roughness [m]
real :: z0_t !< thermal roughness [m]
real :: Rib_conv_lim !< Ri-bulk convection critical value [n/d]
real :: Cm !< transfer coefficient for momentum [n/d]
real :: Ct !< transfer coefficient for heat [n/d]
real :: Km !< eddy viscosity coeff. at h [m^2/s]
real :: Pr_t_inv !< inverse turbulent Prandtl number at h [n/d]

Evgeny Mortikov
committed
end type
!> @brief surface flux output data
!> &details using arrays as output

Evgeny Mortikov
committed
type, public :: sfxDataVecType
real, allocatable :: zeta(:) !< = z/L [n/d]
real, allocatable :: Rib(:) !< bulk Richardson number [n/d]
real, allocatable :: Re(:) !< Reynolds number [n/d]
real, allocatable :: B(:) !< = log(z0_m / z0_h) [n/d]
real, allocatable :: z0_m(:) !< aerodynamic roughness [m]
real, allocatable :: z0_t(:) !< thermal roughness [m]
real, allocatable :: Rib_conv_lim(:) !< Ri-bulk convection critical value [n/d]
real, allocatable :: Cm(:) !< transfer coefficient for momentum [n/d]
real, allocatable :: Ct(:) !< transfer coefficient for heat [n/d]
real, allocatable :: Km(:) !< eddy viscosity coeff. at h [m^2/s]
real, allocatable :: Pr_t_inv(:) !< inverse turbulent Prandtl number at h [n/d]

Evgeny Mortikov
committed
end type
! --------------------------------------------------------------------------------
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
! --------------------------------------------------------------------------------
subroutine allocate_meteo_vec(meteo, n)
!> @brief allocate meteo data vector
! ----------------------------------------------------------------------------
type (meteoDataVecType), intent(inout) :: meteo
integer, intent(in) :: n
! ----------------------------------------------------------------------------
allocate(meteo%h(n))
allocate(meteo%U(n))
allocate(meteo%dT(n))
allocate(meteo%Tsemi(n))
allocate(meteo%dQ(n))
allocate(meteo%z0_m(n))
end subroutine allocate_meteo_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine deallocate_meteo_vec(meteo)
!> @brief deallocate meteo data vector
! ----------------------------------------------------------------------------
type (meteoDataVecType), intent(inout) :: meteo
! ----------------------------------------------------------------------------
deallocate(meteo%h)
deallocate(meteo%U)
deallocate(meteo%dT)
deallocate(meteo%Tsemi)
deallocate(meteo%dQ)
deallocate(meteo%z0_m)
end subroutine deallocate_meteo_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine allocate_sfx_vec(sfx, n)
!> @brief allocate surface fluxes data vector
! ----------------------------------------------------------------------------
type (sfxDataVecType), intent(inout) :: sfx
integer, intent(in) :: n
! ----------------------------------------------------------------------------
allocate(sfx%zeta(n))
allocate(sfx%Rib(n))
allocate(sfx%Re(n))
allocate(sfx%B(n))
allocate(sfx%z0_m(n))
allocate(sfx%z0_t(n))
allocate(sfx%Rib_conv_lim(n))
allocate(sfx%Cm(n))
allocate(sfx%Ct(n))
allocate(sfx%Km(n))
allocate(sfx%Pr_t_inv(n))
end subroutine allocate_sfx_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
subroutine deallocate_sfx_vec(sfx)
! ----------------------------------------------------------------------------
type (sfxDataVecType), intent(inout) :: sfx
! ----------------------------------------------------------------------------
deallocate(sfx%zeta)
deallocate(sfx%Rib)
deallocate(sfx%Re)
deallocate(sfx%B)
deallocate(sfx%z0_m)
deallocate(sfx%z0_t)
deallocate(sfx%Rib_conv_lim)
deallocate(sfx%Cm)
deallocate(sfx%Ct)
deallocate(sfx%Km)
deallocate(sfx%Pr_t_inv)
end subroutine deallocate_sfx_vec
! --------------------------------------------------------------------------------
! --------------------------------------------------------------------------------
!> @brief helper subroutine for copying data in sfxDataVecType
subroutine push_sfx_data(sfx, sfx_cell, idx)
! ----------------------------------------------------------------------------
type (sfxDataVecType), intent(inout) :: sfx
type (sfxDataType), intent(in) :: sfx_cell
integer, intent(in) :: idx
! ----------------------------------------------------------------------------
sfx%zeta(idx) = sfx_cell%zeta
sfx%Rib(idx) = sfx_cell%Rib
sfx%Re(idx) = sfx_cell%Re
sfx%B(idx) = sfx_cell%B
sfx%z0_m(idx) = sfx_cell%z0_m
sfx%z0_t(idx) = sfx_cell%z0_t
sfx%Rib_conv_lim(idx) = sfx_cell%Rib_conv_lim
sfx%Cm(idx) = sfx_cell%Cm
sfx%Ct(idx) = sfx_cell%Ct
sfx%Km(idx) = sfx_cell%Km
sfx%Pr_t_inv(idx) = sfx_cell%Pr_t_inv
end subroutine push_sfx_data
! --------------------------------------------------------------------------------

Evgeny Mortikov
committed
end module sfx_data