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