/[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.23 - (hide annotations) (download)
Wed Aug 8 22:22:42 2012 UTC (11 years, 10 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint63r, checkpoint63s, checkpoint64
Changes since 1.22: +10 -8 lines
- added run time flag mxlSurfFlag to include the code that ensure
  mixing between first and second level (previously included as a comment)

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

  ViewVC Help
Powered by ViewVC 1.1.22