27module monin_obukhov_mod
29use constants_mod,
only: grav, vonkarm
30use mpp_mod,
only: input_nml_file
31use fms_mod,
only: error_mesg, fatal, check_nml_error, &
32 mpp_pe, mpp_root_pe, stdlog, &
36use platform_mod,
only: r4_kind, r8_kind
41 public :: monin_obukhov_init
42 public :: monin_obukhov_end
52 module procedure mo_drag_0d_r4, mo_drag_0d_r8
53 module procedure mo_drag_1d_r4, mo_drag_1d_r8
54 module procedure mo_drag_2d_r4, mo_drag_2d_r8
60 module procedure mo_profile_0d_r4, mo_profile_0d_r8
61 module procedure mo_profile_1d_r4, mo_profile_1d_r8
62 module procedure mo_profile_2d_r4, mo_profile_2d_r8
63 module procedure mo_profile_0d_n_r4, mo_profile_0d_n_r8
64 module procedure mo_profile_1d_n_r4, mo_profile_1d_n_r8
65 module procedure mo_profile_2d_n_r4, mo_profile_2d_n_r8
70 module procedure mo_diff_0d_n_r4, mo_diff_0d_n_r8
71 module procedure mo_diff_0d_1_r4, mo_diff_0d_1_r8
72 module procedure mo_diff_1d_n_r4, mo_diff_1d_n_r8
73 module procedure mo_diff_1d_1_r4, mo_diff_1d_1_r8
74 module procedure mo_diff_2d_n_r4, mo_diff_2d_n_r8
75 module procedure mo_diff_2d_1_r4, mo_diff_2d_1_r8
80 module procedure stable_mix_0d_r4, stable_mix_0d_r8
81 module procedure stable_mix_1d_r4, stable_mix_1d_r8
82 module procedure stable_mix_2d_r4, stable_mix_2d_r8
83 module procedure stable_mix_3d_r4, stable_mix_3d_r8
87 module procedure mo_integral_m_r4, mo_integral_m_r8
91 module procedure mo_integral_tq_r4, mo_integral_tq_r8
95 module procedure mo_derivative_m_r4, mo_derivative_m_r8
99 module procedure mo_derivative_t_r4, mo_derivative_t_r8
107#include<file_version.h>
113real(kind=r8_kind) :: rich_crit = 2.0_r8_kind
114real(kind=r8_kind) :: drag_min_heat = 1.0e-05_r8_kind
115real(kind=r8_kind) :: drag_min_moist = 1.0e-05_r8_kind
116real(kind=r8_kind) :: drag_min_mom = 1.0e-05_r8_kind
117logical :: neutral = .false.
118integer :: stable_option = 1
119real(kind=r8_kind) :: zeta_trans = 0.5_r8_kind
120logical :: new_mo_option = .false.
123namelist /monin_obukhov_nml/ rich_crit, neutral, drag_min_heat, &
124 drag_min_moist, drag_min_mom, &
125 stable_option, zeta_trans, new_mo_option
131real(kind=r8_kind),
parameter :: small = 1.0e-04_r8_kind
132real(kind=r8_kind) :: b_stab, r_crit, lambda, rich_trans
133real(kind=r8_kind) :: sqrt_drag_min_heat, sqrt_drag_min_moist, sqrt_drag_min_mom
134logical :: module_is_initialized = .false.
141subroutine monin_obukhov_init
143integer :: ierr, io, logunit
147 read (input_nml_file, nml=monin_obukhov_nml, iostat=io)
148 ierr = check_nml_error(io,
"monin_obukhov_nml")
152 if ( mpp_pe() == mpp_root_pe() )
then
153 call write_version_number(
'MONIN_OBUKOV_MOD', version)
155 write (logunit, nml=monin_obukhov_nml)
160if(rich_crit.le.0.25_r8_kind)
call error_mesg( &
161 'MONIN_OBUKHOV_INIT in MONIN_OBUKHOV_MOD', &
162 'rich_crit in monin_obukhov_mod must be > 0.25', fatal)
164if(drag_min_heat.le.0.0_r8_kind)
call error_mesg( &
165 'MONIN_OBUKHOV_INIT in MONIN_OBUKHOV_MOD', &
166 'drag_min_heat in monin_obukhov_mod must be >= 0.0', fatal)
168if(drag_min_moist.le.0.0_r8_kind)
call error_mesg( &
169 'MONIN_OBUKHOV_INIT in MONIN_OBUKHOV_MOD', &
170 'drag_min_moist in monin_obukhov_mod must be >= 0.0', fatal)
172if(drag_min_mom.le.0.0_r8_kind)
call error_mesg( &
173 'MONIN_OBUKHOV_INIT in MONIN_OBUKHOV_MOD', &
174 'drag_min_mom in monin_obukhov_mod must be >= 0.0', fatal)
176if(stable_option < 1 .or. stable_option > 2)
call error_mesg( &
177 'MONIN_OBUKHOV_INIT in MONIN_OBUKHOV_MOD', &
178 'the only allowable values of stable_option are 1 and 2', fatal)
180if(stable_option == 2 .and. zeta_trans < 0)
call error_mesg( &
181 'MONIN_OBUKHOV_INIT in MONIN_OBUKHOV_MOD', &
182 'zeta_trans must be positive', fatal)
184b_stab = 1.0_r8_kind/rich_crit
185r_crit = 0.95_r8_kind*rich_crit
188sqrt_drag_min_heat = 0.0_r8_kind
189if(drag_min_heat.ne.0.0_r8_kind) sqrt_drag_min_heat = sqrt(drag_min_heat)
191sqrt_drag_min_moist = 0.0_r8_kind
192if(drag_min_moist.ne.0.0_r8_kind) sqrt_drag_min_moist = sqrt(drag_min_moist)
194sqrt_drag_min_mom = 0.0_r8_kind
195if(drag_min_mom.ne.0.0_r8_kind) sqrt_drag_min_mom = sqrt(drag_min_mom)
197lambda = 1.0_r8_kind + (5.0_r8_kind - b_stab)*zeta_trans
198rich_trans = zeta_trans/(1.0_r8_kind + 5.0_r8_kind*zeta_trans)
200module_is_initialized = .true.
203end subroutine monin_obukhov_init
207subroutine monin_obukhov_end
209module_is_initialized = .false.
211end subroutine monin_obukhov_end
215#include "monin_obukhov_r4.fh"
216#include "monin_obukhov_r8.fh"
218end module monin_obukhov_mod
Compute surface drag coefficients.