/[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.1 - (show annotations) (download)
Wed Jun 21 19:45:52 2000 UTC (23 years, 11 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint29, checkpoint30
Packaged KPP mixing scheme.

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

  ViewVC Help
Powered by ViewVC 1.1.22