/[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.18 - (hide annotations) (download)
Wed Aug 11 03:32:29 2010 UTC (13 years, 10 months ago) by gforget
Branch: MAIN
Changes since 1.17: +2 -2 lines
bug fix. gTKE is a local array, whereas
SOLVE_TRIDIAGONAL expects a global array along with bi,bj.
So 1,1 (not bi,bj) needs to be passed as argument to
SOLVE_TRIDIAGONAL in this case.

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

  ViewVC Help
Powered by ViewVC 1.1.22