/[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.48 - (show annotations) (download)
Wed Jan 20 01:32:50 2010 UTC (14 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63g, checkpoint63, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint62b
Changes since 1.47: +9 -3 lines
avoid unused variables

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

  ViewVC Help
Powered by ViewVC 1.1.22