/[MITgcm]/MITgcm_contrib/atnguyen/code_21Dec2012_saltplume/kpp_routines.F
ViewVC logotype

Annotation of /MITgcm_contrib/atnguyen/code_21Dec2012_saltplume/kpp_routines.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (hide annotations) (download)
Sun Apr 20 04:03:07 2014 UTC (11 years, 3 months ago) by atn
Branch: MAIN
salt plume volume fraction code mod in progress, installment 01

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

  ViewVC Help
Powered by ViewVC 1.1.22