/[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.10 - (hide annotations) (download)
Sat Jul 13 03:12:30 2002 UTC (21 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint46g_pre, checkpoint46f_post, checkpoint46b_post, checkpoint46d_pre, checkpoint46a_post, checkpoint46e_pre, checkpoint46b_pre, checkpoint46c_pre, checkpoint46, checkpoint46a_pre, checkpoint46c_post, checkpoint46e_post, checkpoint46d_post
Changes since 1.9: +14 -4 lines
Merging from release1_p5
o Adjoint-related bug fixes in kpp:
  - kpp_calc: sore of kpphbl avoids recomputation/call to S/R kppmix
  - kpp_routines: store of Rib avoids partial recomputation bug of TAF.

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

  ViewVC Help
Powered by ViewVC 1.1.22