/[MITgcm]/MITgcm_contrib/atnguyen/code_21Dec2012_saltplume/kpp_routines.F
ViewVC logotype

Contents of /MITgcm_contrib/atnguyen/code_21Dec2012_saltplume/kpp_routines.F

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


Revision 1.6 - (show annotations) (download)
Sat May 3 06:10:54 2014 UTC (11 years, 2 months ago) by atn
Branch: MAIN
CVS Tags: HEAD
Changes since 1.5: +3 -3 lines
bug fix

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

  ViewVC Help
Powered by ViewVC 1.1.22