/[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.21 - (hide annotations) (download)
Wed Jun 27 22:39:09 2012 UTC (11 years, 11 months ago) by jmc
Branch: MAIN
Changes since 1.20: +27 -22 lines
pass sigmaR as argument (but not yet used); will save several RHO computation

1 jmc 1.21 C $Header: /u/gcmpack/MITgcm/pkg/ggl90/ggl90_calc.F,v 1.20 2012/03/15 15:23:22 jmc 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
14     C !DESCRIPTION: \bv
15 jmc 1.12 C *==========================================================*
16 mlosch 1.1 C | SUBROUTINE GGL90_CALC |
17     C | o Compute all GGL90 fields defined in GGL90.h |
18 jmc 1.12 C *==========================================================*
19 mlosch 1.1 C | Equation numbers refer to |
20     C | Gaspar et al. (1990), JGR 95 (C9), pp 16,179 |
21     C | Some parts of the implementation follow Blanke and |
22     C | Delecuse (1993), JPO, and OPA code, in particular the |
23     C | computation of the |
24     C | mixing length = max(min(lk,depth),lkmin) |
25 jmc 1.12 C *==========================================================*
26 mlosch 1.1
27     C global parameters updated by ggl90_calc
28 jmc 1.12 C GGL90TKE :: sub-grid turbulent kinetic energy (m^2/s^2)
29     C GGL90viscAz :: GGL90 eddy viscosity coefficient (m^2/s)
30     C GGL90diffKzT :: GGL90 diffusion coefficient for temperature (m^2/s)
31 mlosch 1.1 C \ev
32    
33     C !USES: ============================================================
34 jmc 1.12 IMPLICIT NONE
35 mlosch 1.1 #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 jmc 1.12 C Routine arguments
45 jmc 1.21 C bi, bj :: Current tile indices
46     C sigmaR :: Vertical gradient of iso-neutral density
47 jmc 1.12 C myTime :: Current time in simulation
48 jmc 1.19 C myIter :: Current time-step number
49 jmc 1.12 C myThid :: My Thread Id number
50 mlosch 1.1 INTEGER bi, bj
51 jmc 1.21 _RL sigmaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
52 jmc 1.12 _RL myTime
53 jmc 1.19 INTEGER myIter
54 mlosch 1.1 INTEGER myThid
55 jmc 1.12 CEOP
56 mlosch 1.1
57     #ifdef ALLOW_GGL90
58    
59     C !LOCAL VARIABLES: ====================================================
60 jmc 1.12 C Local constants
61 jmc 1.19 C iMin,iMax,jMin,jMax :: index boundaries of computation domain
62     C i, j, k, kp1,km1 :: array computation indices
63     C kSurf, kBottom :: vertical indices of domain boundaries
64     C explDissFac :: explicit Dissipation Factor (in [0-1])
65     C implDissFac :: implicit Dissipation Factor (in [0-1])
66     C uStarSquare :: square of friction velocity
67     C verticalShear :: (squared) vertical shear of horizontal velocity
68     C Nsquare :: squared buoyancy freqency
69     C RiNumber :: local Richardson number
70     C KappaM :: (local) viscosity parameter (eq.10)
71     C KappaH :: (local) diffusivity parameter for temperature (eq.11)
72     C KappaE :: (local) diffusivity parameter for TKE (eq.15)
73     C TKEdissipation :: dissipation of TKE
74     C GGL90mixingLength:: mixing length of scheme following Banke+Delecuse
75     C rMixingLength:: inverse of mixing length
76     C totalDepth :: thickness of water column (inverse of recip_Rcol)
77     C TKEPrandtlNumber :: here, an empirical function of the Richardson number
78     C rhoK, rhoKm1 :: density at layer k and km1 (relative to k)
79 mlosch 1.1 INTEGER iMin ,iMax ,jMin ,jMax
80 jmc 1.19 INTEGER i, j, k, kp1, km1, kSurf, kBottom
81     _RL explDissFac, implDissFac
82 mlosch 1.1 _RL uStarSquare
83     _RL verticalShear
84     _RL KappaM, KappaH
85 dfer 1.11 c _RL Nsquare
86     _RL Nsquare(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
87 mlosch 1.1 _RL deltaTggl90
88 dfer 1.11 c _RL SQRTTKE
89     _RL SQRTTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90 mlosch 1.1 _RL RiNumber
91     _RL TKEdissipation
92     _RL tempU, tempV, prTemp
93 gforget 1.16 _RL MaxLength, tmpmlx, tmpVisc
94 mlosch 1.1 _RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
95     _RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
96 dfer 1.13 _RL rMixingLength (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
97     _RL mxLength_Dn (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
98 mlosch 1.1 _RL KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
99     _RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100     _RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
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 dfer 1.13 C- tri-diagonal matrix
104 jmc 1.20 _RL a3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
105     _RL b3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
106     _RL c3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
107 dfer 1.15 INTEGER errCode
108 mlosch 1.2 #ifdef ALLOW_GGL90_HORIZDIFF
109 jmc 1.19 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 mlosch 1.2 _RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113     _RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114     _RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115     _RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 jmc 1.19 _RL gTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 mlosch 1.2 #endif /* ALLOW_GGL90_HORIZDIFF */
118 dfer 1.11 #ifdef ALLOW_GGL90_SMOOTH
119 gforget 1.16 _RL p4, p8, p16
120 dfer 1.11 p4=0.25 _d 0
121     p8=0.125 _d 0
122     p16=0.0625 _d 0
123     #endif
124 mlosch 1.1 iMin = 2-OLx
125     iMax = sNx+OLx-1
126     jMin = 2-OLy
127     jMax = sNy+OLy-1
128    
129     C set separate time step (should be deltaTtracer)
130 jmc 1.4 deltaTggl90 = dTtracerLev(1)
131 jmc 1.12
132 mlosch 1.1 kSurf = 1
133 jmc 1.19 C explicit/implicit timestepping weights for dissipation
134     explDissFac = 0. _d 0
135     implDissFac = 1. _d 0 - explDissFac
136 mlosch 1.1
137     C Initialize local fields
138 jmc 1.19 DO k = 1, Nr
139 jmc 1.21 DO j=1-OLy,sNy+OLy
140     DO i=1-OLx,sNx+OLx
141 jmc 1.19 KappaE(i,j,k) = 0. _d 0
142     TKEPrandtlNumber(i,j,k) = 1. _d 0
143     GGL90mixingLength(i,j,k) = GGL90mixingLengthMin
144     GGL90visctmp(i,j,k) = 0. _d 0
145 jmc 1.20 #ifndef SOLVE_DIAGONAL_LOWMEMORY
146     a3d(i,j,k) = 0. _d 0
147     b3d(i,j,k) = 1. _d 0
148     c3d(i,j,k) = 0. _d 0
149     #endif
150 mlosch 1.1 ENDDO
151 jmc 1.8 ENDDO
152 mlosch 1.1 ENDDO
153 jmc 1.21 DO j=1-OLy,sNy+OLy
154     DO i=1-OLx,sNx+OLx
155 jmc 1.19 rhoK(i,j) = 0. _d 0
156     rhoKm1(i,j) = 0. _d 0
157     totalDepth(i,j) = Ro_surf(i,j,bi,bj) - R_low(i,j,bi,bj)
158 jmc 1.14 rMixingLength(i,j,1) = 0. _d 0
159 jmc 1.19 mxLength_Dn(i,j,1) = GGL90mixingLengthMin
160 jmc 1.14 SQRTTKE(i,j,1) = SQRT( GGL90TKE(i,j,1,bi,bj) )
161 mlosch 1.1 ENDDO
162     ENDDO
163    
164     C start k-loop
165 jmc 1.19 DO k = 2, Nr
166     km1 = k-1
167     c kp1 = MIN(Nr,k+1)
168 jmc 1.8 CALL FIND_RHO_2D(
169 jmc 1.19 I iMin, iMax, jMin, jMax, k,
170     I theta(1-OLx,1-OLy,km1,bi,bj), salt(1-OLx,1-OLy,km1,bi,bj),
171 mlosch 1.1 O rhoKm1,
172 jmc 1.19 I km1, bi, bj, myThid )
173 jmc 1.8
174     CALL FIND_RHO_2D(
175 jmc 1.19 I iMin, iMax, jMin, jMax, k,
176     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
177 mlosch 1.1 O rhoK,
178 jmc 1.19 I k, bi, bj, myThid )
179     DO j=jMin,jMax
180     DO i=iMin,iMax
181     SQRTTKE(i,j,k)=SQRT( GGL90TKE(i,j,k,bi,bj) )
182 jmc 1.14
183 mlosch 1.1 C buoyancy frequency
184 jmc 1.21 Nsquare(i,j,k) = gravity*gravitySign*recip_rhoConst
185     & * sigmaR(i,j,k)
186 jmc 1.19 Nsquare(i,j,k) = - gravity*recip_rhoConst*recip_drC(k)
187     & * ( rhoKm1(i,j) - rhoK(i,j) )*maskC(i,j,k,bi,bj)
188 dfer 1.11 cC vertical shear term (dU/dz)^2+(dV/dz)^2
189 jmc 1.19 c tempU= .5 _d 0*( uVel(i,j,km1,bi,bj)+uVel(i+1,j,km1,bi,bj)
190     c & -( uVel(i,j,k ,bi,bj)+uVel(i+1,j,k ,bi,bj)) )
191     c & *recip_drC(k)
192     c tempV= .5 _d 0*( vVel(i,j,km1,bi,bj)+vVel(i,j+1,km1,bi,bj)
193     c & -( vVel(i,j,k ,bi,bj)+vVel(i,j+1,k ,bi,bj)) )
194     c & *recip_drC(k)
195 dfer 1.11 c verticalShear = tempU*tempU + tempV*tempV
196     c RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
197     cC compute Prandtl number (always greater than 0)
198     c prTemp = 1. _d 0
199     c IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
200 jmc 1.19 c TKEPrandtlNumber(i,j,k) = MIN(10. _d 0,prTemp)
201 dfer 1.11 C mixing length
202 jmc 1.19 GGL90mixingLength(i,j,k) = SQRTTWO *
203 dfer 1.11 & SQRTTKE(i,j,k)/SQRT( MAX(Nsquare(i,j,k),GGL90eps) )
204     ENDDO
205     ENDDO
206     ENDDO
207    
208 dfer 1.13 C- Impose upper and lower bound for mixing length
209 dfer 1.11 IF ( mxlMaxFlag .EQ. 0 ) THEN
210 jmc 1.21
211 dfer 1.11 DO k=2,Nr
212 jmc 1.19 DO j=jMin,jMax
213     DO i=iMin,iMax
214     MaxLength=totalDepth(i,j)
215     GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k),
216 dfer 1.13 & MaxLength)
217 dfer 1.11 ENDDO
218     ENDDO
219     ENDDO
220 dfer 1.13
221     DO k=2,Nr
222 jmc 1.19 DO j=jMin,jMax
223     DO i=iMin,iMax
224     GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k),
225 dfer 1.13 & GGL90mixingLengthMin)
226 jmc 1.19 rMixingLength(i,j,k) = 1. _d 0 / GGL90mixingLength(i,j,k)
227 dfer 1.13 ENDDO
228     ENDDO
229     ENDDO
230    
231 dfer 1.11 ELSEIF ( mxlMaxFlag .EQ. 1 ) THEN
232 jmc 1.21
233 dfer 1.11 DO k=2,Nr
234 jmc 1.19 DO j=jMin,jMax
235     DO i=iMin,iMax
236     MaxLength=MIN(Ro_surf(i,j,bi,bj)-rF(k),rF(k)-R_low(i,j,bi,bj))
237 dfer 1.11 c MaxLength=MAX(MaxLength,20. _d 0)
238 jmc 1.19 GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k),
239 dfer 1.13 & MaxLength)
240     ENDDO
241     ENDDO
242     ENDDO
243    
244     DO k=2,Nr
245 jmc 1.19 DO j=jMin,jMax
246     DO i=iMin,iMax
247     GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k),
248 dfer 1.13 & GGL90mixingLengthMin)
249 jmc 1.19 rMixingLength(i,j,k) = 1. _d 0 / GGL90mixingLength(i,j,k)
250 dfer 1.11 ENDDO
251     ENDDO
252     ENDDO
253 dfer 1.13
254 dfer 1.11 ELSEIF ( mxlMaxFlag .EQ. 2 ) THEN
255 jmc 1.21
256 gforget 1.16 cgf ensure mixing between first and second level
257 jmc 1.19 c DO j=jMin,jMax
258     c DO i=iMin,iMax
259     c GGL90mixingLength(i,j,2)=drF(1)
260 gforget 1.16 c ENDDO
261     c ENDDO
262     cgf
263 dfer 1.11 DO k=2,Nr
264 jmc 1.19 DO j=jMin,jMax
265     DO i=iMin,iMax
266     GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k),
267     & GGL90mixingLength(i,j,k-1)+drF(k-1))
268 dfer 1.11 ENDDO
269     ENDDO
270     ENDDO
271 jmc 1.19 DO j=jMin,jMax
272     DO i=iMin,iMax
273     GGL90mixingLength(i,j,Nr) = MIN(GGL90mixingLength(i,j,Nr),
274 dfer 1.11 & GGL90mixingLengthMin+drF(Nr))
275     ENDDO
276     ENDDO
277     DO k=Nr-1,2,-1
278 jmc 1.19 DO j=jMin,jMax
279     DO i=iMin,iMax
280     GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k),
281     & GGL90mixingLength(i,j,k+1)+drF(k))
282 jmc 1.12 ENDDO
283     ENDDO
284 dfer 1.11 ENDDO
285    
286 dfer 1.13 DO k=2,Nr
287 jmc 1.19 DO j=jMin,jMax
288     DO i=iMin,iMax
289     GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k),
290 dfer 1.13 & GGL90mixingLengthMin)
291 jmc 1.19 rMixingLength(i,j,k) = 1. _d 0 / GGL90mixingLength(i,j,k)
292 dfer 1.13 ENDDO
293     ENDDO
294     ENDDO
295    
296     ELSEIF ( mxlMaxFlag .EQ. 3 ) THEN
297 jmc 1.21
298 dfer 1.13 DO k=2,Nr
299 jmc 1.19 DO j=jMin,jMax
300     DO i=iMin,iMax
301     mxLength_Dn(i,j,k) = MIN(GGL90mixingLength(i,j,k),
302     & mxLength_Dn(i,j,k-1)+drF(k-1))
303 dfer 1.13 ENDDO
304     ENDDO
305     ENDDO
306 jmc 1.19 DO j=jMin,jMax
307     DO i=iMin,iMax
308     GGL90mixingLength(i,j,Nr) = MIN(GGL90mixingLength(i,j,Nr),
309 dfer 1.13 & GGL90mixingLengthMin+drF(Nr))
310     ENDDO
311     ENDDO
312     DO k=Nr-1,2,-1
313 jmc 1.19 DO j=jMin,jMax
314     DO i=iMin,iMax
315     GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k),
316     & GGL90mixingLength(i,j,k+1)+drF(k))
317 dfer 1.13 ENDDO
318 dfer 1.11 ENDDO
319     ENDDO
320 dfer 1.13
321     DO k=2,Nr
322 jmc 1.19 DO j=jMin,jMax
323     DO i=iMin,iMax
324     GGL90mixingLength(i,j,k) = MIN(GGL90mixingLength(i,j,k),
325     & mxLength_Dn(i,j,k))
326     tmpmlx = SQRT( GGL90mixingLength(i,j,k)*mxLength_Dn(i,j,k) )
327 dfer 1.13 tmpmlx = MAX( tmpmlx, GGL90mixingLengthMin)
328 jmc 1.19 rMixingLength(i,j,k) = 1. _d 0 / tmpmlx
329 dfer 1.13 ENDDO
330     ENDDO
331     ENDDO
332    
333     ELSE
334     STOP 'GGL90_CALC: Wrong mxlMaxFlag (mixing length limit)'
335     ENDIF
336    
337     C- Impose minimum mixing length (to avoid division by zero)
338     c DO k=2,Nr
339 jmc 1.19 c DO j=jMin,jMax
340     c DO i=iMin,iMax
341     c GGL90mixingLength(i,j,k) = MAX(GGL90mixingLength(i,j,k),
342 dfer 1.13 c & GGL90mixingLengthMin)
343 jmc 1.19 c rMixingLength(i,j,k) = 1. _d 0 /GGL90mixingLength(i,j,k)
344 dfer 1.13 c ENDDO
345     c ENDDO
346     c ENDDO
347 dfer 1.11
348     DO k=2,Nr
349 jmc 1.19 km1 = k-1
350 dfer 1.11
351 jmc 1.19 #ifdef ALLOW_GGL90_HORIZDIFF
352     IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN
353 jmc 1.8 C horizontal diffusion of TKE (requires an exchange in
354     C do_fields_blocking_exchanges)
355 mlosch 1.2 C common factors
356 jmc 1.21 DO j=1-OLy,sNy+OLy
357     DO i=1-OLx,sNx+OLx
358 mlosch 1.2 xA(i,j) = _dyG(i,j,bi,bj)
359     & *drF(k)*_hFacW(i,j,k,bi,bj)
360     yA(i,j) = _dxG(i,j,bi,bj)
361     & *drF(k)*_hFacS(i,j,k,bi,bj)
362     ENDDO
363 jmc 1.8 ENDDO
364 mlosch 1.2 C Compute diffusive fluxes
365     C ... across x-faces
366 jmc 1.21 DO j=1-OLy,sNy+OLy
367     dfx(1-OLx,j)=0. _d 0
368     DO i=1-OLx+1,sNx+OLx
369 mlosch 1.2 dfx(i,j) = -GGL90diffTKEh*xA(i,j)
370     & *_recip_dxC(i,j,bi,bj)
371     & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj))
372     & *CosFacU(j,bi,bj)
373     ENDDO
374     ENDDO
375     C ... across y-faces
376 jmc 1.21 DO i=1-OLx,sNx+OLx
377     dfy(i,1-OLy)=0. _d 0
378 mlosch 1.2 ENDDO
379 jmc 1.21 DO j=1-OLy+1,sNy+OLy
380     DO i=1-OLx,sNx+OLx
381 mlosch 1.2 dfy(i,j) = -GGL90diffTKEh*yA(i,j)
382     & *_recip_dyC(i,j,bi,bj)
383     & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj))
384     #ifdef ISOTROPIC_COS_SCALING
385     & *CosFacV(j,bi,bj)
386     #endif /* ISOTROPIC_COS_SCALING */
387     ENDDO
388 jmc 1.8 ENDDO
389 mlosch 1.2 C Compute divergence of fluxes
390 jmc 1.21 DO j=1-OLy,sNy+OLy-1
391     DO i=1-OLx,sNx+OLx-1
392 jmc 1.19 gTKE(i,j) =
393     & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
394 mlosch 1.2 & *( (dfx(i+1,j)-dfx(i,j))
395 jmc 1.8 & +(dfy(i,j+1)-dfy(i,j))
396 jmc 1.19 & )
397 jmc 1.8 ENDDO
398 mlosch 1.2 ENDDO
399 jmc 1.19 C end if GGL90diffTKEh .eq. 0.
400     ENDIF
401     #endif /* ALLOW_GGL90_HORIZDIFF */
402    
403     DO j=jMin,jMax
404     DO i=iMin,iMax
405     C vertical shear term (dU/dz)^2+(dV/dz)^2
406     tempU= .5 _d 0*( uVel(i,j,km1,bi,bj)+uVel(i+1,j,km1,bi,bj)
407     & -( uVel(i,j,k ,bi,bj)+uVel(i+1,j,k ,bi,bj)) )
408     & *recip_drC(k)
409     tempV= .5 _d 0*( vVel(i,j,km1,bi,bj)+vVel(i,j+1,km1,bi,bj)
410     & -( vVel(i,j,k ,bi,bj)+vVel(i,j+1,k ,bi,bj)) )
411     & *recip_drC(k)
412     verticalShear = tempU*tempU + tempV*tempV
413     RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
414     C compute Prandtl number (always greater than 0)
415     prTemp = 1. _d 0
416     IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
417     TKEPrandtlNumber(i,j,k) = MIN(10. _d 0,prTemp)
418     c TKEPrandtlNumber(i,j,k) = 1. _d 0
419    
420     C viscosity and diffusivity
421     KappaM = GGL90ck*GGL90mixingLength(i,j,k)*SQRTTKE(i,j,k)
422     GGL90visctmp(i,j,k) = MAX(KappaM,diffKrNrT(k))
423     & * maskC(i,j,k,bi,bj)
424     c note: storing GGL90visctmp like this, and using it later to compute
425     c GGL9rdiffKr etc. is robust in case of smoothing (e.g. see OPA)
426     KappaM = MAX(KappaM,viscArNr(k)) * maskC(i,j,k,bi,bj)
427     KappaH = KappaM/TKEPrandtlNumber(i,j,k)
428     KappaE(i,j,k) = GGL90alpha * KappaM * maskC(i,j,k,bi,bj)
429    
430     C dissipation term
431     TKEdissipation = explDissFac*GGL90ceps
432     & *SQRTTKE(i,j,k)*rMixingLength(i,j,k)
433     & *GGL90TKE(i,j,k,bi,bj)
434     C partial update with sum of explicit contributions
435     GGL90TKE(i,j,k,bi,bj) = GGL90TKE(i,j,k,bi,bj)
436     & + deltaTggl90*(
437     & + KappaM*verticalShear
438     & - KappaH*Nsquare(i,j,k)
439     & - TKEdissipation
440     & )
441     ENDDO
442 mlosch 1.2 ENDDO
443 jmc 1.19
444     #ifdef ALLOW_GGL90_HORIZDIFF
445     IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN
446     C-- Add horiz. diffusion tendency
447     DO j=jMin,jMax
448     DO i=iMin,iMax
449     GGL90TKE(i,j,k,bi,bj) = GGL90TKE(i,j,k,bi,bj)
450     & + gTKE(i,j)*deltaTggl90
451     ENDDO
452     ENDDO
453     ENDIF
454 mlosch 1.2 #endif /* ALLOW_GGL90_HORIZDIFF */
455    
456 jmc 1.19 C-- end of k loop
457     ENDDO
458    
459 mlosch 1.2 C ============================================
460     C Implicit time step to update TKE for k=1,Nr;
461     C TKE(Nr+1)=0 by default
462     C ============================================
463 mlosch 1.1 C set up matrix
464 mlosch 1.2 C-- Lower diagonal
465 mlosch 1.1 DO j=jMin,jMax
466     DO i=iMin,iMax
467 jmc 1.20 a3d(i,j,1) = 0. _d 0
468 mlosch 1.1 ENDDO
469     ENDDO
470     DO k=2,Nr
471 dfer 1.15 km1=MAX(2,k-1)
472 mlosch 1.1 DO j=jMin,jMax
473     DO i=iMin,iMax
474 dfer 1.15 C- We keep recip_hFacC in the diffusive flux calculation,
475     C- but no hFacC in TKE volume control
476     C- No need for maskC(k-1) with recip_hFacC(k-1)
477 jmc 1.20 a3d(i,j,k) = -deltaTggl90
478 dfer 1.11 & *recip_drF(k-1)*recip_hFacC(i,j,k-1,bi,bj)
479 dfer 1.10 & *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1))
480 dfer 1.15 & *recip_drC(k)*maskC(i,j,k,bi,bj)
481 mlosch 1.1 ENDDO
482     ENDDO
483     ENDDO
484 mlosch 1.2 C-- Upper diagonal
485 mlosch 1.1 DO j=jMin,jMax
486     DO i=iMin,iMax
487 jmc 1.20 c3d(i,j,1) = 0. _d 0
488 mlosch 1.1 ENDDO
489     ENDDO
490 dfer 1.11 DO k=2,Nr
491 mlosch 1.1 DO j=jMin,jMax
492     DO i=iMin,iMax
493 dfer 1.15 kp1=MAX(1,MIN(klowC(i,j,bi,bj),k+1))
494     C- We keep recip_hFacC in the diffusive flux calculation,
495     C- but no hFacC in TKE volume control
496     C- No need for maskC(k) with recip_hFacC(k)
497 jmc 1.20 c3d(i,j,k) = -deltaTggl90
498 dfer 1.11 & *recip_drF( k ) * recip_hFacC(i,j,k,bi,bj)
499     & *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1))
500 dfer 1.15 & *recip_drC(k)*maskC(i,j,k-1,bi,bj)
501 mlosch 1.1 ENDDO
502     ENDDO
503     ENDDO
504 mlosch 1.2 C-- Center diagonal
505 mlosch 1.1 DO k=1,Nr
506 dfer 1.15 km1 = MAX(k-1,1)
507 mlosch 1.1 DO j=jMin,jMax
508     DO i=iMin,iMax
509 jmc 1.20 b3d(i,j,k) = 1. _d 0 - c3d(i,j,k) - a3d(i,j,k)
510 jmc 1.19 & + implDissFac*deltaTggl90*GGL90ceps*SQRTTKE(i,j,k)
511     & * rMixingLength(i,j,k)
512 dfer 1.15 & * maskC(i,j,k,bi,bj)*maskC(i,j,km1,bi,bj)
513 mlosch 1.1 ENDDO
514     ENDDO
515     ENDDO
516     C end set up matrix
517    
518     C Apply boundary condition
519 dfer 1.15 kp1 = MIN(Nr,kSurf+1)
520 jmc 1.19 DO j=jMin,jMax
521     DO i=iMin,iMax
522 mlosch 1.1 C estimate friction velocity uStar from surface forcing
523 jmc 1.8 uStarSquare = SQRT(
524 jmc 1.19 & ( .5 _d 0*( surfaceForcingU(i, j, bi,bj)
525     & + surfaceForcingU(i+1,j, bi,bj) ) )**2
526     & + ( .5 _d 0*( surfaceForcingV(i, j, bi,bj)
527     & + surfaceForcingV(i, j+1,bi,bj) ) )**2
528 mlosch 1.1 & )
529     C Dirichlet surface boundary condition for TKE
530 jmc 1.19 GGL90TKE(i,j,kSurf,bi,bj) = maskC(i,j,kSurf,bi,bj)
531     & *MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare)
532     GGL90TKE(i,j,kp1,bi,bj) = GGL90TKE(i,j,kp1,bi,bj)
533 jmc 1.20 & - a3d(i,j,kp1)*GGL90TKE(i,j,kSurf,bi,bj)
534     a3d(i,j,kp1) = 0. _d 0
535 mlosch 1.1 C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom
536 jmc 1.19 kBottom = MAX(kLowC(i,j,bi,bj),1)
537     GGL90TKE(i,j,kBottom,bi,bj) = GGL90TKE(i,j,kBottom,bi,bj)
538 jmc 1.20 & - GGL90TKEbottom*c3d(i,j,kBottom)
539     c3d(i,j,kBottom) = 0. _d 0
540 mlosch 1.1 ENDDO
541 jmc 1.8 ENDDO
542 jmc 1.14
543 jmc 1.19 C solve tri-diagonal system
544 dfer 1.15 CALL SOLVE_TRIDIAGONAL( iMin,iMax, jMin,jMax,
545 jmc 1.20 I a3d, b3d, c3d,
546 jmc 1.19 U GGL90TKE,
547 dfer 1.15 O errCode,
548 jmc 1.19 I bi, bj, myThid )
549 jmc 1.14
550 jmc 1.19 DO k=1,Nr
551     DO j=jMin,jMax
552     DO i=iMin,iMax
553 mlosch 1.1 C impose minimum TKE to avoid numerical undershoots below zero
554 jmc 1.19 GGL90TKE(i,j,k,bi,bj) = maskC(i,j,k,bi,bj)
555     & *MAX( GGL90TKE(i,j,k,bi,bj), GGL90TKEmin )
556 mlosch 1.1 ENDDO
557     ENDDO
558 jmc 1.8 ENDDO
559 dfer 1.11
560 mlosch 1.2 C end of time step
561     C ===============================
562 dfer 1.11
563 jmc 1.19 DO k=2,Nr
564     DO j=1,sNy
565     DO i=1,sNx
566 gforget 1.16 #ifdef ALLOW_GGL90_SMOOTH
567     tmpVisc=
568 dfer 1.11 & (
569 gforget 1.16 & p4 * GGL90visctmp(i ,j ,k) * mskCor(i ,j ,bi,bj)
570     & +p8 *( GGL90visctmp(i-1,j ,k) * mskCor(i-1,j ,bi,bj)
571     & + GGL90visctmp(i ,j-1,k) * mskCor(i ,j-1,bi,bj)
572     & + GGL90visctmp(i+1,j ,k) * mskCor(i+1,j ,bi,bj)
573     & + GGL90visctmp(i ,j+1,k) * mskCor(i ,j+1,bi,bj))
574     & +p16*( GGL90visctmp(i+1,j+1,k) * mskCor(i+1,j+1,bi,bj)
575     & + GGL90visctmp(i+1,j-1,k) * mskCor(i+1,j-1,bi,bj)
576     & + GGL90visctmp(i-1,j+1,k) * mskCor(i-1,j+1,bi,bj)
577     & + GGL90visctmp(i-1,j-1,k) * mskCor(i-1,j-1,bi,bj))
578 dfer 1.11 & )
579     & /(p4
580     & +p8 *( maskC(i-1,j ,k,bi,bj) * mskCor(i-1,j ,bi,bj)
581     & + maskC(i ,j-1,k,bi,bj) * mskCor(i ,j-1,bi,bj)
582     & + maskC(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj)
583     & + maskC(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj))
584     & +p16*( maskC(i+1,j+1,k,bi,bj) * mskCor(i+1,j+1,bi,bj)
585     & + maskC(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj)
586     & + maskC(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
587     & + maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj))
588     & )*maskC(i,j,k,bi,bj)*mskCor(i,j,bi,bj)
589 gforget 1.16 #else
590 jmc 1.19 tmpVisc = GGL90visctmp(i,j,k)
591 gforget 1.16 #endif
592     tmpVisc = MIN(tmpVisc/TKEPrandtlNumber(i,j,k),GGL90diffMax)
593 jmc 1.19 GGL90diffKr(i,j,k,bi,bj)= MAX( tmpVisc , diffKrNrT(k) )
594 dfer 1.11 ENDDO
595     ENDDO
596     ENDDO
597 gforget 1.16
598 jmc 1.19 DO k=2,Nr
599     DO j=1,sNy
600     DO i=1,sNx+1
601 gforget 1.16 #ifdef ALLOW_GGL90_SMOOTH
602 jmc 1.19 tmpVisc =
603 gforget 1.16 & (
604     & p4 *(GGL90visctmp(i ,j ,k) * mskCor(i ,j ,bi,bj)
605     & +GGL90visctmp(i-1,j ,k) * mskCor(i-1,j ,bi,bj))
606     & +p8 *(GGL90visctmp(i-1,j-1,k) * mskCor(i-1,j-1,bi,bj)
607     & +GGL90visctmp(i-1,j+1,k) * mskCor(i-1,j+1,bi,bj)
608     & +GGL90visctmp(i ,j-1,k) * mskCor(i ,j-1,bi,bj)
609     & +GGL90visctmp(i ,j+1,k) * mskCor(i ,j+1,bi,bj))
610     & )
611     & /(p4 * 2. _d 0
612     & +p8 *( maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj)
613     & + maskC(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
614     & + maskC(i ,j-1,k,bi,bj) * mskCor(i ,j-1,bi,bj)
615     & + maskC(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj))
616     & )
617     & *maskC(i ,j,k,bi,bj)*mskCor(i ,j,bi,bj)
618     & *maskC(i-1,j,k,bi,bj)*mskCor(i-1,j,bi,bj)
619     #else
620     tmpVisc = _maskW(i,j,k,bi,bj) *
621     & (.5 _d 0*(GGL90visctmp(i,j,k)
622     & +GGL90visctmp(i-1,j,k))
623     & )
624 dfer 1.11 #endif
625 jmc 1.19 tmpVisc = MIN( tmpVisc , GGL90viscMax )
626     GGL90viscArU(i,j,k,bi,bj) = MAX( tmpVisc, viscArNr(k) )
627 gforget 1.16 ENDDO
628     ENDDO
629     ENDDO
630    
631 jmc 1.19 DO k=2,Nr
632     DO j=1,sNy+1
633     DO i=1,sNx
634 gforget 1.16 #ifdef ALLOW_GGL90_SMOOTH
635     tmpVisc =
636     & (
637     & p4 *(GGL90visctmp(i ,j ,k) * mskCor(i ,j ,bi,bj)
638     & +GGL90visctmp(i ,j-1,k) * mskCor(i ,j-1,bi,bj))
639     & +p8 *(GGL90visctmp(i-1,j ,k) * mskCor(i-1,j ,bi,bj)
640     & +GGL90visctmp(i-1,j-1,k) * mskCor(i-1,j-1,bi,bj)
641     & +GGL90visctmp(i+1,j ,k) * mskCor(i+1,j ,bi,bj)
642     & +GGL90visctmp(i+1,j-1,k) * mskCor(i+1,j-1,bi,bj))
643     & )
644     & /(p4 * 2. _d 0
645     & +p8 *( maskC(i-1,j ,k,bi,bj) * mskCor(i-1,j ,bi,bj)
646     & + maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj)
647     & + maskC(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj)
648     & + maskC(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj))
649     & )
650     & *maskC(i,j ,k,bi,bj)*mskCor(i,j ,bi,bj)
651     & *maskC(i,j-1,k,bi,bj)*mskCor(i,j-1,bi,bj)
652     #else
653     tmpVisc = _maskS(i,j,k,bi,bj) *
654     & (.5 _d 0*(GGL90visctmp(i,j,k)
655     & +GGL90visctmp(i,j-1,k))
656     & )
657    
658     #endif
659     tmpVisc = MIN( tmpVisc , GGL90viscMax )
660 jmc 1.19 GGL90viscArV(i,j,k,bi,bj) = MAX( tmpVisc, viscArNr(k) )
661 gforget 1.16 ENDDO
662     ENDDO
663     ENDDO
664 dfer 1.10
665     #ifdef ALLOW_DIAGNOSTICS
666     IF ( useDiagnostics ) THEN
667     CALL DIAGNOSTICS_FILL( GGL90TKE ,'GGL90TKE',
668     & 0,Nr, 1, bi, bj, myThid )
669 gforget 1.16 CALL DIAGNOSTICS_FILL( GGL90viscArU,'GGL90ArU',
670     & 0,Nr, 1, bi, bj, myThid )
671     CALL DIAGNOSTICS_FILL( GGL90viscArV,'GGL90ArV',
672 dfer 1.10 & 0,Nr, 1, bi, bj, myThid )
673     CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ',
674     & 0,Nr, 1, bi, bj, myThid )
675     CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl',
676     & 0,Nr, 2, bi, bj, myThid )
677     CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx',
678     & 0,Nr, 2, bi, bj, myThid )
679     ENDIF
680     #endif
681 mlosch 1.1
682     #endif /* ALLOW_GGL90 */
683    
684     RETURN
685     END

  ViewVC Help
Powered by ViewVC 1.1.22