1 |
mlosch |
1.1 |
C $Header: $ |
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 buoyFreq - buoyancy freqency |
69 |
|
|
C TKEdissipation - dissipation of TKE |
70 |
|
|
C GGL90mixingLength- mixing length of scheme following Banke+Delecuse |
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 KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
90 |
|
|
_RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
91 |
|
|
_RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
92 |
|
|
_RL totalDepth (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
93 |
|
|
_RL gTKE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
94 |
|
|
C tri-diagonal matrix |
95 |
|
|
_RL a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
96 |
|
|
_RL b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
97 |
|
|
_RL c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
98 |
|
|
CEOP |
99 |
|
|
iMin = 2-OLx |
100 |
|
|
iMax = sNx+OLx-1 |
101 |
|
|
jMin = 2-OLy |
102 |
|
|
jMax = sNy+OLy-1 |
103 |
|
|
|
104 |
|
|
C set separate time step (should be deltaTtracer) |
105 |
|
|
deltaTggl90 = deltaTtracer |
106 |
|
|
C |
107 |
|
|
kSurf = 1 |
108 |
|
|
C implicit timestepping weights for dissipation |
109 |
|
|
ab15 = 1.5 _d 0 |
110 |
|
|
ab05 = -0.5 _d 0 |
111 |
|
|
ab15 = 1. _d 0 |
112 |
|
|
ab05 = 0. _d 0 |
113 |
|
|
|
114 |
|
|
C Initialize local fields |
115 |
|
|
DO K = 1, Nr |
116 |
|
|
DO J=1-Oly,sNy+Oly |
117 |
|
|
DO I=1-Olx,sNx+Olx |
118 |
|
|
gTKE(I,J,K) = 0. _d 0 |
119 |
|
|
KappaE(I,J,K) = 0. _d 0 |
120 |
|
|
TKEPrandtlNumber(I,J,K) = 0. _d 0 |
121 |
|
|
GGL90mixingLength(I,J,K) = 0. _d 0 |
122 |
|
|
ENDDO |
123 |
|
|
ENDDO |
124 |
|
|
ENDDO |
125 |
|
|
DO J=1-Oly,sNy+Oly |
126 |
|
|
DO I=1-Olx,sNx+Olx |
127 |
|
|
rhoK (I,J) = 0. _d 0 |
128 |
|
|
rhoKm1 (I,J) = 0. _d 0 |
129 |
|
|
totalDepth(I,J) = 0. _d 0 |
130 |
|
|
IF ( recip_Rcol(I,J,bi,bj) .NE. 0. ) |
131 |
|
|
& totalDepth(I,J) = 1./recip_Rcol(I,J,bi,bj) |
132 |
|
|
ENDDO |
133 |
|
|
ENDDO |
134 |
|
|
|
135 |
|
|
C start k-loop |
136 |
|
|
DO K = 2, Nr |
137 |
|
|
Km1 = K-1 |
138 |
|
|
Kp1 = MIN(Nr,K+1) |
139 |
|
|
CALL FIND_RHO( |
140 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, Km1, K, |
141 |
|
|
I theta, salt, |
142 |
|
|
O rhoKm1, |
143 |
|
|
I myThid ) |
144 |
|
|
CALL FIND_RHO( |
145 |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, K, |
146 |
|
|
I theta, salt, |
147 |
|
|
O rhoK, |
148 |
|
|
I myThid ) |
149 |
|
|
DO J=jMin,jMax |
150 |
|
|
DO I=iMin,iMax |
151 |
|
|
SQRTTKE=SQRT( GGL90TKE(I,J,K,bi,bj) ) |
152 |
|
|
C |
153 |
|
|
C buoyancy frequency |
154 |
|
|
C |
155 |
|
|
Nsquare = - gravity*recip_rhoConst*recip_drC(K) |
156 |
|
|
& * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj) |
157 |
|
|
C vertical shear term (dU/dz)^2+(dV/dz)^2 |
158 |
|
|
tempu= .5*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj) |
159 |
|
|
& - (uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) ) |
160 |
|
|
& *recip_drC(K) |
161 |
|
|
tempv= .5*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj) |
162 |
|
|
& - (vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) ) |
163 |
|
|
& *recip_drC(K) |
164 |
|
|
verticalShear = tempU*tempU + tempV*tempV |
165 |
|
|
RiNumber = MAX(Nsquare,0. _d 0)/(verticalShear+GGL90eps) |
166 |
|
|
C compute Prandtl number (always greater than 0) |
167 |
|
|
prTemp = 1. _d 0 |
168 |
|
|
IF ( RiNumber .GE. 0.2 ) prTemp = 5.0 * RiNumber |
169 |
|
|
TKEPrandtlNumber(I,J,K) = MIN(10.,prTemp) |
170 |
|
|
C mixing length |
171 |
|
|
GGL90mixingLength(I,J,K) = |
172 |
|
|
& SQRTTKE/SQRT( MAX(Nsquare,GGL90eps) ) |
173 |
|
|
C impose upper bound for mixing length (total depth) |
174 |
|
|
GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K), |
175 |
|
|
& totalDepth(I,J)) |
176 |
|
|
C impose minimum mixing length (to avoid division by zero) |
177 |
|
|
GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K), |
178 |
|
|
& GGL90mixingLengthMin) |
179 |
|
|
C viscosity of last timestep |
180 |
|
|
KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE |
181 |
|
|
KappaE(I,J,K) = KappaM*GGL90alpha |
182 |
|
|
C dissipation term |
183 |
|
|
TKEdissipation = ab05*GGL90ceps |
184 |
|
|
& *SQRTTKE/GGL90mixingLength(I,J,K) |
185 |
|
|
& *GGL90TKE(I,J,K,bi,bj) |
186 |
|
|
C sum up contributions to form the right hand side |
187 |
|
|
gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj) |
188 |
|
|
& + deltaTggl90*( |
189 |
|
|
& + KappaM*verticalShear |
190 |
|
|
& - KappaM*Nsquare/TKEPrandtlNumber(I,J,K) |
191 |
|
|
& - TKEdissipation |
192 |
|
|
& ) |
193 |
|
|
ENDDO |
194 |
|
|
ENDDO |
195 |
|
|
ENDDO |
196 |
|
|
C |
197 |
|
|
C Implicit time step to update TKE for k=1,Nr; TKE(Nr+1)=0 by default |
198 |
|
|
C |
199 |
|
|
C set up matrix |
200 |
|
|
C-- Old aLower |
201 |
|
|
DO j=jMin,jMax |
202 |
|
|
DO i=iMin,iMax |
203 |
|
|
a(i,j,1) = 0. _d 0 |
204 |
|
|
ENDDO |
205 |
|
|
ENDDO |
206 |
|
|
DO k=2,Nr |
207 |
|
|
km1=MAX(1,k-1) |
208 |
|
|
DO j=jMin,jMax |
209 |
|
|
DO i=iMin,iMax |
210 |
|
|
a(i,j,k) = -deltaTggl90 |
211 |
|
|
& *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj) |
212 |
|
|
& *.5*(KappaE(i,j, k )+KappaE(i,j,km1)) |
213 |
|
|
& *recip_drC(k) |
214 |
|
|
IF (recip_hFacI(i,j,km1,bi,bj).EQ.0.) a(i,j,k)=0. |
215 |
|
|
ENDDO |
216 |
|
|
ENDDO |
217 |
|
|
ENDDO |
218 |
|
|
|
219 |
|
|
C-- Old aUpper |
220 |
|
|
DO j=jMin,jMax |
221 |
|
|
DO i=iMin,iMax |
222 |
|
|
c(i,j,1) = 0. _d 0 |
223 |
|
|
c(i,j,Nr) = 0. _d 0 |
224 |
|
|
ENDDO |
225 |
|
|
ENDDO |
226 |
|
|
CML DO k=1,Nr-1 |
227 |
|
|
DO k=2,Nr |
228 |
|
|
kp1=min(Nr,k+1) |
229 |
|
|
DO j=jMin,jMax |
230 |
|
|
DO i=iMin,iMax |
231 |
|
|
c(i,j,k) = -deltaTggl90 |
232 |
|
|
& *recip_drF( k )*recip_hFacI(i,j,k,bi,bj) |
233 |
|
|
& *.5*(KappaE(i,j,k)+KappaE(i,j,kp1)) |
234 |
|
|
& *recip_drC(k) |
235 |
|
|
C IF (recip_hFacI(i,j,kp1,bi,bj).EQ.0.) c(i,j,k)=0. |
236 |
|
|
ENDDO |
237 |
|
|
ENDDO |
238 |
|
|
ENDDO |
239 |
|
|
|
240 |
|
|
C-- Old aCenter |
241 |
|
|
DO k=1,Nr |
242 |
|
|
DO j=jMin,jMax |
243 |
|
|
DO i=iMin,iMax |
244 |
|
|
b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k) |
245 |
|
|
& + ab15*deltaTggl90*GGL90ceps*SQRT(GGL90TKE(I,J,K,bi,bj)) |
246 |
|
|
& /GGL90mixingLength(I,J,K) |
247 |
|
|
ENDDO |
248 |
|
|
ENDDO |
249 |
|
|
ENDDO |
250 |
|
|
C end set up matrix |
251 |
|
|
|
252 |
|
|
C |
253 |
|
|
C Apply boundary condition |
254 |
|
|
C |
255 |
|
|
DO J=jMin,jMax |
256 |
|
|
DO I=iMin,iMax |
257 |
|
|
C estimate friction velocity uStar from surface forcing |
258 |
|
|
uStarSquare = SQRT( |
259 |
|
|
& ( .5*( surfaceForcingU(I, J, bi,bj) |
260 |
|
|
& + surfaceForcingU(I+1,J, bi,bj) ) )**2 |
261 |
|
|
& + ( .5*( surfaceForcingV(I, J, bi,bj) |
262 |
|
|
& + surfaceForcingV(I, J+1,bi,bj) ) )**2 |
263 |
|
|
& ) |
264 |
|
|
C Dirichlet surface boundary condition for TKE |
265 |
|
|
gTKE(I,J,kSurf) = MAX(GGL90TKEmin,GGL90m2*uStarSquare) |
266 |
|
|
& *maskC(I,J,kSurf,bi,bj) |
267 |
|
|
C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom |
268 |
|
|
kBottom = MIN(MAX(kLowC(I,J,bi,bj),1),Nr) |
269 |
|
|
gTKE(I,J,kBottom) = gTKE(I,J,kBottom) |
270 |
|
|
& - GGL90TKEbottom*c(I,J,kBottom) |
271 |
|
|
ENDDO |
272 |
|
|
ENDDO |
273 |
|
|
C |
274 |
|
|
C solve tri-diagonal system, and store solution on gTKE (previously rhs) |
275 |
|
|
C |
276 |
|
|
CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax, |
277 |
|
|
I a, b, c, |
278 |
|
|
U gTKE, |
279 |
|
|
I myThid ) |
280 |
|
|
C |
281 |
|
|
C now update TKE |
282 |
|
|
C |
283 |
|
|
DO K=1,Nr |
284 |
|
|
DO J=jMin,jMax |
285 |
|
|
DO I=iMin,iMax |
286 |
|
|
C impose minimum TKE to avoid numerical undershoots below zero |
287 |
|
|
GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin ) |
288 |
|
|
& * maskC(I,J,K,bi,bj) |
289 |
|
|
C |
290 |
|
|
C end of time step |
291 |
|
|
C |
292 |
|
|
ENDDO |
293 |
|
|
ENDDO |
294 |
|
|
ENDDO |
295 |
|
|
C |
296 |
|
|
C compute viscosity coefficients |
297 |
|
|
C |
298 |
|
|
DO K=2,Nr |
299 |
|
|
DO J=jMin,jMax |
300 |
|
|
DO I=iMin,iMax |
301 |
|
|
C Eq. (11), (18) and (21) |
302 |
|
|
KappaM = GGL90ck*GGL90mixingLength(I,J,K)* |
303 |
|
|
& SQRT( GGL90TKE(I,J,K,bi,bj) ) |
304 |
|
|
KappaH = KappaM/TKEPrandtlNumber(I,J,K) |
305 |
|
|
C Set a minium (= background) value |
306 |
|
|
KappaM = MAX(KappaM,viscAr) |
307 |
|
|
KappaH = MAX(KappaH,diffKrT) |
308 |
|
|
C Set a maximum and mask land point |
309 |
|
|
GGL90viscAr(I,J,K,bi,bj) = MIN(KappaM,GGL90viscMax) |
310 |
|
|
& * maskC(I,J,K,bi,bj) |
311 |
|
|
GGL90diffKr(I,J,K,bi,bj) = MIN(KappaH,GGL90diffMax) |
312 |
|
|
& * maskC(I,J,K,bi,bj) |
313 |
|
|
ENDDO |
314 |
|
|
ENDDO |
315 |
|
|
C end third k-loop |
316 |
|
|
ENDDO |
317 |
|
|
|
318 |
|
|
#endif /* ALLOW_GGL90 */ |
319 |
|
|
|
320 |
|
|
RETURN |
321 |
|
|
END |
322 |
|
|
|