/[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.4 - (hide annotations) (download)
Mon Nov 13 16:37:02 2000 UTC (23 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: branch-atmos-merge-start, checkpoint33, checkpoint32, branch-atmos-merge-phase1, branch-atmos-merge-phase3, branch-atmos-merge-phase2
Branch point for: branch-atmos-merge
Changes since 1.3: +246 -370 lines
Modified and fixed version. Tested with verification/natl_box.

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

  ViewVC Help
Powered by ViewVC 1.1.22