1 |
C $Header: /u/gcmpack/MITgcm_contrib/ecco_darwin/v4_llc270/code_darwin/ggl90_calc.F,v 1.2 2019/08/24 13:10:27 dcarroll Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "GGL90_OPTIONS.h" |
5 |
CBOP |
6 |
C !ROUTINE: GGL90_CALC |
7 |
|
8 |
C !INTERFACE: ====================================================== |
9 |
SUBROUTINE GGL90_CALC( |
10 |
I bi, bj, sigmaR, myTime, myIter, myThid ) |
11 |
|
12 |
C !DESCRIPTION: \bv |
13 |
C *==========================================================* |
14 |
C | SUBROUTINE GGL90_CALC | |
15 |
C | o Compute all GGL90 fields defined in GGL90.h | |
16 |
C *==========================================================* |
17 |
C | Equation numbers refer to | |
18 |
C | Gaspar et al. (1990), JGR 95 (C9), pp 16,179 | |
19 |
C | Some parts of the implementation follow Blanke and | |
20 |
C | Delecuse (1993), JPO, and OPA code, in particular the | |
21 |
C | computation of the | |
22 |
C | mixing length = max(min(lk,depth),lkmin) | |
23 |
C *==========================================================* |
24 |
|
25 |
C global parameters updated by ggl90_calc |
26 |
C GGL90TKE :: sub-grid turbulent kinetic energy (m^2/s^2) |
27 |
C GGL90viscAz :: GGL90 eddy viscosity coefficient (m^2/s) |
28 |
C GGL90diffKzT :: GGL90 diffusion coefficient for temperature (m^2/s) |
29 |
C \ev |
30 |
|
31 |
C !USES: ============================================================ |
32 |
IMPLICIT NONE |
33 |
#include "SIZE.h" |
34 |
#include "EEPARAMS.h" |
35 |
#include "PARAMS.h" |
36 |
#include "DYNVARS.h" |
37 |
#include "FFIELDS.h" |
38 |
#include "GRID.h" |
39 |
#include "GGL90.h" |
40 |
|
41 |
C !INPUT PARAMETERS: =================================================== |
42 |
C Routine arguments |
43 |
C bi, bj :: Current tile indices |
44 |
C sigmaR :: Vertical gradient of iso-neutral density |
45 |
C myTime :: Current time in simulation |
46 |
C myIter :: Current time-step number |
47 |
C myThid :: My Thread Id number |
48 |
INTEGER bi, bj |
49 |
_RL sigmaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
50 |
_RL myTime |
51 |
INTEGER myIter |
52 |
INTEGER myThid |
53 |
|
54 |
#ifdef ALLOW_GGL90 |
55 |
|
56 |
C !LOCAL VARIABLES: ==================================================== |
57 |
C Local constants |
58 |
C iMin,iMax,jMin,jMax :: index boundaries of computation domain |
59 |
C i, j, k, kp1,km1 :: array computation indices |
60 |
C kSurf, kBottom :: vertical indices of domain boundaries |
61 |
C hFac/hFacI :: fractional thickness of W-cell |
62 |
C explDissFac :: explicit Dissipation Factor (in [0-1]) |
63 |
C implDissFac :: implicit Dissipation Factor (in [0-1]) |
64 |
C uStarSquare :: square of friction velocity |
65 |
C verticalShear :: (squared) vertical shear of horizontal velocity |
66 |
C Nsquare :: squared buoyancy freqency |
67 |
C RiNumber :: local Richardson number |
68 |
C KappaM :: (local) viscosity parameter (eq.10) |
69 |
C KappaH :: (local) diffusivity parameter for temperature (eq.11) |
70 |
C KappaE :: (local) diffusivity parameter for TKE (eq.15) |
71 |
C TKEdissipation :: dissipation of TKE |
72 |
C GGL90mixingLength:: mixing length of scheme following Banke+Delecuse |
73 |
C rMixingLength:: inverse of mixing length |
74 |
C totalDepth :: thickness of water column (inverse of recip_Rcol) |
75 |
C TKEPrandtlNumber :: here, an empirical function of the Richardson number |
76 |
INTEGER iMin ,iMax ,jMin ,jMax |
77 |
INTEGER i, j, k, kp1, km1, kSurf, kBottom |
78 |
_RL explDissFac, implDissFac |
79 |
_RL uStarSquare |
80 |
_RL verticalShear(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
81 |
_RL KappaM(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
82 |
_RL KappaH |
83 |
c _RL Nsquare |
84 |
_RL Nsquare(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
85 |
_RL deltaTggl90 |
86 |
c _RL SQRTTKE |
87 |
_RL SQRTTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
88 |
_RL RiNumber |
89 |
#ifdef ALLOW_GGL90_IDEMIX |
90 |
_RL IDEMIX_RiNumber |
91 |
#endif |
92 |
_RL TKEdissipation |
93 |
_RL tempU, tempUp, tempV, tempVp, prTemp |
94 |
_RL MaxLength, tmpmlx, tmpVisc |
95 |
_RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
96 |
_RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
97 |
_RL rMixingLength (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
98 |
_RL mxLength_Dn (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
99 |
_RL KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
100 |
_RL totalDepth (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
101 |
_RL GGL90visctmp (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
102 |
#ifdef ALLOW_DIAGNOSTICS |
103 |
_RL surf_flx_tke (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
104 |
#endif /* ALLOW_DIAGNOSTICS */ |
105 |
C hFac(I) :: fractional thickness of W-cell |
106 |
_RL hFac |
107 |
#ifdef ALLOW_GGL90_IDEMIX |
108 |
_RL hFacI(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
109 |
#endif /* ALLOW_GGL90_IDEMIX */ |
110 |
_RL recip_hFacI(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
111 |
C- tri-diagonal matrix |
112 |
_RL a3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
113 |
_RL b3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
114 |
_RL c3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
115 |
INTEGER errCode |
116 |
#ifdef ALLOW_GGL90_HORIZDIFF |
117 |
C xA, yA :: area of lateral faces |
118 |
C dfx, dfy :: diffusive flux across lateral faces |
119 |
C gTKE :: right hand side of diffusion equation |
120 |
_RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
121 |
_RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
122 |
_RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
123 |
_RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
124 |
_RL gTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
125 |
#endif /* ALLOW_GGL90_HORIZDIFF */ |
126 |
#ifdef ALLOW_GGL90_SMOOTH |
127 |
_RL p4, p8, p16 |
128 |
#endif |
129 |
CEOP |
130 |
|
131 |
PARAMETER( iMin = 2-OLx, iMax = sNx+OLx-1 ) |
132 |
PARAMETER( jMin = 2-OLy, jMax = sNy+OLy-1 ) |
133 |
#ifdef ALLOW_GGL90_SMOOTH |
134 |
p4 = 0.25 _d 0 |
135 |
p8 = 0.125 _d 0 |
136 |
p16 = 0.0625 _d 0 |
137 |
#endif |
138 |
|
139 |
C set separate time step (should be deltaTtracer) |
140 |
deltaTggl90 = dTtracerLev(1) |
141 |
|
142 |
kSurf = 1 |
143 |
C explicit/implicit timestepping weights for dissipation |
144 |
explDissFac = 0. _d 0 |
145 |
implDissFac = 1. _d 0 - explDissFac |
146 |
|
147 |
C For nonlinear free surface and especially with r*-coordinates, the |
148 |
C hFacs change every timestep, so we need to update them here in the |
149 |
C case of using IDEMIX. |
150 |
DO K=1,Nr |
151 |
Km1 = MAX(K-1,1) |
152 |
DO j=1-OLy,sNy+OLy |
153 |
DO i=1-OLx,sNx+OLx |
154 |
hFac = |
155 |
& MIN(.5 _d 0,_hFacC(i,j,km1,bi,bj) ) + |
156 |
& MIN(.5 _d 0,_hFacC(i,j,k ,bi,bj) ) |
157 |
recip_hFacI(I,J,K)=0. _d 0 |
158 |
IF ( hFac .NE. 0. _d 0 ) |
159 |
& recip_hFacI(I,J,K)=1. _d 0/hFac |
160 |
#ifdef ALLOW_GGL90_IDEMIX |
161 |
hFacI(i,j,k) = hFac |
162 |
#endif /* ALLOW_GGL90_IDEMIX */ |
163 |
ENDDO |
164 |
ENDDO |
165 |
ENDDO |
166 |
|
167 |
C Initialize local fields |
168 |
DO k = 1, Nr |
169 |
DO j=1-OLy,sNy+OLy |
170 |
DO i=1-OLx,sNx+OLx |
171 |
rMixingLength(i,j,k) = 0. _d 0 |
172 |
mxLength_Dn(i,j,k) = 0. _d 0 |
173 |
GGL90visctmp(i,j,k) = 0. _d 0 |
174 |
KappaE(i,j,k) = 0. _d 0 |
175 |
TKEPrandtlNumber(i,j,k) = 1. _d 0 |
176 |
GGL90mixingLength(i,j,k) = GGL90mixingLengthMin |
177 |
GGL90visctmp(i,j,k) = 0. _d 0 |
178 |
#ifndef SOLVE_DIAGONAL_LOWMEMORY |
179 |
a3d(i,j,k) = 0. _d 0 |
180 |
b3d(i,j,k) = 1. _d 0 |
181 |
c3d(i,j,k) = 0. _d 0 |
182 |
#endif |
183 |
Nsquare(i,j,k) = 0. _d 0 |
184 |
SQRTTKE(i,j,k) = 0. _d 0 |
185 |
ENDDO |
186 |
ENDDO |
187 |
ENDDO |
188 |
DO j=1-OLy,sNy+OLy |
189 |
DO i=1-OLx,sNx+OLx |
190 |
KappaM(i,j) = 0. _d 0 |
191 |
verticalShear(i,j) = 0. _d 0 |
192 |
totalDepth(i,j) = Ro_surf(i,j,bi,bj) - R_low(i,j,bi,bj) |
193 |
rMixingLength(i,j,1) = 0. _d 0 |
194 |
mxLength_Dn(i,j,1) = GGL90mixingLengthMin |
195 |
SQRTTKE(i,j,1) = SQRT( GGL90TKE(i,j,1,bi,bj) ) |
196 |
#ifdef ALLOW_GGL90_HORIZDIFF |
197 |
xA(i,j) = 0. _d 0 |
198 |
yA(i,j) = 0. _d 0 |
199 |
dfx(i,j) = 0. _d 0 |
200 |
dfy(i,j) = 0. _d 0 |
201 |
gTKE(i,j) = 0. _d 0 |
202 |
#endif /* ALLOW_GGL90_HORIZDIFF */ |
203 |
ENDDO |
204 |
ENDDO |
205 |
|
206 |
#ifdef ALLOW_GGL90_IDEMIX |
207 |
IF ( useIDEMIX) CALL GGL90_IDEMIX( |
208 |
& bi, bj, hFacI, recip_hFacI, sigmaR, myTime, myIter, myThid ) |
209 |
#endif /* ALLOW_GGL90_IDEMIX */ |
210 |
|
211 |
DO k = 2, Nr |
212 |
DO j=jMin,jMax |
213 |
DO i=iMin,iMax |
214 |
SQRTTKE(i,j,k)=SQRT( GGL90TKE(i,j,k,bi,bj) ) |
215 |
|
216 |
C buoyancy frequency |
217 |
Nsquare(i,j,k) = gravity*gravitySign*recip_rhoConst |
218 |
& * sigmaR(i,j,k) |
219 |
C vertical shear term (dU/dz)^2+(dV/dz)^2 is computed later |
220 |
C to save some memory |
221 |
C mixing length |
222 |
GGL90mixingLength(i,j,k) = SQRTTWO * |
223 |
& SQRTTKE(i,j,k)/SQRT( MAX(Nsquare(i,j,k),GGL90eps) ) |
224 |
ENDDO |
225 |
ENDDO |
226 |
ENDDO |
227 |
|
228 |
C- ensure mixing between first and second level |
229 |
IF (mxlSurfFlag) THEN |
230 |
DO j=jMin,jMax |
231 |
DO i=iMin,iMax |
232 |
GGL90mixingLength(i,j,2)=drF(1) |
233 |
ENDDO |
234 |
ENDDO |
235 |
ENDIF |
236 |
|
237 |
C-- Impose upper and lower bound for mixing length |
238 |
C-- Impose minimum mixing length to avoid division by zero |
239 |
IF ( mxlMaxFlag .EQ. 0 ) THEN |
240 |
|
241 |
DO k=2,Nr |
242 |
DO j=jMin,jMax |
243 |
DO i=iMin,iMax |
244 |
MaxLength=totalDepth(i,j) |
245 |
GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
246 |
& MaxLength) |
247 |
ENDDO |
248 |
ENDDO |
249 |
ENDDO |
250 |
|
251 |
DO k=2,Nr |
252 |
DO j=jMin,jMax |
253 |
DO i=iMin,iMax |
254 |
GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k), |
255 |
& GGL90mixingLengthMin) |
256 |
rMixingLength(i,j,k) = 1. _d 0 / GGL90mixingLength(i,j,k) |
257 |
ENDDO |
258 |
ENDDO |
259 |
ENDDO |
260 |
|
261 |
ELSEIF ( mxlMaxFlag .EQ. 1 ) THEN |
262 |
|
263 |
DO k=2,Nr |
264 |
DO j=jMin,jMax |
265 |
DO i=iMin,iMax |
266 |
MaxLength=MIN(Ro_surf(i,j,bi,bj)-rF(k),rF(k)-R_low(i,j,bi,bj)) |
267 |
c MaxLength=MAX(MaxLength,20. _d 0) |
268 |
GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
269 |
& MaxLength) |
270 |
ENDDO |
271 |
ENDDO |
272 |
ENDDO |
273 |
|
274 |
DO k=2,Nr |
275 |
DO j=jMin,jMax |
276 |
DO i=iMin,iMax |
277 |
GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k), |
278 |
& GGL90mixingLengthMin) |
279 |
rMixingLength(i,j,k) = 1. _d 0 / GGL90mixingLength(i,j,k) |
280 |
ENDDO |
281 |
ENDDO |
282 |
ENDDO |
283 |
|
284 |
ELSEIF ( mxlMaxFlag .EQ. 2 ) THEN |
285 |
|
286 |
DO k=2,Nr |
287 |
DO j=jMin,jMax |
288 |
DO i=iMin,iMax |
289 |
GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
290 |
& GGL90mixingLength(i,j,k-1)+drF(k-1)) |
291 |
ENDDO |
292 |
ENDDO |
293 |
ENDDO |
294 |
DO j=jMin,jMax |
295 |
DO i=iMin,iMax |
296 |
GGL90mixingLength(i,j,Nr) = MIN(GGL90mixingLength(i,j,Nr), |
297 |
& GGL90mixingLengthMin+drF(Nr)) |
298 |
ENDDO |
299 |
ENDDO |
300 |
DO k=Nr-1,2,-1 |
301 |
DO j=jMin,jMax |
302 |
DO i=iMin,iMax |
303 |
GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
304 |
& GGL90mixingLength(i,j,k+1)+drF(k)) |
305 |
ENDDO |
306 |
ENDDO |
307 |
ENDDO |
308 |
|
309 |
DO k=2,Nr |
310 |
DO j=jMin,jMax |
311 |
DO i=iMin,iMax |
312 |
GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k), |
313 |
& GGL90mixingLengthMin) |
314 |
rMixingLength(i,j,k) = 1. _d 0 / GGL90mixingLength(i,j,k) |
315 |
ENDDO |
316 |
ENDDO |
317 |
ENDDO |
318 |
|
319 |
ELSEIF ( mxlMaxFlag .EQ. 3 ) THEN |
320 |
|
321 |
DO k=2,Nr |
322 |
DO j=jMin,jMax |
323 |
DO i=iMin,iMax |
324 |
mxLength_Dn(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
325 |
& mxLength_Dn(i,j,k-1)+drF(k-1)) |
326 |
ENDDO |
327 |
ENDDO |
328 |
ENDDO |
329 |
DO j=jMin,jMax |
330 |
DO i=iMin,iMax |
331 |
GGL90mixingLength(i,j,Nr) = MIN(GGL90mixingLength(i,j,Nr), |
332 |
& GGL90mixingLengthMin+drF(Nr)) |
333 |
ENDDO |
334 |
ENDDO |
335 |
DO k=Nr-1,2,-1 |
336 |
DO j=jMin,jMax |
337 |
DO i=iMin,iMax |
338 |
GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
339 |
& GGL90mixingLength(i,j,k+1)+drF(k)) |
340 |
ENDDO |
341 |
ENDDO |
342 |
ENDDO |
343 |
|
344 |
DO k=2,Nr |
345 |
DO j=jMin,jMax |
346 |
DO i=iMin,iMax |
347 |
GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k), |
348 |
& mxLength_Dn(i,j,k)) |
349 |
tmpmlx = SQRT( GGL90mixingLength(i,j,k)*mxLength_Dn(i,j,k) ) |
350 |
tmpmlx = MAX( tmpmlx, GGL90mixingLengthMin) |
351 |
rMixingLength(i,j,k) = 1. _d 0 / tmpmlx |
352 |
ENDDO |
353 |
ENDDO |
354 |
ENDDO |
355 |
|
356 |
ELSE |
357 |
STOP 'GGL90_CALC: Wrong mxlMaxFlag (mixing length limit)' |
358 |
ENDIF |
359 |
|
360 |
C start "proper" k-loop (the code above was moved out and up to |
361 |
C implemement various mixing parameters efficiently) |
362 |
DO k=2,Nr |
363 |
km1 = k-1 |
364 |
|
365 |
#ifdef ALLOW_GGL90_HORIZDIFF |
366 |
IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN |
367 |
C horizontal diffusion of TKE (requires an exchange in |
368 |
C do_fields_blocking_exchanges) |
369 |
C common factors |
370 |
DO j=1-OLy,sNy+OLy |
371 |
DO i=1-OLx,sNx+OLx |
372 |
xA(i,j) = _dyG(i,j,bi,bj)*drC(k)* |
373 |
& (min(.5 _d 0,_hFacW(i,j,k-1,bi,bj) ) + |
374 |
& min(.5 _d 0,_hFacW(i,j,k ,bi,bj) ) ) |
375 |
yA(i,j) = _dxG(i,j,bi,bj)*drC(k)* |
376 |
& (min(.5 _d 0,_hFacS(i,j,k-1,bi,bj) ) + |
377 |
& min(.5 _d 0,_hFacS(i,j,k ,bi,bj) ) ) |
378 |
ENDDO |
379 |
ENDDO |
380 |
C Compute diffusive fluxes |
381 |
C ... across x-faces |
382 |
DO j=1-OLy,sNy+OLy |
383 |
dfx(1-OLx,j)=0. _d 0 |
384 |
DO i=1-OLx+1,sNx+OLx |
385 |
dfx(i,j) = -GGL90diffTKEh*xA(i,j) |
386 |
& *_recip_dxC(i,j,bi,bj) |
387 |
& *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj)) |
388 |
#ifdef ISOTROPIC_COS_SCALING |
389 |
& *CosFacU(j,bi,bj) |
390 |
#endif /* ISOTROPIC_COS_SCALING */ |
391 |
ENDDO |
392 |
ENDDO |
393 |
C ... across y-faces |
394 |
DO i=1-OLx,sNx+OLx |
395 |
dfy(i,1-OLy)=0. _d 0 |
396 |
ENDDO |
397 |
DO j=1-OLy+1,sNy+OLy |
398 |
DO i=1-OLx,sNx+OLx |
399 |
dfy(i,j) = -GGL90diffTKEh*yA(i,j) |
400 |
& *_recip_dyC(i,j,bi,bj) |
401 |
& *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj)) |
402 |
#ifdef ISOTROPIC_COS_SCALING |
403 |
& *CosFacV(j,bi,bj) |
404 |
#endif /* ISOTROPIC_COS_SCALING */ |
405 |
ENDDO |
406 |
ENDDO |
407 |
C Compute divergence of fluxes |
408 |
DO j=1-OLy,sNy+OLy-1 |
409 |
DO i=1-OLx,sNx+OLx-1 |
410 |
gTKE(i,j) = -recip_drC(k)*recip_rA(i,j,bi,bj) |
411 |
& *recip_hFacI(i,j,k) |
412 |
& *((dfx(i+1,j)-dfx(i,j)) |
413 |
& + (dfy(i,j+1)-dfy(i,j)) ) |
414 |
ENDDO |
415 |
ENDDO |
416 |
C end if GGL90diffTKEh .eq. 0. |
417 |
ENDIF |
418 |
#endif /* ALLOW_GGL90_HORIZDIFF */ |
419 |
|
420 |
C viscosity and diffusivity |
421 |
DO j=jMin,jMax |
422 |
DO i=iMin,iMax |
423 |
KappaM(i,j) = GGL90ck*GGL90mixingLength(i,j,k)*SQRTTKE(i,j,k) |
424 |
GGL90visctmp(i,j,k) = MAX(KappaM(i,j),diffKrNrS(k)) |
425 |
& * maskC(i,j,k,bi,bj) |
426 |
C note: storing GGL90visctmp like this, and using it later to compute |
427 |
C GGL9rdiffKr etc. is robust in case of smoothing (e.g. see OPA) |
428 |
KappaM(i,j) = MAX(KappaM(i,j),viscArNr(k)) * maskC(i,j,k,bi,bj) |
429 |
ENDDO |
430 |
ENDDO |
431 |
|
432 |
C compute vertical shear (dU/dz)^2+(dV/dz)^2 |
433 |
IF ( calcMeanVertShear ) THEN |
434 |
C by averaging (@ grid-cell center) the 4 vertical shear compon @ U,V pos. |
435 |
DO j=jMin,jMax |
436 |
DO i=iMin,iMax |
437 |
tempU = ( uVel( i ,j,km1,bi,bj) - uVel( i ,j,k,bi,bj) ) |
438 |
tempUp = ( uVel(i+1,j,km1,bi,bj) - uVel(i+1,j,k,bi,bj) ) |
439 |
tempV = ( vVel(i, j ,km1,bi,bj) - vVel(i, j ,k,bi,bj) ) |
440 |
tempVp = ( vVel(i,j+1,km1,bi,bj) - vVel(i,j+1,k,bi,bj) ) |
441 |
verticalShear(i,j) = ( |
442 |
& ( tempU*tempU + tempUp*tempUp )*halfRL |
443 |
& + ( tempV*tempV + tempVp*tempVp )*halfRL |
444 |
& )*recip_drC(k)*recip_drC(k) |
445 |
ENDDO |
446 |
ENDDO |
447 |
ELSE |
448 |
C from the averaged flow at grid-cell center (2 compon x 2 pos.) |
449 |
DO j=jMin,jMax |
450 |
DO i=iMin,iMax |
451 |
tempU = ( ( uVel(i,j,km1,bi,bj) + uVel(i+1,j,km1,bi,bj) ) |
452 |
& -( uVel(i,j,k ,bi,bj) + uVel(i+1,j,k ,bi,bj) ) |
453 |
& )*halfRL*recip_drC(k) |
454 |
tempV = ( ( vVel(i,j,km1,bi,bj) + vVel(i,j+1,km1,bi,bj) ) |
455 |
& -( vVel(i,j,k ,bi,bj) + vVel(i,j+1,k ,bi,bj) ) |
456 |
& )*halfRL*recip_drC(k) |
457 |
verticalShear(i,j) = tempU*tempU + tempV*tempV |
458 |
ENDDO |
459 |
ENDDO |
460 |
ENDIF |
461 |
|
462 |
C compute Prandtl number (always greater than 0) |
463 |
#ifdef ALLOW_GGL90_IDEMIX |
464 |
IF ( useIDEMIX ) THEN |
465 |
DO j=jMin,jMax |
466 |
DO i=iMin,iMax |
467 |
C account for partical cell factor in vertical shear: |
468 |
verticalShear(i,j) = verticalShear(i,j) |
469 |
& * recip_hFacI(i,j,k)*recip_hFacI(i,j,k) |
470 |
RiNumber = MAX(Nsquare(i,j,k),0. _d 0) |
471 |
& /(verticalShear(i,j)+GGL90eps) |
472 |
CML IDEMIX_RiNumber = 1./GGL90eps |
473 |
IDEMIX_RiNumber = MAX( KappaM(i,j)*Nsquare(i,j,k), 0. _d 0)/ |
474 |
& (GGL90eps+IDEMIX_tau_d(i,j,k,bi,bj)*IDEMIX_E(i,j,k,bi,bj)**2) |
475 |
prTemp = MIN(5.*RiNumber, 6.6 _d 0*IDEMIX_RiNumber) |
476 |
TKEPrandtlNumber(i,j,k) = MIN(10. _d 0,prTemp) |
477 |
TKEPrandtlNumber(i,j,k) = MAX( 1. _d 0,TKEPrandtlNumber(i,j,k)) |
478 |
ENDDO |
479 |
ENDDO |
480 |
ELSE |
481 |
#else /* ndef ALLOW_GGL90_IDEMIX */ |
482 |
IF (.TRUE.) THEN |
483 |
#endif /* ALLOW_GGL90_IDEMIX */ |
484 |
DO j=jMin,jMax |
485 |
DO i=iMin,iMax |
486 |
RiNumber = MAX(Nsquare(i,j,k),0. _d 0) |
487 |
& /(verticalShear(i,j)+GGL90eps) |
488 |
prTemp = 1. _d 0 |
489 |
IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber |
490 |
TKEPrandtlNumber(i,j,k) = MIN(10. _d 0,prTemp) |
491 |
ENDDO |
492 |
ENDDO |
493 |
ENDIF |
494 |
|
495 |
DO j=jMin,jMax |
496 |
DO i=iMin,iMax |
497 |
C diffusivity |
498 |
KappaH = KappaM(i,j)/TKEPrandtlNumber(i,j,k) |
499 |
KappaE(i,j,k) = GGL90alpha * KappaM(i,j) * maskC(i,j,k,bi,bj) |
500 |
|
501 |
C dissipation term |
502 |
TKEdissipation = explDissFac*GGL90ceps |
503 |
& *SQRTTKE(i,j,k)*rMixingLength(i,j,k) |
504 |
& *GGL90TKE(i,j,k,bi,bj) |
505 |
C partial update with sum of explicit contributions |
506 |
GGL90TKE(i,j,k,bi,bj) = GGL90TKE(i,j,k,bi,bj) |
507 |
& + deltaTggl90*( |
508 |
& + KappaM(i,j)*verticalShear(i,j) |
509 |
& - KappaH*Nsquare(i,j,k) |
510 |
& - TKEdissipation |
511 |
& ) |
512 |
ENDDO |
513 |
ENDDO |
514 |
|
515 |
#ifdef ALLOW_GGL90_IDEMIX |
516 |
IF ( useIDEMIX ) THEN |
517 |
C add IDEMIX contribution to the turbulent kinetic energy |
518 |
DO j=jMin,jMax |
519 |
DO i=iMin,iMax |
520 |
GGL90TKE(i,j,k,bi,bj) = GGL90TKE(i,j,k,bi,bj) |
521 |
& + deltaTggl90*( |
522 |
& + IDEMIX_tau_d(i,j,k,bi,bj)*IDEMIX_E(i,j,k,bi,bj)**2 |
523 |
& ) |
524 |
ENDDO |
525 |
ENDDO |
526 |
ENDIF |
527 |
#endif /* ALLOW_GGL90_IDEMIX */ |
528 |
|
529 |
#ifdef ALLOW_GGL90_HORIZDIFF |
530 |
IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN |
531 |
C-- Add horiz. diffusion tendency |
532 |
DO j=jMin,jMax |
533 |
DO i=iMin,iMax |
534 |
GGL90TKE(i,j,k,bi,bj) = GGL90TKE(i,j,k,bi,bj) |
535 |
& + gTKE(i,j)*deltaTggl90 |
536 |
ENDDO |
537 |
ENDDO |
538 |
ENDIF |
539 |
#endif /* ALLOW_GGL90_HORIZDIFF */ |
540 |
|
541 |
C-- end of k loop |
542 |
ENDDO |
543 |
|
544 |
C ============================================ |
545 |
C Implicit time step to update TKE for k=1,Nr; |
546 |
C TKE(Nr+1)=0 by default |
547 |
C ============================================ |
548 |
C set up matrix |
549 |
C-- Lower diagonal |
550 |
DO j=jMin,jMax |
551 |
DO i=iMin,iMax |
552 |
a3d(i,j,1) = 0. _d 0 |
553 |
ENDDO |
554 |
ENDDO |
555 |
DO k=2,Nr |
556 |
km1=MAX(2,k-1) |
557 |
DO j=jMin,jMax |
558 |
DO i=iMin,iMax |
559 |
C- We keep recip_hFacC in the diffusive flux calculation, |
560 |
C- but no hFacC in TKE volume control |
561 |
C- No need for maskC(k-1) with recip_hFacC(k-1) |
562 |
a3d(i,j,k) = -deltaTggl90 |
563 |
& *recip_drF(k-1)*recip_hFacC(i,j,k-1,bi,bj) |
564 |
& *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1)) |
565 |
& *recip_drC(k)*maskC(i,j,k,bi,bj) |
566 |
ENDDO |
567 |
ENDDO |
568 |
ENDDO |
569 |
C-- Upper diagonal |
570 |
DO j=jMin,jMax |
571 |
DO i=iMin,iMax |
572 |
c3d(i,j,1) = 0. _d 0 |
573 |
ENDDO |
574 |
ENDDO |
575 |
DO k=2,Nr |
576 |
DO j=jMin,jMax |
577 |
DO i=iMin,iMax |
578 |
kp1=MAX(1,MIN(klowC(i,j,bi,bj),k+1)) |
579 |
C- We keep recip_hFacC in the diffusive flux calculation, |
580 |
C- but no hFacC in TKE volume control |
581 |
C- No need for maskC(k) with recip_hFacC(k) |
582 |
c3d(i,j,k) = -deltaTggl90 |
583 |
& *recip_drF( k ) * recip_hFacC(i,j,k,bi,bj) |
584 |
& *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1)) |
585 |
& *recip_drC(k)*maskC(i,j,k-1,bi,bj) |
586 |
ENDDO |
587 |
ENDDO |
588 |
ENDDO |
589 |
|
590 |
#ifdef ALLOW_GGL90_IDEMIX |
591 |
IF ( useIDEMIX ) THEN |
592 |
DO k=2,Nr |
593 |
DO j=jMin,jMax |
594 |
DO i=iMin,iMax |
595 |
a3d(i,j,k) = a3d(i,j,k)*recip_hFacI(i,j,k) |
596 |
c3d(i,j,k) = c3d(i,j,k)*recip_hFacI(i,j,k) |
597 |
ENDDO |
598 |
ENDDO |
599 |
ENDDO |
600 |
ENDIF |
601 |
#endif /* ALLOW_GGL90_IDEMIX */ |
602 |
|
603 |
IF (.NOT.GGL90_dirichlet) THEN |
604 |
C Neumann bottom boundary condition for TKE: no flux from bottom |
605 |
DO j=jMin,jMax |
606 |
DO i=iMin,iMax |
607 |
kBottom = MAX(kLowC(i,j,bi,bj),1) |
608 |
c3d(i,j,kBottom) = 0. _d 0 |
609 |
ENDDO |
610 |
ENDDO |
611 |
ENDIF |
612 |
|
613 |
C-- Center diagonal |
614 |
DO k=1,Nr |
615 |
km1 = MAX(k-1,1) |
616 |
DO j=jMin,jMax |
617 |
DO i=iMin,iMax |
618 |
b3d(i,j,k) = 1. _d 0 - c3d(i,j,k) - a3d(i,j,k) |
619 |
& + implDissFac*deltaTggl90*GGL90ceps*SQRTTKE(i,j,k) |
620 |
& * rMixingLength(i,j,k) |
621 |
& * maskC(i,j,k,bi,bj)*maskC(i,j,km1,bi,bj) |
622 |
ENDDO |
623 |
ENDDO |
624 |
ENDDO |
625 |
C end set up matrix |
626 |
|
627 |
C Apply boundary condition |
628 |
kp1 = MIN(Nr,kSurf+1) |
629 |
DO j=jMin,jMax |
630 |
DO i=iMin,iMax |
631 |
C estimate friction velocity uStar from surface forcing |
632 |
uStarSquare = SQRT( |
633 |
& ( .5 _d 0*( surfaceForcingU(i, j, bi,bj) |
634 |
& + surfaceForcingU(i+1,j, bi,bj) ) )**2 |
635 |
& + ( .5 _d 0*( surfaceForcingV(i, j, bi,bj) |
636 |
& + surfaceForcingV(i, j+1,bi,bj) ) )**2 |
637 |
& ) |
638 |
C Dirichlet surface boundary condition for TKE |
639 |
GGL90TKE(i,j,kSurf,bi,bj) = maskC(i,j,kSurf,bi,bj) |
640 |
& *MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare) |
641 |
GGL90TKE(i,j,kp1,bi,bj) = GGL90TKE(i,j,kp1,bi,bj) |
642 |
& - a3d(i,j,kp1)*GGL90TKE(i,j,kSurf,bi,bj) |
643 |
a3d(i,j,kp1) = 0. _d 0 |
644 |
ENDDO |
645 |
ENDDO |
646 |
|
647 |
IF (GGL90_dirichlet) THEN |
648 |
C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom |
649 |
DO j=jMin,jMax |
650 |
DO i=iMin,iMax |
651 |
kBottom = MAX(kLowC(i,j,bi,bj),1) |
652 |
GGL90TKE(i,j,kBottom,bi,bj) = GGL90TKE(i,j,kBottom,bi,bj) |
653 |
& - GGL90TKEbottom*c3d(i,j,kBottom) |
654 |
c3d(i,j,kBottom) = 0. _d 0 |
655 |
ENDDO |
656 |
ENDDO |
657 |
ENDIF |
658 |
|
659 |
C solve tri-diagonal system |
660 |
errCode = -1 |
661 |
CALL SOLVE_TRIDIAGONAL( iMin,iMax, jMin,jMax, |
662 |
I a3d, b3d, c3d, |
663 |
U GGL90TKE(1-OLx,1-OLy,1,bi,bj), |
664 |
O errCode, |
665 |
I bi, bj, myThid ) |
666 |
|
667 |
DO k=1,Nr |
668 |
DO j=jMin,jMax |
669 |
DO i=iMin,iMax |
670 |
C impose minimum TKE to avoid numerical undershoots below zero |
671 |
GGL90TKE(i,j,k,bi,bj) = maskC(i,j,k,bi,bj) |
672 |
& *MAX( GGL90TKE(i,j,k,bi,bj), GGL90TKEmin ) |
673 |
ENDDO |
674 |
ENDDO |
675 |
ENDDO |
676 |
|
677 |
C end of time step |
678 |
C =============================== |
679 |
|
680 |
DO k=2,Nr |
681 |
DO j=1,sNy |
682 |
DO i=1,sNx |
683 |
#ifdef ALLOW_GGL90_SMOOTH |
684 |
tmpVisc = ( |
685 |
& p4 * GGL90visctmp(i ,j ,k)*mskCor(i ,j ,bi,bj) |
686 |
& +p8 *( ( GGL90visctmp(i-1,j ,k)*mskCor(i-1,j ,bi,bj) |
687 |
& + GGL90visctmp(i+1,j ,k)*mskCor(i+1,j ,bi,bj) ) |
688 |
& + ( GGL90visctmp(i ,j-1,k)*mskCor(i ,j-1,bi,bj) |
689 |
& + GGL90visctmp(i ,j+1,k)*mskCor(i ,j+1,bi,bj) ) ) |
690 |
& +p16*( ( GGL90visctmp(i+1,j+1,k)*mskCor(i+1,j+1,bi,bj) |
691 |
& + GGL90visctmp(i-1,j-1,k)*mskCor(i-1,j-1,bi,bj) ) |
692 |
& + ( GGL90visctmp(i+1,j-1,k)*mskCor(i+1,j-1,bi,bj) |
693 |
& + GGL90visctmp(i-1,j+1,k)*mskCor(i-1,j+1,bi,bj) ) ) |
694 |
& )/( |
695 |
& p4 |
696 |
& +p8 *( ( maskC(i-1,j ,k,bi,bj)*mskCor(i-1,j ,bi,bj) |
697 |
& + maskC(i+1,j ,k,bi,bj)*mskCor(i+1,j ,bi,bj) ) |
698 |
& + ( maskC(i ,j-1,k,bi,bj)*mskCor(i ,j-1,bi,bj) |
699 |
& + maskC(i ,j+1,k,bi,bj)*mskCor(i ,j+1,bi,bj) ) ) |
700 |
& +p16*( ( maskC(i+1,j+1,k,bi,bj)* mskCor(i+1,j+1,bi,bj) |
701 |
& + maskC(i-1,j-1,k,bi,bj)*mskCor(i-1,j-1,bi,bj) ) |
702 |
& + ( maskC(i+1,j-1,k,bi,bj)*mskCor(i+1,j-1,bi,bj) |
703 |
& + maskC(i-1,j+1,k,bi,bj)*mskCor(i-1,j+1,bi,bj) ) ) |
704 |
& )*maskC(i,j,k,bi,bj)*mskCor(i,j,bi,bj) |
705 |
#else |
706 |
tmpVisc = GGL90visctmp(i,j,k) |
707 |
#endif |
708 |
tmpVisc = MIN(tmpVisc/TKEPrandtlNumber(i,j,k),GGL90diffMax) |
709 |
GGL90diffKr(i,j,k,bi,bj)= MAX( tmpVisc , diffKrNrS(k) ) |
710 |
ENDDO |
711 |
ENDDO |
712 |
ENDDO |
713 |
|
714 |
DO k=2,Nr |
715 |
DO j=1,sNy |
716 |
DO i=1,sNx+1 |
717 |
#ifdef ALLOW_GGL90_SMOOTH |
718 |
tmpVisc = ( |
719 |
& p4 *( GGL90visctmp(i-1,j ,k)*mskCor(i-1,j ,bi,bj) |
720 |
& + GGL90visctmp(i ,j ,k)*mskCor(i ,j ,bi,bj) ) |
721 |
& +p8 *( ( GGL90visctmp(i-1,j-1,k)*mskCor(i-1,j-1,bi,bj) |
722 |
& + GGL90visctmp(i ,j-1,k)*mskCor(i ,j-1,bi,bj) ) |
723 |
& + ( GGL90visctmp(i-1,j+1,k)*mskCor(i-1,j+1,bi,bj) |
724 |
& + GGL90visctmp(i ,j+1,k)*mskCor(i ,j+1,bi,bj) ) ) |
725 |
& )/( |
726 |
& p4 * 2. _d 0 |
727 |
& +p8 *( ( maskC(i-1,j-1,k,bi,bj)*mskCor(i-1,j-1,bi,bj) |
728 |
& + maskC(i ,j-1,k,bi,bj)*mskCor(i ,j-1,bi,bj) ) |
729 |
& + ( maskC(i-1,j+1,k,bi,bj)*mskCor(i-1,j+1,bi,bj) |
730 |
& + maskC(i ,j+1,k,bi,bj)*mskCor(i ,j+1,bi,bj) ) ) |
731 |
& )*maskC(i-1,j,k,bi,bj)*mskCor(i-1,j,bi,bj) |
732 |
& *maskC(i ,j,k,bi,bj)*mskCor(i ,j,bi,bj) |
733 |
#else |
734 |
tmpVisc = _maskW(i,j,k,bi,bj) * halfRL |
735 |
& *( GGL90visctmp(i-1,j,k) |
736 |
& + GGL90visctmp(i,j,k) ) |
737 |
#endif |
738 |
tmpVisc = MIN( tmpVisc , GGL90viscMax ) |
739 |
GGL90viscArU(i,j,k,bi,bj) = MAX( tmpVisc, viscArNr(k) ) |
740 |
ENDDO |
741 |
ENDDO |
742 |
ENDDO |
743 |
|
744 |
DO k=2,Nr |
745 |
DO j=1,sNy+1 |
746 |
DO i=1,sNx |
747 |
#ifdef ALLOW_GGL90_SMOOTH |
748 |
tmpVisc = ( |
749 |
& p4 *( GGL90visctmp(i ,j-1,k)*mskCor(i ,j-1,bi,bj) |
750 |
& + GGL90visctmp(i ,j ,k)*mskCor(i ,j ,bi,bj) ) |
751 |
& +p8 *( ( GGL90visctmp(i-1,j-1,k)*mskCor(i-1,j-1,bi,bj) |
752 |
& + GGL90visctmp(i-1,j ,k)*mskCor(i-1,j ,bi,bj) ) |
753 |
& + ( GGL90visctmp(i+1,j-1,k)*mskCor(i+1,j-1,bi,bj) |
754 |
& + GGL90visctmp(i+1,j ,k)*mskCor(i+1,j ,bi,bj) ) ) |
755 |
& )/( |
756 |
& p4 * 2. _d 0 |
757 |
& +p8 *( ( maskC(i-1,j-1,k,bi,bj)*mskCor(i-1,j-1,bi,bj) |
758 |
& + maskC(i-1,j ,k,bi,bj)*mskCor(i-1,j ,bi,bj) ) |
759 |
& + ( maskC(i+1,j-1,k,bi,bj)*mskCor(i+1,j-1,bi,bj) |
760 |
& + maskC(i+1,j ,k,bi,bj)*mskCor(i+1,j ,bi,bj) ) ) |
761 |
& )*maskC(i,j-1,k,bi,bj)*mskCor(i,j-1,bi,bj) |
762 |
& *maskC(i,j ,k,bi,bj)*mskCor(i,j ,bi,bj) |
763 |
#else |
764 |
tmpVisc = _maskS(i,j,k,bi,bj) * halfRL |
765 |
& *( GGL90visctmp(i,j-1,k) |
766 |
& + GGL90visctmp(i,j,k) ) |
767 |
#endif |
768 |
tmpVisc = MIN( tmpVisc , GGL90viscMax ) |
769 |
GGL90viscArV(i,j,k,bi,bj) = MAX( tmpVisc, viscArNr(k) ) |
770 |
ENDDO |
771 |
ENDDO |
772 |
ENDDO |
773 |
|
774 |
DO k=1,Nr |
775 |
DO j=jMin,jMax |
776 |
DO i=iMin,iMax |
777 |
mixingLength(i,j,k,bi,bj) = GGL90mixingLength(i,j,k) / |
778 |
& drF(k) |
779 |
ENDDO |
780 |
ENDDO |
781 |
ENDDO |
782 |
|
783 |
|
784 |
#ifdef ALLOW_DIAGNOSTICS |
785 |
IF ( useDiagnostics ) THEN |
786 |
CALL DIAGNOSTICS_FILL( GGL90TKE ,'GGL90TKE', |
787 |
& 0,Nr, 1, bi, bj, myThid ) |
788 |
CALL DIAGNOSTICS_FILL( GGL90viscArU,'GGL90ArU', |
789 |
& 0,Nr, 1, bi, bj, myThid ) |
790 |
CALL DIAGNOSTICS_FILL( GGL90viscArV,'GGL90ArV', |
791 |
& 0,Nr, 1, bi, bj, myThid ) |
792 |
CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ', |
793 |
& 0,Nr, 1, bi, bj, myThid ) |
794 |
CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl', |
795 |
& 0,Nr, 2, bi, bj, myThid ) |
796 |
CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx', |
797 |
& 0,Nr, 2, bi, bj, myThid ) |
798 |
|
799 |
kp1 = MIN(Nr,kSurf+1) |
800 |
DO j=jMin,jMax |
801 |
DO i=iMin,iMax |
802 |
C diagnose surface flux of TKE |
803 |
surf_flx_tke(i,j) =(GGL90TKE(i,j,kSurf,bi,bj)- |
804 |
& GGL90TKE(i,j,kp1,bi,bj)) |
805 |
& *recip_drF(kSurf)*recip_hFacC(i,j,kSurf,bi,bj) |
806 |
& *KappaE(i,j,kp1) |
807 |
ENDDO |
808 |
ENDDO |
809 |
CALL DIAGNOSTICS_FILL( surf_flx_tke,'GGL90flx', |
810 |
& 0, 1, 2, bi, bj, myThid ) |
811 |
|
812 |
k=kSurf |
813 |
DO j=jMin,jMax |
814 |
DO i=iMin,iMax |
815 |
C diagnose work done by the wind |
816 |
surf_flx_tke(i,j) = |
817 |
& halfRL*( surfaceForcingU(i, j,bi,bj)*uVel(i ,j,k,bi,bj) |
818 |
& +surfaceForcingU(i+1,j,bi,bj)*uVel(i+1,j,k,bi,bj)) |
819 |
& + halfRL*( surfaceForcingV(i,j, bi,bj)*vVel(i,j ,k,bi,bj) |
820 |
& +surfaceForcingV(i,j+1,bi,bj)*vVel(i,j+1,k,bi,bj)) |
821 |
ENDDO |
822 |
ENDDO |
823 |
CALL DIAGNOSTICS_FILL( surf_flx_tke,'GGL90tau', |
824 |
& 0, 1, 2, bi, bj, myThid ) |
825 |
|
826 |
ENDIF |
827 |
#endif /* ALLOW_DIAGNOSTICS */ |
828 |
|
829 |
#endif /* ALLOW_GGL90 */ |
830 |
|
831 |
RETURN |
832 |
END |