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