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