1 |
! $Header: /u/gcmpack/MITgcm/pkg/atm_phys/vert_turb_driver_mod.F90,v 1.1 2013/05/08 22:14:15 jmc Exp $ |
2 |
! $Name: $ |
3 |
|
4 |
module vert_turb_driver_mod |
5 |
|
6 |
!----------------------------------------------------------------------- |
7 |
! |
8 |
! driver for compuing vertical diffusion coefficients |
9 |
! |
10 |
! choose either: |
11 |
! 1) mellor-yamada 2.5 (with tke) |
12 |
! 2) non-local K scheme |
13 |
! |
14 |
!----------------------------------------------------------------------- |
15 |
!---------------- modules --------------------- |
16 |
|
17 |
use my25_turb_mod, only: my25_turb_init, my25_turb_end, & |
18 |
my25_turb, tke_surf, tke |
19 |
|
20 |
use diffusivity_mod, only: diffusivity, molecular_diff |
21 |
|
22 |
use shallow_conv_mod, only: shallow_conv_init, shallow_conv |
23 |
|
24 |
!use diag_manager_mod, only: register_diag_field, send_data |
25 |
|
26 |
!use time_manager_mod, only: time_type, get_time, operator(-) |
27 |
|
28 |
use gcm_params_mod, only: gcm_LEN_MBUF, gcm_SQZ_R, gcm_stdMsgUnit |
29 |
use constants_mod, only: rdgas, rvgas, kappa |
30 |
|
31 |
!use fms_mod, only: error_mesg, open_namelist_file, file_exist, & |
32 |
! check_nml_error, mpp_root_pe, & |
33 |
! mpp_pe, close_file, FATAL, stdlog, & |
34 |
! write_version_number |
35 |
|
36 |
implicit none |
37 |
private |
38 |
|
39 |
!---------------- interfaces --------------------- |
40 |
|
41 |
public vert_turb_driver_init, vert_turb_driver_end, vert_turb_driver |
42 |
|
43 |
!----------------------------------------------------------------------- |
44 |
!--------------------- version number ---------------------------------- |
45 |
|
46 |
character(len=128) :: version = '$Id: vert_turb_driver_mod.F90,v 1.1 2013/05/08 22:14:15 jmc Exp $' |
47 |
character(len=128) :: tag = '$Name: $' |
48 |
|
49 |
!----------------------------------------------------------------------- |
50 |
real, parameter :: p00 = 1000.0E2 |
51 |
real, parameter :: p00inv = 1./p00 |
52 |
real, parameter :: d622 = rdgas/rvgas |
53 |
real, parameter :: d378 = 1.-d622 |
54 |
real, parameter :: d608 = d378/d622 |
55 |
|
56 |
!---------------- private data ------------------- |
57 |
|
58 |
logical :: do_init = .true. |
59 |
|
60 |
real :: gust_zi = 1000. ! constant for computed gustiness (meters) |
61 |
|
62 |
!----------------------------------------------------------------------- |
63 |
!-------------------- namelist ----------------------------------------- |
64 |
|
65 |
logical :: do_shallow_conv = .false. |
66 |
logical :: do_mellor_yamada = .true. |
67 |
logical :: use_tau = .true. |
68 |
logical :: do_molecular_diffusion = .false. |
69 |
|
70 |
character(len=24) :: gust_scheme = 'constant' ! valid schemes are: |
71 |
! => 'constant' |
72 |
! => 'beljaars' |
73 |
real :: constant_gust = 1.0 |
74 |
|
75 |
namelist /vert_turb_driver_nml/ do_shallow_conv, do_mellor_yamada, & |
76 |
gust_scheme, constant_gust, use_tau, & |
77 |
do_molecular_diffusion |
78 |
|
79 |
!-------------------- diagnostics fields ------------------------------- |
80 |
|
81 |
integer :: id_tke, id_lscale, id_lscale_0, id_z_pbl, id_gust, & |
82 |
id_diff_t, id_diff_m, id_diff_sc, id_z_full, id_z_half,& |
83 |
id_uwnd, id_vwnd |
84 |
|
85 |
real :: missing_value = -999. |
86 |
|
87 |
character(len=9) :: mod_name = 'vert_turb' |
88 |
|
89 |
!----------------------------------------------------------------------- |
90 |
|
91 |
contains |
92 |
|
93 |
!####################################################################### |
94 |
|
95 |
subroutine vert_turb_driver (is, js, Time, Time_next, dt, frac_land, & |
96 |
p_half, p_full, z_half, z_full, & |
97 |
u_star, b_star, rough, & |
98 |
uu, vv, tt, qq, & |
99 |
diff_t, diff_m, gust, & |
100 |
myThid, mask, kbot ) |
101 |
! u, v, t, q, um, vm, tm, qm, & |
102 |
! udt, vdt, tdt, qdt, & |
103 |
|
104 |
!----------------------------------------------------------------------- |
105 |
integer, intent(in) :: is, js |
106 |
!type(time_type), intent(in) :: Time, Time_next |
107 |
real, intent(in) :: Time, Time_next |
108 |
real, intent(in) :: dt |
109 |
real, intent(in), dimension(:,:) :: frac_land, & |
110 |
u_star, b_star, rough |
111 |
real, intent(in), dimension(:,:,:) :: p_half, p_full, & |
112 |
z_half, z_full, & |
113 |
uu, vv, tt, qq |
114 |
! u, v, t, q, um, vm, tm, qm, & |
115 |
! udt, vdt, tdt, qdt |
116 |
real, intent(out), dimension(:,:,:) :: diff_t, diff_m |
117 |
real, intent(out), dimension(:,:) :: gust |
118 |
integer, intent(in) :: myThid |
119 |
real, intent(in),optional, dimension(:,:,:) :: mask |
120 |
integer, intent(in),optional, dimension(:,:) :: kbot |
121 |
!----------------------------------------------------------------------- |
122 |
real , dimension(size(tt,1),size(tt,2),size(tt,3)) :: ape, thv |
123 |
!logical, dimension(size(tt,1),size(tt,2),size(tt,3)+1) :: lmask |
124 |
!real , dimension(size(tt,1),size(tt,2),size(tt,3)+1) :: diag3 |
125 |
real , dimension(size(tt,1),size(tt,2),size(tt,3)+1) :: el |
126 |
real , dimension(size(tt,1),size(tt,2)) :: el0, z_pbl |
127 |
real , dimension(size(diff_t,1),size(diff_t,2), & |
128 |
size(diff_t,3)) :: diff_sc |
129 |
!real , dimension(size(t,1),size(t,2),size(t,3)) :: tt, qq, uu, vv |
130 |
real :: dt_tke |
131 |
integer :: ie, je, nlev, sec, day |
132 |
logical :: used |
133 |
!----------------------------------------------------------------------- |
134 |
!----------------------- vertical turbulence --------------------------- |
135 |
!----------------------------------------------------------------------- |
136 |
|
137 |
if (do_init) CALL PRINT_ERROR( & |
138 |
'vert_turb_driver in vert_turb_driver_mod'// & |
139 |
'initialization has not been called', myThid ) |
140 |
! if (do_init) call error_mesg & |
141 |
! ('vert_turb_driver in vert_turb_driver_mod', & |
142 |
! 'initialization has not been called', FATAL) |
143 |
|
144 |
nlev = size(p_full,3) |
145 |
ie = is + size(p_full,1) - 1 |
146 |
je = js + size(p_full,2) - 1 |
147 |
|
148 |
!----------------------------------------------------------------------- |
149 |
!---- set up state variable used by this module ---- |
150 |
|
151 |
! if (use_tau) then |
152 |
! !-- variables at time tau |
153 |
! uu = u |
154 |
! vv = v |
155 |
! tt = t |
156 |
! qq = q |
157 |
! else |
158 |
! !-- variables at time tau+1 |
159 |
! uu = um + dt*udt |
160 |
! vv = vm + dt*vdt |
161 |
! tt = tm + dt*tdt |
162 |
! qq = qm + dt*qdt |
163 |
! endif |
164 |
|
165 |
!----------------------------------------------------------------------- |
166 |
|
167 |
!--------------------------- |
168 |
if (do_mellor_yamada) then |
169 |
!--------------------------- |
170 |
|
171 |
! ----- time step for prognostic tke calculation ----- |
172 |
! call get_time (Time_next-Time, sec, day) |
173 |
! dt_tke = real(sec+day*86400) |
174 |
dt_tke = dt |
175 |
|
176 |
! ----- virtual temp ---------- |
177 |
ape(:,:,:)=(p_full(:,:,:)*p00inv)**(-kappa) |
178 |
thv(:,:,:)=tt(:,:,:)*(qq(:,:,:)*d608+1.0)*ape(:,:,:) |
179 |
if (present(mask)) where (mask < 0.5) thv = 200. |
180 |
|
181 |
! --------------------- update tke----------------------------------- |
182 |
! ---- compute surface tke -------- |
183 |
! ---- compute tke, master length scale (el0), ------------- |
184 |
! ---- length scale (el), and vert mix coeffs (diff_t,diff_m) ---- |
185 |
|
186 |
call tke_surf (u_star, tke(is:ie,js:je,:), kbot=kbot) |
187 |
|
188 |
if ( id_z_pbl > 0 ) then |
189 |
!------ compute pbl depth from k_profile if diagnostic needed ----- |
190 |
call my25_turb (dt_tke, frac_land, p_half, p_full, thv, uu, vv, & |
191 |
z_half, z_full, rough, tke(is:ie,js:je,:), & |
192 |
el0, el, diff_m, diff_t, & |
193 |
mask=mask, kbot=kbot, & |
194 |
ustar=u_star,bstar=b_star,h=z_pbl) |
195 |
else |
196 |
call my25_turb (dt_tke, frac_land, p_half, p_full, thv, uu, vv, & |
197 |
z_half, z_full, rough, tke(is:ie,js:je,:), & |
198 |
el0, el, diff_m, diff_t, & |
199 |
mask=mask, kbot=kbot) |
200 |
end if |
201 |
|
202 |
!--------------------------- |
203 |
else |
204 |
!-------------------------------------------------------------------- |
205 |
!----------- compute molecular diffusion, if desired --------------- |
206 |
|
207 |
if (do_molecular_diffusion) then |
208 |
call molecular_diff (tt, p_half, diff_m, diff_t, myThid ) |
209 |
else |
210 |
diff_m = 0.0 |
211 |
diff_t = 0.0 |
212 |
endif |
213 |
|
214 |
!--------------------------- |
215 |
!------------------- non-local K scheme -------------- |
216 |
|
217 |
call diffusivity ( tt, qq, uu, vv, p_full, p_half, z_full, z_half, & |
218 |
u_star, b_star, z_pbl, diff_m, diff_t, myThid, & |
219 |
kbot = kbot) |
220 |
|
221 |
!--------------------------- |
222 |
endif |
223 |
!----------------------------------------------------------------------- |
224 |
!------------------ shallow convection ???? ---------------------------- |
225 |
|
226 |
if (do_shallow_conv) then |
227 |
call shallow_conv (tt, qq, p_full, p_half, diff_sc, kbot) |
228 |
diff_t = diff_t + diff_sc |
229 |
endif |
230 |
|
231 |
!----------------------------------------------------------------------- |
232 |
!------------- define gustiness ------------ |
233 |
|
234 |
if ( trim(gust_scheme) == 'constant' ) then |
235 |
gust = constant_gust |
236 |
else if ( trim(gust_scheme) == 'beljaars' ) then |
237 |
! --- from Beljaars (1994) and Beljaars and Viterbo (1999) --- |
238 |
where (b_star > 0.) |
239 |
gust = (u_star*b_star*gust_zi)**(1./3.) |
240 |
elsewhere |
241 |
gust = 0. |
242 |
endwhere |
243 |
endif |
244 |
|
245 |
!----------------------------------------------------------------------- |
246 |
!------------------------ diagnostics section -------------------------- |
247 |
|
248 |
if (do_mellor_yamada) then |
249 |
|
250 |
! --- set up local mask for fields with surface data --- |
251 |
if ( present(mask) ) then |
252 |
! lmask(:,:,1) = .true. |
253 |
! lmask(:,:,2:nlev+1) = mask(:,:,1:nlev) > 0.5 |
254 |
else |
255 |
! lmask = .true. |
256 |
endif |
257 |
|
258 |
!------- tke -------------------------------- |
259 |
if ( id_tke > 0 ) then |
260 |
! used = send_data ( id_tke, tke(is:ie,js:je,:), Time_next, & |
261 |
! is, js, 1, mask=lmask ) |
262 |
endif |
263 |
|
264 |
!------- length scale (at half levels) ------ |
265 |
if ( id_lscale > 0 ) then |
266 |
! used = send_data ( id_lscale, el, Time_next, is, js, 1, & |
267 |
! mask=lmask ) |
268 |
endif |
269 |
|
270 |
!------- master length scale ------- |
271 |
if ( id_lscale_0 > 0 ) then |
272 |
! used = send_data ( id_lscale_0, el0, Time_next, is, js ) |
273 |
endif |
274 |
|
275 |
end if |
276 |
|
277 |
!------- boundary layer depth ------- |
278 |
if ( id_z_pbl > 0 ) then |
279 |
! used = send_data ( id_z_pbl, z_pbl, Time_next, is, js ) |
280 |
endif |
281 |
|
282 |
!------- gustiness ------- |
283 |
if ( id_gust > 0 ) then |
284 |
! used = send_data ( id_gust, gust, Time_next, is, js ) |
285 |
endif |
286 |
|
287 |
!------- output diffusion coefficients --------- |
288 |
|
289 |
if ( id_diff_t > 0 .or. id_diff_m > 0 .or. id_diff_sc > 0 ) then |
290 |
! --- set up local mask for fields without surface data --- |
291 |
if (present(mask)) then |
292 |
! lmask(:,:,1:nlev) = mask(:,:,1:nlev) > 0.5 |
293 |
! lmask(:,:,nlev+1) = .false. |
294 |
else |
295 |
! lmask(:,:,1:nlev) = .true. |
296 |
! lmask(:,:,nlev+1) = .false. |
297 |
endif |
298 |
! -- dummy data at surface -- |
299 |
! diag3(:,:,nlev+1)=0.0 |
300 |
endif |
301 |
|
302 |
!------- diffusion coefficient for heat/moisture ------- |
303 |
if ( id_diff_t > 0 ) then |
304 |
! diag3(:,:,1:nlev) = diff_t(:,:,1:nlev) |
305 |
! used = send_data ( id_diff_t, diag3, Time_next, is, js, 1, mask=lmask ) |
306 |
endif |
307 |
|
308 |
!------- diffusion coefficient for momentum ------- |
309 |
if ( id_diff_m > 0 ) then |
310 |
! diag3(:,:,1:nlev) = diff_m(:,:,1:nlev) |
311 |
! used = send_data ( id_diff_m, diag3, Time_next, is, js, 1, mask=lmask ) |
312 |
endif |
313 |
|
314 |
!------- diffusion coefficient for shallow conv ------- |
315 |
if (do_shallow_conv) then |
316 |
if ( id_diff_sc > 0 ) then |
317 |
! diag3(:,:,1:nlev) = diff_sc(:,:,1:nlev) |
318 |
! used = send_data ( id_diff_sc, diag3, Time_next, is, js, 1, mask=lmask) |
319 |
endif |
320 |
endif |
321 |
|
322 |
!--- geopotential height relative to the surface on full and half levels ---- |
323 |
|
324 |
if ( id_z_half > 0 ) then |
325 |
!--- set up local mask for fields with surface data --- |
326 |
if ( present(mask) ) then |
327 |
! lmask(:,:,1) = .true. |
328 |
! lmask(:,:,2:nlev+1) = mask(:,:,1:nlev) > 0.5 |
329 |
else |
330 |
! lmask = .true. |
331 |
endif |
332 |
! used = send_data ( id_z_half, z_half, Time_next, is, js, 1, mask=lmask ) |
333 |
endif |
334 |
|
335 |
if ( id_z_full > 0 ) then |
336 |
! used = send_data ( id_z_full, z_full, Time_next, is, js, 1, rmask=mask) |
337 |
endif |
338 |
|
339 |
!--- zonal and meridional wind on mass grid ------- |
340 |
|
341 |
if ( id_uwnd > 0 ) then |
342 |
! used = send_data ( id_uwnd, uu, Time_next, is, js, 1, rmask=mask) |
343 |
endif |
344 |
|
345 |
if ( id_vwnd > 0 ) then |
346 |
! used = send_data ( id_vwnd, vv, Time_next, is, js, 1, rmask=mask) |
347 |
endif |
348 |
|
349 |
!----------------------------------------------------------------------- |
350 |
|
351 |
end subroutine vert_turb_driver |
352 |
|
353 |
!####################################################################### |
354 |
|
355 |
subroutine vert_turb_driver_init (id, jd, kd, axes, Time, myThid) |
356 |
|
357 |
!----------------------------------------------------------------------- |
358 |
integer, intent(in) :: id, jd, kd, axes(4) |
359 |
! type(time_type), intent(in) :: Time |
360 |
real, intent(in) :: Time |
361 |
integer, intent(in) :: myThid |
362 |
!----------------------------------------------------------------------- |
363 |
integer, dimension(3) :: full = (/1,2,3/), half = (/1,2,4/) |
364 |
integer :: ierr, unit, io |
365 |
|
366 |
integer :: iUnit |
367 |
CHARACTER*(gcm_LEN_MBUF) :: msgBuf |
368 |
|
369 |
if (.not.do_init) CALL PRINT_ERROR( & |
370 |
'vert_turb_driver_init in vert_turb_driver_mod'// & |
371 |
'attempting to call initialization twice', myThid ) |
372 |
! if (.not.do_init) & |
373 |
! call error_mesg & |
374 |
! ('vert_turb_driver_init in vert_turb_driver_mod', & |
375 |
! 'attempting to call initialization twice', FATAL) |
376 |
|
377 |
!----------------------------------------------------------------------- |
378 |
!--------------- read namelist ------------------ |
379 |
|
380 |
! _BARRIER |
381 |
! _BEGIN_MASTER(myThid) |
382 |
CALL BARRIER(myThid) |
383 |
IF ( myThid.EQ.1 ) THEN |
384 |
|
385 |
WRITE(msgBuf,'(A)') 'VERT_TURB_DRIVER_INIT: opening data.atm_gray' |
386 |
CALL PRINT_MESSAGE( msgBuf, gcm_stdMsgUnit, gcm_SQZ_R, myThid ) |
387 |
CALL OPEN_COPY_DATA_FILE( & |
388 |
'data.atm_gray', 'VERT_TURB_DRIVER_INIT', & |
389 |
iUnit, & |
390 |
myThid ) |
391 |
! Read parameters from open data file |
392 |
READ(UNIT=iUnit,NML=vert_turb_driver_nml) |
393 |
WRITE(msgBuf,'(A)') & |
394 |
'VERT_TURB_DRIVER_INIT: finished reading data.atm_gray' |
395 |
CALL PRINT_MESSAGE( msgBuf, gcm_stdMsgUnit, gcm_SQZ_R, myThid ) |
396 |
! Close the open data file |
397 |
#ifdef SINGLE_DISK_IO |
398 |
CLOSE(iUnit) |
399 |
#else |
400 |
CLOSE(iUnit,STATUS='DELETE') |
401 |
#endif /* SINGLE_DISK_IO */ |
402 |
|
403 |
! if (file_exist('input.nml')) then |
404 |
! unit = open_namelist_file () |
405 |
! ierr=1; do while (ierr /= 0) |
406 |
! read (unit, nml=vert_turb_driver_nml, iostat=io, end=10) |
407 |
! ierr = check_nml_error (io, 'vert_turb_driver_nml') |
408 |
! enddo |
409 |
! 10 call close_file (unit) |
410 |
! endif |
411 |
|
412 |
!---------- output namelist -------------------------------------------- |
413 |
|
414 |
! call write_version_number(version, tag) |
415 |
! if ( mpp_pe() == mpp_root_pe() ) write (stdlog(), nml=vert_turb_driver_nml) |
416 |
|
417 |
! --- check namelist option --- |
418 |
if ( trim(gust_scheme) /= 'constant' .and. & |
419 |
trim(gust_scheme) /= 'beljaars' ) CALL PRINT_ERROR( & |
420 |
'vert_turb_driver_mod'//'invalid value for namelist'// & |
421 |
' variable GUST_SCHEME', myThid ) |
422 |
! trim(gust_scheme) /= 'beljaars' ) call error_mesg & |
423 |
! ('vert_turb_driver_mod', 'invalid value for namelist & |
424 |
! &variable GUST_SCHEME', FATAL) |
425 |
|
426 |
if (do_molecular_diffusion .and. do_mellor_yamada) CALL PRINT_ERROR( & |
427 |
'vert_turb_driver_mod'//' cannot activate'// & |
428 |
' molecular diffusion with mellor_yamada', myThid ) |
429 |
! if (do_molecular_diffusion .and. do_mellor_yamada) & |
430 |
! call error_mesg ( 'vert_turb_driver_mod', 'cannot activate & |
431 |
! &molecular diffusion with mellor_yamada', FATAL) |
432 |
|
433 |
!----------------------------------------------------------------------- |
434 |
|
435 |
if (do_mellor_yamada) call my25_turb_init (id, jd, kd) |
436 |
|
437 |
if (do_shallow_conv) call shallow_conv_init (kd) |
438 |
|
439 |
!----------------------------------------------------------------------- |
440 |
!----- initialize diagnostic fields ----- |
441 |
|
442 |
id_tke = 0 |
443 |
id_lscale = 0 |
444 |
id_lscale_0 = 0 |
445 |
id_z_pbl = 0 |
446 |
id_gust = 0 |
447 |
id_diff_t = 0 |
448 |
id_diff_m = 0 |
449 |
id_diff_sc = 0 |
450 |
id_z_full = 0 |
451 |
id_z_half = 0 |
452 |
id_uwnd = 0 |
453 |
id_vwnd = 0 |
454 |
! id_uwnd = register_diag_field ( mod_name, 'uwnd', axes(full), Time, & |
455 |
! 'zonal wind on mass grid', 'meters/second' , & |
456 |
!missing_value=missing_value ) |
457 |
|
458 |
! id_vwnd = register_diag_field ( mod_name, 'vwnd', axes(full), Time, & |
459 |
! 'meridional wind on mass grid', 'meters/second' , & |
460 |
!missing_value=missing_value ) |
461 |
|
462 |
! id_z_full = & |
463 |
! register_diag_field ( mod_name, 'z_full', axes(full), Time, & |
464 |
! 'geopotential height relative to surface at full levels', & |
465 |
! 'meters' , missing_value=missing_value ) |
466 |
|
467 |
! id_z_half = & |
468 |
! register_diag_field ( mod_name, 'z_half', axes(half), Time, & |
469 |
! 'geopotential height relative to surface at half levels', & |
470 |
! 'meters' , missing_value=missing_value ) |
471 |
|
472 |
if (do_mellor_yamada) then |
473 |
|
474 |
! id_tke = & |
475 |
! register_diag_field ( mod_name, 'tke', axes(half), Time, & |
476 |
! 'turbulent kinetic energy', 'm2/s2' , & |
477 |
! missing_value=missing_value ) |
478 |
|
479 |
! id_lscale = & |
480 |
! register_diag_field ( mod_name, 'lscale', axes(half), Time, & |
481 |
! 'turbulent length scale', 'm' , & |
482 |
! missing_value=missing_value ) |
483 |
|
484 |
! id_lscale_0 = & |
485 |
! register_diag_field ( mod_name, 'lscale_0', axes(1:2), Time, & |
486 |
! 'master length scale', 'm' ) |
487 |
endif |
488 |
|
489 |
! id_z_pbl = & |
490 |
! register_diag_field ( mod_name, 'z_pbl', axes(1:2), Time, & |
491 |
! 'depth of planetary boundary layer', 'm' ) |
492 |
|
493 |
! id_gust = & |
494 |
! register_diag_field ( mod_name, 'gust', axes(1:2), Time, & |
495 |
! 'wind gustiness in surface layer', 'm/s' ) |
496 |
|
497 |
! id_diff_t = & |
498 |
! register_diag_field ( mod_name, 'diff_t', axes(half), Time, & |
499 |
! 'vert diff coeff for temp', 'm2/s' , & |
500 |
! missing_value=missing_value ) |
501 |
|
502 |
! id_diff_m = & |
503 |
! register_diag_field ( mod_name, 'diff_m', axes(half), Time, & |
504 |
! 'vert diff coeff for momentum', 'm2/s' , & |
505 |
! missing_value=missing_value ) |
506 |
|
507 |
if (do_shallow_conv) then |
508 |
|
509 |
! id_diff_sc = & |
510 |
! register_diag_field ( mod_name, 'diff_sc', axes(half), Time, & |
511 |
! 'vert diff coeff for shallow conv', 'm2/s' , & |
512 |
! missing_value=missing_value ) |
513 |
endif |
514 |
|
515 |
!----------------------------------------------------------------------- |
516 |
|
517 |
do_init=.false. |
518 |
|
519 |
ENDIF |
520 |
CALL BARRIER(myThid) |
521 |
|
522 |
!----------------------------------------------------------------------- |
523 |
|
524 |
end subroutine vert_turb_driver_init |
525 |
|
526 |
!####################################################################### |
527 |
|
528 |
subroutine vert_turb_driver_end(myThid) |
529 |
|
530 |
integer, intent(in) :: myThid |
531 |
|
532 |
!----------------------------------------------------------------------- |
533 |
if (do_mellor_yamada) call my25_turb_end |
534 |
!----------------------------------------------------------------------- |
535 |
|
536 |
end subroutine vert_turb_driver_end |
537 |
|
538 |
!####################################################################### |
539 |
|
540 |
end module vert_turb_driver_mod |