/[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.39 - (hide annotations) (download)
Mon Sep 29 15:02:08 2008 UTC (15 years, 8 months ago) by dfer
Branch: MAIN
CVS Tags: checkpoint61e
Changes since 1.38: +2 -2 lines
Add missing comma (didn't compile if ALLOW_KPP_VERTICALLY_SMOOTH defined)

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

  ViewVC Help
Powered by ViewVC 1.1.22