/[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.27 - (show annotations) (download)
Thu Feb 19 15:44:12 2015 UTC (9 years, 7 months ago) by mlosch
Branch: MAIN
Changes since 1.26: +103 -17 lines
add IDEMIX (Olbers and Eden, 2013, Eden and Olbers, 2014):
  - code provided by Carsten Eden as an extension of ggl90
  - so far the code is turned on within ggl90 by setting a CPP-flag at
    compile time; a runtime flag implementation is not yet complete

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

  ViewVC Help
Powered by ViewVC 1.1.22