/[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.54 - (show annotations) (download)
Thu Sep 11 19:23:23 2014 UTC (9 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65l, checkpoint65m, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e
Changes since 1.53: +95 -27 lines
Include explicitly AUTODIFF_OPTIONS.h (in case we don't use ECCO_CPPOPTIONS.h)
and put storage dir within #ifdef ALLOW_AUTODIFF_TAMC

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

  ViewVC Help
Powered by ViewVC 1.1.22