/[MITgcm]/MITgcm/pkg/ggl90/ggl90_calc.F
ViewVC logotype

Contents of /MITgcm/pkg/ggl90/ggl90_calc.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.13 - (show annotations) (download)
Tue Oct 27 15:23:22 2009 UTC (14 years, 7 months ago) by dfer
Branch: MAIN
Changes since 1.12: +99 -23 lines
Add option for mixing length min/max

1 C $Header: /u/gcmpack/MITgcm/pkg/ggl90/ggl90_calc.F,v 1.12 2009/10/08 20:07:18 jmc Exp $
2 C $Name: $
3
4 #include "GGL90_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GGL90_CALC
8
9 C !INTERFACE: ======================================================
10 subroutine GGL90_CALC(
11 I bi, bj, myTime, myThid )
12
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | SUBROUTINE GGL90_CALC |
16 C | o Compute all GGL90 fields defined in GGL90.h |
17 C *==========================================================*
18 C | Equation numbers refer to |
19 C | Gaspar et al. (1990), JGR 95 (C9), pp 16,179 |
20 C | Some parts of the implementation follow Blanke and |
21 C | Delecuse (1993), JPO, and OPA code, in particular the |
22 C | computation of the |
23 C | mixing length = max(min(lk,depth),lkmin) |
24 C *==========================================================*
25
26 C global parameters updated by ggl90_calc
27 C GGL90TKE :: sub-grid turbulent kinetic energy (m^2/s^2)
28 C GGL90viscAz :: GGL90 eddy viscosity coefficient (m^2/s)
29 C GGL90diffKzT :: GGL90 diffusion coefficient for temperature (m^2/s)
30 C
31 C \ev
32
33 C !USES: ============================================================
34 IMPLICIT NONE
35 #include "SIZE.h"
36 #include "EEPARAMS.h"
37 #include "PARAMS.h"
38 #include "DYNVARS.h"
39 #include "GGL90.h"
40 #include "FFIELDS.h"
41 #include "GRID.h"
42
43 C !INPUT PARAMETERS: ===================================================
44 C Routine arguments
45 C bi, bj :: array indices on which to apply calculations
46 C myTime :: Current time in simulation
47 C myThid :: My Thread Id number
48 INTEGER bi, bj
49 _RL myTime
50 INTEGER myThid
51 CEOP
52
53 #ifdef ALLOW_GGL90
54
55 C !LOCAL VARIABLES: ====================================================
56 C Local constants
57 C iMin, iMax, jMin, jMax, I, J - array computation indices
58 C K, Kp1, km1, kSurf, kBottom - vertical loop indices
59 C ab15, ab05 - weights for implicit timestepping
60 C uStarSquare - square of friction velocity
61 C verticalShear - (squared) vertical shear of horizontal velocity
62 C Nsquare - squared buoyancy freqency
63 C RiNumber - local Richardson number
64 C KappaM - (local) viscosity parameter (eq.10)
65 C KappaH - (local) diffusivity parameter for temperature (eq.11)
66 C KappaE - (local) diffusivity parameter for TKE (eq.15)
67 C TKEdissipation - dissipation of TKE
68 C GGL90mixingLength- mixing length of scheme following Banke+Delecuse
69 C rMixingLength- inverse of mixing length
70 C totalDepth - thickness of water column (inverse of recip_Rcol)
71 C TKEPrandtlNumber - here, an empirical function of the Richardson number
72 C rhoK, rhoKm1 - density at layer K and Km1 (relative to K)
73 C gTKE - right hand side of implicit equation
74 INTEGER iMin ,iMax ,jMin ,jMax
75 INTEGER I, J, K, Kp1, Km1, kSurf, kBottom
76 _RL ab15, ab05
77 _RL uStarSquare
78 _RL verticalShear
79 _RL KappaM, KappaH
80 c _RL Nsquare
81 _RL Nsquare(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
82 _RL deltaTggl90
83 c _RL SQRTTKE
84 _RL SQRTTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
85 _RL RiNumber
86 _RL TKEdissipation
87 _RL tempU, tempV, prTemp
88 _RL MaxLength, tmpmlx
89 _RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90 _RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
91 _RL rMixingLength (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
92 _RL mxLength_Dn (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
93 _RL KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
94 _RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL totalDepth (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL gTKE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
98 C- tri-diagonal matrix
99 _RL a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
100 _RL b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
101 _RL c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102 #ifdef ALLOW_GGL90_HORIZDIFF
103 C- xA, yA - area of lateral faces
104 C- dfx, dfy - diffusive flux across lateral faces
105 _RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107 _RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108 _RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109 #endif /* ALLOW_GGL90_HORIZDIFF */
110 #ifdef ALLOW_GGL90_SMOOTH
111 _RL p4, p8, p16, tmpdiffKrS
112 p4=0.25 _d 0
113 p8=0.125 _d 0
114 p16=0.0625 _d 0
115 #endif
116 iMin = 2-OLx
117 iMax = sNx+OLx-1
118 jMin = 2-OLy
119 jMax = sNy+OLy-1
120
121 C set separate time step (should be deltaTtracer)
122 deltaTggl90 = dTtracerLev(1)
123
124 kSurf = 1
125 C implicit timestepping weights for dissipation
126 ab15 = 1.5 _d 0
127 ab05 = -0.5 _d 0
128 ab15 = 1. _d 0
129 ab05 = 0. _d 0
130
131 C Initialize local fields
132 DO K = 1, Nr
133 DO J=1-Oly,sNy+Oly
134 DO I=1-Olx,sNx+Olx
135 gTKE(I,J,K) = 0. _d 0
136 KappaE(I,J,K) = 0. _d 0
137 TKEPrandtlNumber(I,J,K) = 1. _d 0
138 GGL90mixingLength(I,J,K) = GGL90mixingLengthMin
139 rMixingLength(I,J,K) = 0. _d 0
140 ENDDO
141 ENDDO
142 ENDDO
143 DO J=1-Oly,sNy+Oly
144 DO I=1-Olx,sNx+Olx
145 rhoK(I,J) = 0. _d 0
146 rhoKm1(I,J) = 0. _d 0
147 totalDepth(I,J) = Ro_surf(i,j,bi,bj) - R_low(i,j,bi,bj)
148 mxLength_Dn(I,J,1) = GGL90mixingLengthMin
149 ENDDO
150 ENDDO
151
152 C start k-loop
153 DO K = 2, Nr
154 Km1 = K-1
155 c Kp1 = MIN(Nr,K+1)
156 CALL FIND_RHO_2D(
157 I iMin, iMax, jMin, jMax, K,
158 I theta(1-OLx,1-OLy,Km1,bi,bj), salt(1-OLx,1-OLy,Km1,bi,bj),
159 O rhoKm1,
160 I Km1, bi, bj, myThid )
161
162 CALL FIND_RHO_2D(
163 I iMin, iMax, jMin, jMax, K,
164 I theta(1-OLx,1-OLy,K,bi,bj), salt(1-OLx,1-OLy,K,bi,bj),
165 O rhoK,
166 I K, bi, bj, myThid )
167 DO J=jMin,jMax
168 DO I=iMin,iMax
169 SQRTTKE(i,j,k)=SQRT( GGL90TKE(I,J,K,bi,bj) )
170 C
171 C buoyancy frequency
172 C
173 Nsquare(i,j,k) = - gravity*recip_rhoConst*recip_drC(K)
174 & * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj)
175 cC vertical shear term (dU/dz)^2+(dV/dz)^2
176 c tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
177 c & -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) )
178 c & *recip_drC(K)
179 c tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
180 c & -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) )
181 c & *recip_drC(K)
182 c verticalShear = tempU*tempU + tempV*tempV
183 c RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
184 cC compute Prandtl number (always greater than 0)
185 c prTemp = 1. _d 0
186 c IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
187 c TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp)
188 C mixing length
189 GGL90mixingLength(I,J,K) = SQRTTWO *
190 & SQRTTKE(i,j,k)/SQRT( MAX(Nsquare(i,j,k),GGL90eps) )
191 ENDDO
192 ENDDO
193 ENDDO
194
195 C- Impose upper and lower bound for mixing length
196 IF ( mxlMaxFlag .EQ. 0 ) THEN
197 C-
198 DO k=2,Nr
199 DO J=jMin,jMax
200 DO I=iMin,iMax
201 MaxLength=totalDepth(I,J)
202 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
203 & MaxLength)
204 ENDDO
205 ENDDO
206 ENDDO
207
208 DO k=2,Nr
209 DO J=jMin,jMax
210 DO I=iMin,iMax
211 GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
212 & GGL90mixingLengthMin)
213 rMixingLength(I,J,K) = 1. _d 0 / GGL90mixingLength(I,J,K)
214 ENDDO
215 ENDDO
216 ENDDO
217
218 ELSEIF ( mxlMaxFlag .EQ. 1 ) THEN
219 C-
220 DO k=2,Nr
221 DO J=jMin,jMax
222 DO I=iMin,iMax
223 MaxLength=MIN(Ro_surf(I,J,bi,bj)-rF(k),rF(k)-R_low(I,J,bi,bj))
224 c MaxLength=MAX(MaxLength,20. _d 0)
225 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
226 & MaxLength)
227 ENDDO
228 ENDDO
229 ENDDO
230
231 DO k=2,Nr
232 DO J=jMin,jMax
233 DO I=iMin,iMax
234 GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
235 & GGL90mixingLengthMin)
236 rMixingLength(I,J,K) = 1. _d 0 / GGL90mixingLength(I,J,K)
237 ENDDO
238 ENDDO
239 ENDDO
240
241 ELSEIF ( mxlMaxFlag .EQ. 2 ) THEN
242 C-
243 DO k=2,Nr
244 DO J=jMin,jMax
245 DO I=iMin,iMax
246 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
247 & GGL90mixingLength(I,J,K-1)+drF(k-1))
248 ENDDO
249 ENDDO
250 ENDDO
251 DO J=jMin,jMax
252 DO I=iMin,iMax
253 GGL90mixingLength(I,J,Nr) = MIN(GGL90mixingLength(I,J,Nr),
254 & GGL90mixingLengthMin+drF(Nr))
255 ENDDO
256 ENDDO
257 DO k=Nr-1,2,-1
258 DO J=jMin,jMax
259 DO I=iMin,iMax
260 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
261 & GGL90mixingLength(I,J,K+1)+drF(k))
262 ENDDO
263 ENDDO
264 ENDDO
265
266 DO k=2,Nr
267 DO J=jMin,jMax
268 DO I=iMin,iMax
269 GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
270 & GGL90mixingLengthMin)
271 rMixingLength(I,J,K) = 1. _d 0 / GGL90mixingLength(I,J,K)
272 ENDDO
273 ENDDO
274 ENDDO
275
276 ELSEIF ( mxlMaxFlag .EQ. 3 ) THEN
277 C-
278 DO k=2,Nr
279 DO J=jMin,jMax
280 DO I=iMin,iMax
281 mxLength_Dn(I,J,K) = MIN(GGL90mixingLength(I,J,K),
282 & mxLength_Dn(I,J,K-1)+drF(k-1))
283 ENDDO
284 ENDDO
285 ENDDO
286 DO J=jMin,jMax
287 DO I=iMin,iMax
288 GGL90mixingLength(I,J,Nr) = MIN(GGL90mixingLength(I,J,Nr),
289 & GGL90mixingLengthMin+drF(Nr))
290 ENDDO
291 ENDDO
292 DO k=Nr-1,2,-1
293 DO J=jMin,jMax
294 DO I=iMin,iMax
295 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
296 & GGL90mixingLength(I,J,K+1)+drF(k))
297 ENDDO
298 ENDDO
299 ENDDO
300
301 DO k=2,Nr
302 DO J=jMin,jMax
303 DO I=iMin,iMax
304 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
305 & mxLength_Dn(I,J,K))
306 tmpmlx = SQRT(GGL90mixingLength(I,J,K) * mxLength_Dn(I,J,K))
307 tmpmlx = MAX( tmpmlx, GGL90mixingLengthMin)
308 rMixingLength(I,J,K) = 1. _d 0 / tmpmlx
309 ENDDO
310 ENDDO
311 ENDDO
312
313 ELSE
314 STOP 'GGL90_CALC: Wrong mxlMaxFlag (mixing length limit)'
315 ENDIF
316
317 C- Impose minimum mixing length (to avoid division by zero)
318 c DO k=2,Nr
319 c DO J=jMin,jMax
320 c DO I=iMin,iMax
321 c GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
322 c & GGL90mixingLengthMin)
323 c rMixingLength(I,J,K) = 1. _d 0 /GGL90mixingLength(I,J,K)
324 c ENDDO
325 c ENDDO
326 c ENDDO
327
328 DO k=2,Nr
329 Km1 = K-1
330 DO J=jMin,jMax
331 DO I=iMin,iMax
332 C vertical shear term (dU/dz)^2+(dV/dz)^2
333 tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
334 & -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) )
335 & *recip_drC(K)
336 tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
337 & -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) )
338 & *recip_drC(K)
339 verticalShear = tempU*tempU + tempV*tempV
340 RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
341 C compute Prandtl number (always greater than 0)
342 prTemp = 1. _d 0
343 IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
344 TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp)
345 c TKEPrandtlNumber(I,J,K) = 1. _d 0
346
347 C viscosity and diffusivity
348 KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE(i,j,k)
349 KappaH = KappaM/TKEPrandtlNumber(I,J,K)
350
351 C Set a minium (= background) and maximum value
352 KappaM = MAX(KappaM,viscArNr(k))
353 KappaH = MAX(KappaH,diffKrNrT(k))
354 KappaM = MIN(KappaM,GGL90viscMax)
355 KappaH = MIN(KappaH,GGL90diffMax)
356
357 C Mask land points and storage
358 GGL90viscAr(I,J,K,bi,bj) = KappaM * maskC(I,J,K,bi,bj)
359 GGL90diffKr(I,J,K,bi,bj) = KappaH * maskC(I,J,K,bi,bj)
360 KappaE(I,J,K) = GGL90alpha * GGL90viscAr(I,J,K,bi,bj)
361
362 C dissipation term
363 TKEdissipation = ab05*GGL90ceps
364 & *SQRTTKE(i,j,k)*rMixingLength(I,J,K)
365 & *GGL90TKE(I,J,K,bi,bj)
366 C sum up contributions to form the right hand side
367 gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj)
368 & + deltaTggl90*(
369 & + KappaM*verticalShear
370 & - KappaH*Nsquare(i,j,k)
371 & - TKEdissipation
372 & )
373 ENDDO
374 ENDDO
375 ENDDO
376
377 C horizontal diffusion of TKE (requires an exchange in
378 C do_fields_blocking_exchanges)
379 #ifdef ALLOW_GGL90_HORIZDIFF
380 IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN
381 DO K=2,Nr
382 C common factors
383 DO j=1-Oly,sNy+Oly
384 DO i=1-Olx,sNx+Olx
385 xA(i,j) = _dyG(i,j,bi,bj)
386 & *drF(k)*_hFacW(i,j,k,bi,bj)
387 yA(i,j) = _dxG(i,j,bi,bj)
388 & *drF(k)*_hFacS(i,j,k,bi,bj)
389 ENDDO
390 ENDDO
391 C Compute diffusive fluxes
392 C ... across x-faces
393 DO j=1-Oly,sNy+Oly
394 dfx(1-Olx,j)=0. _d 0
395 DO i=1-Olx+1,sNx+Olx
396 dfx(i,j) = -GGL90diffTKEh*xA(i,j)
397 & *_recip_dxC(i,j,bi,bj)
398 & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj))
399 & *CosFacU(j,bi,bj)
400 ENDDO
401 ENDDO
402 C ... across y-faces
403 DO i=1-Olx,sNx+Olx
404 dfy(i,1-Oly)=0. _d 0
405 ENDDO
406 DO j=1-Oly+1,sNy+Oly
407 DO i=1-Olx,sNx+Olx
408 dfy(i,j) = -GGL90diffTKEh*yA(i,j)
409 & *_recip_dyC(i,j,bi,bj)
410 & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj))
411 #ifdef ISOTROPIC_COS_SCALING
412 & *CosFacV(j,bi,bj)
413 #endif /* ISOTROPIC_COS_SCALING */
414 ENDDO
415 ENDDO
416 C Compute divergence of fluxes
417 DO j=1-Oly,sNy+Oly-1
418 DO i=1-Olx,sNx+Olx-1
419 gTKE(i,j,k)=gTKE(i,j,k)
420 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
421 & *( (dfx(i+1,j)-dfx(i,j))
422 & +(dfy(i,j+1)-dfy(i,j))
423 & )
424 ENDDO
425 ENDDO
426 C end of k-loop
427 ENDDO
428 C end if GGL90diffTKEh .eq. 0.
429 ENDIF
430 #endif /* ALLOW_GGL90_HORIZDIFF */
431
432 C ============================================
433 C Implicit time step to update TKE for k=1,Nr;
434 C TKE(Nr+1)=0 by default
435 C ============================================
436 C set up matrix
437 C-- Lower diagonal
438 DO j=jMin,jMax
439 DO i=iMin,iMax
440 a(i,j,1) = 0. _d 0
441 ENDDO
442 ENDDO
443 DO k=2,Nr
444 km1=max(2,k-1)
445 DO j=jMin,jMax
446 DO i=iMin,iMax
447 a(i,j,k) = -deltaTggl90
448 c & *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj)
449 & *recip_drF(k-1)*recip_hFacC(i,j,k-1,bi,bj)
450 & *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1))
451 & *recip_drC(k)*maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj)
452 ENDDO
453 ENDDO
454 ENDDO
455 C-- Upper diagonal
456 DO j=jMin,jMax
457 DO i=iMin,iMax
458 c(i,j,1) = 0. _d 0
459 ENDDO
460 ENDDO
461 DO k=2,Nr
462 DO j=jMin,jMax
463 DO i=iMin,iMax
464 kp1=min(klowC(i,j,bi,bj),k+1)
465 c(i,j,k) = -deltaTggl90
466 c & *recip_drF( k )*recip_hFacI(i,j,k,bi,bj)
467 & *recip_drF( k ) * recip_hFacC(i,j,k,bi,bj)
468 & *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1))
469 & *recip_drC(k)*maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj)
470 ENDDO
471 ENDDO
472 ENDDO
473 C-- Center diagonal
474 DO k=1,Nr
475 DO j=jMin,jMax
476 DO i=iMin,iMax
477 b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)
478 & + ab15*deltaTggl90*GGL90ceps*SQRTTKE(I,J,K)
479 & *rMixingLength(I,J,K)*maskC(i,j,k,bi,bj)
480 ENDDO
481 ENDDO
482 ENDDO
483 C end set up matrix
484
485 C
486 C Apply boundary condition
487 C
488 DO J=jMin,jMax
489 DO I=iMin,iMax
490 C estimate friction velocity uStar from surface forcing
491 uStarSquare = SQRT(
492 & ( .5 _d 0*( surfaceForcingU(I, J, bi,bj)
493 & + surfaceForcingU(I+1,J, bi,bj) ) )**2
494 & + ( .5 _d 0*( surfaceForcingV(I, J, bi,bj)
495 & + surfaceForcingV(I, J+1,bi,bj) ) )**2
496 & )
497 C Dirichlet surface boundary condition for TKE
498 gTKE(I,J,kSurf) = MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare)
499 & *maskC(I,J,kSurf,bi,bj)
500 C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom
501 kBottom = MAX(kLowC(I,J,bi,bj),1)
502 gTKE(I,J,kBottom) = gTKE(I,J,kBottom)
503 & - GGL90TKEbottom*c(I,J,kBottom)
504 c(I,J,kBottom) = 0. _d 0
505 ENDDO
506 ENDDO
507 C
508 C solve tri-diagonal system, and store solution on gTKE (previously rhs)
509 C
510 CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax,
511 I a, b, c,
512 U gTKE,
513 I myThid )
514 C
515 C now update TKE
516 C
517 DO K=1,Nr
518 DO J=jMin,jMax
519 DO I=iMin,iMax
520 C impose minimum TKE to avoid numerical undershoots below zero
521 GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin )
522 & * maskC(I,J,K,bi,bj)
523 ENDDO
524 ENDDO
525 ENDDO
526
527 C end of time step
528 C ===============================
529
530 #ifdef ALLOW_GGL90_SMOOTH
531 DO K=1,Nr
532 DO J=jMin,jMax
533 DO I=iMin,iMax
534 tmpdiffKrS=
535 & (
536 & p4 * GGL90viscAr(i ,j ,k,bi,bj) * mskCor(i ,j ,bi,bj)
537 & +p8 *( GGL90viscAr(i-1,j ,k,bi,bj) * mskCor(i-1,j ,bi,bj)
538 & + GGL90viscAr(i ,j-1,k,bi,bj) * mskCor(i ,j-1,bi,bj)
539 & + GGL90viscAr(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj)
540 & + GGL90viscAr(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj))
541 & +p16*( GGL90viscAr(i+1,j+1,k,bi,bj) * mskCor(i+1,j+1,bi,bj)
542 & + GGL90viscAr(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj)
543 & + GGL90viscAr(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
544 & + GGL90viscAr(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj))
545 & )
546 & /(p4
547 & +p8 *( maskC(i-1,j ,k,bi,bj) * mskCor(i-1,j ,bi,bj)
548 & + maskC(i ,j-1,k,bi,bj) * mskCor(i ,j-1,bi,bj)
549 & + maskC(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj)
550 & + maskC(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj))
551 & +p16*( maskC(i+1,j+1,k,bi,bj) * mskCor(i+1,j+1,bi,bj)
552 & + maskC(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj)
553 & + maskC(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
554 & + maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj))
555 & )*maskC(i,j,k,bi,bj)*mskCor(i,j,bi,bj)
556 & /TKEPrandtlNumber(i,j,k)
557 GGL90diffKrS(I,J,K,bi,bj)= MAX( tmpdiffKrS , diffKrNrT(k) )
558 ENDDO
559 ENDDO
560 ENDDO
561 #endif
562
563 #ifdef ALLOW_DIAGNOSTICS
564 IF ( useDiagnostics ) THEN
565 CALL DIAGNOSTICS_FILL( GGL90TKE ,'GGL90TKE',
566 & 0,Nr, 1, bi, bj, myThid )
567 CALL DIAGNOSTICS_FILL( GGL90viscAr,'GGL90Ar ',
568 & 0,Nr, 1, bi, bj, myThid )
569 CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ',
570 & 0,Nr, 1, bi, bj, myThid )
571 CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl',
572 & 0,Nr, 2, bi, bj, myThid )
573 CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx',
574 & 0,Nr, 2, bi, bj, myThid )
575 ENDIF
576 #endif
577
578 #endif /* ALLOW_GGL90 */
579
580 RETURN
581 END

  ViewVC Help
Powered by ViewVC 1.1.22