/[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.40 - (hide annotations) (download)
Thu Oct 23 18:11:00 2008 UTC (15 years, 8 months ago) by dfer
Branch: MAIN
Changes since 1.39: +14 -8 lines
Add Bulk Richardson number diagnostics

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

  ViewVC Help
Powered by ViewVC 1.1.22