/[MITgcm]/MITgcm/pkg/kpp/kpp_routines.F
ViewVC logotype

Contents of /MITgcm/pkg/kpp/kpp_routines.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.23 - (show annotations) (download)
Fri Apr 29 18:47:02 2005 UTC (19 years ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57s_post, checkpoint58b_post, checkpoint57y_post, checkpoint57r_post, checkpoint57i_post, checkpoint58, checkpoint58f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint58y_post, checkpoint58t_post, checkpoint58m_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57h_pre, checkpoint58w_post, checkpoint57h_post, checkpoint57y_pre, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, checkpoint58r_post, checkpoint58n_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, checkpoint58k_post, checkpoint58v_post, checkpoint58l_post, checkpoint57h_done, checkpoint57j_post, checkpoint58g_post, checkpoint58x_post, checkpoint58h_post, checkpoint58j_post, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post, checkpoint58i_post, checkpoint57x_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.22: +57 -2 lines
o added diagnostics to pkg/kpp, including computation of mixed layer
  depth based on a temperature/density criterion
o updated verification/natl_box to test the new pkg/kpp diagnostics

1 C $Header: /u/gcmpack/MITgcm/pkg/kpp/kpp_routines.F,v 1.22 2005/02/10 05:41:44 heimbach Exp $
2 C $Name: $
3
4 #include "KPP_OPTIONS.h"
5
6 C-- File kpp_routines.F: subroutines needed to implement
7 C-- KPP vertical mixing scheme
8 C-- Contents
9 C-- o KPPMIX - Main driver and interface routine.
10 C-- o BLDEPTH - Determine oceanic planetary boundary layer depth.
11 C-- o WSCALE - Compute turbulent velocity scales.
12 C-- o RI_IWMIX - Compute interior viscosity diffusivity coefficients.
13 C-- o Z121 - Apply 121 vertical smoothing.
14 C-- o KPP_SMOOTH_HORIZ - Apply horizontal smoothing to KPP array.
15 C-- o SMOOTH_HORIZ - Apply horizontal smoothing to global array.
16 C-- o BLMIX - Boundary layer mixing coefficients.
17 C-- o ENHANCE - Enhance diffusivity at boundary layer interface.
18 C-- o STATEKPP - Compute buoyancy-related input arrays.
19
20 c***********************************************************************
21
22 SUBROUTINE KPPMIX (
23 I mytime, mythid
24 I , kmtj, shsq, dvsq, ustar
25 I , bo, bosol, dbloc, Ritop, coriol
26 I , ikppkey
27 O , diffus
28 U , ghat
29 O , hbl )
30
31 c-----------------------------------------------------------------------
32 c
33 c Main driver subroutine for kpp vertical mixing scheme and
34 c interface to greater ocean model
35 c
36 c written by: bill large, june 6, 1994
37 c modified by: jan morzel, june 30, 1994
38 c bill large, august 11, 1994
39 c bill large, january 25, 1995 : "dVsq" and 1d code
40 c detlef stammer, august 1997 : for use with MIT GCM Classic
41 c d. menemenlis, june 1998 : for use with MIT GCM UV
42 c
43 c-----------------------------------------------------------------------
44
45 IMPLICIT NONE
46
47 #include "SIZE.h"
48 #include "EEPARAMS.h"
49 #include "PARAMS.h"
50 #include "DYNVARS.h"
51 #include "FFIELDS.h"
52 #include "KPP_PARAMS.h"
53
54 c input
55 c myTime - current time in simulation
56 c myThid - thread number for this instance of the routine
57 c kmtj (imt) - number of vertical layers on this row
58 c shsq (imt,Nr) - (local velocity shear)^2 ((m/s)^2)
59 c dvsq (imt,Nr) - (velocity shear re sfc)^2 ((m/s)^2)
60 c ustar (imt) - surface friction velocity (m/s)
61 c bo (imt) - surface turbulent buoy. forcing (m^2/s^3)
62 c bosol (imt) - radiative buoyancy forcing (m^2/s^3)
63 c dbloc (imt,Nr) - local delta buoyancy across interfaces (m/s^2)
64 c dblocSm(imt,Nr) - horizontally smoothed dbloc (m/s^2)
65 c stored in ghat to save space
66 c Ritop (imt,Nr) - numerator of bulk Richardson Number
67 c (zref-z) * delta buoyancy w.r.t. surface ((m/s)^2)
68 c coriol (imt) - Coriolis parameter (1/s)
69 c note: there is a conversion from 2-D to 1-D for input output variables,
70 c e.g., hbl(sNx,sNy) -> hbl(imt),
71 c where hbl(i,j) -> hbl((j-1)*sNx+i)
72
73 _RL mytime
74 integer mythid
75 integer kmtj (imt )
76 _KPP_RL shsq (imt,Nr)
77 _KPP_RL dvsq (imt,Nr)
78 _KPP_RL ustar (imt )
79 _KPP_RL bo (imt )
80 _KPP_RL bosol (imt )
81 _KPP_RL dbloc (imt,Nr)
82 _KPP_RL Ritop (imt,Nr)
83 _KPP_RL coriol(imt )
84
85 integer ikppkey
86
87 c output
88 c diffus (imt,1) - vertical viscosity coefficient (m^2/s)
89 c diffus (imt,2) - vertical scalar diffusivity (m^2/s)
90 c diffus (imt,3) - vertical temperature diffusivity (m^2/s)
91 c ghat (imt) - nonlocal transport coefficient (s/m^2)
92 c hbl (imt) - mixing layer depth (m)
93
94 _KPP_RL diffus(imt,0:Nrp1,mdiff)
95 _KPP_RL ghat (imt,Nr)
96 _KPP_RL hbl (imt)
97
98 #ifdef ALLOW_KPP
99
100 c local
101 c kbl (imt ) - index of first grid level below hbl
102 c bfsfc (imt ) - surface buoyancy forcing (m^2/s^3)
103 c casea (imt ) - 1 in case A; 0 in case B
104 c stable (imt ) - 1 in stable forcing; 0 if unstable
105 c dkm1 (imt, mdiff) - boundary layer diffusivity at kbl-1 level
106 c blmc (imt,Nr,mdiff) - boundary layer mixing coefficients
107 c sigma (imt ) - normalized depth (d / hbl)
108 c Rib (imt,Nr ) - bulk Richardson number
109
110 integer kbl (imt )
111 _KPP_RL bfsfc (imt )
112 _KPP_RL casea (imt )
113 _KPP_RL stable(imt )
114 _KPP_RL dkm1 (imt, mdiff)
115 _KPP_RL blmc (imt,Nr,mdiff)
116 _KPP_RL sigma (imt )
117 _KPP_RL Rib (imt,Nr )
118
119 integer i, k, md
120
121 c-----------------------------------------------------------------------
122 c compute interior mixing coefficients everywhere, due to constant
123 c internal wave activity, static instability, and local shear
124 c instability.
125 c (ghat is temporary storage for horizontally smoothed dbloc)
126 c-----------------------------------------------------------------------
127
128 cph(
129 cph these storings avoid recomp. of Ri_iwmix
130 CADJ STORE ghat = comlev1_kpp, key = ikppkey
131 CADJ STORE dbloc = comlev1_kpp, key = ikppkey
132 cph)
133 call Ri_iwmix (
134 I kmtj, shsq, dbloc, ghat
135 I , ikppkey
136 O , diffus )
137
138 cph(
139 cph these storings avoid recomp. of Ri_iwmix
140 cph DESPITE TAFs 'not necessary' warning!
141 CADJ STORE dbloc = comlev1_kpp, key = ikppkey
142 CADJ STORE shsq = comlev1_kpp, key = ikppkey
143 CADJ STORE ghat = comlev1_kpp, key = ikppkey
144 CADJ STORE diffus = comlev1_kpp, key = ikppkey
145 cph)
146
147 c-----------------------------------------------------------------------
148 c set seafloor values to zero and fill extra "Nrp1" coefficients
149 c for blmix
150 c-----------------------------------------------------------------------
151
152 do md = 1, mdiff
153 do i = 1,imt
154 do k=kmtj(i),Nrp1
155 diffus(i,k,md) = 0.0
156 end do
157 end do
158 end do
159
160 c-----------------------------------------------------------------------
161 c compute boundary layer mixing coefficients:
162 c
163 c diagnose the new boundary layer depth
164 c-----------------------------------------------------------------------
165
166 call bldepth (
167 I mytime, mythid
168 I , kmtj
169 I , dvsq, dbloc, Ritop, ustar, bo, bosol, coriol
170 I , ikppkey
171 O , hbl, bfsfc, stable, casea, kbl, Rib, sigma
172 & )
173
174 CADJ STORE hbl,bfsfc,stable,casea,kbl = comlev1_kpp, key = ikppkey
175
176 c-----------------------------------------------------------------------
177 c compute boundary layer diffusivities
178 c-----------------------------------------------------------------------
179
180 call blmix (
181 I ustar, bfsfc, hbl, stable, casea, diffus, kbl
182 O , dkm1, blmc, ghat, sigma, ikppkey
183 & )
184 cph(
185 CADJ STORE dkm1,blmc,ghat = comlev1_kpp, key = ikppkey
186 CADJ STORE hbl, kbl, diffus, casea = comlev1_kpp, key = ikppkey
187 cph)
188
189 c-----------------------------------------------------------------------
190 c enhance diffusivity at interface kbl - 1
191 c-----------------------------------------------------------------------
192
193 call enhance (
194 I dkm1, hbl, kbl, diffus, casea
195 U , ghat
196 O , blmc )
197
198 cph(
199 cph avoids recomp. of enhance
200 CADJ STORE blmc = comlev1_kpp, key = ikppkey
201 cph)
202
203 c-----------------------------------------------------------------------
204 c combine interior and boundary layer coefficients and nonlocal term
205 c !!!NOTE!!! In shallow (2-level) regions and for shallow mixed layers
206 c (< 1 level), diffusivity blmc can become negative. The max-s below
207 c are a hack until this problem is properly diagnosed and fixed.
208 c-----------------------------------------------------------------------
209 do k = 1, Nr
210 do i = 1, imt
211 if (k .lt. kbl(i)) then
212 diffus(i,k,1) = max ( blmc(i,k,1), viscAr )
213 diffus(i,k,2) = max ( blmc(i,k,2), diffKrNrS(k) )
214 diffus(i,k,3) = max ( blmc(i,k,3), diffKrNrT(k) )
215 else
216 ghat(i,k) = 0.
217 endif
218 end do
219 end do
220
221 #endif /* ALLOW_KPP */
222
223 return
224 end
225
226 c*************************************************************************
227
228 subroutine bldepth (
229 I mytime, mythid
230 I , kmtj
231 I , dvsq, dbloc, Ritop, ustar, bo, bosol, coriol
232 I , ikppkey
233 O , hbl, bfsfc, stable, casea, kbl, Rib, sigma
234 & )
235
236 c the oceanic planetary boundary layer depth, hbl, is determined as
237 c the shallowest depth where the bulk Richardson number is
238 c equal to the critical value, Ricr.
239 c
240 c bulk Richardson numbers are evaluated by computing velocity and
241 c buoyancy differences between values at zgrid(kl) < 0 and surface
242 c reference values.
243 c in this configuration, the reference values are equal to the
244 c values in the surface layer.
245 c when using a very fine vertical grid, these values should be
246 c computed as the vertical average of velocity and buoyancy from
247 c the surface down to epsilon*zgrid(kl).
248 c
249 c when the bulk Richardson number at k exceeds Ricr, hbl is
250 c linearly interpolated between grid levels zgrid(k) and zgrid(k-1).
251 c
252 c The water column and the surface forcing are diagnosed for
253 c stable/ustable forcing conditions, and where hbl is relative
254 c to grid points (caseA), so that conditional branches can be
255 c avoided in later subroutines.
256 c
257 IMPLICIT NONE
258
259 #include "SIZE.h"
260 #include "EEPARAMS.h"
261 #include "PARAMS.h"
262 #include "KPP_PARAMS.h"
263 #include "FFIELDS.h"
264
265 c input
266 c------
267 c myTime : current time in simulation
268 c myThid : thread number for this instance of the routine
269 c kmtj : number of vertical layers
270 c dvsq : (velocity shear re sfc)^2 ((m/s)^2)
271 c dbloc : local delta buoyancy across interfaces (m/s^2)
272 c Ritop : numerator of bulk Richardson Number
273 c =(z-zref)*dbsfc, where dbsfc=delta
274 c buoyancy with respect to surface ((m/s)^2)
275 c ustar : surface friction velocity (m/s)
276 c bo : surface turbulent buoyancy forcing (m^2/s^3)
277 c bosol : radiative buoyancy forcing (m^2/s^3)
278 c coriol : Coriolis parameter (1/s)
279 _RL mytime
280 integer mythid
281 integer kmtj(imt)
282 _KPP_RL dvsq (imt,Nr)
283 _KPP_RL dbloc (imt,Nr)
284 _KPP_RL Ritop (imt,Nr)
285 _KPP_RL ustar (imt)
286 _KPP_RL bo (imt)
287 _KPP_RL bosol (imt)
288 _KPP_RL coriol(imt)
289 integer ikppkey, kkppkey
290
291 c output
292 c--------
293 c hbl : boundary layer depth (m)
294 c bfsfc : Bo+radiation absorbed to d=hbf*hbl (m^2/s^3)
295 c stable : =1 in stable forcing; =0 unstable
296 c casea : =1 in case A, =0 in case B
297 c kbl : -1 of first grid level below hbl
298 c Rib : Bulk Richardson number
299 c sigma : normalized depth (d/hbl)
300 _KPP_RL hbl (imt)
301 _KPP_RL bfsfc (imt)
302 _KPP_RL stable(imt)
303 _KPP_RL casea (imt)
304 integer kbl (imt)
305 _KPP_RL Rib (imt,Nr)
306 _KPP_RL sigma (imt)
307
308 #ifdef ALLOW_KPP
309
310 c local
311 c-------
312 c wm, ws : turbulent velocity scales (m/s)
313 _KPP_RL wm(imt), ws(imt)
314 _RL worka(imt)
315
316 _KPP_RL bvsq, vtsq, hekman, hmonob, hlimit, tempVar1, tempVar2
317 integer i, kl
318
319 _KPP_RL p5 , eins
320 parameter ( p5=0.5, eins=1.0 )
321 _RL minusone
322 parameter ( minusone=-1.0 )
323
324 c find bulk Richardson number at every grid level until > Ricr
325 c
326 c note: the reference depth is -epsilon/2.*zgrid(k), but the reference
327 c u,v,t,s values are simply the surface layer values,
328 c and not the averaged values from 0 to 2*ref.depth,
329 c which is necessary for very fine grids(top layer < 2m thickness)
330 c note: max values when Ricr never satisfied are
331 c kbl(i)=kmtj(i) and hbl(i)=-zgrid(kmtj(i))
332
333 c initialize hbl and kbl to bottomed out values
334
335 do i = 1, imt
336 Rib(i,1) = 0.0
337 kbl(i) = max(kmtj(i),1)
338 hbl(i) = -zgrid(kbl(i))
339 end do
340
341 do kl = 2, Nr
342
343 #ifdef ALLOW_AUTODIFF_TAMC
344 kkppkey = (ikppkey-1)*Nr + kl
345 #endif
346
347 c compute bfsfc = sw fraction at hbf * zgrid
348
349 do i = 1, imt
350 worka(i) = zgrid(kl)
351 end do
352 CADJ store worka = comlev1_kpp_k, key = kkppkey
353 call SWFRAC(
354 I imt, hbf,
355 I mytime, mythid,
356 U worka )
357 CADJ store worka = comlev1_kpp_k, key = kkppkey
358
359
360 do i = 1, imt
361
362 c use caseA as temporary array
363
364 casea(i) = -zgrid(kl)
365
366 c compute bfsfc= Bo + radiative contribution down to hbf * hbl
367
368 bfsfc(i) = bo(i) + bosol(i)*(1. - worka(i))
369 stable(i) = p5 + sign(p5,bfsfc(i))
370 sigma(i) = stable(i) + (1. - stable(i)) * epsilon
371
372 end do
373
374 c-----------------------------------------------------------------------
375 c compute velocity scales at sigma, for hbl= caseA = -zgrid(kl)
376 c-----------------------------------------------------------------------
377
378 call wscale (
379 I sigma, casea, ustar, bfsfc,
380 O wm, ws )
381 CADJ store ws = comlev1_kpp_k, key = kkppkey
382
383 do i = 1, imt
384
385 c-----------------------------------------------------------------------
386 c compute the turbulent shear contribution to Rib
387 c-----------------------------------------------------------------------
388
389 bvsq = p5 *
390 1 ( dbloc(i,kl-1) / (zgrid(kl-1)-zgrid(kl ))+
391 2 dbloc(i,kl ) / (zgrid(kl )-zgrid(kl+1)))
392
393 if (bvsq .eq. 0.) then
394 vtsq = 0.0
395 else
396 vtsq = -zgrid(kl) * ws(i) * sqrt(abs(bvsq)) * Vtc
397 endif
398
399 c compute bulk Richardson number at new level
400 c note: Ritop needs to be zero on land and ocean bottom
401 c points so that the following if statement gets triggered
402 c correctly; otherwise, hbl might get set to (big) negative
403 c values, that might exceed the limit for the "exp" function
404 c in "SWFRAC"
405
406 c
407 c rg: assignment to double precision variable to avoid overflow
408 c ph: test for zero nominator
409 c
410
411 tempVar1 = dvsq(i,kl) + vtsq
412 tempVar2 = max(tempVar1, phepsi)
413 Rib(i,kl) = Ritop(i,kl) / tempVar2
414
415 end do
416 end do
417
418 cph(
419 cph without this store, there is a recomputation error for
420 cph rib in adbldepth (probably partial recomputation problem)
421 CADJ store Rib = comlev1_kpp
422 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy),Nr /)
423 cph)
424
425 do kl = 2, Nr
426 do i = 1, imt
427 if (kbl(i).eq.kmtj(i) .and. Rib(i,kl).gt.Ricr) kbl(i) = kl
428 end do
429 end do
430
431 CADJ store kbl = comlev1_kpp
432 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
433
434 do i = 1, imt
435 kl = kbl(i)
436 c linearly interpolate to find hbl where Rib = Ricr
437 if (kl.gt.1 .and. kl.lt.kmtj(i)) then
438 tempVar1 = (Rib(i,kl)-Rib(i,kl-1))
439 hbl(i) = -zgrid(kl-1) + (zgrid(kl-1)-zgrid(kl)) *
440 1 (Ricr - Rib(i,kl-1)) / tempVar1
441 endif
442 end do
443
444 CADJ store hbl = comlev1_kpp
445 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
446
447 c-----------------------------------------------------------------------
448 c find stability and buoyancy forcing for boundary layer
449 c-----------------------------------------------------------------------
450
451 do i = 1, imt
452 worka(i) = hbl(i)
453 end do
454 CADJ store worka = comlev1_kpp
455 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
456 call SWFRAC(
457 I imt, minusone,
458 I mytime, mythid,
459 U worka )
460 CADJ store worka = comlev1_kpp
461 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
462
463 do i = 1, imt
464 bfsfc(i) = bo(i) + bosol(i) * (1. - worka(i))
465 end do
466 CADJ store bfsfc = comlev1_kpp
467 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
468
469 c-- ensure bfsfc is never 0
470 do i = 1, imt
471 stable(i) = p5 + sign( p5, bfsfc(i) )
472 bfsfc(i) = sign(eins,bfsfc(i))*max(phepsi,abs(bfsfc(i)))
473 end do
474
475 cph(
476 cph added stable to store list to avoid extensive recomp.
477 CADJ store bfsfc, stable = comlev1_kpp
478 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
479 cph)
480
481 c-----------------------------------------------------------------------
482 c check hbl limits for hekman or hmonob
483 c ph: test for zero nominator
484 c-----------------------------------------------------------------------
485
486 do i = 1, imt
487 if (bfsfc(i) .gt. 0.0) then
488 hekman = cekman * ustar(i) / max(abs(Coriol(i)),phepsi)
489 hmonob = cmonob * ustar(i)*ustar(i)*ustar(i)
490 & / vonk / bfsfc(i)
491 hlimit = stable(i) * min(hekman,hmonob)
492 & + (stable(i)-1.) * zgrid(Nr)
493 hbl(i) = min(hbl(i),hlimit)
494 end if
495 end do
496 CADJ store hbl = comlev1_kpp
497 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
498
499 do i = 1, imt
500 hbl(i) = max(hbl(i),minKPPhbl)
501 kbl(i) = kmtj(i)
502 end do
503
504 CADJ store hbl = comlev1_kpp
505 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
506
507 c-----------------------------------------------------------------------
508 c find new kbl
509 c-----------------------------------------------------------------------
510
511 do kl = 2, Nr
512 do i = 1, imt
513 if ( kbl(i).eq.kmtj(i) .and. (-zgrid(kl)).gt.hbl(i) ) then
514 kbl(i) = kl
515 endif
516 end do
517 end do
518
519 c-----------------------------------------------------------------------
520 c find stability and buoyancy forcing for final hbl values
521 c-----------------------------------------------------------------------
522
523 do i = 1, imt
524 worka(i) = hbl(i)
525 end do
526 CADJ store worka = comlev1_kpp
527 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
528 call SWFRAC(
529 I imt, minusone,
530 I mytime, mythid,
531 U worka )
532 CADJ store worka = comlev1_kpp
533 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
534
535 do i = 1, imt
536 bfsfc(i) = bo(i) + bosol(i) * (1. - worka(i))
537 end do
538 CADJ store bfsfc = comlev1_kpp
539 CADJ & , key = ikppkey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
540
541 c-- ensures bfsfc is never 0
542 do i = 1, imt
543 stable(i) = p5 + sign( p5, bfsfc(i) )
544 bfsfc(i) = sign(eins,bfsfc(i))*max(phepsi,abs(bfsfc(i)))
545 end do
546
547 c-----------------------------------------------------------------------
548 c determine caseA and caseB
549 c-----------------------------------------------------------------------
550
551 do i = 1, imt
552 casea(i) = p5 +
553 1 sign(p5, -zgrid(kbl(i)) - p5*hwide(kbl(i)) - hbl(i))
554 end do
555
556 #endif /* ALLOW_KPP */
557
558 return
559 end
560
561 c*************************************************************************
562
563 subroutine wscale (
564 I sigma, hbl, ustar, bfsfc,
565 O wm, ws )
566
567 c compute turbulent velocity scales.
568 c use a 2D-lookup table for wm and ws as functions of ustar and
569 c zetahat (=vonk*sigma*hbl*bfsfc).
570 c
571 c note: the lookup table is only used for unstable conditions
572 c (zehat.le.0), in the stable domain wm (=ws) gets computed
573 c directly.
574 c
575 IMPLICIT NONE
576
577 #include "SIZE.h"
578 #include "KPP_PARAMS.h"
579
580 c input
581 c------
582 c sigma : normalized depth (d/hbl)
583 c hbl : boundary layer depth (m)
584 c ustar : surface friction velocity (m/s)
585 c bfsfc : total surface buoyancy flux (m^2/s^3)
586 _KPP_RL sigma(imt)
587 _KPP_RL hbl (imt)
588 _KPP_RL ustar(imt)
589 _KPP_RL bfsfc(imt)
590
591 c output
592 c--------
593 c wm, ws : turbulent velocity scales at sigma
594 _KPP_RL wm(imt), ws(imt)
595
596 #ifdef ALLOW_KPP
597
598 c local
599 c------
600 c zehat : = zeta * ustar**3
601 _KPP_RL zehat
602
603 integer iz, izp1, ju, i, jup1
604 _KPP_RL udiff, zdiff, zfrac, ufrac, fzfrac, wam
605 _KPP_RL wbm, was, wbs, u3, tempVar
606
607 c-----------------------------------------------------------------------
608 c use lookup table for zehat < zmax only; otherwise use
609 c stable formulae
610 c-----------------------------------------------------------------------
611
612 do i = 1, imt
613 zehat = vonk*sigma(i)*hbl(i)*bfsfc(i)
614
615 if (zehat .le. zmax) then
616
617 zdiff = zehat - zmin
618 iz = int( zdiff / deltaz )
619 iz = min( iz, nni )
620 iz = max( iz, 0 )
621 izp1 = iz + 1
622
623 udiff = ustar(i) - umin
624 ju = int( udiff / deltau )
625 ju = min( ju, nnj )
626 ju = max( ju, 0 )
627 jup1 = ju + 1
628
629 zfrac = zdiff / deltaz - float(iz)
630 ufrac = udiff / deltau - float(ju)
631
632 fzfrac= 1. - zfrac
633 wam = fzfrac * wmt(iz,jup1) + zfrac * wmt(izp1,jup1)
634 wbm = fzfrac * wmt(iz,ju ) + zfrac * wmt(izp1,ju )
635 wm(i) = (1.-ufrac) * wbm + ufrac * wam
636
637 was = fzfrac * wst(iz,jup1) + zfrac * wst(izp1,jup1)
638 wbs = fzfrac * wst(iz,ju ) + zfrac * wst(izp1,ju )
639 ws(i) = (1.-ufrac) * wbs + ufrac * was
640
641 else
642
643 u3 = ustar(i) * ustar(i) * ustar(i)
644 tempVar = u3 + conc1 * zehat
645 wm(i) = vonk * ustar(i) * u3 / tempVar
646 ws(i) = wm(i)
647
648 endif
649
650 end do
651
652 #endif /* ALLOW_KPP */
653
654 return
655 end
656
657 c*************************************************************************
658
659 subroutine Ri_iwmix (
660 I kmtj, shsq, dbloc, dblocSm
661 I , ikppkey
662 O , diffus )
663
664 c compute interior viscosity diffusivity coefficients due
665 c to shear instability (dependent on a local Richardson number),
666 c to background internal wave activity, and
667 c to static instability (local Richardson number < 0).
668
669 IMPLICIT NONE
670
671 #include "SIZE.h"
672 #include "EEPARAMS.h"
673 #include "PARAMS.h"
674 #include "KPP_PARAMS.h"
675
676 c input
677 c kmtj (imt) number of vertical layers on this row
678 c shsq (imt,Nr) (local velocity shear)^2 ((m/s)^2)
679 c dbloc (imt,Nr) local delta buoyancy (m/s^2)
680 c dblocSm(imt,Nr) horizontally smoothed dbloc (m/s^2)
681 integer kmtj (imt)
682 _KPP_RL shsq (imt,Nr)
683 _KPP_RL dbloc (imt,Nr)
684 _KPP_RL dblocSm(imt,Nr)
685 integer ikppkey
686
687 c output
688 c diffus(imt,0:Nrp1,1) vertical viscosivity coefficient (m^2/s)
689 c diffus(imt,0:Nrp1,2) vertical scalar diffusivity (m^2/s)
690 c diffus(imt,0:Nrp1,3) vertical temperature diffusivity (m^2/s)
691 _KPP_RL diffus(imt,0:Nrp1,3)
692
693 #ifdef ALLOW_KPP
694
695 c local variables
696 c Rig local Richardson number
697 c fRi, fcon function of Rig
698 _KPP_RL Rig
699 _KPP_RL fRi, fcon
700 _KPP_RL ratio
701 integer i, ki
702 _KPP_RL c1, c0
703
704 #ifdef ALLOW_KPP_VERTICALLY_SMOOTH
705 integer mr
706 CADJ INIT kpp_ri_tape_mr = common, 1
707 #endif
708
709 c constants
710 c1 = 1.0
711 c0 = 0.0
712
713 c-----------------------------------------------------------------------
714 c compute interior gradient Ri at all interfaces ki=1,Nr, (not surface)
715 c use diffus(*,*,1) as temporary storage of Ri to be smoothed
716 c use diffus(*,*,2) as temporary storage for Brunt-Vaisala squared
717 c set values at bottom and below to nearest value above bottom
718 #ifdef ALLOW_AUTODIFF_TAMC
719 C break data flow dependence on diffus
720 diffus(1,1,1) = 0.0
721
722 do ki = 1, Nr
723 do i = 1, imt
724 diffus(i,ki,1) = 0.
725 diffus(i,ki,2) = 0.
726 diffus(i,ki,3) = 0.
727 enddo
728 enddo
729 #endif
730
731
732 do ki = 1, Nr
733 do i = 1, imt
734 if (kmtj(i) .LE. 1 ) then
735 diffus(i,ki,1) = 0.
736 diffus(i,ki,2) = 0.
737 elseif (ki .GE. kmtj(i)) then
738 diffus(i,ki,1) = diffus(i,ki-1,1)
739 diffus(i,ki,2) = diffus(i,ki-1,2)
740 else
741 diffus(i,ki,1) = dblocSm(i,ki) * (zgrid(ki)-zgrid(ki+1))
742 & / max( Shsq(i,ki), phepsi )
743 diffus(i,ki,2) = dbloc(i,ki) / (zgrid(ki)-zgrid(ki+1))
744 endif
745 end do
746 end do
747 CADJ store diffus = comlev1_kpp, key = ikppkey
748
749 c-----------------------------------------------------------------------
750 c vertically smooth Ri
751 #ifdef ALLOW_KPP_VERTICALLY_SMOOTH
752 do mr = 1, num_v_smooth_Ri
753
754 CADJ store diffus(:,:,1) = kpp_ri_tape_mr
755 CADJ & , key=mr, shape=(/ (sNx+2*OLx)*(sNy+2*OLy),Nr+2 /)
756
757 call z121 (
758 U diffus(1,0,1))
759 end do
760 #endif
761
762 c-----------------------------------------------------------------------
763 c after smoothing loop
764
765 do ki = 1, Nr
766 do i = 1, imt
767
768 c evaluate f of Brunt-Vaisala squared for convection, store in fcon
769
770 Rig = max ( diffus(i,ki,2) , BVSQcon )
771 ratio = min ( (BVSQcon - Rig) / BVSQcon, c1 )
772 fcon = c1 - ratio * ratio
773 fcon = fcon * fcon * fcon
774
775 c evaluate f of smooth Ri for shear instability, store in fRi
776
777 Rig = max ( diffus(i,ki,1), c0 )
778 ratio = min ( Rig / Riinfty , c1 )
779 fRi = c1 - ratio * ratio
780 fRi = fRi * fRi * fRi
781
782 c ----------------------------------------------------------------------
783 c evaluate diffusivities and viscosity
784 c mixing due to internal waves, and shear and static instability
785
786 #ifndef EXCLUDE_KPP_SHEAR_MIX
787 if ( .NOT. inAdMode ) then
788 diffus(i,ki,1) = viscAr + fcon * difmcon + fRi * difm0
789 diffus(i,ki,2) = diffKrNrS(ki)+ fcon*difscon + fRi*difs0
790 diffus(i,ki,3) = diffKrNrT(ki)+ fcon*difscon + fRi*difs0
791 else
792 diffus(i,ki,1) = viscAr
793 diffus(i,ki,2) = diffKrNrS(ki)
794 diffus(i,ki,3) = diffKrNrT(ki)
795 endif
796 #else
797 diffus(i,ki,1) = viscAr
798 diffus(i,ki,2) = diffKrNrS(ki)
799 diffus(i,ki,3) = diffKrNrT(ki)
800 #endif
801
802 end do
803 end do
804
805 c ------------------------------------------------------------------------
806 c set surface values to 0.0
807
808 do i = 1, imt
809 diffus(i,0,1) = c0
810 diffus(i,0,2) = c0
811 diffus(i,0,3) = c0
812 end do
813
814 #endif /* ALLOW_KPP */
815
816 return
817 end
818
819 c*************************************************************************
820
821 subroutine z121 (
822 U v )
823
824 c Apply 121 smoothing in k to 2-d array V(i,k=1,Nr)
825 c top (0) value is used as a dummy
826 c bottom (Nrp1) value is set to input value from above.
827
828 c Note that it is important to exclude from the smoothing any points
829 c that are outside the range of the K(Ri) scheme, ie. >0.8, or <0.0.
830 c Otherwise, there is interference with other physics, especially
831 c penetrative convection.
832
833 IMPLICIT NONE
834 #include "SIZE.h"
835 #include "KPP_PARAMS.h"
836
837 c input/output
838 c-------------
839 c v : 2-D array to be smoothed in Nrp1 direction
840 _KPP_RL v(imt,0:Nrp1)
841
842 #ifdef ALLOW_KPP
843
844 c local
845 _KPP_RL zwork, zflag
846 _KPP_RL KRi_range(1:Nrp1)
847 integer i, k, km1, kp1
848
849 _KPP_RL p0 , p25 , p5 , p2
850 parameter ( p0 = 0.0, p25 = 0.25, p5 = 0.5, p2 = 2.0 )
851
852 KRi_range(Nrp1) = p0
853
854 #ifdef ALLOW_AUTODIFF_TAMC
855 C-- dummy assignment to end declaration part for TAMC
856 i = 0
857
858 C-- HPF directive to help TAMC
859 CHPF$ INDEPENDENT
860 CADJ INIT z121tape = common, Nr
861 #endif /* ALLOW_AUTODIFF_TAMC */
862
863 do i = 1, imt
864
865 k = 1
866 CADJ STORE v(i,k) = z121tape
867 v(i,Nrp1) = v(i,Nr)
868
869 do k = 1, Nr
870 KRi_range(k) = p5 + SIGN(p5,v(i,k))
871 KRi_range(k) = KRi_range(k) *
872 & ( p5 + SIGN(p5,(Riinfty-v(i,k))) )
873 end do
874
875 zwork = KRi_range(1) * v(i,1)
876 v(i,1) = p2 * v(i,1) +
877 & KRi_range(1) * KRi_range(2) * v(i,2)
878 zflag = p2 + KRi_range(1) * KRi_range(2)
879 v(i,1) = v(i,1) / zflag
880
881 do k = 2, Nr
882 CADJ STORE v(i,k), zwork = z121tape
883 km1 = k - 1
884 kp1 = k + 1
885 zflag = v(i,k)
886 v(i,k) = p2 * v(i,k) +
887 & KRi_range(k) * KRi_range(kp1) * v(i,kp1) +
888 & KRi_range(k) * zwork
889 zwork = KRi_range(k) * zflag
890 zflag = p2 + KRi_range(k)*(KRi_range(kp1)+KRi_range(km1))
891 v(i,k) = v(i,k) / zflag
892 end do
893
894 end do
895
896 #endif /* ALLOW_KPP */
897
898 return
899 end
900
901 c*************************************************************************
902
903 subroutine kpp_smooth_horiz (
904 I k, bi, bj,
905 U fld )
906
907 c Apply horizontal smoothing to KPP array
908
909 IMPLICIT NONE
910 #include "SIZE.h"
911 #include "GRID.h"
912 #include "KPP_PARAMS.h"
913
914 c input
915 c bi, bj : array indices
916 c k : vertical index used for masking
917 integer k, bi, bj
918
919 c input/output
920 c fld : 2-D array to be smoothed
921 _KPP_RL fld( ibot:itop, jbot:jtop )
922
923 #ifdef ALLOW_KPP
924
925 c local
926 integer i, j, im1, ip1, jm1, jp1
927 _KPP_RL tempVar
928 _KPP_RL fld_tmp( ibot:itop, jbot:jtop )
929
930 integer imin , imax , jmin , jmax
931 parameter( imin=ibot+1, imax=itop-1, jmin=jbot+1, jmax=jtop-1 )
932
933 _KPP_RL p0 , p5 , p25 , p125 , p0625
934 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
935
936 DO j = jmin, jmax
937 jm1 = j-1
938 jp1 = j+1
939 DO i = imin, imax
940 im1 = i-1
941 ip1 = i+1
942 tempVar =
943 & p25 * maskC(i ,j ,k,bi,bj) +
944 & p125 * ( maskC(im1,j ,k,bi,bj) +
945 & maskC(ip1,j ,k,bi,bj) +
946 & maskC(i ,jm1,k,bi,bj) +
947 & maskC(i ,jp1,k,bi,bj) ) +
948 & p0625 * ( maskC(im1,jm1,k,bi,bj) +
949 & maskC(im1,jp1,k,bi,bj) +
950 & maskC(ip1,jm1,k,bi,bj) +
951 & maskC(ip1,jp1,k,bi,bj) )
952 IF ( tempVar .GE. p25 ) THEN
953 fld_tmp(i,j) = (
954 & p25 * fld(i ,j )*maskC(i ,j ,k,bi,bj) +
955 & p125 *(fld(im1,j )*maskC(im1,j ,k,bi,bj) +
956 & fld(ip1,j )*maskC(ip1,j ,k,bi,bj) +
957 & fld(i ,jm1)*maskC(i ,jm1,k,bi,bj) +
958 & fld(i ,jp1)*maskC(i ,jp1,k,bi,bj))+
959 & p0625*(fld(im1,jm1)*maskC(im1,jm1,k,bi,bj) +
960 & fld(im1,jp1)*maskC(im1,jp1,k,bi,bj) +
961 & fld(ip1,jm1)*maskC(ip1,jm1,k,bi,bj) +
962 & fld(ip1,jp1)*maskC(ip1,jp1,k,bi,bj)))
963 & / tempVar
964 ELSE
965 fld_tmp(i,j) = fld(i,j)
966 ENDIF
967 ENDDO
968 ENDDO
969
970 c transfer smoothed field to output array
971 DO j = jmin, jmax
972 DO i = imin, imax
973 fld(i,j) = fld_tmp(i,j)
974 ENDDO
975 ENDDO
976
977 #endif /* ALLOW_KPP */
978
979 return
980 end
981
982 c*************************************************************************
983
984 subroutine smooth_horiz (
985 I k, bi, bj,
986 U fld )
987
988 c Apply horizontal smoothing to global _RL 2-D array
989
990 IMPLICIT NONE
991 #include "SIZE.h"
992 #include "GRID.h"
993 #include "KPP_PARAMS.h"
994
995 c input
996 c bi, bj : array indices
997 c k : vertical index used for masking
998 integer k, bi, bj
999
1000 c input/output
1001 c fld : 2-D array to be smoothed
1002 _RL fld( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
1003
1004 #ifdef ALLOW_KPP
1005
1006 c local
1007 integer i, j, im1, ip1, jm1, jp1
1008 _RL tempVar
1009 _RL fld_tmp( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
1010
1011 integer imin , imax , jmin , jmax
1012 parameter(imin=2-OLx, imax=sNx+OLx-1, jmin=2-OLy, jmax=sNy+OLy-1)
1013
1014 _RL p0 , p5 , p25 , p125 , p0625
1015 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
1016
1017 DO j = jmin, jmax
1018 jm1 = j-1
1019 jp1 = j+1
1020 DO i = imin, imax
1021 im1 = i-1
1022 ip1 = i+1
1023 tempVar =
1024 & p25 * maskC(i ,j ,k,bi,bj) +
1025 & p125 * ( maskC(im1,j ,k,bi,bj) +
1026 & maskC(ip1,j ,k,bi,bj) +
1027 & maskC(i ,jm1,k,bi,bj) +
1028 & maskC(i ,jp1,k,bi,bj) ) +
1029 & p0625 * ( maskC(im1,jm1,k,bi,bj) +
1030 & maskC(im1,jp1,k,bi,bj) +
1031 & maskC(ip1,jm1,k,bi,bj) +
1032 & maskC(ip1,jp1,k,bi,bj) )
1033 IF ( tempVar .GE. p25 ) THEN
1034 fld_tmp(i,j) = (
1035 & p25 * fld(i ,j )*maskC(i ,j ,k,bi,bj) +
1036 & p125 *(fld(im1,j )*maskC(im1,j ,k,bi,bj) +
1037 & fld(ip1,j )*maskC(ip1,j ,k,bi,bj) +
1038 & fld(i ,jm1)*maskC(i ,jm1,k,bi,bj) +
1039 & fld(i ,jp1)*maskC(i ,jp1,k,bi,bj))+
1040 & p0625*(fld(im1,jm1)*maskC(im1,jm1,k,bi,bj) +
1041 & fld(im1,jp1)*maskC(im1,jp1,k,bi,bj) +
1042 & fld(ip1,jm1)*maskC(ip1,jm1,k,bi,bj) +
1043 & fld(ip1,jp1)*maskC(ip1,jp1,k,bi,bj)))
1044 & / tempVar
1045 ELSE
1046 fld_tmp(i,j) = fld(i,j)
1047 ENDIF
1048 ENDDO
1049 ENDDO
1050
1051 c transfer smoothed field to output array
1052 DO j = jmin, jmax
1053 DO i = imin, imax
1054 fld(i,j) = fld_tmp(i,j)
1055 ENDDO
1056 ENDDO
1057
1058 #endif /* ALLOW_KPP */
1059
1060 return
1061 end
1062
1063 c*************************************************************************
1064
1065 subroutine blmix (
1066 I ustar, bfsfc, hbl, stable, casea, diffus, kbl
1067 O , dkm1, blmc, ghat, sigma, ikppkey
1068 & )
1069
1070 c mixing coefficients within boundary layer depend on surface
1071 c forcing and the magnitude and gradient of interior mixing below
1072 c the boundary layer ("matching").
1073 c
1074 c caution: if mixing bottoms out at hbl = -zgrid(Nr) then
1075 c fictitious layer at Nrp1 is needed with small but finite width
1076 c hwide(Nrp1) (eg. epsln = 1.e-20).
1077 c
1078 IMPLICIT NONE
1079
1080 #include "SIZE.h"
1081 #include "KPP_PARAMS.h"
1082
1083 c input
1084 c ustar (imt) surface friction velocity (m/s)
1085 c bfsfc (imt) surface buoyancy forcing (m^2/s^3)
1086 c hbl (imt) boundary layer depth (m)
1087 c stable(imt) = 1 in stable forcing
1088 c casea (imt) = 1 in case A
1089 c diffus(imt,0:Nrp1,mdiff) vertical diffusivities (m^2/s)
1090 c kbl(imt) -1 of first grid level below hbl
1091 _KPP_RL ustar (imt)
1092 _KPP_RL bfsfc (imt)
1093 _KPP_RL hbl (imt)
1094 _KPP_RL stable(imt)
1095 _KPP_RL casea (imt)
1096 _KPP_RL diffus(imt,0:Nrp1,mdiff)
1097 integer kbl(imt)
1098
1099 c output
1100 c dkm1 (imt,mdiff) boundary layer difs at kbl-1 level
1101 c blmc (imt,Nr,mdiff) boundary layer mixing coefficients (m^2/s)
1102 c ghat (imt,Nr) nonlocal scalar transport
1103 c sigma(imt) normalized depth (d / hbl)
1104 _KPP_RL dkm1 (imt,mdiff)
1105 _KPP_RL blmc (imt,Nr,mdiff)
1106 _KPP_RL ghat (imt,Nr)
1107 _KPP_RL sigma(imt)
1108 integer ikppkey, kkppkey
1109
1110 #ifdef ALLOW_KPP
1111
1112 c local
1113 c gat1*(imt) shape function at sigma = 1
1114 c dat1*(imt) derivative of shape function at sigma = 1
1115 c ws(imt), wm(imt) turbulent velocity scales (m/s)
1116 _KPP_RL gat1m(imt), gat1s(imt), gat1t(imt)
1117 _KPP_RL dat1m(imt), dat1s(imt), dat1t(imt)
1118 _KPP_RL ws(imt), wm(imt)
1119 integer i, kn, ki
1120 _KPP_RL R, dvdzup, dvdzdn, viscp
1121 _KPP_RL difsp, diftp, visch, difsh, difth
1122 _KPP_RL f1, sig, a1, a2, a3, delhat
1123 _KPP_RL Gm, Gs, Gt
1124 _KPP_RL tempVar
1125
1126 _KPP_RL p0 , eins
1127 parameter (p0=0.0, eins=1.0)
1128
1129 c-----------------------------------------------------------------------
1130 c compute velocity scales at hbl
1131 c-----------------------------------------------------------------------
1132
1133 do i = 1, imt
1134 sigma(i) = stable(i) * 1.0 + (1. - stable(i)) * epsilon
1135 end do
1136
1137 CADJ STORE sigma = comlev1_kpp, key = ikppkey
1138 call wscale (
1139 I sigma, hbl, ustar, bfsfc,
1140 O wm, ws )
1141 CADJ STORE wm = comlev1_kpp, key = ikppkey
1142 CADJ STORE ws = comlev1_kpp, key = ikppkey
1143
1144 do i = 1, imt
1145 wm(i) = sign(eins,wm(i))*max(phepsi,abs(wm(i)))
1146 ws(i) = sign(eins,ws(i))*max(phepsi,abs(ws(i)))
1147 end do
1148 CADJ STORE wm = comlev1_kpp, key = ikppkey
1149 CADJ STORE ws = comlev1_kpp, key = ikppkey
1150
1151 do i = 1, imt
1152
1153 kn = int(caseA(i)+phepsi) *(kbl(i) -1) +
1154 $ (1 - int(caseA(i)+phepsi)) * kbl(i)
1155
1156 c-----------------------------------------------------------------------
1157 c find the interior viscosities and derivatives at hbl(i)
1158 c-----------------------------------------------------------------------
1159
1160 delhat = 0.5*hwide(kn) - zgrid(kn) - hbl(i)
1161 R = 1.0 - delhat / hwide(kn)
1162 dvdzup = (diffus(i,kn-1,1) - diffus(i,kn ,1)) / hwide(kn)
1163 dvdzdn = (diffus(i,kn ,1) - diffus(i,kn+1,1)) / hwide(kn+1)
1164 viscp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1165 1 R * (dvdzdn + abs(dvdzdn)) )
1166
1167 dvdzup = (diffus(i,kn-1,2) - diffus(i,kn ,2)) / hwide(kn)
1168 dvdzdn = (diffus(i,kn ,2) - diffus(i,kn+1,2)) / hwide(kn+1)
1169 difsp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1170 1 R * (dvdzdn + abs(dvdzdn)) )
1171
1172 dvdzup = (diffus(i,kn-1,3) - diffus(i,kn ,3)) / hwide(kn)
1173 dvdzdn = (diffus(i,kn ,3) - diffus(i,kn+1,3)) / hwide(kn+1)
1174 diftp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1175 1 R * (dvdzdn + abs(dvdzdn)) )
1176
1177 visch = diffus(i,kn,1) + viscp * delhat
1178 difsh = diffus(i,kn,2) + difsp * delhat
1179 difth = diffus(i,kn,3) + diftp * delhat
1180
1181 f1 = stable(i) * conc1 * bfsfc(i) /
1182 & max(ustar(i)**4,phepsi)
1183 gat1m(i) = visch / hbl(i) / wm(i)
1184 dat1m(i) = -viscp / wm(i) + f1 * visch
1185
1186 gat1s(i) = difsh / hbl(i) / ws(i)
1187 dat1s(i) = -difsp / ws(i) + f1 * difsh
1188
1189 gat1t(i) = difth / hbl(i) / ws(i)
1190 dat1t(i) = -diftp / ws(i) + f1 * difth
1191
1192 end do
1193 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1194 CADJ STORE gat1m = comlev1_kpp, key = ikppkey
1195 CADJ STORE gat1s = comlev1_kpp, key = ikppkey
1196 CADJ STORE gat1t = comlev1_kpp, key = ikppkey
1197 CADJ STORE dat1m = comlev1_kpp, key = ikppkey
1198 CADJ STORE dat1s = comlev1_kpp, key = ikppkey
1199 CADJ STORE dat1t = comlev1_kpp, key = ikppkey
1200 #endif
1201 do i = 1, imt
1202 dat1m(i) = min(dat1m(i),p0)
1203 dat1s(i) = min(dat1s(i),p0)
1204 dat1t(i) = min(dat1t(i),p0)
1205 end do
1206 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1207 CADJ STORE dat1m = comlev1_kpp, key = ikppkey
1208 CADJ STORE dat1s = comlev1_kpp, key = ikppkey
1209 CADJ STORE dat1t = comlev1_kpp, key = ikppkey
1210 #endif
1211
1212 do ki = 1, Nr
1213
1214 #ifdef ALLOW_AUTODIFF_TAMC
1215 kkppkey = (ikppkey-1)*Nr + ki
1216 #endif
1217
1218 c-----------------------------------------------------------------------
1219 c compute turbulent velocity scales on the interfaces
1220 c-----------------------------------------------------------------------
1221
1222 do i = 1, imt
1223 sig = (-zgrid(ki) + 0.5 * hwide(ki)) / hbl(i)
1224 sigma(i) = stable(i)*sig + (1.-stable(i))*min(sig,epsilon)
1225 end do
1226 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1227 CADJ STORE wm = comlev1_kpp_k, key = kkppkey
1228 CADJ STORE ws = comlev1_kpp_k, key = kkppkey
1229 #endif
1230 CADJ STORE sigma = comlev1_kpp_k, key = kkppkey
1231 call wscale (
1232 I sigma, hbl, ustar, bfsfc,
1233 O wm, ws )
1234 CADJ STORE wm = comlev1_kpp_k, key = kkppkey
1235 CADJ STORE ws = comlev1_kpp_k, key = kkppkey
1236
1237 c-----------------------------------------------------------------------
1238 c compute the dimensionless shape functions at the interfaces
1239 c-----------------------------------------------------------------------
1240
1241 do i = 1, imt
1242 sig = (-zgrid(ki) + 0.5 * hwide(ki)) / hbl(i)
1243 a1 = sig - 2.
1244 a2 = 3. - 2. * sig
1245 a3 = sig - 1.
1246
1247 Gm = a1 + a2 * gat1m(i) + a3 * dat1m(i)
1248 Gs = a1 + a2 * gat1s(i) + a3 * dat1s(i)
1249 Gt = a1 + a2 * gat1t(i) + a3 * dat1t(i)
1250
1251 c-----------------------------------------------------------------------
1252 c compute boundary layer diffusivities at the interfaces
1253 c-----------------------------------------------------------------------
1254
1255 blmc(i,ki,1) = hbl(i) * wm(i) * sig * (1. + sig * Gm)
1256 blmc(i,ki,2) = hbl(i) * ws(i) * sig * (1. + sig * Gs)
1257 blmc(i,ki,3) = hbl(i) * ws(i) * sig * (1. + sig * Gt)
1258
1259 c-----------------------------------------------------------------------
1260 c nonlocal transport term = ghat * <ws>o
1261 c-----------------------------------------------------------------------
1262
1263 tempVar = ws(i) * hbl(i)
1264 ghat(i,ki) = (1.-stable(i)) * cg / max(phepsi,tempVar)
1265
1266 end do
1267 end do
1268
1269 c-----------------------------------------------------------------------
1270 c find diffusivities at kbl-1 grid level
1271 c-----------------------------------------------------------------------
1272
1273 do i = 1, imt
1274 sig = -zgrid(kbl(i)-1) / hbl(i)
1275 sigma(i) = stable(i) * sig
1276 & + (1. - stable(i)) * min(sig,epsilon)
1277 end do
1278
1279 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1280 CADJ STORE wm = comlev1_kpp, key = ikppkey
1281 CADJ STORE ws = comlev1_kpp, key = ikppkey
1282 #endif
1283 CADJ STORE sigma = comlev1_kpp, key = ikppkey
1284 call wscale (
1285 I sigma, hbl, ustar, bfsfc,
1286 O wm, ws )
1287 CADJ STORE wm = comlev1_kpp, key = ikppkey
1288 CADJ STORE ws = comlev1_kpp, key = ikppkey
1289
1290 do i = 1, imt
1291 sig = -zgrid(kbl(i)-1) / hbl(i)
1292 a1 = sig - 2.
1293 a2 = 3. - 2. * sig
1294 a3 = sig - 1.
1295 Gm = a1 + a2 * gat1m(i) + a3 * dat1m(i)
1296 Gs = a1 + a2 * gat1s(i) + a3 * dat1s(i)
1297 Gt = a1 + a2 * gat1t(i) + a3 * dat1t(i)
1298 dkm1(i,1) = hbl(i) * wm(i) * sig * (1. + sig * Gm)
1299 dkm1(i,2) = hbl(i) * ws(i) * sig * (1. + sig * Gs)
1300 dkm1(i,3) = hbl(i) * ws(i) * sig * (1. + sig * Gt)
1301 end do
1302
1303 #endif /* ALLOW_KPP */
1304
1305 return
1306 end
1307
1308 c*************************************************************************
1309
1310 subroutine enhance (
1311 I dkm1, hbl, kbl, diffus, casea
1312 U , ghat
1313 O , blmc
1314 & )
1315
1316 c enhance the diffusivity at the kbl-.5 interface
1317
1318 IMPLICIT NONE
1319
1320 #include "SIZE.h"
1321 #include "KPP_PARAMS.h"
1322
1323 c input
1324 c dkm1(imt,mdiff) bl diffusivity at kbl-1 grid level
1325 c hbl(imt) boundary layer depth (m)
1326 c kbl(imt) grid above hbl
1327 c diffus(imt,0:Nrp1,mdiff) vertical diffusivities (m^2/s)
1328 c casea(imt) = 1 in caseA, = 0 in case B
1329 _KPP_RL dkm1 (imt,mdiff)
1330 _KPP_RL hbl (imt)
1331 integer kbl (imt)
1332 _KPP_RL diffus(imt,0:Nrp1,mdiff)
1333 _KPP_RL casea (imt)
1334
1335 c input/output
1336 c nonlocal transport, modified ghat at kbl(i)-1 interface (s/m**2)
1337 _KPP_RL ghat (imt,Nr)
1338
1339 c output
1340 c enhanced bound. layer mixing coeff.
1341 _KPP_RL blmc (imt,Nr,mdiff)
1342
1343 #ifdef ALLOW_KPP
1344
1345 c local
1346 c fraction hbl lies beteen zgrid neighbors
1347 _KPP_RL delta
1348 integer ki, i, md
1349 _KPP_RL dkmp5, dstar
1350
1351 do i = 1, imt
1352 ki = kbl(i)-1
1353 if ((ki .ge. 1) .and. (ki .lt. Nr)) then
1354 delta = (hbl(i) + zgrid(ki)) / (zgrid(ki) - zgrid(ki+1))
1355 do md = 1, mdiff
1356 dkmp5 = casea(i) * diffus(i,ki,md) +
1357 1 (1.- casea(i)) * blmc (i,ki,md)
1358 dstar = (1.- delta)**2 * dkm1(i,md)
1359 & + delta**2 * dkmp5
1360 blmc(i,ki,md) = (1.- delta)*diffus(i,ki,md)
1361 & + delta*dstar
1362 end do
1363 ghat(i,ki) = (1.- casea(i)) * ghat(i,ki)
1364 endif
1365 end do
1366
1367 #endif /* ALLOW_KPP */
1368
1369 return
1370 end
1371
1372 c*************************************************************************
1373
1374 SUBROUTINE STATEKPP (
1375 I ikppkey, bi, bj, myThid,
1376 O RHO1, DBLOC, DBSFC, TTALPHA, SSBETA)
1377 c
1378 c-----------------------------------------------------------------------
1379 c "statekpp" computes all necessary input arrays
1380 c for the kpp mixing scheme
1381 c
1382 c input:
1383 c bi, bj = array indices on which to apply calculations
1384 c
1385 c output:
1386 c rho1 = potential density of surface layer (kg/m^3)
1387 c dbloc = local buoyancy gradient at Nr interfaces
1388 c g/rho{k+1,k+1} * [ drho{k,k+1}-drho{k+1,k+1} ] (m/s^2)
1389 c dbsfc = buoyancy difference with respect to the surface
1390 c g * [ drho{1,k}/rho{1,k} - drho{k,k}/rho{k,k} ] (m/s^2)
1391 c ttalpha= thermal expansion coefficient without 1/rho factor
1392 c d(rho) / d(potential temperature) (kg/m^3/C)
1393 c ssbeta = salt expansion coefficient without 1/rho factor
1394 c d(rho) / d(salinity) (kg/m^3/PSU)
1395 c
1396 c see also subroutines find_rho.F find_alpha.F find_beta.F
1397 c
1398 c written by: jan morzel, feb. 10, 1995 (converted from "sigma" version)
1399 c modified by: d. menemenlis, june 1998 : for use with MIT GCM UV
1400 c
1401 c 28 april 05: added computation of mixed layer depth KPPmld
1402 c for a deltaRho equivalent to a deltaTheta of 0.8 deg C
1403
1404 c-----------------------------------------------------------------------
1405
1406 IMPLICIT NONE
1407
1408 #include "SIZE.h"
1409 #include "EEPARAMS.h"
1410 #include "PARAMS.h"
1411 #include "KPP_PARAMS.h"
1412 #include "DYNVARS.h"
1413 #include "GRID.h"
1414
1415 c-------------- Routine arguments -----------------------------------------
1416 INTEGER bi, bj, myThid
1417 _KPP_RL RHO1 ( ibot:itop, jbot:jtop )
1418 _KPP_RL DBLOC ( ibot:itop, jbot:jtop, Nr )
1419 _KPP_RL DBSFC ( ibot:itop, jbot:jtop, Nr )
1420 _KPP_RL TTALPHA( ibot:itop, jbot:jtop, Nrp1 )
1421 _KPP_RL SSBETA ( ibot:itop, jbot:jtop, Nrp1 )
1422
1423 #ifdef ALLOW_KPP
1424
1425 c--------------------------------------------------------------------------
1426 c
1427 c local arrays:
1428 c
1429 c rhok - density of t(k ) & s(k ) at depth k
1430 c rhokm1 - density of t(k-1) & s(k-1) at depth k
1431 c rho1k - density of t(1 ) & s(1 ) at depth k
1432 c work1,2,3 - work arrays for holding horizontal slabs
1433
1434 _RL RHOK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1435 _RL RHOKM1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1436 _RL RHO1K (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1437 _RL WORK1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1438 _RL WORK2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1439 _RL WORK3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1440 #ifdef ALLOW_DIAGNOSTICS
1441 c KPPMLD - mixed layer depth based on density criterion
1442 _RL KPPMLD(1 :sNx ,1 :sNy )
1443 #endif /* ALLOW_DIAGNOSTICS */
1444
1445 INTEGER I, J, K
1446 INTEGER ikppkey, kkppkey
1447
1448 c calculate density, alpha, beta in surface layer, and set dbsfc to zero
1449
1450 kkppkey = (ikppkey-1)*Nr + 1
1451
1452 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1453 CADJ STORE theta(:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1454 CADJ STORE salt (:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1455 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1456 call FIND_RHO(
1457 I bi, bj, ibot, itop, jbot, jtop, 1, 1,
1458 I theta, salt,
1459 O WORK1,
1460 I myThid )
1461 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1462 CADJ STORE theta(:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1463 CADJ STORE salt (:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1464 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1465
1466 call FIND_ALPHA(
1467 I bi, bj, ibot, itop, jbot, jtop, 1, 1,
1468 O WORK2 )
1469
1470 call FIND_BETA(
1471 I bi, bj, ibot, itop, jbot, jtop, 1, 1,
1472 O WORK3 )
1473
1474 DO J = jbot, jtop
1475 DO I = ibot, itop
1476 RHO1(I,J) = WORK1(I,J) + rhoConst
1477 TTALPHA(I,J,1) = WORK2(I,J)
1478 SSBETA(I,J,1) = WORK3(I,J)
1479 DBSFC(I,J,1) = 0.
1480 END DO
1481 END DO
1482
1483 #ifdef ALLOW_DIAGNOSTICS
1484 c work3 - density of t(1)-.8 & s(1 ) at depth 1
1485 IF ( useDiagnostics ) THEN
1486 DO J = 1, sNy
1487 DO I = 1, sNx
1488 KPPMLD(I,J) = ABS(R_low(I,J,bi,bj))
1489 WORK3 (I,J) = WORK1(I,J) - 0.8 * WORK2(I,J)
1490 END DO
1491 END DO
1492 ENDIF
1493 #endif /* ALLOW_DIAGNOSTICS */
1494
1495 c calculate alpha, beta, and gradients in interior layers
1496
1497 CHPF$ INDEPENDENT, NEW (RHOK,RHOKM1,RHO1K,WORK1,WORK2)
1498 DO K = 2, Nr
1499
1500 kkppkey = (ikppkey-1)*Nr + k
1501
1502 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1503 CADJ STORE theta(:,:,k,bi,bj) = comlev1_kpp_k, key=kkppkey
1504 CADJ STORE salt (:,:,k,bi,bj) = comlev1_kpp_k, key=kkppkey
1505 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1506 call FIND_RHO(
1507 I bi, bj, ibot, itop, jbot, jtop, K, K,
1508 I theta, salt,
1509 O RHOK,
1510 I myThid )
1511
1512 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1513 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_kpp_k, key=kkppkey
1514 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_kpp_k, key=kkppkey
1515 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1516 call FIND_RHO(
1517 I bi, bj, ibot, itop, jbot, jtop, K-1, K,
1518 I theta, salt,
1519 O RHOKM1,
1520 I myThid )
1521
1522 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1523 CADJ STORE theta(:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1524 CADJ STORE salt (:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1525 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1526 call FIND_RHO(
1527 I bi, bj, ibot, itop, jbot, jtop, 1, K,
1528 I theta, salt,
1529 O RHO1K,
1530 I myThid )
1531
1532 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1533 CADJ STORE rhok (:,:) = comlev1_kpp_k, key=kkppkey
1534 CADJ STORE rhokm1(:,:) = comlev1_kpp_k, key=kkppkey
1535 CADJ STORE rho1k (:,:) = comlev1_kpp_k, key=kkppkey
1536 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1537
1538 call FIND_ALPHA(
1539 I bi, bj, ibot, itop, jbot, jtop, K, K,
1540 O WORK1 )
1541
1542 call FIND_BETA(
1543 I bi, bj, ibot, itop, jbot, jtop, K, K,
1544 O WORK2 )
1545
1546 DO J = jbot, jtop
1547 DO I = ibot, itop
1548 TTALPHA(I,J,K) = WORK1 (I,J)
1549 SSBETA(I,J,K) = WORK2 (I,J)
1550 DBLOC(I,J,K-1) = gravity * (RHOK(I,J) - RHOKM1(I,J)) /
1551 & (RHOK(I,J) + rhoConst)
1552 DBSFC(I,J,K) = gravity * (RHOK(I,J) - RHO1K (I,J)) /
1553 & (RHOK(I,J) + rhoConst)
1554 END DO
1555 END DO
1556
1557 #ifdef ALLOW_DIAGNOSTICS
1558 IF ( useDiagnostics ) THEN
1559 c work1 - density of t(k-1) & s(k-1) at depth 1
1560 c work2 - density of t(k ) & s(k ) at depth 1
1561 c work3 - density of t(1)-.8 & s(1 ) at depth 1
1562 call FIND_RHO(
1563 I bi, bj, 1, sNx, 1, sNy, K-1, 1, theta, salt,
1564 O WORK1,
1565 I myThid )
1566 call FIND_RHO(
1567 I bi, bj, 1, sNx, 1, sNy, K , 1, theta, salt,
1568 O WORK2,
1569 I myThid )
1570 DO J = 1, sNy
1571 DO I = 1, sNx
1572 IF ( k .LE. klowC(I,J,bi,bj) .AND.
1573 & WORK1(I,J) .LT. WORK3(I,J) .AND.
1574 & WORK2(I,J) .GE. WORK3(I,J) ) THEN
1575 KPPMLD(I,J) = MIN ( KPPMLD(I,J),
1576 & ABS(((WORK3(I,J)-WORK1(I,J))*rC(k)+
1577 & (WORK2(I,J)-WORK3(I,J))*rC(k-1))/
1578 & (WORK2(I,J)-WORK1(I,J))))
1579 ENDIF
1580 END DO
1581 END DO
1582 ENDIF
1583 #endif /* ALLOW_DIAGNOSTICS */
1584
1585 END DO
1586
1587 c compute arrays for K = Nrp1
1588 DO J = jbot, jtop
1589 DO I = ibot, itop
1590 TTALPHA(I,J,Nrp1) = TTALPHA(I,J,Nr)
1591 SSBETA(I,J,Nrp1) = SSBETA(I,J,Nr)
1592 DBLOC(I,J,Nr) = 0.
1593 END DO
1594 END DO
1595
1596 #ifdef ALLOW_DIAGNOSTICS
1597 IF ( useDiagnostics ) THEN
1598 CALL DIAGNOSTICS_FILL(KPPmld,'KPPmld ',0,1,3,bi,bj,myThid)
1599 ENDIF
1600 #endif /* ALLOW_DIAGNOSTICS */
1601
1602 #endif /* ALLOW_KPP */
1603
1604 RETURN
1605 END

  ViewVC Help
Powered by ViewVC 1.1.22