/[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.5 - (show annotations) (download)
Mon Jan 29 20:09:23 2001 UTC (23 years, 5 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint34
Changes since 1.4: +13 -7 lines
Replaced some storage for TAMC by cheap  recomputation.

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

  ViewVC Help
Powered by ViewVC 1.1.22