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