/[MITgcm]/MITgcm/pkg/kpp/kpp_routines.F
ViewVC logotype

Annotation of /MITgcm/pkg/kpp/kpp_routines.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.8 - (hide annotations) (download)
Sun Mar 25 22:33:55 2001 UTC (23 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint38, checkpoint40pre2, checkpoint40pre4, c37_adj, checkpoint39, checkpoint40pre5
Changes since 1.7: +2 -2 lines
Modifications and additions to enable automatic differentiation.
Detailed info's in doc/notes_c37_adj.txt

1 heimbach 1.8 C $Header: /u/gcmpack/models/MITgcmUV/pkg/kpp/kpp_routines.F,v 1.7 2001/02/04 14:38:50 cnh Exp $
2     C $Name: checkpoint37 $
3 adcroft 1.1
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 heimbach 1.4 C-- o KPP_SMOOTH_HORIZ - Apply horizontal smoothing to KPP array.
15     C-- o SMOOTH_HORIZ - Apply horizontal smoothing to global array.
16 adcroft 1.1 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 heimbach 1.4 c***********************************************************************
21 adcroft 1.1
22     SUBROUTINE KPPMIX (
23 heimbach 1.4 I mytime, mythid
24     I , kmtj, shsq, dvsq, ustar
25     I , bo, bosol, dbloc, Ritop, coriol
26 adcroft 1.1 I , ikey
27     O , diffus
28     U , ghat
29     O , hbl )
30    
31 heimbach 1.4 c-----------------------------------------------------------------------
32 adcroft 1.1 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 heimbach 1.4 c myTime - current time in simulation
56     c myThid - thread number for this instance of the routine
57 adcroft 1.1 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 heimbach 1.4 _RL mytime
74     integer mythid
75 adcroft 1.1 integer kmtj (imt )
76 heimbach 1.4 _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 adcroft 1.1
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 heimbach 1.4 _KPP_RL diffus(imt,0:Nrp1,mdiff)
95     _KPP_RL ghat (imt,Nr)
96     _KPP_RL hbl (imt)
97 adcroft 1.1
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 heimbach 1.4 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 adcroft 1.1
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 heimbach 1.4 CADJ STORE ghat = comlev1_kpp, key = ikey
129    
130 adcroft 1.1 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 heimbach 1.4 I mytime, mythid
156     I , kmtj
157 adcroft 1.1 I , dvsq, dbloc, Ritop, ustar, bo, bosol, coriol
158     I , ikey
159     O , hbl, bfsfc, stable, casea, kbl, Rib, sigma
160     & )
161    
162 heimbach 1.2 CADJ STORE hbl,bfsfc,stable,casea,kbl = comlev1_kpp, key = ikey
163 adcroft 1.1
164     c-----------------------------------------------------------------------
165     c compute boundary layer diffusivities
166     c-----------------------------------------------------------------------
167    
168     call blmix (
169     I ustar, bfsfc, hbl, stable, casea, diffus, kbl
170 heimbach 1.4 O , dkm1, blmc, ghat, sigma, ikey
171 adcroft 1.1 & )
172    
173 heimbach 1.2 CADJ STORE dkm1,blmc,ghat = comlev1_kpp, key = ikey
174 adcroft 1.1
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 heimbach 1.4 I mytime, mythid
209     I , kmtj
210 adcroft 1.1 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 heimbach 1.4 c myTime : current time in simulation
247     c myThid : thread number for this instance of the routine
248 adcroft 1.1 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 heimbach 1.4 _RL mytime
259     integer mythid
260 adcroft 1.1 integer kmtj(imt)
261 heimbach 1.4 _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 adcroft 1.1 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 heimbach 1.4 _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 adcroft 1.1
287     #ifdef ALLOW_KPP
288    
289     c local
290     c-------
291     c wm, ws : turbulent velocity scales (m/s)
292 heimbach 1.4 _KPP_RL wm(imt), ws(imt)
293     _RL worka(imt)
294 adcroft 1.1
295 heimbach 1.4 _KPP_RL bvsq, vtsq, hekman, hmonob, hlimit, tempVar1, tempVar2
296 adcroft 1.1 integer i, kl
297    
298 heimbach 1.4 _KPP_RL p5 , eins
299     parameter ( p5=0.5, eins=1.0 )
300     _RL minusone
301     parameter ( minusone=-1.0 )
302 adcroft 1.1
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 heimbach 1.4 worka(i) = zgrid(kl)
326 adcroft 1.1 end do
327     call SWFRAC(
328 heimbach 1.4 I imt, hbf,
329     I mytime, mythid,
330     U worka )
331 adcroft 1.1
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 heimbach 1.4 bfsfc(i) = bo(i) + bosol(i)*(1. - worka(i))
341 adcroft 1.1 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 heimbach 1.4
382     tempVar1 = dvsq(i,kl) + vtsq
383     tempVar2 = max(tempVar1, phepsi)
384     Rib(i,kl) = Ritop(i,kl) / tempVar2
385    
386 adcroft 1.1 end do
387     end do
388    
389     do kl = 2, Nr
390     do i = 1, imt
391     if (kbl(i).eq.kmtj(i) .and. Rib(i,kl).gt.Ricr) kbl(i) = kl
392     end do
393     end do
394    
395 heimbach 1.2 CADJ store kbl = comlev1_kpp
396 adcroft 1.1 CADJ & , key = ikey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
397    
398     do i = 1, imt
399     kl = kbl(i)
400     c linearly interpolate to find hbl where Rib = Ricr
401     if (kl.gt.1 .and. kl.lt.kmtj(i)) then
402 heimbach 1.4 tempVar1 = (Rib(i,kl)-Rib(i,kl-1))
403 adcroft 1.1 hbl(i) = -zgrid(kl-1) + (zgrid(kl-1)-zgrid(kl)) *
404 heimbach 1.4 1 (Ricr - Rib(i,kl-1)) / tempVar1
405 adcroft 1.1 endif
406     end do
407    
408 heimbach 1.2 CADJ store hbl = comlev1_kpp
409 adcroft 1.1 CADJ & , key = ikey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
410    
411     c-----------------------------------------------------------------------
412     c find stability and buoyancy forcing for boundary layer
413     c-----------------------------------------------------------------------
414    
415 heimbach 1.4 do i = 1, imt
416     worka(i) = hbl(i)
417     end do
418 adcroft 1.1 call SWFRAC(
419 heimbach 1.4 I imt, minusone,
420     I mytime, mythid,
421     U worka )
422 adcroft 1.1
423 heimbach 1.4 do i = 1, imt
424     bfsfc(i) = bo(i) + bosol(i) * (1. - worka(i))
425     end do
426 heimbach 1.2 CADJ store bfsfc = comlev1_kpp
427     CADJ & , key = ikey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
428    
429 heimbach 1.4 c-- ensure bfsfc is never 0
430 adcroft 1.1 do i = 1, imt
431     stable(i) = p5 + sign( p5, bfsfc(i) )
432 heimbach 1.4 bfsfc(i) = sign(eins,bfsfc(i))*max(phepsi,abs(bfsfc(i)))
433     end do
434 adcroft 1.1
435 heimbach 1.4 CADJ store bfsfc = comlev1_kpp
436     CADJ & , key = ikey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
437 adcroft 1.1
438     c-----------------------------------------------------------------------
439     c check hbl limits for hekman or hmonob
440     c ph: test for zero nominator
441     c-----------------------------------------------------------------------
442    
443     do i = 1, imt
444     if (bfsfc(i) .gt. 0.0) then
445     hekman = cekman * ustar(i) / max(abs(Coriol(i)),phepsi)
446     hmonob = cmonob * ustar(i)*ustar(i)*ustar(i)
447     & / vonk / bfsfc(i)
448     hlimit = stable(i) * min(hekman,hmonob)
449     & + (stable(i)-1.) * zgrid(Nr)
450     hbl(i) = min(hbl(i),hlimit)
451 heimbach 1.4 end if
452     end do
453     CADJ store hbl = comlev1_kpp
454     CADJ & , key = ikey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
455    
456     do i = 1, imt
457     hbl(i) = max(hbl(i),minKPPhbl)
458 adcroft 1.1 kbl(i) = kmtj(i)
459     end do
460    
461 heimbach 1.2 CADJ store hbl = comlev1_kpp
462 adcroft 1.1 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 heimbach 1.4 do i = 1, imt
481     worka(i) = hbl(i)
482     end do
483 adcroft 1.1 call SWFRAC(
484 heimbach 1.4 I imt, minusone,
485     I mytime, mythid,
486     U worka )
487 adcroft 1.1
488 heimbach 1.4 do i = 1, imt
489     bfsfc(i) = bo(i) + bosol(i) * (1. - worka(i))
490     end do
491 heimbach 1.2 CADJ store bfsfc = comlev1_kpp
492     CADJ & , key = ikey, shape = (/ (sNx+2*OLx)*(sNy+2*OLy) /)
493    
494 heimbach 1.4 c-- ensures bfsfc is never 0
495 adcroft 1.1 do i = 1, imt
496     stable(i) = p5 + sign( p5, bfsfc(i) )
497 heimbach 1.4 bfsfc(i) = sign(eins,bfsfc(i))*max(phepsi,abs(bfsfc(i)))
498 adcroft 1.1 end do
499    
500     c-----------------------------------------------------------------------
501     c determine caseA and caseB
502     c-----------------------------------------------------------------------
503    
504     do i = 1, imt
505     casea(i) = p5 +
506     1 sign(p5, -zgrid(kbl(i)) - p5*hwide(kbl(i)) - hbl(i))
507     end do
508    
509     #endif /* ALLOW_KPP */
510    
511     return
512     end
513    
514     c*************************************************************************
515    
516     subroutine wscale (
517     I sigma, hbl, ustar, bfsfc,
518     O wm, ws )
519    
520     c compute turbulent velocity scales.
521     c use a 2D-lookup table for wm and ws as functions of ustar and
522     c zetahat (=vonk*sigma*hbl*bfsfc).
523     c
524     c note: the lookup table is only used for unstable conditions
525     c (zehat.le.0), in the stable domain wm (=ws) gets computed
526     c directly.
527     c
528     IMPLICIT NONE
529    
530     #include "SIZE.h"
531     #include "KPP_PARAMS.h"
532    
533     c input
534     c------
535     c sigma : normalized depth (d/hbl)
536     c hbl : boundary layer depth (m)
537     c ustar : surface friction velocity (m/s)
538     c bfsfc : total surface buoyancy flux (m^2/s^3)
539 heimbach 1.4 _KPP_RL sigma(imt)
540     _KPP_RL hbl (imt)
541     _KPP_RL ustar(imt)
542     _KPP_RL bfsfc(imt)
543 adcroft 1.1
544     c output
545     c--------
546     c wm, ws : turbulent velocity scales at sigma
547 heimbach 1.4 _KPP_RL wm(imt), ws(imt)
548 adcroft 1.1
549     #ifdef ALLOW_KPP
550    
551     c local
552     c------
553     c zehat : = zeta * ustar**3
554 heimbach 1.4 _KPP_RL zehat
555 adcroft 1.1
556     integer iz, izp1, ju, i, jup1
557 heimbach 1.4 _KPP_RL udiff, zdiff, zfrac, ufrac, fzfrac, wam
558     _KPP_RL wbm, was, wbs, u3, tempVar
559 adcroft 1.1
560     c-----------------------------------------------------------------------
561     c use lookup table for zehat < zmax only; otherwise use
562     c stable formulae
563     c-----------------------------------------------------------------------
564    
565     do i = 1, imt
566     zehat = vonk*sigma(i)*hbl(i)*bfsfc(i)
567    
568     if (zehat .le. zmax) then
569    
570     zdiff = zehat - zmin
571     iz = int( zdiff / deltaz )
572     iz = min( iz, nni )
573     iz = max( iz, 0 )
574     izp1 = iz + 1
575 heimbach 1.4
576 adcroft 1.1 udiff = ustar(i) - umin
577     ju = int( udiff / deltau )
578     ju = min( ju, nnj )
579     ju = max( ju, 0 )
580     jup1 = ju + 1
581 heimbach 1.4
582 adcroft 1.1 zfrac = zdiff / deltaz - float(iz)
583     ufrac = udiff / deltau - float(ju)
584 heimbach 1.4
585 adcroft 1.1 fzfrac= 1. - zfrac
586     wam = fzfrac * wmt(iz,jup1) + zfrac * wmt(izp1,jup1)
587     wbm = fzfrac * wmt(iz,ju ) + zfrac * wmt(izp1,ju )
588     wm(i) = (1.-ufrac) * wbm + ufrac * wam
589 heimbach 1.4
590 adcroft 1.1 was = fzfrac * wst(iz,jup1) + zfrac * wst(izp1,jup1)
591     wbs = fzfrac * wst(iz,ju ) + zfrac * wst(izp1,ju )
592     ws(i) = (1.-ufrac) * wbs + ufrac * was
593 heimbach 1.4
594 adcroft 1.1 else
595 heimbach 1.4
596     u3 = ustar(i) * ustar(i) * ustar(i)
597     tempVar = u3 + conc1 * zehat
598     wm(i) = vonk * ustar(i) * u3 / tempVar
599 adcroft 1.1 ws(i) = wm(i)
600 heimbach 1.4
601 adcroft 1.1 endif
602    
603     end do
604    
605     #endif /* ALLOW_KPP */
606    
607     return
608     end
609    
610     c*************************************************************************
611    
612     subroutine Ri_iwmix (
613     I kmtj, shsq, dbloc, dblocSm
614     I , ikey
615     O , diffus )
616    
617     c compute interior viscosity diffusivity coefficients due
618     c to shear instability (dependent on a local Richardson number),
619     c to background internal wave activity, and
620     c to static instability (local Richardson number < 0).
621    
622     IMPLICIT NONE
623    
624     #include "SIZE.h"
625     #include "EEPARAMS.h"
626     #include "PARAMS.h"
627     #include "KPP_PARAMS.h"
628    
629     c input
630     c kmtj (imt) number of vertical layers on this row
631     c shsq (imt,Nr) (local velocity shear)^2 ((m/s)^2)
632     c dbloc (imt,Nr) local delta buoyancy (m/s^2)
633     c dblocSm(imt,Nr) horizontally smoothed dbloc (m/s^2)
634     integer kmtj (imt)
635 heimbach 1.4 _KPP_RL shsq (imt,Nr)
636     _KPP_RL dbloc (imt,Nr)
637     _KPP_RL dblocSm(imt,Nr)
638 adcroft 1.1 integer ikey
639    
640     c output
641     c diffus(imt,0:Nrp1,1) vertical viscosivity coefficient (m^2/s)
642     c diffus(imt,0:Nrp1,2) vertical scalar diffusivity (m^2/s)
643     c diffus(imt,0:Nrp1,3) vertical temperature diffusivity (m^2/s)
644 heimbach 1.4 _KPP_RL diffus(imt,0:Nrp1,3)
645 adcroft 1.1
646     #ifdef ALLOW_KPP
647    
648     c local variables
649     c Rig local Richardson number
650     c fRi, fcon function of Rig
651 heimbach 1.4 _KPP_RL Rig
652     _KPP_RL fRi, fcon
653     _KPP_RL ratio
654 adcroft 1.1 integer i, ki, mr
655 heimbach 1.4 _KPP_RL c1, c0
656 adcroft 1.1
657 heimbach 1.5 #ifdef ALLOW_KPP_VERTICALLY_SMOOTH
658     CADJ INIT kpp_ri_tape_mr = common, 1
659     #endif
660    
661 adcroft 1.1 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 heimbach 1.5 #ifdef ALLOW_AUTODIFF_TAMC
671     C break data flow dependence on diffus
672     diffus(1,1,1) = 0.0
673     #endif
674 adcroft 1.1
675     do ki = 1, Nr
676     do i = 1, imt
677     if (kmtj(i) .EQ. 0 ) then
678     diffus(i,ki,1) = 0.
679     diffus(i,ki,2) = 0.
680     elseif (ki .GE. kmtj(i)) then
681     diffus(i,ki,1) = diffus(i,ki-1,1)
682     diffus(i,ki,2) = diffus(i,ki-1,2)
683     else
684     diffus(i,ki,1) = dblocSm(i,ki) * (zgrid(ki)-zgrid(ki+1))
685     & / max( Shsq(i,ki), phepsi )
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 heimbach 1.5 #ifdef ALLOW_KPP_VERTICALLY_SMOOTH
694 adcroft 1.1 do mr = 1, num_v_smooth_Ri
695 heimbach 1.2
696 heimbach 1.5 CADJ store diffus(:,:,1) = kpp_ri_tape_mr
697     CADJ & , key=mr, shape=(/ (sNx+2*OLx)*(sNy+2*OLy),Nr+2 /)
698 heimbach 1.2
699 adcroft 1.1 call z121 (
700     U diffus(1,0,1))
701     end do
702 heimbach 1.5 #endif
703 adcroft 1.1
704     c-----------------------------------------------------------------------
705     c after smoothing loop
706    
707     do ki = 1, Nr
708     do i = 1, imt
709    
710     c evaluate f of Brunt-Vaisala squared for convection, store in fcon
711    
712     Rig = max ( diffus(i,ki,2) , BVSQcon )
713     ratio = min ( (BVSQcon - Rig) / BVSQcon, c1 )
714     fcon = c1 - ratio * ratio
715     fcon = fcon * fcon * fcon
716    
717     c evaluate f of smooth Ri for shear instability, store in fRi
718    
719     Rig = max ( diffus(i,ki,1), c0 )
720     ratio = min ( Rig / Riinfty , c1 )
721     fRi = c1 - ratio * ratio
722     fRi = fRi * fRi * fRi
723    
724     c ----------------------------------------------------------------------
725     c evaluate diffusivities and viscosity
726     c mixing due to internal waves, and shear and static instability
727    
728     diffus(i,ki,1) = viscAr + fcon * difmcon + fRi * difm0
729     diffus(i,ki,2) = diffKrS + fcon * difscon + fRi * difs0
730 heimbach 1.2 diffus(i,ki,3) = diffKrT + fcon * difscon + fRi * difs0
731 adcroft 1.1
732     end do
733     end do
734    
735     c ------------------------------------------------------------------------
736     c set surface values to 0.0
737    
738     do i = 1, imt
739     diffus(i,0,1) = c0
740     diffus(i,0,2) = c0
741     diffus(i,0,3) = c0
742     end do
743    
744     #endif /* ALLOW_KPP */
745    
746     return
747     end
748    
749     c*************************************************************************
750    
751     subroutine z121 (
752     U v )
753    
754     c Apply 121 smoothing in k to 2-d array V(i,k=1,Nr)
755     c top (0) value is used as a dummy
756     c bottom (Nrp1) value is set to input value from above.
757    
758 heimbach 1.2 c Note that it is important to exclude from the smoothing any points
759 adcroft 1.1 c that are outside the range of the K(Ri) scheme, ie. >0.8, or <0.0.
760 heimbach 1.2 c Otherwise, there is interference with other physics, especially
761 adcroft 1.1 c penetrative convection.
762    
763     IMPLICIT NONE
764     #include "SIZE.h"
765     #include "KPP_PARAMS.h"
766    
767     c input/output
768     c-------------
769     c v : 2-D array to be smoothed in Nrp1 direction
770 heimbach 1.4 _KPP_RL v(imt,0:Nrp1)
771 adcroft 1.1
772     #ifdef ALLOW_KPP
773    
774     c local
775 heimbach 1.4 _KPP_RL zwork, zflag
776     _KPP_RL KRi_range(1:Nrp1)
777 adcroft 1.1 integer i, k, km1, kp1
778    
779 heimbach 1.4 _KPP_RL p0 , p25 , p5 , p2
780 heimbach 1.2 parameter ( p0 = 0.0, p25 = 0.25, p5 = 0.5, p2 = 2.0 )
781    
782     KRi_range(Nrp1) = p0
783 adcroft 1.1
784 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
785     C-- dummy assignment to end declaration part for TAMC
786     i = 0
787    
788     C-- HPF directive to help TAMC
789     CHPF$ INDEPENDENT
790 heimbach 1.4 CADJ INIT z121tape = common, Nr
791 heimbach 1.2 #endif /* ALLOW_AUTODIFF_TAMC */
792 heimbach 1.4
793 adcroft 1.1 do i = 1, imt
794    
795 heimbach 1.4 k = 1
796     CADJ STORE v(i,k) = z121tape
797 adcroft 1.1 v(i,Nrp1) = v(i,Nr)
798    
799     do k = 1, Nr
800     KRi_range(k) = p5 + SIGN(p5,v(i,k))
801     KRi_range(k) = KRi_range(k) *
802     & ( p5 + SIGN(p5,(Riinfty-v(i,k))) )
803     end do
804    
805     zwork = KRi_range(1) * v(i,1)
806     v(i,1) = p2 * v(i,1) +
807     & KRi_range(1) * KRi_range(2) * v(i,2)
808     zflag = p2 + KRi_range(1) * KRi_range(2)
809     v(i,1) = v(i,1) / zflag
810    
811     do k = 2, Nr
812 heimbach 1.2 CADJ STORE v(i,k), zwork = z121tape
813 adcroft 1.1 km1 = k - 1
814     kp1 = k + 1
815     zflag = v(i,k)
816     v(i,k) = p2 * v(i,k) +
817     & KRi_range(k) * KRi_range(kp1) * v(i,kp1) +
818     & KRi_range(k) * zwork
819     zwork = KRi_range(k) * zflag
820     zflag = p2 + KRi_range(k)*(KRi_range(kp1)+KRi_range(km1))
821     v(i,k) = v(i,k) / zflag
822     end do
823    
824     end do
825    
826     #endif /* ALLOW_KPP */
827    
828     return
829     end
830    
831     c*************************************************************************
832    
833 heimbach 1.4 subroutine kpp_smooth_horiz (
834 adcroft 1.1 I k, bi, bj,
835 heimbach 1.4 U fld )
836 adcroft 1.1
837 heimbach 1.4 c Apply horizontal smoothing to KPP array
838 adcroft 1.1
839     IMPLICIT NONE
840     #include "SIZE.h"
841     #include "KPP_PARAMS.h"
842    
843     c input
844 heimbach 1.4 c bi, bj : array indices
845     c k : vertical index used for masking
846 adcroft 1.1 integer k, bi, bj
847    
848 heimbach 1.4 c input/output
849     c fld : 2-D array to be smoothed
850     _KPP_RL fld( ibot:itop, jbot:jtop )
851 adcroft 1.1
852     #ifdef ALLOW_KPP
853    
854     c local
855     integer i, j, im1, ip1, jm1, jp1
856 heimbach 1.4 _KPP_RL tempVar
857     _KPP_RL fld_tmp( ibot:itop, jbot:jtop )
858 adcroft 1.1
859 heimbach 1.4 integer imin , imax , jmin , jmax
860     parameter( imin=ibot+1, imax=itop-1, jmin=jbot+1, jmax=jtop-1 )
861 adcroft 1.1
862 heimbach 1.4 _KPP_RL p0 , p5 , p25 , p125 , p0625
863 adcroft 1.1 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
864    
865     DO j = jmin, jmax
866     jm1 = j-1
867     jp1 = j+1
868     DO i = imin, imax
869     im1 = i-1
870     ip1 = i+1
871     tempVar =
872     & p25 * pMask(i ,j ,k,bi,bj) +
873     & p125 * ( pMask(im1,j ,k,bi,bj) +
874     & pMask(ip1,j ,k,bi,bj) +
875     & pMask(i ,jm1,k,bi,bj) +
876     & pMask(i ,jp1,k,bi,bj) ) +
877     & p0625 * ( pMask(im1,jm1,k,bi,bj) +
878     & pMask(im1,jp1,k,bi,bj) +
879     & pMask(ip1,jm1,k,bi,bj) +
880     & pMask(ip1,jp1,k,bi,bj) )
881     IF ( tempVar .GE. p25 ) THEN
882     fld_tmp(i,j) = (
883 heimbach 1.4 & p25 * fld(i ,j )*pMask(i ,j ,k,bi,bj) +
884     & p125 *(fld(im1,j )*pMask(im1,j ,k,bi,bj) +
885     & fld(ip1,j )*pMask(ip1,j ,k,bi,bj) +
886     & fld(i ,jm1)*pMask(i ,jm1,k,bi,bj) +
887     & fld(i ,jp1)*pMask(i ,jp1,k,bi,bj))+
888     & p0625*(fld(im1,jm1)*pMask(im1,jm1,k,bi,bj) +
889     & fld(im1,jp1)*pMask(im1,jp1,k,bi,bj) +
890     & fld(ip1,jm1)*pMask(ip1,jm1,k,bi,bj) +
891     & fld(ip1,jp1)*pMask(ip1,jp1,k,bi,bj)))
892 adcroft 1.1 & / tempVar
893     ELSE
894 heimbach 1.4 fld_tmp(i,j) = fld(i,j)
895 adcroft 1.1 ENDIF
896     ENDDO
897     ENDDO
898    
899     c transfer smoothed field to output array
900     DO j = jmin, jmax
901     DO i = imin, imax
902 heimbach 1.4 fld(i,j) = fld_tmp(i,j)
903 adcroft 1.1 ENDDO
904     ENDDO
905    
906     #endif /* ALLOW_KPP */
907    
908     return
909     end
910    
911     c*************************************************************************
912    
913 heimbach 1.4 subroutine smooth_horiz (
914 adcroft 1.1 I k, bi, bj,
915 heimbach 1.4 U fld )
916 adcroft 1.1
917 heimbach 1.4 c Apply horizontal smoothing to global _RL 2-D array
918 adcroft 1.1
919     IMPLICIT NONE
920     #include "SIZE.h"
921     #include "KPP_PARAMS.h"
922    
923     c input
924 heimbach 1.4 c bi, bj : array indices
925     c k : vertical index used for masking
926 adcroft 1.1 integer k, bi, bj
927    
928 heimbach 1.4 c input/output
929     c fld : 2-D array to be smoothed
930     _RL fld( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
931 adcroft 1.1
932     #ifdef ALLOW_KPP
933    
934     c local
935     integer i, j, im1, ip1, jm1, jp1
936     _RL tempVar
937 heimbach 1.4 _RL fld_tmp( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
938 adcroft 1.1
939 heimbach 1.4 integer imin , imax , jmin , jmax
940     parameter(imin=2-OLx, imax=sNx+OLx-1, jmin=2-OLy, jmax=sNy+OLy-1)
941 adcroft 1.1
942 heimbach 1.4 _RL p0 , p5 , p25 , p125 , p0625
943 adcroft 1.1 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
944    
945     DO j = jmin, jmax
946     jm1 = j-1
947     jp1 = j+1
948     DO i = imin, imax
949     im1 = i-1
950     ip1 = i+1
951     tempVar =
952     & p25 * pMask(i ,j ,k,bi,bj) +
953     & p125 * ( pMask(im1,j ,k,bi,bj) +
954     & pMask(ip1,j ,k,bi,bj) +
955     & pMask(i ,jm1,k,bi,bj) +
956     & pMask(i ,jp1,k,bi,bj) ) +
957     & p0625 * ( pMask(im1,jm1,k,bi,bj) +
958     & pMask(im1,jp1,k,bi,bj) +
959     & pMask(ip1,jm1,k,bi,bj) +
960     & pMask(ip1,jp1,k,bi,bj) )
961     IF ( tempVar .GE. p25 ) THEN
962     fld_tmp(i,j) = (
963 heimbach 1.4 & p25 * fld(i ,j )*pMask(i ,j ,k,bi,bj) +
964     & p125 *(fld(im1,j )*pMask(im1,j ,k,bi,bj) +
965     & fld(ip1,j )*pMask(ip1,j ,k,bi,bj) +
966     & fld(i ,jm1)*pMask(i ,jm1,k,bi,bj) +
967     & fld(i ,jp1)*pMask(i ,jp1,k,bi,bj))+
968     & p0625*(fld(im1,jm1)*pMask(im1,jm1,k,bi,bj) +
969     & fld(im1,jp1)*pMask(im1,jp1,k,bi,bj) +
970     & fld(ip1,jm1)*pMask(ip1,jm1,k,bi,bj) +
971     & fld(ip1,jp1)*pMask(ip1,jp1,k,bi,bj)))
972 adcroft 1.1 & / tempVar
973     ELSE
974 heimbach 1.4 fld_tmp(i,j) = fld(i,j)
975 adcroft 1.1 ENDIF
976     ENDDO
977     ENDDO
978    
979     c transfer smoothed field to output array
980     DO j = jmin, jmax
981     DO i = imin, imax
982 heimbach 1.4 fld(i,j) = fld_tmp(i,j)
983 adcroft 1.1 ENDDO
984     ENDDO
985    
986     #endif /* ALLOW_KPP */
987    
988     return
989     end
990    
991     c*************************************************************************
992    
993     subroutine blmix (
994     I ustar, bfsfc, hbl, stable, casea, diffus, kbl
995 heimbach 1.4 O , dkm1, blmc, ghat, sigma, ikey
996 adcroft 1.1 & )
997    
998     c mixing coefficients within boundary layer depend on surface
999     c forcing and the magnitude and gradient of interior mixing below
1000     c the boundary layer ("matching").
1001     c
1002     c caution: if mixing bottoms out at hbl = -zgrid(Nr) then
1003     c fictitious layer at Nrp1 is needed with small but finite width
1004     c hwide(Nrp1) (eg. epsln = 1.e-20).
1005     c
1006     IMPLICIT NONE
1007    
1008     #include "SIZE.h"
1009     #include "KPP_PARAMS.h"
1010    
1011     c input
1012     c ustar (imt) surface friction velocity (m/s)
1013     c bfsfc (imt) surface buoyancy forcing (m^2/s^3)
1014     c hbl (imt) boundary layer depth (m)
1015     c stable(imt) = 1 in stable forcing
1016     c casea (imt) = 1 in case A
1017     c diffus(imt,0:Nrp1,mdiff) vertical diffusivities (m^2/s)
1018     c kbl(imt) -1 of first grid level below hbl
1019 heimbach 1.4 _KPP_RL ustar (imt)
1020     _KPP_RL bfsfc (imt)
1021     _KPP_RL hbl (imt)
1022     _KPP_RL stable(imt)
1023     _KPP_RL casea (imt)
1024     _KPP_RL diffus(imt,0:Nrp1,mdiff)
1025 adcroft 1.1 integer kbl(imt)
1026    
1027     c output
1028     c dkm1 (imt,mdiff) boundary layer difs at kbl-1 level
1029     c blmc (imt,Nr,mdiff) boundary layer mixing coefficients (m^2/s)
1030     c ghat (imt,Nr) nonlocal scalar transport
1031     c sigma(imt) normalized depth (d / hbl)
1032 heimbach 1.4 _KPP_RL dkm1 (imt,mdiff)
1033     _KPP_RL blmc (imt,Nr,mdiff)
1034     _KPP_RL ghat (imt,Nr)
1035     _KPP_RL sigma(imt)
1036     integer ikey
1037 adcroft 1.1
1038     #ifdef ALLOW_KPP
1039    
1040     c local
1041 heimbach 1.2 c gat1*(imt) shape function at sigma = 1
1042     c dat1*(imt) derivative of shape function at sigma = 1
1043 adcroft 1.1 c ws(imt), wm(imt) turbulent velocity scales (m/s)
1044 heimbach 1.4 _KPP_RL gat1m(imt), gat1s(imt), gat1t(imt)
1045     _KPP_RL dat1m(imt), dat1s(imt), dat1t(imt)
1046     _KPP_RL ws(imt), wm(imt)
1047 adcroft 1.1 integer i, kn, ki
1048 heimbach 1.4 _KPP_RL R, dvdzup, dvdzdn, viscp
1049     _KPP_RL difsp, diftp, visch, difsh, difth
1050     _KPP_RL f1, sig, a1, a2, a3, delhat
1051     _KPP_RL Gm, Gs, Gt
1052     _KPP_RL tempVar
1053 adcroft 1.1
1054 heimbach 1.4 _KPP_RL p0 , eins
1055 heimbach 1.2 parameter (p0=0.0, eins=1.0)
1056 adcroft 1.1
1057     c-----------------------------------------------------------------------
1058     c compute velocity scales at hbl
1059     c-----------------------------------------------------------------------
1060    
1061     do i = 1, imt
1062     sigma(i) = stable(i) * 1.0 + (1. - stable(i)) * epsilon
1063     end do
1064    
1065     call wscale (
1066     I sigma, hbl, ustar, bfsfc,
1067     O wm, ws )
1068    
1069     do i = 1, imt
1070 heimbach 1.4 wm(i) = sign(eins,wm(i))*max(phepsi,abs(wm(i)))
1071     ws(i) = sign(eins,ws(i))*max(phepsi,abs(ws(i)))
1072     end do
1073     CADJ STORE wm = comlev1_kpp, key = ikey
1074     CADJ STORE ws = comlev1_kpp, key = ikey
1075    
1076     do i = 1, imt
1077 adcroft 1.1
1078     kn = int(caseA(i)+phepsi) *(kbl(i) -1) +
1079     $ (1 - int(caseA(i)+phepsi)) * kbl(i)
1080    
1081     c-----------------------------------------------------------------------
1082     c find the interior viscosities and derivatives at hbl(i)
1083     c-----------------------------------------------------------------------
1084    
1085     delhat = 0.5*hwide(kn) - zgrid(kn) - hbl(i)
1086     R = 1.0 - delhat / hwide(kn)
1087     dvdzup = (diffus(i,kn-1,1) - diffus(i,kn ,1)) / hwide(kn)
1088     dvdzdn = (diffus(i,kn ,1) - diffus(i,kn+1,1)) / hwide(kn+1)
1089     viscp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1090     1 R * (dvdzdn + abs(dvdzdn)) )
1091    
1092     dvdzup = (diffus(i,kn-1,2) - diffus(i,kn ,2)) / hwide(kn)
1093     dvdzdn = (diffus(i,kn ,2) - diffus(i,kn+1,2)) / hwide(kn+1)
1094     difsp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1095     1 R * (dvdzdn + abs(dvdzdn)) )
1096    
1097     dvdzup = (diffus(i,kn-1,3) - diffus(i,kn ,3)) / hwide(kn)
1098     dvdzdn = (diffus(i,kn ,3) - diffus(i,kn+1,3)) / hwide(kn+1)
1099     diftp = 0.5 * ( (1.-R) * (dvdzup + abs(dvdzup)) +
1100     1 R * (dvdzdn + abs(dvdzdn)) )
1101    
1102     visch = diffus(i,kn,1) + viscp * delhat
1103     difsh = diffus(i,kn,2) + difsp * delhat
1104     difth = diffus(i,kn,3) + diftp * delhat
1105    
1106     f1 = stable(i) * conc1 * bfsfc(i) /
1107     & max(ustar(i)**4,phepsi)
1108 heimbach 1.2 gat1m(i) = visch / hbl(i) / wm(i)
1109     dat1m(i) = -viscp / wm(i) + f1 * visch
1110     dat1m(i) = min(dat1m(i),p0)
1111 adcroft 1.1
1112 heimbach 1.2 gat1s(i) = difsh / hbl(i) / ws(i)
1113     dat1s(i) = -difsp / ws(i) + f1 * difsh
1114     dat1s(i) = min(dat1s(i),p0)
1115 adcroft 1.1
1116 heimbach 1.2 gat1t(i) = difth / hbl(i) / ws(i)
1117     dat1t(i) = -diftp / ws(i) + f1 * difth
1118     dat1t(i) = min(dat1t(i),p0)
1119 adcroft 1.1
1120     end do
1121    
1122     do ki = 1, Nr
1123    
1124     c-----------------------------------------------------------------------
1125     c compute turbulent velocity scales on the interfaces
1126     c-----------------------------------------------------------------------
1127    
1128     do i = 1, imt
1129     sig = (-zgrid(ki) + 0.5 * hwide(ki)) / hbl(i)
1130     sigma(i) = stable(i)*sig + (1.-stable(i))*min(sig,epsilon)
1131     end do
1132     call wscale (
1133     I sigma, hbl, ustar, bfsfc,
1134     O wm, ws )
1135    
1136     c-----------------------------------------------------------------------
1137     c compute the dimensionless shape functions at the interfaces
1138     c-----------------------------------------------------------------------
1139    
1140     do i = 1, imt
1141     sig = (-zgrid(ki) + 0.5 * hwide(ki)) / hbl(i)
1142     a1 = sig - 2.
1143     a2 = 3. - 2. * sig
1144     a3 = sig - 1.
1145    
1146 heimbach 1.2 Gm = a1 + a2 * gat1m(i) + a3 * dat1m(i)
1147     Gs = a1 + a2 * gat1s(i) + a3 * dat1s(i)
1148     Gt = a1 + a2 * gat1t(i) + a3 * dat1t(i)
1149 adcroft 1.1
1150     c-----------------------------------------------------------------------
1151     c compute boundary layer diffusivities at the interfaces
1152     c-----------------------------------------------------------------------
1153    
1154 heimbach 1.2 blmc(i,ki,1) = hbl(i) * wm(i) * sig * (1. + sig * Gm)
1155 adcroft 1.1 blmc(i,ki,2) = hbl(i) * ws(i) * sig * (1. + sig * Gs)
1156     blmc(i,ki,3) = hbl(i) * ws(i) * sig * (1. + sig * Gt)
1157    
1158     c-----------------------------------------------------------------------
1159     c nonlocal transport term = ghat * <ws>o
1160     c-----------------------------------------------------------------------
1161 heimbach 1.4
1162     tempVar = ws(i) * hbl(i)
1163     ghat(i,ki) = (1.-stable(i)) * cg / max(phepsi,tempVar)
1164 adcroft 1.1
1165     end do
1166     end do
1167    
1168     c-----------------------------------------------------------------------
1169     c find diffusivities at kbl-1 grid level
1170     c-----------------------------------------------------------------------
1171    
1172     do i = 1, imt
1173     sig = -zgrid(kbl(i)-1) / hbl(i)
1174     sigma(i) = stable(i) * sig
1175     & + (1. - stable(i)) * min(sig,epsilon)
1176     end do
1177    
1178     call wscale (
1179     I sigma, hbl, ustar, bfsfc,
1180     O wm, ws )
1181    
1182     do i = 1, imt
1183     sig = -zgrid(kbl(i)-1) / hbl(i)
1184     a1 = sig - 2.
1185     a2 = 3. - 2. * sig
1186     a3 = sig - 1.
1187 heimbach 1.2 Gm = a1 + a2 * gat1m(i) + a3 * dat1m(i)
1188     Gs = a1 + a2 * gat1s(i) + a3 * dat1s(i)
1189     Gt = a1 + a2 * gat1t(i) + a3 * dat1t(i)
1190     dkm1(i,1) = hbl(i) * wm(i) * sig * (1. + sig * Gm)
1191 adcroft 1.1 dkm1(i,2) = hbl(i) * ws(i) * sig * (1. + sig * Gs)
1192     dkm1(i,3) = hbl(i) * ws(i) * sig * (1. + sig * Gt)
1193     end do
1194    
1195     #endif /* ALLOW_KPP */
1196    
1197     return
1198     end
1199    
1200     c*************************************************************************
1201    
1202     subroutine enhance (
1203     I dkm1, hbl, kbl, diffus, casea
1204     U , ghat
1205     O , blmc
1206     & )
1207    
1208     c enhance the diffusivity at the kbl-.5 interface
1209    
1210     IMPLICIT NONE
1211    
1212     #include "SIZE.h"
1213     #include "KPP_PARAMS.h"
1214    
1215     c input
1216     c dkm1(imt,mdiff) bl diffusivity at kbl-1 grid level
1217     c hbl(imt) boundary layer depth (m)
1218     c kbl(imt) grid above hbl
1219     c diffus(imt,0:Nrp1,mdiff) vertical diffusivities (m^2/s)
1220     c casea(imt) = 1 in caseA, = 0 in case B
1221 heimbach 1.4 _KPP_RL dkm1 (imt,mdiff)
1222     _KPP_RL hbl (imt)
1223 adcroft 1.1 integer kbl (imt)
1224 heimbach 1.4 _KPP_RL diffus(imt,0:Nrp1,mdiff)
1225     _KPP_RL casea (imt)
1226 adcroft 1.1
1227     c input/output
1228     c nonlocal transport, modified ghat at kbl(i)-1 interface (s/m**2)
1229 heimbach 1.4 _KPP_RL ghat (imt,Nr)
1230 adcroft 1.1
1231     c output
1232     c enhanced bound. layer mixing coeff.
1233 heimbach 1.4 _KPP_RL blmc (imt,Nr,mdiff)
1234 adcroft 1.1
1235     #ifdef ALLOW_KPP
1236    
1237     c local
1238     c fraction hbl lies beteen zgrid neighbors
1239 heimbach 1.4 _KPP_RL delta
1240 adcroft 1.1 integer ki, i, md
1241 heimbach 1.4 _KPP_RL dkmp5, dstar
1242 adcroft 1.1
1243     do i = 1, imt
1244     ki = kbl(i)-1
1245     if ((ki .ge. 1) .and. (ki .lt. Nr)) then
1246     delta = (hbl(i) + zgrid(ki)) / (zgrid(ki) - zgrid(ki+1))
1247     do md = 1, mdiff
1248     dkmp5 = casea(i) * diffus(i,ki,md) +
1249     1 (1.- casea(i)) * blmc (i,ki,md)
1250     dstar = (1.- delta)**2 * dkm1(i,md)
1251     & + delta**2 * dkmp5
1252     blmc(i,ki,md) = (1.- delta)*diffus(i,ki,md)
1253     & + delta*dstar
1254     end do
1255     ghat(i,ki) = (1.- casea(i)) * ghat(i,ki)
1256     endif
1257     end do
1258    
1259     #endif /* ALLOW_KPP */
1260    
1261     return
1262     end
1263    
1264     c*************************************************************************
1265    
1266     SUBROUTINE STATEKPP (
1267     I bi, bj, myThid,
1268     O RHO1, DBLOC, DBSFC, TTALPHA, SSBETA)
1269     c
1270     c-----------------------------------------------------------------------
1271     c "statekpp" computes all necessary input arrays
1272     c for the kpp mixing scheme
1273     c
1274     c input:
1275     c bi, bj = array indices on which to apply calculations
1276     c
1277     c output:
1278     c rho1 = potential density of surface layer (kg/m^3)
1279     c dbloc = local buoyancy gradient at Nr interfaces
1280     c g/rho{k+1,k+1} * [ drho{k,k+1}-drho{k+1,k+1} ] (m/s^2)
1281     c dbsfc = buoyancy difference with respect to the surface
1282     c g * [ drho{1,k}/rho{1,k} - drho{k,k}/rho{k,k} ] (m/s^2)
1283     c ttalpha= thermal expansion coefficient without 1/rho factor
1284     c d(rho) / d(potential temperature) (kg/m^3/C)
1285     c ssbeta = salt expansion coefficient without 1/rho factor
1286     c d(rho) / d(salinity) (kg/m^3/PSU)
1287     c
1288     c see also subroutines find_rho.F find_alpha.F find_beta.F
1289     c
1290     c written by: jan morzel, feb. 10, 1995 (converted from "sigma" version)
1291     c modified by: d. menemenlis, june 1998 : for use with MIT GCM UV
1292     c
1293     c-----------------------------------------------------------------------
1294    
1295     IMPLICIT NONE
1296    
1297     #include "SIZE.h"
1298     #include "EEPARAMS.h"
1299     #include "PARAMS.h"
1300     #include "KPP_PARAMS.h"
1301 adcroft 1.6 #include "DYNVARS.h"
1302 adcroft 1.1
1303     c-------------- Routine arguments -----------------------------------------
1304     INTEGER bi, bj, myThid
1305 heimbach 1.4 _KPP_RL RHO1 ( ibot:itop, jbot:jtop )
1306     _KPP_RL DBLOC ( ibot:itop, jbot:jtop, Nr )
1307     _KPP_RL DBSFC ( ibot:itop, jbot:jtop, Nr )
1308     _KPP_RL TTALPHA( ibot:itop, jbot:jtop, Nrp1 )
1309     _KPP_RL SSBETA ( ibot:itop, jbot:jtop, Nrp1 )
1310 adcroft 1.1
1311     #ifdef ALLOW_KPP
1312    
1313     c--------------------------------------------------------------------------
1314     c
1315     c local arrays:
1316     c
1317     c rhok - density of t(k ) & s(k ) at depth k
1318     c rhokm1 - density of t(k-1) & s(k-1) at depth k
1319     c rho1k - density of t(1 ) & s(1 ) at depth k
1320     c work1, work2 - work arrays for holding horizontal slabs
1321    
1322     _RL RHOK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1323     _RL RHOKM1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1324     _RL RHO1K (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1325     _RL WORK1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1326     _RL WORK2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1327     _RL WORK3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
1328     INTEGER I, J, K
1329    
1330     c calculate density, alpha, beta in surface layer, and set dbsfc to zero
1331    
1332     call FIND_RHO(
1333 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, 1, 1, eosType,
1334 adcroft 1.6 I theta, salt,
1335 adcroft 1.1 O WORK1,
1336     I myThid )
1337    
1338     call FIND_ALPHA(
1339 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, 1, 1, eosType,
1340 adcroft 1.1 O WORK2 )
1341    
1342     call FIND_BETA(
1343 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, 1, 1, eosType,
1344 adcroft 1.1 O WORK3 )
1345    
1346 heimbach 1.4 DO J = jbot, jtop
1347     DO I = ibot, itop
1348 adcroft 1.1 RHO1(I,J) = WORK1(I,J) + rhonil
1349     TTALPHA(I,J,1) = WORK2(I,J)
1350     SSBETA(I,J,1) = WORK3(I,J)
1351     DBSFC(I,J,1) = 0.
1352     END DO
1353     END DO
1354    
1355     c calculate alpha, beta, and gradients in interior layers
1356    
1357 heimbach 1.2 CHPF$ INDEPENDENT, NEW (RHOK,RHOKM1,RHO1K,WORK1,WORK2)
1358 adcroft 1.1 DO K = 2, Nr
1359    
1360     call FIND_RHO(
1361 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, K, K, eosType,
1362 adcroft 1.6 I theta, salt,
1363 adcroft 1.1 O RHOK,
1364     I myThid )
1365    
1366     call FIND_RHO(
1367 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, K-1, K, eosType,
1368 adcroft 1.6 I theta, salt,
1369 adcroft 1.1 O RHOKM1,
1370     I myThid )
1371    
1372     call FIND_RHO(
1373 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, 1, K, eosType,
1374 adcroft 1.6 I theta, salt,
1375 adcroft 1.1 O RHO1K,
1376     I myThid )
1377    
1378     call FIND_ALPHA(
1379 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, K, K, eosType,
1380 adcroft 1.1 O WORK1 )
1381    
1382     call FIND_BETA(
1383 heimbach 1.4 I bi, bj, ibot, itop, jbot, jtop, K, K, eosType,
1384 adcroft 1.1 O WORK2 )
1385    
1386 heimbach 1.4 DO J = jbot, jtop
1387     DO I = ibot, itop
1388 adcroft 1.1 TTALPHA(I,J,K) = WORK1 (I,J)
1389     SSBETA(I,J,K) = WORK2 (I,J)
1390     DBLOC(I,J,K-1) = gravity * (RHOK(I,J) - RHOKM1(I,J)) /
1391     & (RHOK(I,J) + rhonil)
1392     DBSFC(I,J,K) = gravity * (RHOK(I,J) - RHO1K (I,J)) /
1393     & (RHOK(I,J) + rhonil)
1394     END DO
1395     END DO
1396    
1397     END DO
1398    
1399     c compute arrays for K = Nrp1
1400 heimbach 1.4 DO J = jbot, jtop
1401     DO I = ibot, itop
1402 adcroft 1.1 TTALPHA(I,J,Nrp1) = TTALPHA(I,J,Nr)
1403     SSBETA(I,J,Nrp1) = SSBETA(I,J,Nr)
1404     DBLOC(I,J,Nr) = 0.
1405     END DO
1406     END DO
1407    
1408     #endif /* ALLOW_KPP */
1409    
1410     RETURN
1411     END

  ViewVC Help
Powered by ViewVC 1.1.22