/[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.28 - (show annotations) (download)
Sun Apr 29 19:30:57 2007 UTC (17 years, 1 month ago) by molod
Branch: MAIN
Changes since 1.27: +4 -7 lines
Begin modifications for heterogeneous bcs in kpp - zero diff arg list changes 1

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

  ViewVC Help
Powered by ViewVC 1.1.22