/[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.44 - (hide annotations) (download)
Thu May 28 22:59:23 2009 UTC (15 years, 1 month ago) by dfer
Branch: MAIN
CVS Tags: checkpoint61q, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p
Changes since 1.43: +4 -1 lines
Add run-time flag to remove limit on hbl under stable conditions
(default behavior unchanged)

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

  ViewVC Help
Powered by ViewVC 1.1.22