/[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.13 - (show annotations) (download)
Sat Dec 28 10:11:11 2002 UTC (21 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint48e_post, checkpoint48b_post, checkpoint48c_pre, checkpoint48d_pre, checkpoint47i_post, checkpoint48d_post, checkpoint47g_post, checkpoint48a_post, checkpoint47j_post, checkpoint48c_post, checkpoint47f_post, checkpoint48, checkpoint47h_post
Changes since 1.12: +3 -2 lines
checkpoint47f_post
Merging from release1_p10:
o modifications for using pkg/exf with pkg/seaice
  - pkg/seaice CPP options SEAICE_EXTERNAL_FORCING
    and SEAICE_EXTERNAL_FLUXES
  - pkg/exf CPP options EXF_READ_EVAP and
    EXF_NO_BULK_COMPUTATIONS
  - usage examples are Experiments 8 and 9 in
    verification/lab_sea/README
  - verification/lab_sea default experiment now uses
    pkg/gmredi, pkg/kpp, pkg/seaice, and pkg/exf

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

  ViewVC Help
Powered by ViewVC 1.1.22