1 |
C $Header: /u/gcmpack/MITgcm/pkg/ggl90/ggl90_calc.F,v 1.9 2008/09/29 13:47:23 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 |
IMPLICIT NONE |
26 |
C |
27 |
C-------------------------------------------------------------------- |
28 |
|
29 |
C global parameters updated by ggl90_calc |
30 |
C GGL90TKE - sub-grid turbulent kinetic energy (m^2/s^2) |
31 |
C GGL90viscAz - GGL90 eddy viscosity coefficient (m^2/s) |
32 |
C GGL90diffKzT - GGL90 diffusion coefficient for temperature (m^2/s) |
33 |
C |
34 |
C \ev |
35 |
|
36 |
C !USES: ============================================================ |
37 |
#include "SIZE.h" |
38 |
#include "EEPARAMS.h" |
39 |
#include "PARAMS.h" |
40 |
#include "DYNVARS.h" |
41 |
#include "GGL90.h" |
42 |
#include "FFIELDS.h" |
43 |
#include "GRID.h" |
44 |
|
45 |
C !INPUT PARAMETERS: =================================================== |
46 |
c Routine arguments |
47 |
c bi, bj - array indices on which to apply calculations |
48 |
c myTime - Current time in simulation |
49 |
|
50 |
INTEGER bi, bj |
51 |
INTEGER myThid |
52 |
_RL myTime |
53 |
|
54 |
#ifdef ALLOW_GGL90 |
55 |
|
56 |
C !LOCAL VARIABLES: ==================================================== |
57 |
c Local constants |
58 |
C iMin, iMax, jMin, jMax, I, J - array computation indices |
59 |
C K, Kp1, km1, kSurf, kBottom - vertical loop indices |
60 |
C ab15, ab05 - weights for implicit timestepping |
61 |
C uStarSquare - square of friction velocity |
62 |
C verticalShear - (squared) vertical shear of horizontal velocity |
63 |
C Nsquare - squared buoyancy freqency |
64 |
C RiNumber - local Richardson number |
65 |
C KappaM - (local) viscosity parameter (eq.10) |
66 |
C KappaH - (local) diffusivity parameter for temperature (eq.11) |
67 |
C KappaE - (local) diffusivity parameter for TKE (eq.15) |
68 |
C TKEdissipation - dissipation of TKE |
69 |
C GGL90mixingLength- mixing length of scheme following Banke+Delecuse |
70 |
C rMixingLength- inverse of mixing length |
71 |
C totalDepth - thickness of water column (inverse of recip_Rcol) |
72 |
C TKEPrandtlNumber - here, an empirical function of the Richardson number |
73 |
C rhoK, rhoKm1 - density at layer K and Km1 (relative to K) |
74 |
C gTKE - right hand side of implicit equation |
75 |
INTEGER iMin ,iMax ,jMin ,jMax |
76 |
INTEGER I, J, K, Kp1, Km1, kSurf, kBottom |
77 |
_RL ab15, ab05 |
78 |
_RL uStarSquare |
79 |
_RL verticalShear |
80 |
_RL KappaM, KappaH |
81 |
_RL Nsquare |
82 |
_RL deltaTggl90 |
83 |
_RL SQRTTKE |
84 |
_RL RiNumber |
85 |
_RL TKEdissipation |
86 |
_RL tempU, tempV, prTemp |
87 |
_RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
88 |
_RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
89 |
_RL rMixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
90 |
_RL KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
91 |
_RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
92 |
_RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
93 |
_RL totalDepth (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
94 |
_RL gTKE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
95 |
C tri-diagonal matrix |
96 |
_RL a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
97 |
_RL b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
98 |
_RL c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
99 |
#ifdef ALLOW_GGL90_HORIZDIFF |
100 |
C xA, yA - area of lateral faces |
101 |
C dfx, dfy - diffusive flux across lateral faces |
102 |
_RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
103 |
_RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
104 |
_RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
105 |
_RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
106 |
#endif /* ALLOW_GGL90_HORIZDIFF */ |
107 |
CEOP |
108 |
iMin = 2-OLx |
109 |
iMax = sNx+OLx-1 |
110 |
jMin = 2-OLy |
111 |
jMax = sNy+OLy-1 |
112 |
|
113 |
C set separate time step (should be deltaTtracer) |
114 |
deltaTggl90 = dTtracerLev(1) |
115 |
C |
116 |
kSurf = 1 |
117 |
C implicit timestepping weights for dissipation |
118 |
ab15 = 1.5 _d 0 |
119 |
ab05 = -0.5 _d 0 |
120 |
ab15 = 1. _d 0 |
121 |
ab05 = 0. _d 0 |
122 |
|
123 |
C Initialize local fields |
124 |
DO K = 1, Nr |
125 |
DO J=1-Oly,sNy+Oly |
126 |
DO I=1-Olx,sNx+Olx |
127 |
gTKE(I,J,K) = 0. _d 0 |
128 |
KappaE(I,J,K) = 0. _d 0 |
129 |
TKEPrandtlNumber(I,J,K) = 0. _d 0 |
130 |
GGL90mixingLength(I,J,K) = GGL90mixingLengthMin |
131 |
rMixingLength(I,J,K) = 0. _d 0 |
132 |
ENDDO |
133 |
ENDDO |
134 |
ENDDO |
135 |
DO J=1-Oly,sNy+Oly |
136 |
DO I=1-Olx,sNx+Olx |
137 |
rhoK (I,J) = 0. _d 0 |
138 |
rhoKm1 (I,J) = 0. _d 0 |
139 |
totalDepth(I,J) = 0. _d 0 |
140 |
IF ( recip_Rcol(I,J,bi,bj) .NE. 0. _d 0 ) |
141 |
& totalDepth(I,J) = 1./recip_Rcol(I,J,bi,bj) |
142 |
ENDDO |
143 |
ENDDO |
144 |
|
145 |
C start k-loop |
146 |
DO K = 2, Nr |
147 |
Km1 = K-1 |
148 |
Kp1 = MIN(Nr,K+1) |
149 |
CALL FIND_RHO_2D( |
150 |
I iMin, iMax, jMin, jMax, K, |
151 |
I theta(1-OLx,1-OLy,Km1,bi,bj), salt(1-OLx,1-OLy,Km1,bi,bj), |
152 |
O rhoKm1, |
153 |
I Km1, bi, bj, myThid ) |
154 |
|
155 |
CALL FIND_RHO_2D( |
156 |
I iMin, iMax, jMin, jMax, K, |
157 |
I theta(1-OLx,1-OLy,K,bi,bj), salt(1-OLx,1-OLy,K,bi,bj), |
158 |
O rhoK, |
159 |
I K, bi, bj, myThid ) |
160 |
DO J=jMin,jMax |
161 |
DO I=iMin,iMax |
162 |
SQRTTKE=SQRT( GGL90TKE(I,J,K,bi,bj) ) |
163 |
C |
164 |
C buoyancy frequency |
165 |
C |
166 |
Nsquare = - gravity*recip_rhoConst*recip_drC(K) |
167 |
& * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj) |
168 |
C vertical shear term (dU/dz)^2+(dV/dz)^2 |
169 |
tempu= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj) |
170 |
& -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) ) |
171 |
& *recip_drC(K) |
172 |
tempv= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj) |
173 |
& -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) ) |
174 |
& *recip_drC(K) |
175 |
verticalShear = tempU*tempU + tempV*tempV |
176 |
RiNumber = MAX(Nsquare,0. _d 0)/(verticalShear+GGL90eps) |
177 |
C compute Prandtl number (always greater than 0) |
178 |
prTemp = 1. _d 0 |
179 |
IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber |
180 |
TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp) |
181 |
C mixing length |
182 |
GGL90mixingLength(I,J,K) = SQRTTWO * |
183 |
& SQRTTKE/SQRT( MAX(Nsquare,GGL90eps) ) |
184 |
C impose upper bound for mixing length (total depth) |
185 |
GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K), |
186 |
& totalDepth(I,J)) |
187 |
C impose minimum mixing length (to avoid division by zero) |
188 |
GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K), |
189 |
& GGL90mixingLengthMin) |
190 |
rMixingLength(I,J,K) = 1. _d 0 /GGL90mixingLength(I,J,K) |
191 |
C viscosity of last timestep |
192 |
KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE |
193 |
KappaH = KappaM/TKEPrandtlNumber(I,J,K) |
194 |
|
195 |
C Set a minium (= background) and maximum value |
196 |
KappaM = MAX(KappaM,viscAr) |
197 |
KappaH = MAX(KappaH,diffKrNrT(k)) |
198 |
KappaM = MIN(KappaM,GGL90viscMax) |
199 |
KappaH = MIN(KappaH,GGL90diffMax) |
200 |
|
201 |
C Mask land points and save |
202 |
KappaE(I,J,K) = KappaM*GGL90alpha |
203 |
GGL90viscAr(I,J,K,bi,bj) = KappaM * maskC(I,J,K,bi,bj) |
204 |
GGL90diffKr(I,J,K,bi,bj) = KappaH * maskC(I,J,K,bi,bj) |
205 |
|
206 |
C dissipation term |
207 |
TKEdissipation = ab05*GGL90ceps |
208 |
& *SQRTTKE*rMixingLength(I,J,K) |
209 |
& *GGL90TKE(I,J,K,bi,bj) |
210 |
C sum up contributions to form the right hand side |
211 |
gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj) |
212 |
& + deltaTggl90*( |
213 |
& + KappaM*verticalShear |
214 |
& - KappaH*Nsquare |
215 |
c & - KappaM*Nsquare/TKEPrandtlNumber(I,J,K) |
216 |
& - TKEdissipation |
217 |
& ) |
218 |
ENDDO |
219 |
ENDDO |
220 |
ENDDO |
221 |
C horizontal diffusion of TKE (requires an exchange in |
222 |
C do_fields_blocking_exchanges) |
223 |
#ifdef ALLOW_GGL90_HORIZDIFF |
224 |
IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN |
225 |
DO K=2,Nr |
226 |
C common factors |
227 |
DO j=1-Oly,sNy+Oly |
228 |
DO i=1-Olx,sNx+Olx |
229 |
xA(i,j) = _dyG(i,j,bi,bj) |
230 |
& *drF(k)*_hFacW(i,j,k,bi,bj) |
231 |
yA(i,j) = _dxG(i,j,bi,bj) |
232 |
& *drF(k)*_hFacS(i,j,k,bi,bj) |
233 |
ENDDO |
234 |
ENDDO |
235 |
C Compute diffusive fluxes |
236 |
C ... across x-faces |
237 |
DO j=1-Oly,sNy+Oly |
238 |
dfx(1-Olx,j)=0. _d 0 |
239 |
DO i=1-Olx+1,sNx+Olx |
240 |
dfx(i,j) = -GGL90diffTKEh*xA(i,j) |
241 |
& *_recip_dxC(i,j,bi,bj) |
242 |
& *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj)) |
243 |
& *CosFacU(j,bi,bj) |
244 |
ENDDO |
245 |
ENDDO |
246 |
C ... across y-faces |
247 |
DO i=1-Olx,sNx+Olx |
248 |
dfy(i,1-Oly)=0. _d 0 |
249 |
ENDDO |
250 |
DO j=1-Oly+1,sNy+Oly |
251 |
DO i=1-Olx,sNx+Olx |
252 |
dfy(i,j) = -GGL90diffTKEh*yA(i,j) |
253 |
& *_recip_dyC(i,j,bi,bj) |
254 |
& *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj)) |
255 |
#ifdef ISOTROPIC_COS_SCALING |
256 |
& *CosFacV(j,bi,bj) |
257 |
#endif /* ISOTROPIC_COS_SCALING */ |
258 |
ENDDO |
259 |
ENDDO |
260 |
C Compute divergence of fluxes |
261 |
DO j=1-Oly,sNy+Oly-1 |
262 |
DO i=1-Olx,sNx+Olx-1 |
263 |
gTKE(i,j,k)=gTKE(i,j,k) |
264 |
& -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj) |
265 |
& *( (dfx(i+1,j)-dfx(i,j)) |
266 |
& +(dfy(i,j+1)-dfy(i,j)) |
267 |
& ) |
268 |
ENDDO |
269 |
ENDDO |
270 |
C end of k-loop |
271 |
ENDDO |
272 |
C end if GGL90diffTKEh .eq. 0. |
273 |
ENDIF |
274 |
#endif /* ALLOW_GGL90_HORIZDIFF */ |
275 |
|
276 |
C ============================================ |
277 |
C Implicit time step to update TKE for k=1,Nr; |
278 |
C TKE(Nr+1)=0 by default |
279 |
C ============================================ |
280 |
C set up matrix |
281 |
C-- Lower diagonal |
282 |
DO j=jMin,jMax |
283 |
DO i=iMin,iMax |
284 |
a(i,j,1) = 0. _d 0 |
285 |
ENDDO |
286 |
ENDDO |
287 |
DO k=2,Nr |
288 |
km1=MAX(1,k-1) |
289 |
DO j=jMin,jMax |
290 |
DO i=iMin,iMax |
291 |
a(i,j,k) = -deltaTggl90 |
292 |
& *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj) |
293 |
& *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1)) |
294 |
& *recip_drC(k) |
295 |
IF (recip_hFacI(i,j,km1,bi,bj).EQ.0. _d 0) a(i,j,k)=0. _d 0 |
296 |
ENDDO |
297 |
ENDDO |
298 |
ENDDO |
299 |
C-- Upper diagonal |
300 |
DO j=jMin,jMax |
301 |
DO i=iMin,iMax |
302 |
c(i,j,1) = 0. _d 0 |
303 |
c(i,j,Nr) = 0. _d 0 |
304 |
ENDDO |
305 |
ENDDO |
306 |
CML DO k=1,Nr-1 |
307 |
DO k=2,Nr-1 |
308 |
kp1=min(Nr,k+1) |
309 |
DO j=jMin,jMax |
310 |
DO i=iMin,iMax |
311 |
c(i,j,k) = -deltaTggl90 |
312 |
& *recip_drF( k )*recip_hFacI(i,j,k,bi,bj) |
313 |
& *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1)) |
314 |
& *recip_drC(k) |
315 |
IF (recip_hFacI(i,j,kp1,bi,bj).EQ.0. _d 0) c(i,j,k)=0. _d 0 |
316 |
ENDDO |
317 |
ENDDO |
318 |
ENDDO |
319 |
C-- Center diagonal |
320 |
DO k=1,Nr |
321 |
DO j=jMin,jMax |
322 |
DO i=iMin,iMax |
323 |
b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k) |
324 |
& + ab15*deltaTggl90*GGL90ceps*SQRT(GGL90TKE(I,J,K,bi,bj)) |
325 |
& *rMixingLength(I,J,K) |
326 |
ENDDO |
327 |
ENDDO |
328 |
ENDDO |
329 |
C end set up matrix |
330 |
|
331 |
C |
332 |
C Apply boundary condition |
333 |
C |
334 |
DO J=jMin,jMax |
335 |
DO I=iMin,iMax |
336 |
C estimate friction velocity uStar from surface forcing |
337 |
uStarSquare = SQRT( |
338 |
& ( .5 _d 0*( surfaceForcingU(I, J, bi,bj) |
339 |
& + surfaceForcingU(I+1,J, bi,bj) ) )**2 |
340 |
& + ( .5 _d 0*( surfaceForcingV(I, J, bi,bj) |
341 |
& + surfaceForcingV(I, J+1,bi,bj) ) )**2 |
342 |
& ) |
343 |
C Dirichlet surface boundary condition for TKE |
344 |
gTKE(I,J,kSurf) = MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare) |
345 |
& *maskC(I,J,kSurf,bi,bj) |
346 |
C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom |
347 |
kBottom = MIN(MAX(kLowC(I,J,bi,bj),1),Nr) |
348 |
gTKE(I,J,kBottom) = gTKE(I,J,kBottom) |
349 |
& - GGL90TKEbottom*c(I,J,kBottom) |
350 |
ENDDO |
351 |
ENDDO |
352 |
C |
353 |
C solve tri-diagonal system, and store solution on gTKE (previously rhs) |
354 |
C |
355 |
CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax, |
356 |
I a, b, c, |
357 |
U gTKE, |
358 |
I myThid ) |
359 |
C |
360 |
C now update TKE |
361 |
C |
362 |
DO K=1,Nr |
363 |
DO J=jMin,jMax |
364 |
DO I=iMin,iMax |
365 |
C impose minimum TKE to avoid numerical undershoots below zero |
366 |
GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin ) |
367 |
& * maskC(I,J,K,bi,bj) |
368 |
ENDDO |
369 |
ENDDO |
370 |
ENDDO |
371 |
C |
372 |
C end of time step |
373 |
C =============================== |
374 |
C compute viscosity coefficients |
375 |
C |
376 |
c DO K=2,Nr |
377 |
c DO J=jMin,jMax |
378 |
c DO I=iMin,iMax |
379 |
C Eq. (11), (18) and (21) |
380 |
c KappaM = GGL90ck*GGL90mixingLength(I,J,K)* |
381 |
c & SQRT( GGL90TKE(I,J,K,bi,bj) ) |
382 |
c KappaH = KappaM/TKEPrandtlNumber(I,J,K) |
383 |
C Set a minium (= background) value |
384 |
c KappaM = MAX(KappaM,viscAr) |
385 |
c KappaH = MAX(KappaH,diffKrNrT(k)) |
386 |
C Set a maximum and mask land point |
387 |
c GGL90viscAr(I,J,K,bi,bj) = MIN(KappaM,GGL90viscMax) |
388 |
c & * maskC(I,J,K,bi,bj) |
389 |
c GGL90diffKr(I,J,K,bi,bj) = MIN(KappaH,GGL90diffMax) |
390 |
c & * maskC(I,J,K,bi,bj) |
391 |
c ENDDO |
392 |
c ENDDO |
393 |
C end third k-loop |
394 |
c ENDDO |
395 |
|
396 |
#ifdef ALLOW_DIAGNOSTICS |
397 |
IF ( useDiagnostics ) THEN |
398 |
CALL DIAGNOSTICS_FILL( GGL90TKE ,'GGL90TKE', |
399 |
& 0,Nr, 1, bi, bj, myThid ) |
400 |
CALL DIAGNOSTICS_FILL( GGL90viscAr,'GGL90Ar ', |
401 |
& 0,Nr, 1, bi, bj, myThid ) |
402 |
CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ', |
403 |
& 0,Nr, 1, bi, bj, myThid ) |
404 |
CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl', |
405 |
& 0,Nr, 2, bi, bj, myThid ) |
406 |
CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx', |
407 |
& 0,Nr, 2, bi, bj, myThid ) |
408 |
ENDIF |
409 |
#endif |
410 |
|
411 |
#endif /* ALLOW_GGL90 */ |
412 |
|
413 |
RETURN |
414 |
END |
415 |
|