1 |
dfer |
1.11 |
C $Header: /u/gcmpack/MITgcm/pkg/ggl90/ggl90_calc.F,v 1.10 2008/10/07 19:35:10 dfer 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 TKEdissipation - dissipation of TKE |
69 |
|
|
C GGL90mixingLength- mixing length of scheme following Banke+Delecuse |
70 |
mlosch |
1.7 |
C rMixingLength- inverse of mixing length |
71 |
mlosch |
1.1 |
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 |
dfer |
1.11 |
c _RL Nsquare |
82 |
|
|
_RL Nsquare(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
83 |
mlosch |
1.1 |
_RL deltaTggl90 |
84 |
dfer |
1.11 |
c _RL SQRTTKE |
85 |
|
|
_RL SQRTTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
86 |
mlosch |
1.1 |
_RL RiNumber |
87 |
|
|
_RL TKEdissipation |
88 |
|
|
_RL tempU, tempV, prTemp |
89 |
dfer |
1.11 |
_RL MaxLength |
90 |
mlosch |
1.1 |
_RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
91 |
|
|
_RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
92 |
mlosch |
1.7 |
_RL rMixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
93 |
mlosch |
1.1 |
_RL KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
94 |
|
|
_RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
95 |
|
|
_RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
96 |
|
|
_RL totalDepth (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
97 |
|
|
_RL gTKE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
98 |
|
|
C tri-diagonal matrix |
99 |
|
|
_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 |
mlosch |
1.2 |
#ifdef ALLOW_GGL90_HORIZDIFF |
103 |
|
|
C xA, yA - area of lateral faces |
104 |
|
|
C dfx, dfy - diffusive flux across lateral faces |
105 |
|
|
_RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
106 |
|
|
_RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
107 |
|
|
_RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
108 |
|
|
_RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
109 |
|
|
#endif /* ALLOW_GGL90_HORIZDIFF */ |
110 |
dfer |
1.11 |
#ifdef ALLOW_GGL90_SMOOTH |
111 |
|
|
_RL p4, p8, p16, tmpdiffKrS |
112 |
|
|
p4=0.25 _d 0 |
113 |
|
|
p8=0.125 _d 0 |
114 |
|
|
p16=0.0625 _d 0 |
115 |
|
|
#endif |
116 |
mlosch |
1.1 |
CEOP |
117 |
|
|
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.8 |
C |
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 |
|
|
TKEPrandtlNumber(I,J,K) = 0. _d 0 |
139 |
mlosch |
1.7 |
GGL90mixingLength(I,J,K) = GGL90mixingLengthMin |
140 |
|
|
rMixingLength(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 |
jmc |
1.8 |
rhoK (I,J) = 0. _d 0 |
147 |
|
|
rhoKm1 (I,J) = 0. _d 0 |
148 |
dfer |
1.11 |
totalDepth(I,J) = Ro_surf(i,j,bi,bj) - R_low(i,j,bi,bj) |
149 |
mlosch |
1.1 |
ENDDO |
150 |
|
|
ENDDO |
151 |
|
|
|
152 |
|
|
C start k-loop |
153 |
|
|
DO K = 2, Nr |
154 |
|
|
Km1 = K-1 |
155 |
dfer |
1.11 |
c Kp1 = MIN(Nr,K+1) |
156 |
jmc |
1.8 |
CALL FIND_RHO_2D( |
157 |
|
|
I iMin, iMax, jMin, jMax, K, |
158 |
|
|
I theta(1-OLx,1-OLy,Km1,bi,bj), salt(1-OLx,1-OLy,Km1,bi,bj), |
159 |
mlosch |
1.1 |
O rhoKm1, |
160 |
jmc |
1.8 |
I Km1, bi, bj, myThid ) |
161 |
|
|
|
162 |
|
|
CALL FIND_RHO_2D( |
163 |
|
|
I iMin, iMax, jMin, jMax, K, |
164 |
|
|
I theta(1-OLx,1-OLy,K,bi,bj), salt(1-OLx,1-OLy,K,bi,bj), |
165 |
mlosch |
1.1 |
O rhoK, |
166 |
jmc |
1.8 |
I K, bi, bj, myThid ) |
167 |
mlosch |
1.1 |
DO J=jMin,jMax |
168 |
|
|
DO I=iMin,iMax |
169 |
dfer |
1.11 |
SQRTTKE(i,j,k)=SQRT( GGL90TKE(I,J,K,bi,bj) ) |
170 |
mlosch |
1.1 |
C |
171 |
|
|
C buoyancy frequency |
172 |
|
|
C |
173 |
dfer |
1.11 |
Nsquare(i,j,k) = - gravity*recip_rhoConst*recip_drC(K) |
174 |
mlosch |
1.1 |
& * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj) |
175 |
dfer |
1.11 |
cC vertical shear term (dU/dz)^2+(dV/dz)^2 |
176 |
|
|
c tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj) |
177 |
|
|
c & -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) ) |
178 |
|
|
c & *recip_drC(K) |
179 |
|
|
c tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj) |
180 |
|
|
c & -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) ) |
181 |
|
|
c & *recip_drC(K) |
182 |
|
|
c verticalShear = tempU*tempU + tempV*tempV |
183 |
|
|
c RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps) |
184 |
|
|
cC compute Prandtl number (always greater than 0) |
185 |
|
|
c prTemp = 1. _d 0 |
186 |
|
|
c IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber |
187 |
|
|
c TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp) |
188 |
|
|
C mixing length |
189 |
|
|
GGL90mixingLength(I,J,K) = SQRTTWO * |
190 |
|
|
& SQRTTKE(i,j,k)/SQRT( MAX(Nsquare(i,j,k),GGL90eps) ) |
191 |
|
|
ENDDO |
192 |
|
|
ENDDO |
193 |
|
|
ENDDO |
194 |
|
|
|
195 |
|
|
C- Impose upper bound for mixing length (total depth) |
196 |
|
|
IF ( mxlMaxFlag .EQ. 0 ) THEN |
197 |
|
|
DO k=2,Nr |
198 |
|
|
DO J=jMin,jMax |
199 |
|
|
DO I=iMin,iMax |
200 |
|
|
MaxLength=totalDepth(I,J) |
201 |
|
|
GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K), |
202 |
|
|
& MaxLength) |
203 |
|
|
ENDDO |
204 |
|
|
ENDDO |
205 |
|
|
ENDDO |
206 |
|
|
ELSEIF ( mxlMaxFlag .EQ. 1 ) THEN |
207 |
|
|
DO k=2,Nr |
208 |
|
|
DO J=jMin,jMax |
209 |
|
|
DO I=iMin,iMax |
210 |
|
|
MaxLength=MIN(Ro_surf(I,J,bi,bj)-rF(k),rF(k)-R_low(I,J,bi,bj)) |
211 |
|
|
c MaxLength=MAX(MaxLength,20. _d 0) |
212 |
|
|
GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K), |
213 |
|
|
& MaxLength) |
214 |
|
|
ENDDO |
215 |
|
|
ENDDO |
216 |
|
|
ENDDO |
217 |
|
|
ELSEIF ( mxlMaxFlag .EQ. 2 ) THEN |
218 |
|
|
DO k=2,Nr |
219 |
|
|
DO J=jMin,jMax |
220 |
|
|
DO I=iMin,iMax |
221 |
|
|
GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K), |
222 |
|
|
& GGL90mixingLength(I,J,K-1)+drF(k-1)) |
223 |
|
|
ENDDO |
224 |
|
|
ENDDO |
225 |
|
|
ENDDO |
226 |
|
|
DO J=jMin,jMax |
227 |
|
|
DO I=iMin,iMax |
228 |
|
|
GGL90mixingLength(I,J,Nr) = MIN(GGL90mixingLength(I,J,Nr), |
229 |
|
|
& GGL90mixingLengthMin+drF(Nr)) |
230 |
|
|
ENDDO |
231 |
|
|
ENDDO |
232 |
|
|
DO k=Nr-1,2,-1 |
233 |
|
|
DO J=jMin,jMax |
234 |
|
|
DO I=iMin,iMax |
235 |
|
|
GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K), |
236 |
|
|
& GGL90mixingLength(I,J,K+1)+drF(k)) |
237 |
|
|
ENDDO |
238 |
|
|
ENDDO |
239 |
|
|
ENDDO |
240 |
|
|
ELSE |
241 |
|
|
STOP 'GGL90_CALC: Wrong mxlMaxFlag (mixing lenght limit)' |
242 |
|
|
ENDIF |
243 |
|
|
|
244 |
|
|
C- Impose minimum mixing length (to avoid division by zero) |
245 |
|
|
DO k=2,Nr |
246 |
|
|
DO J=jMin,jMax |
247 |
|
|
DO I=iMin,iMax |
248 |
|
|
GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K), |
249 |
|
|
& GGL90mixingLengthMin) |
250 |
|
|
rMixingLength(I,J,K) = 1. _d 0 /GGL90mixingLength(I,J,K) |
251 |
|
|
ENDDO |
252 |
|
|
ENDDO |
253 |
|
|
ENDDO |
254 |
|
|
|
255 |
|
|
DO k=2,Nr |
256 |
|
|
Km1 = K-1 |
257 |
|
|
DO J=jMin,jMax |
258 |
|
|
DO I=iMin,iMax |
259 |
mlosch |
1.1 |
C vertical shear term (dU/dz)^2+(dV/dz)^2 |
260 |
dfer |
1.11 |
tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj) |
261 |
dfer |
1.10 |
& -( uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) ) |
262 |
mlosch |
1.1 |
& *recip_drC(K) |
263 |
dfer |
1.11 |
tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj) |
264 |
dfer |
1.10 |
& -( vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) ) |
265 |
mlosch |
1.1 |
& *recip_drC(K) |
266 |
|
|
verticalShear = tempU*tempU + tempV*tempV |
267 |
dfer |
1.11 |
RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps) |
268 |
mlosch |
1.1 |
C compute Prandtl number (always greater than 0) |
269 |
|
|
prTemp = 1. _d 0 |
270 |
dfer |
1.10 |
IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber |
271 |
|
|
TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp) |
272 |
dfer |
1.11 |
|
273 |
|
|
C viscosity and diffusivity |
274 |
|
|
KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE(i,j,k) |
275 |
dfer |
1.10 |
KappaH = KappaM/TKEPrandtlNumber(I,J,K) |
276 |
|
|
|
277 |
|
|
C Set a minium (= background) and maximum value |
278 |
|
|
KappaM = MAX(KappaM,viscAr) |
279 |
|
|
KappaH = MAX(KappaH,diffKrNrT(k)) |
280 |
|
|
KappaM = MIN(KappaM,GGL90viscMax) |
281 |
|
|
KappaH = MIN(KappaH,GGL90diffMax) |
282 |
|
|
|
283 |
dfer |
1.11 |
C Mask land points and storage |
284 |
dfer |
1.10 |
GGL90viscAr(I,J,K,bi,bj) = KappaM * maskC(I,J,K,bi,bj) |
285 |
|
|
GGL90diffKr(I,J,K,bi,bj) = KappaH * maskC(I,J,K,bi,bj) |
286 |
dfer |
1.11 |
KappaE(I,J,K) = GGL90alpha * GGL90viscAr(I,J,K,bi,bj) |
287 |
dfer |
1.10 |
|
288 |
mlosch |
1.1 |
C dissipation term |
289 |
|
|
TKEdissipation = ab05*GGL90ceps |
290 |
dfer |
1.11 |
& *SQRTTKE(i,j,k)*rMixingLength(I,J,K) |
291 |
jmc |
1.8 |
& *GGL90TKE(I,J,K,bi,bj) |
292 |
mlosch |
1.1 |
C sum up contributions to form the right hand side |
293 |
jmc |
1.8 |
gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj) |
294 |
mlosch |
1.1 |
& + deltaTggl90*( |
295 |
|
|
& + KappaM*verticalShear |
296 |
dfer |
1.11 |
& - KappaH*Nsquare(i,j,k) |
297 |
jmc |
1.8 |
& - TKEdissipation |
298 |
mlosch |
1.1 |
& ) |
299 |
jmc |
1.8 |
ENDDO |
300 |
mlosch |
1.1 |
ENDDO |
301 |
mlosch |
1.2 |
ENDDO |
302 |
dfer |
1.11 |
|
303 |
jmc |
1.8 |
C horizontal diffusion of TKE (requires an exchange in |
304 |
|
|
C do_fields_blocking_exchanges) |
305 |
mlosch |
1.2 |
#ifdef ALLOW_GGL90_HORIZDIFF |
306 |
|
|
IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN |
307 |
|
|
DO K=2,Nr |
308 |
|
|
C common factors |
309 |
|
|
DO j=1-Oly,sNy+Oly |
310 |
|
|
DO i=1-Olx,sNx+Olx |
311 |
|
|
xA(i,j) = _dyG(i,j,bi,bj) |
312 |
|
|
& *drF(k)*_hFacW(i,j,k,bi,bj) |
313 |
|
|
yA(i,j) = _dxG(i,j,bi,bj) |
314 |
|
|
& *drF(k)*_hFacS(i,j,k,bi,bj) |
315 |
|
|
ENDDO |
316 |
jmc |
1.8 |
ENDDO |
317 |
mlosch |
1.2 |
C Compute diffusive fluxes |
318 |
|
|
C ... across x-faces |
319 |
|
|
DO j=1-Oly,sNy+Oly |
320 |
dfer |
1.10 |
dfx(1-Olx,j)=0. _d 0 |
321 |
mlosch |
1.2 |
DO i=1-Olx+1,sNx+Olx |
322 |
|
|
dfx(i,j) = -GGL90diffTKEh*xA(i,j) |
323 |
|
|
& *_recip_dxC(i,j,bi,bj) |
324 |
|
|
& *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj)) |
325 |
|
|
& *CosFacU(j,bi,bj) |
326 |
|
|
ENDDO |
327 |
|
|
ENDDO |
328 |
|
|
C ... across y-faces |
329 |
|
|
DO i=1-Olx,sNx+Olx |
330 |
dfer |
1.10 |
dfy(i,1-Oly)=0. _d 0 |
331 |
mlosch |
1.2 |
ENDDO |
332 |
|
|
DO j=1-Oly+1,sNy+Oly |
333 |
|
|
DO i=1-Olx,sNx+Olx |
334 |
|
|
dfy(i,j) = -GGL90diffTKEh*yA(i,j) |
335 |
|
|
& *_recip_dyC(i,j,bi,bj) |
336 |
|
|
& *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj)) |
337 |
|
|
#ifdef ISOTROPIC_COS_SCALING |
338 |
|
|
& *CosFacV(j,bi,bj) |
339 |
|
|
#endif /* ISOTROPIC_COS_SCALING */ |
340 |
|
|
ENDDO |
341 |
jmc |
1.8 |
ENDDO |
342 |
mlosch |
1.2 |
C Compute divergence of fluxes |
343 |
|
|
DO j=1-Oly,sNy+Oly-1 |
344 |
|
|
DO i=1-Olx,sNx+Olx-1 |
345 |
|
|
gTKE(i,j,k)=gTKE(i,j,k) |
346 |
|
|
& -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj) |
347 |
|
|
& *( (dfx(i+1,j)-dfx(i,j)) |
348 |
jmc |
1.8 |
& +(dfy(i,j+1)-dfy(i,j)) |
349 |
mlosch |
1.2 |
& ) |
350 |
jmc |
1.8 |
ENDDO |
351 |
mlosch |
1.2 |
ENDDO |
352 |
jmc |
1.8 |
C end of k-loop |
353 |
mlosch |
1.2 |
ENDDO |
354 |
|
|
C end if GGL90diffTKEh .eq. 0. |
355 |
|
|
ENDIF |
356 |
|
|
#endif /* ALLOW_GGL90_HORIZDIFF */ |
357 |
|
|
|
358 |
|
|
C ============================================ |
359 |
|
|
C Implicit time step to update TKE for k=1,Nr; |
360 |
|
|
C TKE(Nr+1)=0 by default |
361 |
|
|
C ============================================ |
362 |
mlosch |
1.1 |
C set up matrix |
363 |
mlosch |
1.2 |
C-- Lower diagonal |
364 |
mlosch |
1.1 |
DO j=jMin,jMax |
365 |
|
|
DO i=iMin,iMax |
366 |
jmc |
1.8 |
a(i,j,1) = 0. _d 0 |
367 |
mlosch |
1.1 |
ENDDO |
368 |
|
|
ENDDO |
369 |
|
|
DO k=2,Nr |
370 |
dfer |
1.11 |
km1=max(2,k-1) |
371 |
mlosch |
1.1 |
DO j=jMin,jMax |
372 |
|
|
DO i=iMin,iMax |
373 |
|
|
a(i,j,k) = -deltaTggl90 |
374 |
dfer |
1.11 |
c & *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj) |
375 |
|
|
& *recip_drF(k-1)*recip_hFacC(i,j,k-1,bi,bj) |
376 |
dfer |
1.10 |
& *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1)) |
377 |
dfer |
1.11 |
& *recip_drC(k)*maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj) |
378 |
mlosch |
1.1 |
ENDDO |
379 |
|
|
ENDDO |
380 |
|
|
ENDDO |
381 |
mlosch |
1.2 |
C-- Upper diagonal |
382 |
mlosch |
1.1 |
DO j=jMin,jMax |
383 |
|
|
DO i=iMin,iMax |
384 |
|
|
c(i,j,1) = 0. _d 0 |
385 |
|
|
ENDDO |
386 |
|
|
ENDDO |
387 |
dfer |
1.11 |
DO k=2,Nr |
388 |
mlosch |
1.1 |
DO j=jMin,jMax |
389 |
|
|
DO i=iMin,iMax |
390 |
dfer |
1.11 |
kp1=min(klowC(i,j,bi,bj),k+1) |
391 |
mlosch |
1.1 |
c(i,j,k) = -deltaTggl90 |
392 |
dfer |
1.11 |
c & *recip_drF( k )*recip_hFacI(i,j,k,bi,bj) |
393 |
|
|
& *recip_drF( k ) * recip_hFacC(i,j,k,bi,bj) |
394 |
|
|
& *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1)) |
395 |
|
|
& *recip_drC(k)*maskC(i,j,k,bi,bj)*maskC(i,j,k-1,bi,bj) |
396 |
mlosch |
1.1 |
ENDDO |
397 |
|
|
ENDDO |
398 |
|
|
ENDDO |
399 |
mlosch |
1.2 |
C-- Center diagonal |
400 |
mlosch |
1.1 |
DO k=1,Nr |
401 |
|
|
DO j=jMin,jMax |
402 |
|
|
DO i=iMin,iMax |
403 |
|
|
b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k) |
404 |
|
|
& + ab15*deltaTggl90*GGL90ceps*SQRT(GGL90TKE(I,J,K,bi,bj)) |
405 |
dfer |
1.11 |
& *rMixingLength(I,J,K)*maskC(i,j,k,bi,bj) |
406 |
mlosch |
1.1 |
ENDDO |
407 |
|
|
ENDDO |
408 |
|
|
ENDDO |
409 |
|
|
C end set up matrix |
410 |
|
|
|
411 |
|
|
C |
412 |
|
|
C Apply boundary condition |
413 |
jmc |
1.8 |
C |
414 |
mlosch |
1.1 |
DO J=jMin,jMax |
415 |
|
|
DO I=iMin,iMax |
416 |
|
|
C estimate friction velocity uStar from surface forcing |
417 |
jmc |
1.8 |
uStarSquare = SQRT( |
418 |
dfer |
1.10 |
& ( .5 _d 0*( surfaceForcingU(I, J, bi,bj) |
419 |
mlosch |
1.1 |
& + surfaceForcingU(I+1,J, bi,bj) ) )**2 |
420 |
dfer |
1.10 |
& + ( .5 _d 0*( surfaceForcingV(I, J, bi,bj) |
421 |
mlosch |
1.1 |
& + surfaceForcingV(I, J+1,bi,bj) ) )**2 |
422 |
|
|
& ) |
423 |
|
|
C Dirichlet surface boundary condition for TKE |
424 |
mlosch |
1.6 |
gTKE(I,J,kSurf) = MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare) |
425 |
mlosch |
1.1 |
& *maskC(I,J,kSurf,bi,bj) |
426 |
|
|
C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom |
427 |
dfer |
1.11 |
kBottom = MAX(kLowC(I,J,bi,bj),1) |
428 |
jmc |
1.8 |
gTKE(I,J,kBottom) = gTKE(I,J,kBottom) |
429 |
dfer |
1.11 |
& - GGL90TKEbottom*c(I,J,kBottom) |
430 |
|
|
c(I,J,kBottom) = 0. _d 0 |
431 |
mlosch |
1.1 |
ENDDO |
432 |
jmc |
1.8 |
ENDDO |
433 |
mlosch |
1.1 |
C |
434 |
|
|
C solve tri-diagonal system, and store solution on gTKE (previously rhs) |
435 |
|
|
C |
436 |
|
|
CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax, |
437 |
|
|
I a, b, c, |
438 |
|
|
U gTKE, |
439 |
|
|
I myThid ) |
440 |
|
|
C |
441 |
|
|
C now update TKE |
442 |
jmc |
1.8 |
C |
443 |
mlosch |
1.1 |
DO K=1,Nr |
444 |
|
|
DO J=jMin,jMax |
445 |
|
|
DO I=iMin,iMax |
446 |
|
|
C impose minimum TKE to avoid numerical undershoots below zero |
447 |
jmc |
1.8 |
GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin ) |
448 |
mlosch |
1.1 |
& * maskC(I,J,K,bi,bj) |
449 |
|
|
ENDDO |
450 |
|
|
ENDDO |
451 |
jmc |
1.8 |
ENDDO |
452 |
dfer |
1.11 |
|
453 |
mlosch |
1.2 |
C end of time step |
454 |
|
|
C =============================== |
455 |
dfer |
1.11 |
|
456 |
|
|
#ifdef ALLOW_GGL90_SMOOTH |
457 |
|
|
DO K=1,Nr |
458 |
|
|
DO J=jMin,jMax |
459 |
|
|
DO I=iMin,iMax |
460 |
|
|
tmpdiffKrS= |
461 |
|
|
& ( |
462 |
|
|
& p4 * GGL90viscAr(i ,j ,k,bi,bj) * mskCor(i ,j ,bi,bj) |
463 |
|
|
& +p8 *( 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 |
|
|
& + GGL90viscAr(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj) |
466 |
|
|
& + GGL90viscAr(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj)) |
467 |
|
|
& +p16*( 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 |
|
|
& + GGL90viscAr(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj) |
470 |
|
|
& + GGL90viscAr(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj)) |
471 |
|
|
& ) |
472 |
|
|
& /(p4 |
473 |
|
|
& +p8 *( 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 |
|
|
& + maskC(i+1,j ,k,bi,bj) * mskCor(i+1,j ,bi,bj) |
476 |
|
|
& + maskC(i ,j+1,k,bi,bj) * mskCor(i ,j+1,bi,bj)) |
477 |
|
|
& +p16*( 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-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj) |
480 |
|
|
& + maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj)) |
481 |
|
|
& )*maskC(i,j,k,bi,bj)*mskCor(i,j,bi,bj) |
482 |
|
|
& /TKEPrandtlNumber(i,j,k) |
483 |
|
|
GGL90diffKrS(I,J,K,bi,bj)= MAX( tmpdiffKrS , diffKrNrT(k) ) |
484 |
|
|
ENDDO |
485 |
|
|
ENDDO |
486 |
|
|
ENDDO |
487 |
|
|
#endif |
488 |
dfer |
1.10 |
|
489 |
|
|
#ifdef ALLOW_DIAGNOSTICS |
490 |
|
|
IF ( useDiagnostics ) THEN |
491 |
|
|
CALL DIAGNOSTICS_FILL( GGL90TKE ,'GGL90TKE', |
492 |
|
|
& 0,Nr, 1, bi, bj, myThid ) |
493 |
|
|
CALL DIAGNOSTICS_FILL( GGL90viscAr,'GGL90Ar ', |
494 |
|
|
& 0,Nr, 1, bi, bj, myThid ) |
495 |
|
|
CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ', |
496 |
|
|
& 0,Nr, 1, bi, bj, myThid ) |
497 |
|
|
CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl', |
498 |
|
|
& 0,Nr, 2, bi, bj, myThid ) |
499 |
|
|
CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx', |
500 |
|
|
& 0,Nr, 2, bi, bj, myThid ) |
501 |
|
|
ENDIF |
502 |
|
|
#endif |
503 |
mlosch |
1.1 |
|
504 |
|
|
#endif /* ALLOW_GGL90 */ |
505 |
|
|
|
506 |
|
|
RETURN |
507 |
|
|
END |
508 |
|
|
|