24 implicit none ;
private
26 #include <MOM_memory.h>
28 public mixedlayer_restrat
29 public mixedlayer_restrat_init
30 public mixedlayer_restrat_register_restarts
39 real :: ml_restrat_coef
43 real :: ml_restrat_coef2
48 logical :: mle_use_pbl_mld
51 real :: mle_mld_decay_time
52 real :: mle_mld_decay_time2
53 real :: mle_density_diff
57 real :: mle_mld_stretch
59 logical :: debug = .false.
63 real,
dimension(:,:),
pointer :: &
64 mld_filtered => null(), &
65 mld_filtered_slow => null()
69 integer :: id_urestrat_time = -1
70 integer :: id_vrestrat_time = -1
71 integer :: id_uhml = -1
72 integer :: id_vhml = -1
73 integer :: id_mld = -1
74 integer :: id_rml = -1
75 integer :: id_udml = -1
76 integer :: id_vdml = -1
77 integer :: id_uml = -1
78 integer :: id_vml = -1
83 character(len=40) :: mdl =
"MOM_mixed_layer_restrat"
90 subroutine mixedlayer_restrat(h, uhtr, vhtr, tv, forces, dt, MLD, VarMix, G, GV, US, CS)
94 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
95 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uhtr
97 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vhtr
101 real,
intent(in) :: dt
102 real,
dimension(:,:),
pointer :: mld
107 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
108 "Module must be initialized before it is used.")
111 call mixedlayer_restrat_bml(h, uhtr, vhtr, tv, forces, dt, g, gv, us, cs)
113 call mixedlayer_restrat_general(h, uhtr, vhtr, tv, forces, dt, mld, varmix, g, gv, us, cs)
116 end subroutine mixedlayer_restrat
119 subroutine mixedlayer_restrat_general(h, uhtr, vhtr, tv, forces, dt, MLD_in, VarMix, G, GV, US, CS)
124 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
125 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uhtr
127 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vhtr
131 real,
intent(in) :: dt
132 real,
dimension(:,:),
pointer :: MLD_in
137 real :: uhml(SZIB_(G),SZJ_(G),SZK_(G))
138 real :: vhml(SZI_(G),SZJB_(G),SZK_(G))
139 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
142 real,
dimension(SZI_(G),SZJ_(G)) :: &
150 real :: rho_ml(SZI_(G))
161 real :: Ihtot,Ihtot_slow
167 real :: uDml(SZIB_(G))
168 real :: vDml(SZI_(G))
169 real :: uDml_slow(SZIB_(G))
170 real :: vDml_slow(SZI_(G))
171 real :: utimescale_diag(SZIB_(G),SZJ_(G))
172 real :: vtimescale_diag(SZI_(G),SZJB_(G))
174 real :: uDml_diag(SZIB_(G),SZJ_(G)), vDml_diag(SZI_(G),SZJB_(G))
175 real,
dimension(SZI_(G)) :: rhoSurf, deltaRhoAtKm1, deltaRhoAtK
176 real,
dimension(SZI_(G)) :: dK, dKm1
177 real,
dimension(SZI_(G)) :: pRef_MLD
179 real,
dimension(SZI_(G)) :: rhoAtK, rho1, d1, pRef_N2
182 real :: hAtVel, zpa, zpb, dh, res_scaling_fac
184 logical :: line_is_empty, keep_going, res_upscale
185 integer,
dimension(2) :: EOSdom
186 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
188 real :: PSI, PSI1, z, BOTTOP, XP, DD
191 psi1(z) = max(0., (1. - (2.*z+1.)**2 ) * (1. + (5./21.)*(2.*z+1.)**2) )
192 bottop(z) = 0.5*(1.-sign(1.,z+0.5))
193 xp(z) = max(0., min(1., (-z-0.5)*2./(1.+2.*cs%MLE_tail_dh) ) )
194 dd(z) = (1.-3.*(xp(z)**2)+2.*(xp(z)**3))**(1.+2.*cs%MLE_tail_dh)
195 psi(z) = max( psi1(z), dd(z)*bottop(z) )
197 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
198 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
200 if (.not.
associated(tv%eqn_of_state))
call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
201 "An equation of state must be used with this module.")
202 if (.not.
associated(varmix) .and. cs%front_length>0.)
call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
203 "The resolution argument, Rd/dx, was not associated.")
205 if (cs%MLE_density_diff > 0.)
then
208 eosdom(:) = eos_domain(g%HI, halo=1)
210 dk(:) = 0.5 * h(:,j,1)
211 call calculate_density(tv%T(:,j,1), tv%S(:,j,1), pref_mld, rhosurf, tv%eqn_of_state, eosdom)
216 dk(:) = dk(:) + 0.5 * ( h(:,j,k) + h(:,j,k-1) )
218 deltarhoatkm1(:) = deltarhoatk(:)
219 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pref_mld, deltarhoatk, tv%eqn_of_state, eosdom)
221 deltarhoatk(i) = deltarhoatk(i) - rhosurf(i)
224 ddrho = deltarhoatk(i) - deltarhoatkm1(i)
225 if ((mld_fast(i,j)==0.) .and. (ddrho>0.) .and. &
226 (deltarhoatkm1(i)<cs%MLE_density_diff) .and. (deltarhoatk(i)>=cs%MLE_density_diff))
then
227 afac = ( cs%MLE_density_diff - deltarhoatkm1(i) ) / ddrho
228 mld_fast(i,j) = dk(i) * afac + dkm1(i) * (1. - afac)
233 mld_fast(i,j) = cs%MLE_MLD_stretch * mld_fast(i,j)
234 if ((mld_fast(i,j)==0.) .and. (deltarhoatk(i)<cs%MLE_density_diff)) &
235 mld_fast(i,j) = dk(i)
238 elseif (cs%MLE_use_PBL_MLD)
then
239 if (.not.
associated(mld_in))
call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
240 "Argument MLD_in was not associated!")
241 do j = js-1, je+1 ;
do i = is-1, ie+1
242 mld_fast(i,j) = (cs%MLE_MLD_stretch * gv%Z_to_H) * mld_in(i,j)
245 call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
246 "No MLD to use for MLE parameterization.")
250 if (cs%MLE_MLD_decay_time>0.)
then
252 call hchksum(cs%MLD_filtered,
'mixed_layer_restrat: MLD_filtered', g%HI, haloshift=1, scale=gv%H_to_m)
253 call hchksum(mld_in,
'mixed_layer_restrat: MLD in', g%HI, haloshift=1, scale=us%Z_to_m)
255 afac = cs%MLE_MLD_decay_time / ( dt + cs%MLE_MLD_decay_time )
256 bfac = dt / ( dt + cs%MLE_MLD_decay_time )
257 do j = js-1, je+1 ;
do i = is-1, ie+1
261 cs%MLD_filtered(i,j) = max( mld_fast(i,j), bfac*mld_fast(i,j) + afac*cs%MLD_filtered(i,j) )
262 mld_fast(i,j) = cs%MLD_filtered(i,j)
267 if (cs%MLE_MLD_decay_time2>0.)
then
269 call hchksum(cs%MLD_filtered_slow,
'mixed_layer_restrat: MLD_filtered_slow',g%HI,haloshift=1,scale=gv%H_to_m)
270 call hchksum(mld_fast,
'mixed_layer_restrat: MLD fast',g%HI,haloshift=1,scale=gv%H_to_m)
272 afac = cs%MLE_MLD_decay_time2 / ( dt + cs%MLE_MLD_decay_time2 )
273 bfac = dt / ( dt + cs%MLE_MLD_decay_time2 )
274 do j = js-1, je+1 ;
do i = is-1, ie+1
278 cs%MLD_filtered_slow(i,j) = max( mld_fast(i,j), bfac*mld_fast(i,j) + afac*cs%MLD_filtered_slow(i,j) )
279 mld_slow(i,j) = cs%MLD_filtered_slow(i,j)
282 do j = js-1, je+1 ;
do i = is-1, ie+1
283 mld_slow(i,j) = mld_fast(i,j)
287 udml(:) = 0.0 ; vdml(:) = 0.0
288 udml_slow(:) = 0.0 ; vdml_slow(:) = 0.0
290 g_rho0 = gv%g_Earth / gv%Rho0
291 h_neglect = gv%H_subroundoff
292 dz_neglect = gv%H_subroundoff*gv%H_to_Z
293 if (cs%front_length>0.)
then
295 i_lfront = 1. / cs%front_length
297 res_upscale = .false.
301 eosdom(:) = eos_domain(g%HI, halo=1)
314 htot_fast(i,j) = 0.0 ; rml_av_fast(i,j) = 0.0
315 htot_slow(i,j) = 0.0 ; rml_av_slow(i,j) = 0.0
320 h_avail(i,j,k) = max(i4dt*g%areaT(i,j)*(h(i,j,k)-gv%Angstrom_H),0.0)
323 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p0, rho_ml(:), tv%eqn_of_state, eosdom)
324 line_is_empty = .true.
326 if (htot_fast(i,j) < mld_fast(i,j))
then
327 dh = min( h(i,j,k), mld_fast(i,j)-htot_fast(i,j) )
328 rml_av_fast(i,j) = rml_av_fast(i,j) + dh*rho_ml(i)
329 htot_fast(i,j) = htot_fast(i,j) + dh
330 line_is_empty = .false.
332 if (htot_slow(i,j) < mld_slow(i,j))
then
333 dh = min( h(i,j,k), mld_slow(i,j)-htot_slow(i,j) )
334 rml_av_slow(i,j) = rml_av_slow(i,j) + dh*rho_ml(i)
335 htot_slow(i,j) = htot_slow(i,j) + dh
336 line_is_empty = .false.
339 if (line_is_empty) keep_going=.false.
344 rml_av_fast(i,j) = -(g_rho0*rml_av_fast(i,j)) / (htot_fast(i,j) + h_neglect)
345 rml_av_slow(i,j) = -(g_rho0*rml_av_slow(i,j)) / (htot_slow(i,j) + h_neglect)
350 call hchksum(h,
'mixed_layer_restrat: h', g%HI, haloshift=1, scale=gv%H_to_m)
351 call hchksum(forces%ustar,
'mixed_layer_restrat: u*', g%HI, haloshift=1, scale=us%Z_to_m*us%s_to_T)
352 call hchksum(mld_fast,
'mixed_layer_restrat: MLD', g%HI, haloshift=1, scale=gv%H_to_m)
353 call hchksum(rml_av_fast,
'mixed_layer_restrat: rml', g%HI, haloshift=1, &
354 scale=us%m_to_Z*us%L_T_to_m_s**2)
363 do j=js,je ;
do i=is-1,ie
364 u_star = 0.5*(forces%ustar(i,j) + forces%ustar(i+1,j))
365 absf = 0.5*(abs(g%CoriolisBu(i,j-1)) + abs(g%CoriolisBu(i,j)))
367 if (res_upscale) res_scaling_fac = &
368 ( sqrt( 0.5 * ( g%dxCu(i,j)**2 + g%dyCu(i,j)**2 ) ) * i_lfront ) &
369 * min( 1., 0.5*( varmix%Rd_dx_h(i,j) + varmix%Rd_dx_h(i+1,j) ) )
374 h_vel = 0.5*((htot_fast(i,j) + htot_fast(i+1,j)) + h_neglect) * gv%H_to_Z
375 mom_mixrate = (0.41*9.8696)*u_star**2 / &
376 (absf*h_vel**2 + 4.0*(h_vel+dz_neglect)*u_star)
377 timescale = 0.0625 * (absf + 2.0*mom_mixrate) / (absf**2 + mom_mixrate**2)
378 timescale = timescale * cs%ml_restrat_coef
379 if (res_upscale) timescale = timescale * res_scaling_fac
380 udml(i) = timescale * g%mask2dCu(i,j)*g%dyCu(i,j)*g%IdxCu(i,j) * &
381 (rml_av_fast(i+1,j)-rml_av_fast(i,j)) * (h_vel**2 * gv%Z_to_H)
383 h_vel = 0.5*((htot_slow(i,j) + htot_slow(i+1,j)) + h_neglect) * gv%H_to_Z
384 mom_mixrate = (0.41*9.8696)*u_star**2 / &
385 (absf*h_vel**2 + 4.0*(h_vel+dz_neglect)*u_star)
386 timescale = 0.0625 * (absf + 2.0*mom_mixrate) / (absf**2 + mom_mixrate**2)
387 timescale = timescale * cs%ml_restrat_coef2
388 if (res_upscale) timescale = timescale * res_scaling_fac
389 udml_slow(i) = timescale * g%mask2dCu(i,j)*g%dyCu(i,j)*g%IdxCu(i,j) * &
390 (rml_av_slow(i+1,j)-rml_av_slow(i,j)) * (h_vel**2 * gv%Z_to_H)
392 if (udml(i) + udml_slow(i) == 0.)
then
393 do k=1,nz ; uhml(i,j,k) = 0.0 ;
enddo
395 ihtot = 2.0 / ((htot_fast(i,j) + htot_fast(i+1,j)) + h_neglect)
396 ihtot_slow = 2.0 / ((htot_slow(i,j) + htot_slow(i+1,j)) + h_neglect)
397 zpa = 0.0 ; zpb = 0.0
401 hatvel = 0.5*(h(i,j,k) + h(i+1,j,k))
403 zpa = zpa - (hatvel * ihtot)
404 a(k) = a(k) - psi(zpa)
406 if (a(k)*udml(i) > 0.0)
then
407 if (a(k)*udml(i) > h_avail(i,j,k)) udml(i) = h_avail(i,j,k) / a(k)
408 elseif (a(k)*udml(i) < 0.0)
then
409 if (-a(k)*udml(i) > h_avail(i+1,j,k)) udml(i) = -h_avail(i+1,j,k) / a(k)
414 hatvel = 0.5*(h(i,j,k) + h(i+1,j,k))
416 zpb = zpb - (hatvel * ihtot_slow)
417 b(k) = b(k) - psi(zpb)
419 if (b(k)*udml_slow(i) > 0.0)
then
420 if (b(k)*udml_slow(i) > h_avail(i,j,k) - a(k)*udml(i)) &
421 udml_slow(i) = max( 0., h_avail(i,j,k) - a(k)*udml(i) ) / b(k)
422 elseif (b(k)*udml_slow(i) < 0.0)
then
423 if (-b(k)*udml_slow(i) > h_avail(i+1,j,k) + a(k)*udml(i)) &
424 udml_slow(i) = -max( 0., h_avail(i+1,j,k) + a(k)*udml(i) ) / b(k)
428 uhml(i,j,k) = a(k)*udml(i) + b(k)*udml_slow(i)
429 uhtr(i,j,k) = uhtr(i,j,k) + uhml(i,j,k)*dt
433 utimescale_diag(i,j) = timescale
434 udml_diag(i,j) = udml(i)
439 do j=js-1,je ;
do i=is,ie
440 u_star = 0.5*(forces%ustar(i,j) + forces%ustar(i,j+1))
441 absf = 0.5*(abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j)))
443 if (res_upscale) res_scaling_fac = &
444 ( sqrt( 0.5 * ( (g%dxCv(i,j))**2 + (g%dyCv(i,j))**2 ) ) * i_lfront ) &
445 * min( 1., 0.5*( varmix%Rd_dx_h(i,j) + varmix%Rd_dx_h(i,j+1) ) )
450 h_vel = 0.5*((htot_fast(i,j) + htot_fast(i,j+1)) + h_neglect) * gv%H_to_Z
451 mom_mixrate = (0.41*9.8696)*u_star**2 / &
452 (absf*h_vel**2 + 4.0*(h_vel+dz_neglect)*u_star)
453 timescale = 0.0625 * (absf + 2.0*mom_mixrate) / (absf**2 + mom_mixrate**2)
454 timescale = timescale * cs%ml_restrat_coef
455 if (res_upscale) timescale = timescale * res_scaling_fac
456 vdml(i) = timescale * g%mask2dCv(i,j)*g%dxCv(i,j)*g%IdyCv(i,j) * &
457 (rml_av_fast(i,j+1)-rml_av_fast(i,j)) * (h_vel**2 * gv%Z_to_H)
459 h_vel = 0.5*((htot_slow(i,j) + htot_slow(i,j+1)) + h_neglect) * gv%H_to_Z
460 mom_mixrate = (0.41*9.8696)*u_star**2 / &
461 (absf*h_vel**2 + 4.0*(h_vel+dz_neglect)*u_star)
462 timescale = 0.0625 * (absf + 2.0*mom_mixrate) / (absf**2 + mom_mixrate**2)
463 timescale = timescale * cs%ml_restrat_coef2
464 if (res_upscale) timescale = timescale * res_scaling_fac
465 vdml_slow(i) = timescale * g%mask2dCv(i,j)*g%dxCv(i,j)*g%IdyCv(i,j) * &
466 (rml_av_slow(i,j+1)-rml_av_slow(i,j)) * (h_vel**2 * gv%Z_to_H)
468 if (vdml(i) + vdml_slow(i) == 0.)
then
469 do k=1,nz ; vhml(i,j,k) = 0.0 ;
enddo
471 ihtot = 2.0 / ((htot_fast(i,j) + htot_fast(i,j+1)) + h_neglect)
472 ihtot_slow = 2.0 / ((htot_slow(i,j) + htot_slow(i,j+1)) + h_neglect)
473 zpa = 0.0 ; zpb = 0.0
477 hatvel = 0.5*(h(i,j,k) + h(i,j+1,k))
479 zpa = zpa - (hatvel * ihtot)
480 a(k) = a(k) - psi( zpa )
482 if (a(k)*vdml(i) > 0.0)
then
483 if (a(k)*vdml(i) > h_avail(i,j,k)) vdml(i) = h_avail(i,j,k) / a(k)
484 elseif (a(k)*vdml(i) < 0.0)
then
485 if (-a(k)*vdml(i) > h_avail(i,j+1,k)) vdml(i) = -h_avail(i,j+1,k) / a(k)
490 hatvel = 0.5*(h(i,j,k) + h(i,j+1,k))
492 zpb = zpb - (hatvel * ihtot_slow)
493 b(k) = b(k) - psi(zpb)
495 if (b(k)*vdml_slow(i) > 0.0)
then
496 if (b(k)*vdml_slow(i) > h_avail(i,j,k) - a(k)*vdml(i)) &
497 vdml_slow(i) = max( 0., h_avail(i,j,k) - a(k)*vdml(i) ) / b(k)
498 elseif (b(k)*vdml_slow(i) < 0.0)
then
499 if (-b(k)*vdml_slow(i) > h_avail(i,j+1,k) + a(k)*vdml(i)) &
500 vdml_slow(i) = -max( 0., h_avail(i,j+1,k) + a(k)*vdml(i) ) / b(k)
504 vhml(i,j,k) = a(k)*vdml(i) + b(k)*vdml_slow(i)
505 vhtr(i,j,k) = vhtr(i,j,k) + vhml(i,j,k)*dt
509 vtimescale_diag(i,j) = timescale
510 vdml_diag(i,j) = vdml(i)
514 do j=js,je ;
do k=1,nz ;
do i=is,ie
515 h(i,j,k) = h(i,j,k) - dt*g%IareaT(i,j) * &
516 ((uhml(i,j,k) - uhml(i-1,j,k)) + (vhml(i,j,k) - vhml(i,j-1,k)))
517 enddo ;
enddo ;
enddo
522 if (cs%id_uhml > 0 .or. cs%id_vhml > 0) &
524 call pass_var(h, g%domain, to_west+to_south+omit_corners, halo=1)
525 call diag_update_remap_grids(cs%diag)
528 if (query_averaging_enabled(cs%diag))
then
529 if (cs%id_urestrat_time > 0)
call post_data(cs%id_urestrat_time, utimescale_diag, cs%diag)
530 if (cs%id_vrestrat_time > 0)
call post_data(cs%id_vrestrat_time, vtimescale_diag, cs%diag)
531 if (cs%id_uhml > 0)
call post_data(cs%id_uhml, uhml, cs%diag)
532 if (cs%id_vhml > 0)
call post_data(cs%id_vhml, vhml, cs%diag)
533 if (cs%id_MLD > 0)
call post_data(cs%id_MLD, mld_fast, cs%diag)
534 if (cs%id_Rml > 0)
call post_data(cs%id_Rml, rml_av_fast, cs%diag)
535 if (cs%id_uDml > 0)
call post_data(cs%id_uDml, udml_diag, cs%diag)
536 if (cs%id_vDml > 0)
call post_data(cs%id_vDml, vdml_diag, cs%diag)
538 if (cs%id_uml > 0)
then
539 do j=js,je ;
do i=is-1,ie
540 h_vel = 0.5*((htot_fast(i,j) + htot_fast(i+1,j)) + h_neglect)
541 udml_diag(i,j) = udml_diag(i,j) / (0.01*h_vel) * g%IdyCu(i,j) * (psi(0.)-psi(-.01))
543 call post_data(cs%id_uml, udml_diag, cs%diag)
545 if (cs%id_vml > 0)
then
546 do j=js-1,je ;
do i=is,ie
547 h_vel = 0.5*((htot_fast(i,j) + htot_fast(i,j+1)) + h_neglect)
548 vdml_diag(i,j) = vdml_diag(i,j) / (0.01*h_vel) * g%IdxCv(i,j) * (psi(0.)-psi(-.01))
550 call post_data(cs%id_vml, vdml_diag, cs%diag)
556 call diag_update_remap_grids(cs%diag)
558 end subroutine mixedlayer_restrat_general
562 subroutine mixedlayer_restrat_bml(h, uhtr, vhtr, tv, forces, dt, G, GV, US, CS)
566 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
567 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(inout) :: uhtr
569 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(inout) :: vhtr
573 real,
intent(in) :: dt
576 real :: uhml(SZIB_(G),SZJ_(G),SZK_(G))
577 real :: vhml(SZI_(G),SZJB_(G),SZK_(G))
578 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
581 real,
dimension(SZI_(G),SZJ_(G)) :: &
585 real :: Rho0(SZI_(G))
603 real :: uDml(SZIB_(G))
604 real :: vDml(SZI_(G))
605 real :: utimescale_diag(SZIB_(G),SZJ_(G))
606 real :: vtimescale_diag(SZI_(G),SZJB_(G))
609 real :: uDml_diag(SZIB_(G),SZJ_(G)), vDml_diag(SZI_(G),SZJB_(G))
612 integer,
dimension(2) :: EOSdom
613 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, nkml
614 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
615 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nkml = gv%nkml
617 if (.not.
associated(cs))
call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
618 "Module must be initialized before it is used.")
619 if ((nkml<2) .or. (cs%ml_restrat_coef<=0.0))
return
621 udml(:) = 0.0 ; vdml(:) = 0.0
623 g_rho0 = gv%g_Earth / gv%Rho0
624 use_eos =
associated(tv%eqn_of_state)
625 h_neglect = gv%H_subroundoff
626 dz_neglect = gv%H_subroundoff*gv%H_to_Z
628 if (.not.use_eos)
call mom_error(fatal,
"MOM_mixedlayer_restrat: "// &
629 "An equation of state must be used with this module.")
634 eosdom(:) = eos_domain(g%HI, halo=1)
645 htot(i,j) = 0.0 ; rml_av(i,j) = 0.0
648 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p0, rho0(:), tv%eqn_of_state, eosdom)
650 rml_av(i,j) = rml_av(i,j) + h(i,j,k)*rho0(i)
651 htot(i,j) = htot(i,j) + h(i,j,k)
652 h_avail(i,j,k) = max(i4dt*g%areaT(i,j)*(h(i,j,k)-gv%Angstrom_H),0.0)
657 rml_av(i,j) = (g_rho0*rml_av(i,j)) / (htot(i,j) + h_neglect)
667 do j=js,je;
do i=is-1,ie
668 h_vel = 0.5*(htot(i,j) + htot(i+1,j)) * gv%H_to_Z
670 u_star = 0.5*(forces%ustar(i,j) + forces%ustar(i+1,j))
671 absf = 0.5*(abs(g%CoriolisBu(i,j-1)) + abs(g%CoriolisBu(i,j)))
675 mom_mixrate = (0.41*9.8696)*u_star**2 / &
676 (absf*h_vel**2 + 4.0*(h_vel+dz_neglect)*u_star)
677 timescale = 0.0625 * (absf + 2.0*mom_mixrate) / (absf**2 + mom_mixrate**2)
679 timescale = timescale * cs%ml_restrat_coef
682 udml(i) = timescale * g%mask2dCu(i,j)*g%dyCu(i,j)*g%IdxCu(i,j) * &
683 (rml_av(i+1,j)-rml_av(i,j)) * (h_vel**2 * gv%Z_to_H)
685 if (udml(i) == 0)
then
686 do k=1,nkml ; uhml(i,j,k) = 0.0 ;
enddo
688 i2htot = 1.0 / (htot(i,j) + htot(i+1,j) + h_neglect)
693 hx2 = (h(i,j,k) + h(i+1,j,k) + h_neglect)
694 a(k) = (hx2 * i2htot) * (2.0 - 4.0*(z_topx2+0.5*hx2)*i2htot)
695 z_topx2 = z_topx2 + hx2
696 if (a(k)*udml(i) > 0.0)
then
697 if (a(k)*udml(i) > h_avail(i,j,k)) udml(i) = h_avail(i,j,k) / a(k)
699 if (-a(k)*udml(i) > h_avail(i+1,j,k)) udml(i) = -h_avail(i+1,j,k)/a(k)
703 uhml(i,j,k) = a(k)*udml(i)
704 uhtr(i,j,k) = uhtr(i,j,k) + uhml(i,j,k)*dt
708 udml_diag(i,j) = udml(i)
709 utimescale_diag(i,j) = timescale
714 do j=js-1,je ;
do i=is,ie
715 h_vel = 0.5*(htot(i,j) + htot(i,j+1)) * gv%H_to_Z
717 u_star = 0.5*(forces%ustar(i,j) + forces%ustar(i,j+1))
718 absf = 0.5*(abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j)))
722 mom_mixrate = (0.41*9.8696)*u_star**2 / &
723 (absf*h_vel**2 + 4.0*(h_vel+dz_neglect)*u_star)
724 timescale = 0.0625 * (absf + 2.0*mom_mixrate) / (absf**2 + mom_mixrate**2)
726 timescale = timescale * cs%ml_restrat_coef
729 vdml(i) = timescale * g%mask2dCv(i,j)*g%dxCv(i,j)*g%IdyCv(i,j) * &
730 (rml_av(i,j+1)-rml_av(i,j)) * (h_vel**2 * gv%Z_to_H)
731 if (vdml(i) == 0)
then
732 do k=1,nkml ; vhml(i,j,k) = 0.0 ;
enddo
734 i2htot = 1.0 / (htot(i,j) + htot(i,j+1) + h_neglect)
739 hx2 = (h(i,j,k) + h(i,j+1,k) + h_neglect)
740 a(k) = (hx2 * i2htot) * (2.0 - 4.0*(z_topx2+0.5*hx2)*i2htot)
741 z_topx2 = z_topx2 + hx2
742 if (a(k)*vdml(i) > 0.0)
then
743 if (a(k)*vdml(i) > h_avail(i,j,k)) vdml(i) = h_avail(i,j,k) / a(k)
745 if (-a(k)*vdml(i) > h_avail(i,j+1,k)) vdml(i) = -h_avail(i,j+1,k)/a(k)
749 vhml(i,j,k) = a(k)*vdml(i)
750 vhtr(i,j,k) = vhtr(i,j,k) + vhml(i,j,k)*dt
754 vtimescale_diag(i,j) = timescale
755 vdml_diag(i,j) = vdml(i)
759 do j=js,je ;
do k=1,nkml ;
do i=is,ie
760 h(i,j,k) = h(i,j,k) - dt*g%IareaT(i,j) * &
761 ((uhml(i,j,k) - uhml(i-1,j,k)) + (vhml(i,j,k) - vhml(i,j-1,k)))
762 enddo ;
enddo ;
enddo
767 if (cs%id_uhml > 0 .or. cs%id_vhml > 0) &
769 call pass_var(h, g%domain, to_west+to_south+omit_corners, halo=1)
770 call diag_update_remap_grids(cs%diag)
773 if (query_averaging_enabled(cs%diag) .and. &
774 ((cs%id_urestrat_time > 0) .or. (cs%id_vrestrat_time > 0)))
then
775 call post_data(cs%id_urestrat_time, utimescale_diag, cs%diag)
776 call post_data(cs%id_vrestrat_time, vtimescale_diag, cs%diag)
778 if (query_averaging_enabled(cs%diag) .and. &
779 ((cs%id_uhml>0) .or. (cs%id_vhml>0)))
then
781 do j=js,je ;
do i=isq,ieq ; uhml(i,j,k) = 0.0 ;
enddo ;
enddo
782 do j=jsq,jeq ;
do i=is,ie ; vhml(i,j,k) = 0.0 ;
enddo ;
enddo
784 if (cs%id_uhml > 0)
call post_data(cs%id_uhml, uhml, cs%diag)
785 if (cs%id_vhml > 0)
call post_data(cs%id_vhml, vhml, cs%diag)
786 if (cs%id_MLD > 0)
call post_data(cs%id_MLD, htot, cs%diag)
787 if (cs%id_Rml > 0)
call post_data(cs%id_Rml, rml_av, cs%diag)
788 if (cs%id_uDml > 0)
call post_data(cs%id_uDml, udml_diag, cs%diag)
789 if (cs%id_vDml > 0)
call post_data(cs%id_vDml, vdml_diag, cs%diag)
792 end subroutine mixedlayer_restrat_bml
796 logical function mixedlayer_restrat_init(Time, G, GV, US, param_file, diag, CS, restart_CS)
797 type(time_type),
intent(in) :: time
802 type(
diag_ctrl),
target,
intent(inout) :: diag
809 real :: flux_to_kg_per_s
811 # include "version_variable.h"
815 call get_param(param_file, mdl,
"MIXEDLAYER_RESTRAT", mixedlayer_restrat_init, &
816 default=.false., do_not_log=.true.)
817 call log_version(param_file, mdl, version,
"", all_default=.not.mixedlayer_restrat_init)
818 call get_param(param_file, mdl,
"MIXEDLAYER_RESTRAT", mixedlayer_restrat_init, &
819 "If true, a density-gradient dependent re-stratifying "//&
820 "flow is imposed in the mixed layer. Can be used in ALE mode "//&
821 "without restriction but in layer mode can only be used if "//&
822 "BULKMIXEDLAYER is true.", default=.false.)
823 if (.not. mixedlayer_restrat_init)
return
825 if (.not.
associated(cs))
then
826 call mom_error(fatal,
"mixedlayer_restrat_init called without an associated control structure.")
830 cs%MLE_MLD_decay_time = -9.e9*us%s_to_T
831 cs%MLE_density_diff = -9.e9*us%kg_m3_to_R
832 cs%MLE_tail_dh = -9.e9
833 cs%MLE_use_PBL_MLD = .false.
834 cs%MLE_MLD_stretch = -9.e9
836 call get_param(param_file, mdl,
"DEBUG", cs%debug, default=.false., do_not_log=.true.)
837 call get_param(param_file, mdl,
"FOX_KEMPER_ML_RESTRAT_COEF", cs%ml_restrat_coef, &
838 "A nondimensional coefficient that is proportional to "//&
839 "the ratio of the deformation radius to the dominant "//&
840 "lengthscale of the submesoscale mixed layer "//&
841 "instabilities, times the minimum of the ratio of the "//&
842 "mesoscale eddy kinetic energy to the large-scale "//&
843 "geostrophic kinetic energy or 1 plus the square of the "//&
844 "grid spacing over the deformation radius, as detailed "//&
845 "by Fox-Kemper et al. (2010)", units=
"nondim", default=0.0)
849 call get_param(param_file, mdl,
"FOX_KEMPER_ML_RESTRAT_COEF2", cs%ml_restrat_coef2, &
850 "As for FOX_KEMPER_ML_RESTRAT_COEF but used in a second application "//&
851 "of the MLE restratification parameterization.", units=
"nondim", default=0.0)
852 call get_param(param_file, mdl,
"MLE_FRONT_LENGTH", cs%front_length, &
853 "If non-zero, is the frontal-length scale used to calculate the "//&
854 "upscaling of buoyancy gradients that is otherwise represented "//&
855 "by the parameter FOX_KEMPER_ML_RESTRAT_COEF. If MLE_FRONT_LENGTH is "//&
856 "non-zero, it is recommended to set FOX_KEMPER_ML_RESTRAT_COEF=1.0.",&
857 units=
"m", default=0.0, scale=us%m_to_L)
858 call get_param(param_file, mdl,
"MLE_USE_PBL_MLD", cs%MLE_use_PBL_MLD, &
859 "If true, the MLE parameterization will use the mixed-layer "//&
860 "depth provided by the active PBL parameterization. If false, "//&
861 "MLE will estimate a MLD based on a density difference with the "//&
862 "surface using the parameter MLE_DENSITY_DIFF.", default=.false.)
863 call get_param(param_file, mdl,
"MLE_MLD_DECAY_TIME", cs%MLE_MLD_decay_time, &
864 "The time-scale for a running-mean filter applied to the mixed-layer "//&
865 "depth used in the MLE restratification parameterization. When "//&
866 "the MLD deepens below the current running-mean the running-mean "//&
867 "is instantaneously set to the current MLD.", units=
"s", default=0., scale=us%s_to_T)
868 call get_param(param_file, mdl,
"MLE_MLD_DECAY_TIME2", cs%MLE_MLD_decay_time2, &
869 "The time-scale for a running-mean filter applied to the filtered "//&
870 "mixed-layer depth used in a second MLE restratification parameterization. "//&
871 "When the MLD deepens below the current running-mean the running-mean "//&
872 "is instantaneously set to the current MLD.", units=
"s", default=0., scale=us%s_to_T)
873 if (.not. cs%MLE_use_PBL_MLD)
then
874 call get_param(param_file, mdl,
"MLE_DENSITY_DIFF", cs%MLE_density_diff, &
875 "Density difference used to detect the mixed-layer "//&
876 "depth used for the mixed-layer eddy parameterization "//&
877 "by Fox-Kemper et al. (2010)", units=
"kg/m3", default=0.03, scale=us%kg_m3_to_R)
879 call get_param(param_file, mdl,
"MLE_TAIL_DH", cs%MLE_tail_dh, &
880 "Fraction by which to extend the mixed-layer restratification "//&
881 "depth used for a smoother stream function at the base of "//&
882 "the mixed-layer.", units=
"nondim", default=0.0)
883 call get_param(param_file, mdl,
"MLE_MLD_STRETCH", cs%MLE_MLD_stretch, &
884 "A scaling coefficient for stretching/shrinking the MLD "//&
885 "used in the MLE scheme. This simply multiplies MLD wherever used.",&
886 units=
"nondim", default=1.0)
891 flux_to_kg_per_s = gv%H_to_kg_m2 * us%L_to_m**2 * us%s_to_T
893 cs%id_uhml = register_diag_field(
'ocean_model',
'uhml', diag%axesCuL, time, &
894 'Zonal Thickness Flux to Restratify Mixed Layer',
'kg s-1', &
895 conversion=flux_to_kg_per_s, y_cell_method=
'sum', v_extensive=.true.)
896 cs%id_vhml = register_diag_field(
'ocean_model',
'vhml', diag%axesCvL, time, &
897 'Meridional Thickness Flux to Restratify Mixed Layer',
'kg s-1', &
898 conversion=flux_to_kg_per_s, x_cell_method=
'sum', v_extensive=.true.)
899 cs%id_urestrat_time = register_diag_field(
'ocean_model',
'MLu_restrat_time', diag%axesCu1, time, &
900 'Mixed Layer Zonal Restratification Timescale',
's', conversion=us%T_to_s)
901 cs%id_vrestrat_time = register_diag_field(
'ocean_model',
'MLv_restrat_time', diag%axesCv1, time, &
902 'Mixed Layer Meridional Restratification Timescale',
's', conversion=us%T_to_s)
903 cs%id_MLD = register_diag_field(
'ocean_model',
'MLD_restrat', diag%axesT1, time, &
904 'Mixed Layer Depth as used in the mixed-layer restratification parameterization',
'm', &
905 conversion=gv%H_to_m)
906 cs%id_Rml = register_diag_field(
'ocean_model',
'ML_buoy_restrat', diag%axesT1, time, &
907 'Mixed Layer Buoyancy as used in the mixed-layer restratification parameterization', &
908 'm s2', conversion=us%m_to_Z*(us%L_to_m**2)*(us%s_to_T**2))
909 cs%id_uDml = register_diag_field(
'ocean_model',
'udml_restrat', diag%axesCu1, time, &
910 'Transport stream function amplitude for zonal restratification of mixed layer', &
911 'm3 s-1', conversion=gv%H_to_m*(us%L_to_m**2)*us%s_to_T)
912 cs%id_vDml = register_diag_field(
'ocean_model',
'vdml_restrat', diag%axesCv1, time, &
913 'Transport stream function amplitude for meridional restratification of mixed layer', &
914 'm3 s-1', conversion=gv%H_to_m*(us%L_to_m**2)*us%s_to_T)
915 cs%id_uml = register_diag_field(
'ocean_model',
'uml_restrat', diag%axesCu1, time, &
916 'Surface zonal velocity component of mixed layer restratification', &
917 'm s-1', conversion=us%L_T_to_m_s)
918 cs%id_vml = register_diag_field(
'ocean_model',
'vml_restrat', diag%axesCv1, time, &
919 'Surface meridional velocity component of mixed layer restratification', &
920 'm s-1', conversion=us%L_T_to_m_s)
923 if (cs%MLE_MLD_decay_time>0. .or. cs%MLE_MLD_decay_time2>0.)
then
925 (gv%m_to_H_restart /= 0.0) .and. (gv%m_to_H_restart /= gv%m_to_H))
then
926 h_rescale = gv%m_to_H / gv%m_to_H_restart
927 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
928 cs%MLD_filtered(i,j) = h_rescale * cs%MLD_filtered(i,j)
932 if (cs%MLE_MLD_decay_time2>0.)
then
933 if (
query_initialized(cs%MLD_filtered_slow,
"MLD_MLE_filtered_slow", restart_cs) .and. &
934 (gv%m_to_H_restart /= 0.0) .and. (gv%m_to_H_restart /= gv%m_to_H))
then
935 h_rescale = gv%m_to_H / gv%m_to_H_restart
936 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
937 cs%MLD_filtered_slow(i,j) = h_rescale * cs%MLD_filtered_slow(i,j)
943 if (
associated(cs%MLD_filtered))
call pass_var(cs%MLD_filtered, g%domain)
945 end function mixedlayer_restrat_init
948 subroutine mixedlayer_restrat_register_restarts(HI, param_file, CS, restart_CS)
956 logical :: mixedlayer_restrat_init
959 call get_param(param_file, mdl,
"MIXEDLAYER_RESTRAT", mixedlayer_restrat_init, &
960 default=.false., do_not_log=.true.)
961 if (.not. mixedlayer_restrat_init)
return
964 if (
associated(cs))
call mom_error(fatal, &
965 "mixedlayer_restrat_register_restarts called with an associated control structure.")
968 call get_param(param_file, mdl,
"MLE_MLD_DECAY_TIME", cs%MLE_MLD_decay_time, &
969 default=0., do_not_log=.true.)
970 call get_param(param_file, mdl,
"MLE_MLD_DECAY_TIME2", cs%MLE_MLD_decay_time2, &
971 default=0., do_not_log=.true.)
972 if (cs%MLE_MLD_decay_time>0. .or. cs%MLE_MLD_decay_time2>0.)
then
974 allocate(cs%MLD_filtered(hi%isd:hi%ied,hi%jsd:hi%jed)) ; cs%MLD_filtered(:,:) = 0.
975 vd = var_desc(
"MLD_MLE_filtered",
"m",
"Time-filtered MLD for use in MLE", &
976 hor_grid=
'h', z_grid=
'1')
979 if (cs%MLE_MLD_decay_time2>0.)
then
981 allocate(cs%MLD_filtered_slow(hi%isd:hi%ied,hi%jsd:hi%jed)) ; cs%MLD_filtered_slow(:,:) = 0.
982 vd = var_desc(
"MLD_MLE_filtered_slow",
"m",
"c Slower time-filtered MLD for use in MLE", &
983 hor_grid=
'h', z_grid=
'1')
987 end subroutine mixedlayer_restrat_register_restarts