/[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.3 - (show annotations) (download)
Wed Sep 13 17:07:11 2000 UTC (23 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint31
Changes since 1.2: +0 -0 lines
KPP package without flag kppPackageIsOn.

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

  ViewVC Help
Powered by ViewVC 1.1.22