/[MITgcm]/MITgcm/pkg/kpp/kpp_routines.F
ViewVC logotype

Contents of /MITgcm/pkg/kpp/kpp_routines.F

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


Revision 1.49 - (show annotations) (download)
Wed Jul 4 20:23:23 2012 UTC (11 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64l, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f
Changes since 1.48: +18 -15 lines
move flag "inAdMode" from PARAMS.h to AUTODIFF_PARAMS.h

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

  ViewVC Help
Powered by ViewVC 1.1.22