/[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.21 - (show annotations) (download)
Sun Oct 17 23:05:09 2004 UTC (19 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57b_post, checkpoint56b_post, checkpoint57, checkpoint56, checkpoint55i_post, checkpoint57a_post, checkpoint57c_post, checkpoint57c_pre, checkpoint55j_post, checkpoint55h_post, checkpoint56a_post, checkpoint56c_post, checkpoint57a_pre
Changes since 1.20: +7 -7 lines
allow to set a vertical profile of vertical diffusivity for T & S

1 C $Header: /u/gcmpack/MITgcm/pkg/kpp/kpp_routines.F,v 1.20 2004/07/18 01:18:55 jmc 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 diffus(i,ki,1) = viscAr + fcon * difmcon + fRi * difm0
788 diffus(i,ki,2) = diffKrNrS(ki)+ fcon*difscon + fRi*difs0
789 diffus(i,ki,3) = diffKrNrT(ki)+ fcon*difscon + fRi*difs0
790 #else
791 diffus(i,ki,1) = viscAr
792 diffus(i,ki,2) = diffKrNrS(ki)
793 diffus(i,ki,3) = diffKrNrT(ki)
794 #endif
795
796 end do
797 end do
798
799 c ------------------------------------------------------------------------
800 c set surface values to 0.0
801
802 do i = 1, imt
803 diffus(i,0,1) = c0
804 diffus(i,0,2) = c0
805 diffus(i,0,3) = c0
806 end do
807
808 #endif /* ALLOW_KPP */
809
810 return
811 end
812
813 c*************************************************************************
814
815 subroutine z121 (
816 U v )
817
818 c Apply 121 smoothing in k to 2-d array V(i,k=1,Nr)
819 c top (0) value is used as a dummy
820 c bottom (Nrp1) value is set to input value from above.
821
822 c Note that it is important to exclude from the smoothing any points
823 c that are outside the range of the K(Ri) scheme, ie. >0.8, or <0.0.
824 c Otherwise, there is interference with other physics, especially
825 c penetrative convection.
826
827 IMPLICIT NONE
828 #include "SIZE.h"
829 #include "KPP_PARAMS.h"
830
831 c input/output
832 c-------------
833 c v : 2-D array to be smoothed in Nrp1 direction
834 _KPP_RL v(imt,0:Nrp1)
835
836 #ifdef ALLOW_KPP
837
838 c local
839 _KPP_RL zwork, zflag
840 _KPP_RL KRi_range(1:Nrp1)
841 integer i, k, km1, kp1
842
843 _KPP_RL p0 , p25 , p5 , p2
844 parameter ( p0 = 0.0, p25 = 0.25, p5 = 0.5, p2 = 2.0 )
845
846 KRi_range(Nrp1) = p0
847
848 #ifdef ALLOW_AUTODIFF_TAMC
849 C-- dummy assignment to end declaration part for TAMC
850 i = 0
851
852 C-- HPF directive to help TAMC
853 CHPF$ INDEPENDENT
854 CADJ INIT z121tape = common, Nr
855 #endif /* ALLOW_AUTODIFF_TAMC */
856
857 do i = 1, imt
858
859 k = 1
860 CADJ STORE v(i,k) = z121tape
861 v(i,Nrp1) = v(i,Nr)
862
863 do k = 1, Nr
864 KRi_range(k) = p5 + SIGN(p5,v(i,k))
865 KRi_range(k) = KRi_range(k) *
866 & ( p5 + SIGN(p5,(Riinfty-v(i,k))) )
867 end do
868
869 zwork = KRi_range(1) * v(i,1)
870 v(i,1) = p2 * v(i,1) +
871 & KRi_range(1) * KRi_range(2) * v(i,2)
872 zflag = p2 + KRi_range(1) * KRi_range(2)
873 v(i,1) = v(i,1) / zflag
874
875 do k = 2, Nr
876 CADJ STORE v(i,k), zwork = z121tape
877 km1 = k - 1
878 kp1 = k + 1
879 zflag = v(i,k)
880 v(i,k) = p2 * v(i,k) +
881 & KRi_range(k) * KRi_range(kp1) * v(i,kp1) +
882 & KRi_range(k) * zwork
883 zwork = KRi_range(k) * zflag
884 zflag = p2 + KRi_range(k)*(KRi_range(kp1)+KRi_range(km1))
885 v(i,k) = v(i,k) / zflag
886 end do
887
888 end do
889
890 #endif /* ALLOW_KPP */
891
892 return
893 end
894
895 c*************************************************************************
896
897 subroutine kpp_smooth_horiz (
898 I k, bi, bj,
899 U fld )
900
901 c Apply horizontal smoothing to KPP array
902
903 IMPLICIT NONE
904 #include "SIZE.h"
905 #include "GRID.h"
906 #include "KPP_PARAMS.h"
907
908 c input
909 c bi, bj : array indices
910 c k : vertical index used for masking
911 integer k, bi, bj
912
913 c input/output
914 c fld : 2-D array to be smoothed
915 _KPP_RL fld( ibot:itop, jbot:jtop )
916
917 #ifdef ALLOW_KPP
918
919 c local
920 integer i, j, im1, ip1, jm1, jp1
921 _KPP_RL tempVar
922 _KPP_RL fld_tmp( ibot:itop, jbot:jtop )
923
924 integer imin , imax , jmin , jmax
925 parameter( imin=ibot+1, imax=itop-1, jmin=jbot+1, jmax=jtop-1 )
926
927 _KPP_RL p0 , p5 , p25 , p125 , p0625
928 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
929
930 DO j = jmin, jmax
931 jm1 = j-1
932 jp1 = j+1
933 DO i = imin, imax
934 im1 = i-1
935 ip1 = i+1
936 tempVar =
937 & p25 * maskC(i ,j ,k,bi,bj) +
938 & p125 * ( maskC(im1,j ,k,bi,bj) +
939 & maskC(ip1,j ,k,bi,bj) +
940 & maskC(i ,jm1,k,bi,bj) +
941 & maskC(i ,jp1,k,bi,bj) ) +
942 & p0625 * ( maskC(im1,jm1,k,bi,bj) +
943 & maskC(im1,jp1,k,bi,bj) +
944 & maskC(ip1,jm1,k,bi,bj) +
945 & maskC(ip1,jp1,k,bi,bj) )
946 IF ( tempVar .GE. p25 ) THEN
947 fld_tmp(i,j) = (
948 & p25 * fld(i ,j )*maskC(i ,j ,k,bi,bj) +
949 & p125 *(fld(im1,j )*maskC(im1,j ,k,bi,bj) +
950 & fld(ip1,j )*maskC(ip1,j ,k,bi,bj) +
951 & fld(i ,jm1)*maskC(i ,jm1,k,bi,bj) +
952 & fld(i ,jp1)*maskC(i ,jp1,k,bi,bj))+
953 & p0625*(fld(im1,jm1)*maskC(im1,jm1,k,bi,bj) +
954 & fld(im1,jp1)*maskC(im1,jp1,k,bi,bj) +
955 & fld(ip1,jm1)*maskC(ip1,jm1,k,bi,bj) +
956 & fld(ip1,jp1)*maskC(ip1,jp1,k,bi,bj)))
957 & / tempVar
958 ELSE
959 fld_tmp(i,j) = fld(i,j)
960 ENDIF
961 ENDDO
962 ENDDO
963
964 c transfer smoothed field to output array
965 DO j = jmin, jmax
966 DO i = imin, imax
967 fld(i,j) = fld_tmp(i,j)
968 ENDDO
969 ENDDO
970
971 #endif /* ALLOW_KPP */
972
973 return
974 end
975
976 c*************************************************************************
977
978 subroutine smooth_horiz (
979 I k, bi, bj,
980 U fld )
981
982 c Apply horizontal smoothing to global _RL 2-D array
983
984 IMPLICIT NONE
985 #include "SIZE.h"
986 #include "GRID.h"
987 #include "KPP_PARAMS.h"
988
989 c input
990 c bi, bj : array indices
991 c k : vertical index used for masking
992 integer k, bi, bj
993
994 c input/output
995 c fld : 2-D array to be smoothed
996 _RL fld( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
997
998 #ifdef ALLOW_KPP
999
1000 c local
1001 integer i, j, im1, ip1, jm1, jp1
1002 _RL tempVar
1003 _RL fld_tmp( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
1004
1005 integer imin , imax , jmin , jmax
1006 parameter(imin=2-OLx, imax=sNx+OLx-1, jmin=2-OLy, jmax=sNy+OLy-1)
1007
1008 _RL p0 , p5 , p25 , p125 , p0625
1009 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
1010
1011 DO j = jmin, jmax
1012 jm1 = j-1
1013 jp1 = j+1
1014 DO i = imin, imax
1015 im1 = i-1
1016 ip1 = i+1
1017 tempVar =
1018 & p25 * maskC(i ,j ,k,bi,bj) +
1019 & p125 * ( maskC(im1,j ,k,bi,bj) +
1020 & maskC(ip1,j ,k,bi,bj) +
1021 & maskC(i ,jm1,k,bi,bj) +
1022 & maskC(i ,jp1,k,bi,bj) ) +
1023 & p0625 * ( maskC(im1,jm1,k,bi,bj) +
1024 & maskC(im1,jp1,k,bi,bj) +
1025 & maskC(ip1,jm1,k,bi,bj) +
1026 & maskC(ip1,jp1,k,bi,bj) )
1027 IF ( tempVar .GE. p25 ) THEN
1028 fld_tmp(i,j) = (
1029 & p25 * fld(i ,j )*maskC(i ,j ,k,bi,bj) +
1030 & p125 *(fld(im1,j )*maskC(im1,j ,k,bi,bj) +
1031 & fld(ip1,j )*maskC(ip1,j ,k,bi,bj) +
1032 & fld(i ,jm1)*maskC(i ,jm1,k,bi,bj) +
1033 & fld(i ,jp1)*maskC(i ,jp1,k,bi,bj))+
1034 & p0625*(fld(im1,jm1)*maskC(im1,jm1,k,bi,bj) +
1035 & fld(im1,jp1)*maskC(im1,jp1,k,bi,bj) +
1036 & fld(ip1,jm1)*maskC(ip1,jm1,k,bi,bj) +
1037 & fld(ip1,jp1)*maskC(ip1,jp1,k,bi,bj)))
1038 & / tempVar
1039 ELSE
1040 fld_tmp(i,j) = fld(i,j)
1041 ENDIF
1042 ENDDO
1043 ENDDO
1044
1045 c transfer smoothed field to output array
1046 DO j = jmin, jmax
1047 DO i = imin, imax
1048 fld(i,j) = fld_tmp(i,j)
1049 ENDDO
1050 ENDDO
1051
1052 #endif /* ALLOW_KPP */
1053
1054 return
1055 end
1056
1057 c*************************************************************************
1058
1059 subroutine blmix (
1060 I ustar, bfsfc, hbl, stable, casea, diffus, kbl
1061 O , dkm1, blmc, ghat, sigma, ikppkey
1062 & )
1063
1064 c mixing coefficients within boundary layer depend on surface
1065 c forcing and the magnitude and gradient of interior mixing below
1066 c the boundary layer ("matching").
1067 c
1068 c caution: if mixing bottoms out at hbl = -zgrid(Nr) then
1069 c fictitious layer at Nrp1 is needed with small but finite width
1070 c hwide(Nrp1) (eg. epsln = 1.e-20).
1071 c
1072 IMPLICIT NONE
1073
1074 #include "SIZE.h"
1075 #include "KPP_PARAMS.h"
1076
1077 c input
1078 c ustar (imt) surface friction velocity (m/s)
1079 c bfsfc (imt) surface buoyancy forcing (m^2/s^3)
1080 c hbl (imt) boundary layer depth (m)
1081 c stable(imt) = 1 in stable forcing
1082 c casea (imt) = 1 in case A
1083 c diffus(imt,0:Nrp1,mdiff) vertical diffusivities (m^2/s)
1084 c kbl(imt) -1 of first grid level below hbl
1085 _KPP_RL ustar (imt)
1086 _KPP_RL bfsfc (imt)
1087 _KPP_RL hbl (imt)
1088 _KPP_RL stable(imt)
1089 _KPP_RL casea (imt)
1090 _KPP_RL diffus(imt,0:Nrp1,mdiff)
1091 integer kbl(imt)
1092
1093 c output
1094 c dkm1 (imt,mdiff) boundary layer difs at kbl-1 level
1095 c blmc (imt,Nr,mdiff) boundary layer mixing coefficients (m^2/s)
1096 c ghat (imt,Nr) nonlocal scalar transport
1097 c sigma(imt) normalized depth (d / hbl)
1098 _KPP_RL dkm1 (imt,mdiff)
1099 _KPP_RL blmc (imt,Nr,mdiff)
1100 _KPP_RL ghat (imt,Nr)
1101 _KPP_RL sigma(imt)
1102 integer ikppkey, kkppkey
1103
1104 #ifdef ALLOW_KPP
1105
1106 c local
1107 c gat1*(imt) shape function at sigma = 1
1108 c dat1*(imt) derivative of shape function at sigma = 1
1109 c ws(imt), wm(imt) turbulent velocity scales (m/s)
1110 _KPP_RL gat1m(imt), gat1s(imt), gat1t(imt)
1111 _KPP_RL dat1m(imt), dat1s(imt), dat1t(imt)
1112 _KPP_RL ws(imt), wm(imt)
1113 integer i, kn, ki
1114 _KPP_RL R, dvdzup, dvdzdn, viscp
1115 _KPP_RL difsp, diftp, visch, difsh, difth
1116 _KPP_RL f1, sig, a1, a2, a3, delhat
1117 _KPP_RL Gm, Gs, Gt
1118 _KPP_RL tempVar
1119
1120 _KPP_RL p0 , eins
1121 parameter (p0=0.0, eins=1.0)
1122
1123 c-----------------------------------------------------------------------
1124 c compute velocity scales at hbl
1125 c-----------------------------------------------------------------------
1126
1127 do i = 1, imt
1128 sigma(i) = stable(i) * 1.0 + (1. - stable(i)) * epsilon
1129 end do
1130
1131 CADJ STORE sigma = comlev1_kpp, key = ikppkey
1132 call wscale (
1133 I sigma, hbl, ustar, bfsfc,
1134 O wm, ws )
1135 CADJ STORE wm = comlev1_kpp, key = ikppkey
1136 CADJ STORE ws = comlev1_kpp, key = ikppkey
1137
1138 do i = 1, imt
1139 wm(i) = sign(eins,wm(i))*max(phepsi,abs(wm(i)))
1140 ws(i) = sign(eins,ws(i))*max(phepsi,abs(ws(i)))
1141 end do
1142 CADJ STORE wm = comlev1_kpp, key = ikppkey
1143 CADJ STORE ws = comlev1_kpp, key = ikppkey
1144
1145 do i = 1, imt
1146
1147 kn = int(caseA(i)+phepsi) *(kbl(i) -1) +
1148 $ (1 - int(caseA(i)+phepsi)) * kbl(i)
1149
1150 c-----------------------------------------------------------------------
1151 c find the interior viscosities and derivatives at hbl(i)
1152 c-----------------------------------------------------------------------
1153
1154 delhat = 0.5*hwide(kn) - zgrid(kn) - hbl(i)
1155 R = 1.0 - delhat / hwide(kn)
1156 dvdzup = (diffus(i,kn-1,1) - diffus(i,kn ,1)) / hwide(kn)
1157 dvdzdn = (diffus(i,kn ,1) - diffus(i,kn+1,1)) / hwide(kn+1)
1158 viscp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1159 1 R * (dvdzdn + abs(dvdzdn)) )
1160
1161 dvdzup = (diffus(i,kn-1,2) - diffus(i,kn ,2)) / hwide(kn)
1162 dvdzdn = (diffus(i,kn ,2) - diffus(i,kn+1,2)) / hwide(kn+1)
1163 difsp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1164 1 R * (dvdzdn + abs(dvdzdn)) )
1165
1166 dvdzup = (diffus(i,kn-1,3) - diffus(i,kn ,3)) / hwide(kn)
1167 dvdzdn = (diffus(i,kn ,3) - diffus(i,kn+1,3)) / hwide(kn+1)
1168 diftp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1169 1 R * (dvdzdn + abs(dvdzdn)) )
1170
1171 visch = diffus(i,kn,1) + viscp * delhat
1172 difsh = diffus(i,kn,2) + difsp * delhat
1173 difth = diffus(i,kn,3) + diftp * delhat
1174
1175 f1 = stable(i) * conc1 * bfsfc(i) /
1176 & max(ustar(i)**4,phepsi)
1177 gat1m(i) = visch / hbl(i) / wm(i)
1178 dat1m(i) = -viscp / wm(i) + f1 * visch
1179
1180 gat1s(i) = difsh / hbl(i) / ws(i)
1181 dat1s(i) = -difsp / ws(i) + f1 * difsh
1182
1183 gat1t(i) = difth / hbl(i) / ws(i)
1184 dat1t(i) = -diftp / ws(i) + f1 * difth
1185
1186 end do
1187 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1188 CADJ STORE gat1m = comlev1_kpp, key = ikppkey
1189 CADJ STORE gat1s = comlev1_kpp, key = ikppkey
1190 CADJ STORE gat1t = comlev1_kpp, key = ikppkey
1191 CADJ STORE dat1m = comlev1_kpp, key = ikppkey
1192 CADJ STORE dat1s = comlev1_kpp, key = ikppkey
1193 CADJ STORE dat1t = comlev1_kpp, key = ikppkey
1194 #endif
1195 do i = 1, imt
1196 dat1m(i) = min(dat1m(i),p0)
1197 dat1s(i) = min(dat1s(i),p0)
1198 dat1t(i) = min(dat1t(i),p0)
1199 end do
1200 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1201 CADJ STORE dat1m = comlev1_kpp, key = ikppkey
1202 CADJ STORE dat1s = comlev1_kpp, key = ikppkey
1203 CADJ STORE dat1t = comlev1_kpp, key = ikppkey
1204 #endif
1205
1206 do ki = 1, Nr
1207
1208 #ifdef ALLOW_AUTODIFF_TAMC
1209 kkppkey = (ikppkey-1)*Nr + ki
1210 #endif
1211
1212 c-----------------------------------------------------------------------
1213 c compute turbulent velocity scales on the interfaces
1214 c-----------------------------------------------------------------------
1215
1216 do i = 1, imt
1217 sig = (-zgrid(ki) + 0.5 * hwide(ki)) / hbl(i)
1218 sigma(i) = stable(i)*sig + (1.-stable(i))*min(sig,epsilon)
1219 end do
1220 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1221 CADJ STORE wm = comlev1_kpp_k, key = kkppkey
1222 CADJ STORE ws = comlev1_kpp_k, key = kkppkey
1223 #endif
1224 CADJ STORE sigma = comlev1_kpp_k, key = kkppkey
1225 call wscale (
1226 I sigma, hbl, ustar, bfsfc,
1227 O wm, ws )
1228 CADJ STORE wm = comlev1_kpp_k, key = kkppkey
1229 CADJ STORE ws = comlev1_kpp_k, key = kkppkey
1230
1231 c-----------------------------------------------------------------------
1232 c compute the dimensionless shape functions at the interfaces
1233 c-----------------------------------------------------------------------
1234
1235 do i = 1, imt
1236 sig = (-zgrid(ki) + 0.5 * hwide(ki)) / hbl(i)
1237 a1 = sig - 2.
1238 a2 = 3. - 2. * sig
1239 a3 = sig - 1.
1240
1241 Gm = a1 + a2 * gat1m(i) + a3 * dat1m(i)
1242 Gs = a1 + a2 * gat1s(i) + a3 * dat1s(i)
1243 Gt = a1 + a2 * gat1t(i) + a3 * dat1t(i)
1244
1245 c-----------------------------------------------------------------------
1246 c compute boundary layer diffusivities at the interfaces
1247 c-----------------------------------------------------------------------
1248
1249 blmc(i,ki,1) = hbl(i) * wm(i) * sig * (1. + sig * Gm)
1250 blmc(i,ki,2) = hbl(i) * ws(i) * sig * (1. + sig * Gs)
1251 blmc(i,ki,3) = hbl(i) * ws(i) * sig * (1. + sig * Gt)
1252
1253 c-----------------------------------------------------------------------
1254 c nonlocal transport term = ghat * <ws>o
1255 c-----------------------------------------------------------------------
1256
1257 tempVar = ws(i) * hbl(i)
1258 ghat(i,ki) = (1.-stable(i)) * cg / max(phepsi,tempVar)
1259
1260 end do
1261 end do
1262
1263 c-----------------------------------------------------------------------
1264 c find diffusivities at kbl-1 grid level
1265 c-----------------------------------------------------------------------
1266
1267 do i = 1, imt
1268 sig = -zgrid(kbl(i)-1) / hbl(i)
1269 sigma(i) = stable(i) * sig
1270 & + (1. - stable(i)) * min(sig,epsilon)
1271 end do
1272
1273 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1274 CADJ STORE wm = comlev1_kpp, key = ikppkey
1275 CADJ STORE ws = comlev1_kpp, key = ikppkey
1276 #endif
1277 CADJ STORE sigma = comlev1_kpp, key = ikppkey
1278 call wscale (
1279 I sigma, hbl, ustar, bfsfc,
1280 O wm, ws )
1281 CADJ STORE wm = comlev1_kpp, key = ikppkey
1282 CADJ STORE ws = comlev1_kpp, key = ikppkey
1283
1284 do i = 1, imt
1285 sig = -zgrid(kbl(i)-1) / hbl(i)
1286 a1 = sig - 2.
1287 a2 = 3. - 2. * sig
1288 a3 = sig - 1.
1289 Gm = a1 + a2 * gat1m(i) + a3 * dat1m(i)
1290 Gs = a1 + a2 * gat1s(i) + a3 * dat1s(i)
1291 Gt = a1 + a2 * gat1t(i) + a3 * dat1t(i)
1292 dkm1(i,1) = hbl(i) * wm(i) * sig * (1. + sig * Gm)
1293 dkm1(i,2) = hbl(i) * ws(i) * sig * (1. + sig * Gs)
1294 dkm1(i,3) = hbl(i) * ws(i) * sig * (1. + sig * Gt)
1295 end do
1296
1297 #endif /* ALLOW_KPP */
1298
1299 return
1300 end
1301
1302 c*************************************************************************
1303
1304 subroutine enhance (
1305 I dkm1, hbl, kbl, diffus, casea
1306 U , ghat
1307 O , blmc
1308 & )
1309
1310 c enhance the diffusivity at the kbl-.5 interface
1311
1312 IMPLICIT NONE
1313
1314 #include "SIZE.h"
1315 #include "KPP_PARAMS.h"
1316
1317 c input
1318 c dkm1(imt,mdiff) bl diffusivity at kbl-1 grid level
1319 c hbl(imt) boundary layer depth (m)
1320 c kbl(imt) grid above hbl
1321 c diffus(imt,0:Nrp1,mdiff) vertical diffusivities (m^2/s)
1322 c casea(imt) = 1 in caseA, = 0 in case B
1323 _KPP_RL dkm1 (imt,mdiff)
1324 _KPP_RL hbl (imt)
1325 integer kbl (imt)
1326 _KPP_RL diffus(imt,0:Nrp1,mdiff)
1327 _KPP_RL casea (imt)
1328
1329 c input/output
1330 c nonlocal transport, modified ghat at kbl(i)-1 interface (s/m**2)
1331 _KPP_RL ghat (imt,Nr)
1332
1333 c output
1334 c enhanced bound. layer mixing coeff.
1335 _KPP_RL blmc (imt,Nr,mdiff)
1336
1337 #ifdef ALLOW_KPP
1338
1339 c local
1340 c fraction hbl lies beteen zgrid neighbors
1341 _KPP_RL delta
1342 integer ki, i, md
1343 _KPP_RL dkmp5, dstar
1344
1345 do i = 1, imt
1346 ki = kbl(i)-1
1347 if ((ki .ge. 1) .and. (ki .lt. Nr)) then
1348 delta = (hbl(i) + zgrid(ki)) / (zgrid(ki) - zgrid(ki+1))
1349 do md = 1, mdiff
1350 dkmp5 = casea(i) * diffus(i,ki,md) +
1351 1 (1.- casea(i)) * blmc (i,ki,md)
1352 dstar = (1.- delta)**2 * dkm1(i,md)
1353 & + delta**2 * dkmp5
1354 blmc(i,ki,md) = (1.- delta)*diffus(i,ki,md)
1355 & + delta*dstar
1356 end do
1357 ghat(i,ki) = (1.- casea(i)) * ghat(i,ki)
1358 endif
1359 end do
1360
1361 #endif /* ALLOW_KPP */
1362
1363 return
1364 end
1365
1366 c*************************************************************************
1367
1368 SUBROUTINE STATEKPP (
1369 I ikppkey, bi, bj, myThid,
1370 O RHO1, DBLOC, DBSFC, TTALPHA, SSBETA)
1371 c
1372 c-----------------------------------------------------------------------
1373 c "statekpp" computes all necessary input arrays
1374 c for the kpp mixing scheme
1375 c
1376 c input:
1377 c bi, bj = array indices on which to apply calculations
1378 c
1379 c output:
1380 c rho1 = potential density of surface layer (kg/m^3)
1381 c dbloc = local buoyancy gradient at Nr interfaces
1382 c g/rho{k+1,k+1} * [ drho{k,k+1}-drho{k+1,k+1} ] (m/s^2)
1383 c dbsfc = buoyancy difference with respect to the surface
1384 c g * [ drho{1,k}/rho{1,k} - drho{k,k}/rho{k,k} ] (m/s^2)
1385 c ttalpha= thermal expansion coefficient without 1/rho factor
1386 c d(rho) / d(potential temperature) (kg/m^3/C)
1387 c ssbeta = salt expansion coefficient without 1/rho factor
1388 c d(rho) / d(salinity) (kg/m^3/PSU)
1389 c
1390 c see also subroutines find_rho.F find_alpha.F find_beta.F
1391 c
1392 c written by: jan morzel, feb. 10, 1995 (converted from "sigma" version)
1393 c modified by: d. menemenlis, june 1998 : for use with MIT GCM UV
1394 c
1395 c-----------------------------------------------------------------------
1396
1397 IMPLICIT NONE
1398
1399 #include "SIZE.h"
1400 #include "EEPARAMS.h"
1401 #include "PARAMS.h"
1402 #include "KPP_PARAMS.h"
1403 #include "DYNVARS.h"
1404
1405 c-------------- Routine arguments -----------------------------------------
1406 INTEGER bi, bj, myThid
1407 _KPP_RL RHO1 ( ibot:itop, jbot:jtop )
1408 _KPP_RL DBLOC ( ibot:itop, jbot:jtop, Nr )
1409 _KPP_RL DBSFC ( ibot:itop, jbot:jtop, Nr )
1410 _KPP_RL TTALPHA( ibot:itop, jbot:jtop, Nrp1 )
1411 _KPP_RL SSBETA ( ibot:itop, jbot:jtop, Nrp1 )
1412
1413 #ifdef ALLOW_KPP
1414
1415 c--------------------------------------------------------------------------
1416 c
1417 c local arrays:
1418 c
1419 c rhok - density of t(k ) & s(k ) at depth k
1420 c rhokm1 - density of t(k-1) & s(k-1) at depth k
1421 c rho1k - density of t(1 ) & s(1 ) at depth k
1422 c work1, work2 - work arrays for holding horizontal slabs
1423
1424 _RL RHOK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1425 _RL RHOKM1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1426 _RL RHO1K (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1427 _RL WORK1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1428 _RL WORK2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1429 _RL WORK3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1430 INTEGER I, J, K
1431 INTEGER ikppkey, kkppkey
1432
1433 c calculate density, alpha, beta in surface layer, and set dbsfc to zero
1434
1435 kkppkey = (ikppkey-1)*Nr + 1
1436
1437 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1438 CADJ STORE theta(:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1439 CADJ STORE salt (:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1440 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1441 call FIND_RHO(
1442 I bi, bj, ibot, itop, jbot, jtop, 1, 1,
1443 I theta, salt,
1444 O WORK1,
1445 I myThid )
1446 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1447 CADJ STORE theta(:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1448 CADJ STORE salt (:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1449 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1450
1451 call FIND_ALPHA(
1452 I bi, bj, ibot, itop, jbot, jtop, 1, 1,
1453 O WORK2 )
1454
1455 call FIND_BETA(
1456 I bi, bj, ibot, itop, jbot, jtop, 1, 1,
1457 O WORK3 )
1458
1459 DO J = jbot, jtop
1460 DO I = ibot, itop
1461 RHO1(I,J) = WORK1(I,J) + rhoConst
1462 TTALPHA(I,J,1) = WORK2(I,J)
1463 SSBETA(I,J,1) = WORK3(I,J)
1464 DBSFC(I,J,1) = 0.
1465 END DO
1466 END DO
1467
1468 c calculate alpha, beta, and gradients in interior layers
1469
1470 CHPF$ INDEPENDENT, NEW (RHOK,RHOKM1,RHO1K,WORK1,WORK2)
1471 DO K = 2, Nr
1472
1473 kkppkey = (ikppkey-1)*Nr + k
1474
1475 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1476 CADJ STORE theta(:,:,k,bi,bj) = comlev1_kpp_k, key=kkppkey
1477 CADJ STORE salt (:,:,k,bi,bj) = comlev1_kpp_k, key=kkppkey
1478 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1479 call FIND_RHO(
1480 I bi, bj, ibot, itop, jbot, jtop, K, K,
1481 I theta, salt,
1482 O RHOK,
1483 I myThid )
1484
1485 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1486 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_kpp_k, key=kkppkey
1487 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_kpp_k, key=kkppkey
1488 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1489 call FIND_RHO(
1490 I bi, bj, ibot, itop, jbot, jtop, K-1, K,
1491 I theta, salt,
1492 O RHOKM1,
1493 I myThid )
1494
1495 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1496 CADJ STORE theta(:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1497 CADJ STORE salt (:,:,1,bi,bj) = comlev1_kpp_k, key=kkppkey
1498 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1499 call FIND_RHO(
1500 I bi, bj, ibot, itop, jbot, jtop, 1, K,
1501 I theta, salt,
1502 O RHO1K,
1503 I myThid )
1504
1505 #ifdef KPP_AUTODIFF_EXCESSIVE_STORE
1506 CADJ STORE rhok (:,:) = comlev1_kpp_k, key=kkppkey
1507 CADJ STORE rhokm1(:,:) = comlev1_kpp_k, key=kkppkey
1508 CADJ STORE rho1k (:,:) = comlev1_kpp_k, key=kkppkey
1509 #endif /* KPP_AUTODIFF_EXCESSIVE_STORE */
1510
1511 call FIND_ALPHA(
1512 I bi, bj, ibot, itop, jbot, jtop, K, K,
1513 O WORK1 )
1514
1515 call FIND_BETA(
1516 I bi, bj, ibot, itop, jbot, jtop, K, K,
1517 O WORK2 )
1518
1519 DO J = jbot, jtop
1520 DO I = ibot, itop
1521 TTALPHA(I,J,K) = WORK1 (I,J)
1522 SSBETA(I,J,K) = WORK2 (I,J)
1523 DBLOC(I,J,K-1) = gravity * (RHOK(I,J) - RHOKM1(I,J)) /
1524 & (RHOK(I,J) + rhoConst)
1525 DBSFC(I,J,K) = gravity * (RHOK(I,J) - RHO1K (I,J)) /
1526 & (RHOK(I,J) + rhoConst)
1527 END DO
1528 END DO
1529
1530 END DO
1531
1532 c compute arrays for K = Nrp1
1533 DO J = jbot, jtop
1534 DO I = ibot, itop
1535 TTALPHA(I,J,Nrp1) = TTALPHA(I,J,Nr)
1536 SSBETA(I,J,Nrp1) = SSBETA(I,J,Nr)
1537 DBLOC(I,J,Nr) = 0.
1538 END DO
1539 END DO
1540
1541 #endif /* ALLOW_KPP */
1542
1543 RETURN
1544 END

  ViewVC Help
Powered by ViewVC 1.1.22