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

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

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


Revision 1.12 - (show annotations) (download)
Thu Oct 8 20:07:18 2009 UTC (14 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61w, checkpoint61x
Changes since 1.11: +20 -23 lines
modif for vertical profile of background viscosity

1 C $Header: /u/gcmpack/MITgcm/pkg/ggl90/ggl90_calc.F,v 1.11 2009/01/30 02:23:56 dfer Exp $
2 C $Name: $
3
4 #include "GGL90_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GGL90_CALC
8
9 C !INTERFACE: ======================================================
10 subroutine GGL90_CALC(
11 I bi, bj, myTime, myThid )
12
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | SUBROUTINE GGL90_CALC |
16 C | o Compute all GGL90 fields defined in GGL90.h |
17 C *==========================================================*
18 C | Equation numbers refer to |
19 C | Gaspar et al. (1990), JGR 95 (C9), pp 16,179 |
20 C | Some parts of the implementation follow Blanke and |
21 C | Delecuse (1993), JPO, and OPA code, in particular the |
22 C | computation of the |
23 C | mixing length = max(min(lk,depth),lkmin) |
24 C *==========================================================*
25
26 C global parameters updated by ggl90_calc
27 C GGL90TKE :: sub-grid turbulent kinetic energy (m^2/s^2)
28 C GGL90viscAz :: GGL90 eddy viscosity coefficient (m^2/s)
29 C GGL90diffKzT :: GGL90 diffusion coefficient for temperature (m^2/s)
30 C
31 C \ev
32
33 C !USES: ============================================================
34 IMPLICIT NONE
35 #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 C Routine arguments
45 C bi, bj :: array indices on which to apply calculations
46 C myTime :: Current time in simulation
47 C myThid :: My Thread Id number
48 INTEGER bi, bj
49 _RL myTime
50 INTEGER myThid
51 CEOP
52
53 #ifdef ALLOW_GGL90
54
55 C !LOCAL VARIABLES: ====================================================
56 C Local constants
57 C iMin, iMax, jMin, jMax, I, J - array computation indices
58 C K, Kp1, km1, kSurf, kBottom - vertical loop indices
59 C ab15, ab05 - weights for implicit timestepping
60 C uStarSquare - square of friction velocity
61 C verticalShear - (squared) vertical shear of horizontal velocity
62 C Nsquare - squared buoyancy freqency
63 C RiNumber - local Richardson number
64 C KappaM - (local) viscosity parameter (eq.10)
65 C KappaH - (local) diffusivity parameter for temperature (eq.11)
66 C KappaE - (local) diffusivity parameter for TKE (eq.15)
67 C TKEdissipation - dissipation of TKE
68 C GGL90mixingLength- mixing length of scheme following Banke+Delecuse
69 C rMixingLength- inverse of mixing length
70 C totalDepth - thickness of water column (inverse of recip_Rcol)
71 C TKEPrandtlNumber - here, an empirical function of the Richardson number
72 C rhoK, rhoKm1 - density at layer K and Km1 (relative to K)
73 C gTKE - right hand side of implicit equation
74 INTEGER iMin ,iMax ,jMin ,jMax
75 INTEGER I, J, K, Kp1, Km1, kSurf, kBottom
76 _RL ab15, ab05
77 _RL uStarSquare
78 _RL verticalShear
79 _RL KappaM, KappaH
80 c _RL Nsquare
81 _RL Nsquare(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
82 _RL deltaTggl90
83 c _RL SQRTTKE
84 _RL SQRTTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
85 _RL RiNumber
86 _RL TKEdissipation
87 _RL tempU, tempV, prTemp
88 _RL MaxLength
89 _RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90 _RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
91 _RL rMixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
92 _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 C tri-diagonal matrix
98 _RL a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
99 _RL b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
100 _RL c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
101 #ifdef ALLOW_GGL90_HORIZDIFF
102 C xA, yA - area of lateral faces
103 C dfx, dfy - diffusive flux across lateral faces
104 _RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107 _RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108 #endif /* ALLOW_GGL90_HORIZDIFF */
109 #ifdef ALLOW_GGL90_SMOOTH
110 _RL p4, p8, p16, tmpdiffKrS
111 p4=0.25 _d 0
112 p8=0.125 _d 0
113 p16=0.0625 _d 0
114 #endif
115 iMin = 2-OLx
116 iMax = sNx+OLx-1
117 jMin = 2-OLy
118 jMax = sNy+OLy-1
119
120 C set separate time step (should be deltaTtracer)
121 deltaTggl90 = dTtracerLev(1)
122
123 kSurf = 1
124 C implicit timestepping weights for dissipation
125 ab15 = 1.5 _d 0
126 ab05 = -0.5 _d 0
127 ab15 = 1. _d 0
128 ab05 = 0. _d 0
129
130 C Initialize local fields
131 DO K = 1, Nr
132 DO J=1-Oly,sNy+Oly
133 DO I=1-Olx,sNx+Olx
134 gTKE(I,J,K) = 0. _d 0
135 KappaE(I,J,K) = 0. _d 0
136 TKEPrandtlNumber(I,J,K) = 0. _d 0
137 GGL90mixingLength(I,J,K) = GGL90mixingLengthMin
138 rMixingLength(I,J,K) = 0. _d 0
139 ENDDO
140 ENDDO
141 ENDDO
142 DO J=1-Oly,sNy+Oly
143 DO I=1-Olx,sNx+Olx
144 rhoK (I,J) = 0. _d 0
145 rhoKm1 (I,J) = 0. _d 0
146 totalDepth(I,J) = Ro_surf(i,j,bi,bj) - R_low(i,j,bi,bj)
147 ENDDO
148 ENDDO
149
150 C start k-loop
151 DO K = 2, Nr
152 Km1 = K-1
153 c Kp1 = MIN(Nr,K+1)
154 CALL FIND_RHO_2D(
155 I iMin, iMax, jMin, jMax, K,
156 I theta(1-OLx,1-OLy,Km1,bi,bj), salt(1-OLx,1-OLy,Km1,bi,bj),
157 O rhoKm1,
158 I Km1, bi, bj, myThid )
159
160 CALL FIND_RHO_2D(
161 I iMin, iMax, jMin, jMax, K,
162 I theta(1-OLx,1-OLy,K,bi,bj), salt(1-OLx,1-OLy,K,bi,bj),
163 O rhoK,
164 I K, bi, bj, myThid )
165 DO J=jMin,jMax
166 DO I=iMin,iMax
167 SQRTTKE(i,j,k)=SQRT( GGL90TKE(I,J,K,bi,bj) )
168 C
169 C buoyancy frequency
170 C
171 Nsquare(i,j,k) = - gravity*recip_rhoConst*recip_drC(K)
172 & * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj)
173 cC vertical shear term (dU/dz)^2+(dV/dz)^2
174 c tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
175 c & -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) )
176 c & *recip_drC(K)
177 c tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
178 c & -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) )
179 c & *recip_drC(K)
180 c verticalShear = tempU*tempU + tempV*tempV
181 c RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
182 cC compute Prandtl number (always greater than 0)
183 c prTemp = 1. _d 0
184 c IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
185 c TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp)
186 C mixing length
187 GGL90mixingLength(I,J,K) = SQRTTWO *
188 & SQRTTKE(i,j,k)/SQRT( MAX(Nsquare(i,j,k),GGL90eps) )
189 ENDDO
190 ENDDO
191 ENDDO
192
193 C- Impose upper bound for mixing length (total depth)
194 IF ( mxlMaxFlag .EQ. 0 ) THEN
195 DO k=2,Nr
196 DO J=jMin,jMax
197 DO I=iMin,iMax
198 MaxLength=totalDepth(I,J)
199 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
200 & MaxLength)
201 ENDDO
202 ENDDO
203 ENDDO
204 ELSEIF ( mxlMaxFlag .EQ. 1 ) THEN
205 DO k=2,Nr
206 DO J=jMin,jMax
207 DO I=iMin,iMax
208 MaxLength=MIN(Ro_surf(I,J,bi,bj)-rF(k),rF(k)-R_low(I,J,bi,bj))
209 c MaxLength=MAX(MaxLength,20. _d 0)
210 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
211 & MaxLength)
212 ENDDO
213 ENDDO
214 ENDDO
215 ELSEIF ( mxlMaxFlag .EQ. 2 ) THEN
216 DO k=2,Nr
217 DO J=jMin,jMax
218 DO I=iMin,iMax
219 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
220 & GGL90mixingLength(I,J,K-1)+drF(k-1))
221 ENDDO
222 ENDDO
223 ENDDO
224 DO J=jMin,jMax
225 DO I=iMin,iMax
226 GGL90mixingLength(I,J,Nr) = MIN(GGL90mixingLength(I,J,Nr),
227 & GGL90mixingLengthMin+drF(Nr))
228 ENDDO
229 ENDDO
230 DO k=Nr-1,2,-1
231 DO J=jMin,jMax
232 DO I=iMin,iMax
233 GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
234 & GGL90mixingLength(I,J,K+1)+drF(k))
235 ENDDO
236 ENDDO
237 ENDDO
238 ELSE
239 STOP 'GGL90_CALC: Wrong mxlMaxFlag (mixing lenght limit)'
240 ENDIF
241
242 C- Impose minimum mixing length (to avoid division by zero)
243 DO k=2,Nr
244 DO J=jMin,jMax
245 DO I=iMin,iMax
246 GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
247 & GGL90mixingLengthMin)
248 rMixingLength(I,J,K) = 1. _d 0 /GGL90mixingLength(I,J,K)
249 ENDDO
250 ENDDO
251 ENDDO
252
253 DO k=2,Nr
254 Km1 = K-1
255 DO J=jMin,jMax
256 DO I=iMin,iMax
257 C vertical shear term (dU/dz)^2+(dV/dz)^2
258 tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
259 & -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) )
260 & *recip_drC(K)
261 tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
262 & -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) )
263 & *recip_drC(K)
264 verticalShear = tempU*tempU + tempV*tempV
265 RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
266 C compute Prandtl number (always greater than 0)
267 prTemp = 1. _d 0
268 IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
269 TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp)
270
271 C viscosity and diffusivity
272 KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE(i,j,k)
273 KappaH = KappaM/TKEPrandtlNumber(I,J,K)
274
275 C Set a minium (= background) and maximum value
276 KappaM = MAX(KappaM,viscArNr(k))
277 KappaH = MAX(KappaH,diffKrNrT(k))
278 KappaM = MIN(KappaM,GGL90viscMax)
279 KappaH = MIN(KappaH,GGL90diffMax)
280
281 C Mask land points and storage
282 GGL90viscAr(I,J,K,bi,bj) = KappaM * maskC(I,J,K,bi,bj)
283 GGL90diffKr(I,J,K,bi,bj) = KappaH * maskC(I,J,K,bi,bj)
284 KappaE(I,J,K) = GGL90alpha * GGL90viscAr(I,J,K,bi,bj)
285
286 C dissipation term
287 TKEdissipation = ab05*GGL90ceps
288 & *SQRTTKE(i,j,k)*rMixingLength(I,J,K)
289 & *GGL90TKE(I,J,K,bi,bj)
290 C sum up contributions to form the right hand side
291 gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj)
292 & + deltaTggl90*(
293 & + KappaM*verticalShear
294 & - KappaH*Nsquare(i,j,k)
295 & - TKEdissipation
296 & )
297 ENDDO
298 ENDDO
299 ENDDO
300
301 C horizontal diffusion of TKE (requires an exchange in
302 C do_fields_blocking_exchanges)
303 #ifdef ALLOW_GGL90_HORIZDIFF
304 IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN
305 DO K=2,Nr
306 C common factors
307 DO j=1-Oly,sNy+Oly
308 DO i=1-Olx,sNx+Olx
309 xA(i,j) = _dyG(i,j,bi,bj)
310 & *drF(k)*_hFacW(i,j,k,bi,bj)
311 yA(i,j) = _dxG(i,j,bi,bj)
312 & *drF(k)*_hFacS(i,j,k,bi,bj)
313 ENDDO
314 ENDDO
315 C Compute diffusive fluxes
316 C ... across x-faces
317 DO j=1-Oly,sNy+Oly
318 dfx(1-Olx,j)=0. _d 0
319 DO i=1-Olx+1,sNx+Olx
320 dfx(i,j) = -GGL90diffTKEh*xA(i,j)
321 & *_recip_dxC(i,j,bi,bj)
322 & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj))
323 & *CosFacU(j,bi,bj)
324 ENDDO
325 ENDDO
326 C ... across y-faces
327 DO i=1-Olx,sNx+Olx
328 dfy(i,1-Oly)=0. _d 0
329 ENDDO
330 DO j=1-Oly+1,sNy+Oly
331 DO i=1-Olx,sNx+Olx
332 dfy(i,j) = -GGL90diffTKEh*yA(i,j)
333 & *_recip_dyC(i,j,bi,bj)
334 & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj))
335 #ifdef ISOTROPIC_COS_SCALING
336 & *CosFacV(j,bi,bj)
337 #endif /* ISOTROPIC_COS_SCALING */
338 ENDDO
339 ENDDO
340 C Compute divergence of fluxes
341 DO j=1-Oly,sNy+Oly-1
342 DO i=1-Olx,sNx+Olx-1
343 gTKE(i,j,k)=gTKE(i,j,k)
344 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
345 & *( (dfx(i+1,j)-dfx(i,j))
346 & +(dfy(i,j+1)-dfy(i,j))
347 & )
348 ENDDO
349 ENDDO
350 C end of k-loop
351 ENDDO
352 C end if GGL90diffTKEh .eq. 0.
353 ENDIF
354 #endif /* ALLOW_GGL90_HORIZDIFF */
355
356 C ============================================
357 C Implicit time step to update TKE for k=1,Nr;
358 C TKE(Nr+1)=0 by default
359 C ============================================
360 C set up matrix
361 C-- Lower diagonal
362 DO j=jMin,jMax
363 DO i=iMin,iMax
364 a(i,j,1) = 0. _d 0
365 ENDDO
366 ENDDO
367 DO k=2,Nr
368 km1=max(2,k-1)
369 DO j=jMin,jMax
370 DO i=iMin,iMax
371 a(i,j,k) = -deltaTggl90
372 c & *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj)
373 & *recip_drF(k-1)*recip_hFacC(i,j,k-1,bi,bj)
374 & *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1))
375 & *recip_drC(k)*maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj)
376 ENDDO
377 ENDDO
378 ENDDO
379 C-- Upper diagonal
380 DO j=jMin,jMax
381 DO i=iMin,iMax
382 c(i,j,1) = 0. _d 0
383 ENDDO
384 ENDDO
385 DO k=2,Nr
386 DO j=jMin,jMax
387 DO i=iMin,iMax
388 kp1=min(klowC(i,j,bi,bj),k+1)
389 c(i,j,k) = -deltaTggl90
390 c & *recip_drF( k )*recip_hFacI(i,j,k,bi,bj)
391 & *recip_drF( k ) * recip_hFacC(i,j,k,bi,bj)
392 & *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1))
393 & *recip_drC(k)*maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj)
394 ENDDO
395 ENDDO
396 ENDDO
397 C-- Center diagonal
398 DO k=1,Nr
399 DO j=jMin,jMax
400 DO i=iMin,iMax
401 b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)
402 & + ab15*deltaTggl90*GGL90ceps*SQRT(GGL90TKE(I,J,K,bi,bj))
403 & *rMixingLength(I,J,K)*maskC(i,j,k,bi,bj)
404 ENDDO
405 ENDDO
406 ENDDO
407 C end set up matrix
408
409 C
410 C Apply boundary condition
411 C
412 DO J=jMin,jMax
413 DO I=iMin,iMax
414 C estimate friction velocity uStar from surface forcing
415 uStarSquare = SQRT(
416 & ( .5 _d 0*( surfaceForcingU(I, J, bi,bj)
417 & + surfaceForcingU(I+1,J, bi,bj) ) )**2
418 & + ( .5 _d 0*( surfaceForcingV(I, J, bi,bj)
419 & + surfaceForcingV(I, J+1,bi,bj) ) )**2
420 & )
421 C Dirichlet surface boundary condition for TKE
422 gTKE(I,J,kSurf) = MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare)
423 & *maskC(I,J,kSurf,bi,bj)
424 C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom
425 kBottom = MAX(kLowC(I,J,bi,bj),1)
426 gTKE(I,J,kBottom) = gTKE(I,J,kBottom)
427 & - GGL90TKEbottom*c(I,J,kBottom)
428 c(I,J,kBottom) = 0. _d 0
429 ENDDO
430 ENDDO
431 C
432 C solve tri-diagonal system, and store solution on gTKE (previously rhs)
433 C
434 CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax,
435 I a, b, c,
436 U gTKE,
437 I myThid )
438 C
439 C now update TKE
440 C
441 DO K=1,Nr
442 DO J=jMin,jMax
443 DO I=iMin,iMax
444 C impose minimum TKE to avoid numerical undershoots below zero
445 GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin )
446 & * maskC(I,J,K,bi,bj)
447 ENDDO
448 ENDDO
449 ENDDO
450
451 C end of time step
452 C ===============================
453
454 #ifdef ALLOW_GGL90_SMOOTH
455 DO K=1,Nr
456 DO J=jMin,jMax
457 DO I=iMin,iMax
458 tmpdiffKrS=
459 & (
460 & p4 * GGL90viscAr(i ,j ,k,bi,bj) * mskCor(i ,j ,bi,bj)
461 & +p8 *( GGL90viscAr(i-1,j ,k,bi,bj) * mskCor(i-1,j ,bi,bj)
462 & + GGL90viscAr(i ,j-1,k,bi,bj) * mskCor(i ,j-1,bi,bj)
463 & + GGL90viscAr(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj)
464 & + GGL90viscAr(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj))
465 & +p16*( GGL90viscAr(i+1,j+1,k,bi,bj) * mskCor(i+1,j+1,bi,bj)
466 & + GGL90viscAr(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj)
467 & + GGL90viscAr(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
468 & + GGL90viscAr(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj))
469 & )
470 & /(p4
471 & +p8 *( maskC(i-1,j ,k,bi,bj) * mskCor(i-1,j ,bi,bj)
472 & + maskC(i ,j-1,k,bi,bj) * mskCor(i ,j-1,bi,bj)
473 & + maskC(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj)
474 & + maskC(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj))
475 & +p16*( maskC(i+1,j+1,k,bi,bj) * mskCor(i+1,j+1,bi,bj)
476 & + maskC(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj)
477 & + maskC(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
478 & + maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj))
479 & )*maskC(i,j,k,bi,bj)*mskCor(i,j,bi,bj)
480 & /TKEPrandtlNumber(i,j,k)
481 GGL90diffKrS(I,J,K,bi,bj)= MAX( tmpdiffKrS , diffKrNrT(k) )
482 ENDDO
483 ENDDO
484 ENDDO
485 #endif
486
487 #ifdef ALLOW_DIAGNOSTICS
488 IF ( useDiagnostics ) THEN
489 CALL DIAGNOSTICS_FILL( GGL90TKE ,'GGL90TKE',
490 & 0,Nr, 1, bi, bj, myThid )
491 CALL DIAGNOSTICS_FILL( GGL90viscAr,'GGL90Ar ',
492 & 0,Nr, 1, bi, bj, myThid )
493 CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ',
494 & 0,Nr, 1, bi, bj, myThid )
495 CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl',
496 & 0,Nr, 2, bi, bj, myThid )
497 CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx',
498 & 0,Nr, 2, bi, bj, myThid )
499 ENDIF
500 #endif
501
502 #endif /* ALLOW_GGL90 */
503
504 RETURN
505 END

  ViewVC Help
Powered by ViewVC 1.1.22