/[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.30 - (hide annotations) (download)
Tue May 1 04:09:25 2007 UTC (17 years, 2 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint59a
Changes since 1.29: +53 -28 lines
 add more code to have only Ri-number based mixing in shelf ice caverns
 add myThid to all kpp routines (long overdue)

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

  ViewVC Help
Powered by ViewVC 1.1.22