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

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

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


Revision 1.33 - (hide annotations) (download)
Thu Jul 2 04:42:43 2015 UTC (8 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65n, checkpoint65o
Changes since 1.32: +2 -2 lines
Some minor rearrangement to support AD

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

  ViewVC Help
Powered by ViewVC 1.1.22