/[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.45 - (hide annotations) (download)
Fri Sep 18 11:40:22 2009 UTC (14 years, 9 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint61v
Changes since 1.44: +167 -1 lines
add double diffusive contributions as a hack, turned off by default
for now, and the code can be excluded with a CPP-flag
EXCLUDE_KPP_DOUBLEDIFF just as EXCLUDE_KPP_SHEAR_MIX

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

  ViewVC Help
Powered by ViewVC 1.1.22