/[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.33 - (hide annotations) (download)
Sat Sep 22 17:55:32 2007 UTC (16 years, 9 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint59h
Changes since 1.32: +70 -8 lines
added the salt plume scheme (ALLOW_SALT_PLUME) to pkg/kpp
(code conttributed by An Nguyen)

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

  ViewVC Help
Powered by ViewVC 1.1.22