/[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.26 - (hide annotations) (download)
Mon Apr 23 20:46:49 2007 UTC (17 years, 2 months ago) by dimitri
Branch: MAIN
Changes since 1.25: +43 -32 lines
o bug fixes for vertical diffusivity computations when both KPP and
    3D diffusivity arrays are used.

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

  ViewVC Help
Powered by ViewVC 1.1.22