/[MITgcm]/MITgcm/pkg/gmredi/gmredi_k3d.F
ViewVC logotype

Contents of /MITgcm/pkg/gmredi/gmredi_k3d.F

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


Revision 1.5 - (show annotations) (download)
Thu Aug 8 22:39:36 2013 UTC (10 years, 10 months ago) by m_bates
Branch: MAIN
CVS Tags: checkpoint64m
Changes since 1.4: +3 -3 lines
Namelist option to ignore beta in the calculation of grad(q) for GM_K3D

1 C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_k3d.F,v 1.4 2013/07/11 14:33:23 m_bates Exp $
2 C $Name: $
3 #include "GMREDI_OPTIONS.h"
4
5 C !ROUTINE: GMREDI_K3D
6 C !INTERFACE:
7 SUBROUTINE GMREDI_K3D(
8 I iMin, iMax, jMin, jMax,
9 I sigmaX, sigmaY, sigmaR,
10 I bi, bj, myTime, myThid )
11
12 C !DESCRIPTION: \bv
13 C *==========================================================*
14 C | SUBROUTINE GMREDI_K3D
15 C | o Calculates the 3D diffusivity as per Bates et al. (2013)
16 C *==========================================================*
17 C \ev
18
19 IMPLICIT NONE
20
21 C == Global variables ==
22 #include "SIZE.h"
23 #include "GRID.h"
24 #include "DYNVARS.h"
25 #include "EEPARAMS.h"
26 #include "PARAMS.h"
27 #include "GMREDI.h"
28
29 C !INPUT/OUTPUT PARAMETERS:
30 C == Routine arguments ==
31 C bi, bj :: tile indices
32 C myThid :: My Thread Id. number
33
34 INTEGER iMin,iMax,jMin,jMax
35 _RL sigmaX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
36 _RL sigmaY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
37 _RL sigmaR(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
38 INTEGER bi, bj
39 _RL myTime
40 INTEGER myThid
41
42 #ifdef GM_K3D
43
44 C === Functions ====
45 LOGICAL DIFFERENT_MULTIPLE
46 EXTERNAL DIFFERENT_MULTIPLE
47
48 C !LOCAL VARIABLES:
49 C == Local variables ==
50 INTEGER i,j,k,kk,m
51
52 C update_modes :: Whether to update the eigenmodes
53 LOGICAL update_modes
54
55 C surfk :: index of the depth of the surface layer
56 C kLow_C :: Local version of the index of deepest wet grid cell on tracer grid
57 C kLow_U :: Local version of the index of deepest wet grid cell on U grid
58 C kLow_V :: Local version of the index of deepest wet grid cell on V grid
59 INTEGER surfk(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
60 INTEGER kLow_C(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
61 INTEGER kLow_U(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
62 INTEGER kLow_V(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
63
64 C N2loc :: local N**2
65 C slope :: local slope
66 C Req :: local equatorial deformation radius (m)
67 C Rurms :: a local mixing length used in calculation of urms (m)
68 C deltaH :: local thickness of Eady integration (m)
69 C g_reciprho_sq :: (gravity*recip_rhoConst)**2
70 C M4loc :: local M**4
71 C maxDRhoDz :: maximum value of d(rho)/dz (derived from GM_K3D_minN2)
72 C sigx :: local d(rho)/dx
73 C sigy :: local d(rho)/dy
74 C sigz :: local d(rho)/dz
75 C hsurf :: local surface layer depth
76 C small :: a small number (to avoid floating point exceptions)
77 _RL N2loc
78 _RL slope
79 _RL Req
80 _RL Rurms
81 _RL deltaH
82 _RL g_reciprho_sq
83 _RL M4loc
84 _RL maxDRhoDz
85 _RL sigx, sigy, sigz
86 _RL hsurf
87 _RL small
88
89 C dfdy :: gradient of the Coriolis paramter, df/dy, 1/(m*s)
90 C dfdx :: gradient of the Coriolis paramter, df/dx, 1/(m*s)
91 C gradf :: total gradient of the Coriolis paramter, SQRT(df/dx**2+df/dy**2), 1/(m*s)
92 C RRhines :: The Rhines scale (m)
93 C Rmix :: Mixing length
94 C N2 :: Square of the buoyancy frequency (1/s**2)
95 C N2W :: Square of the buoyancy frequency (1/s**2) averaged to west of grid cell
96 C N2S :: Square of the buoyancy frequency (1/s**2) averaged to south of grid cell
97 C N :: Buoyancy frequency, SQRT(N2)
98 C BVint :: The vertical integral of N (m/s)
99 C ubar :: Zonal velocity on a tracer point (m/s)
100 C vbar :: Meridional velocity on a tracer point (m/s)
101 C Ubaro :: Barotropic velocity (m/s)
102 _RL dfdy( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
103 _RL dfdx( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
104 _RL gradf( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
105 _RL dummy( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
106 _RL RRhines(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
107 _RL Rmix( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
108 _RL N2( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
109 _RL N2W( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
110 _RL N2S( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
111 _RL N( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
112 _RL BVint( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
113 _RL Ubaro( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
114 _RL ubar( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
115 _RL vbar( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
116
117 C Rmid :: Rossby radius (m)
118 C KPV :: Diffusivity (m**2/s)
119 C SlopeX :: isopycnal slope in x direction
120 C SlopeY :: isopycnal slope in y direction
121 C dSigmaDx :: sigmaX averaged onto tracer grid
122 C dSigmaDy :: sigmaY averaged onto tracer grid
123 C tfluxX :: thickness flux in x direction
124 C tfluxY :: thickness flux in y direction
125 C fCoriU :: Coriolis parameter averaged to U points
126 C fCoriV :: Coriolis parameter averaged to V points
127 C cori :: Coriolis parameter forced to be finite near the equator
128 C coriU :: As for cori, but, at U point
129 C coriV :: As for cori, but, at V point
130 C surfkz :: Depth of surface layer (in r units)
131 _RL Rmid(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
132 _RL KPV(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
133 _RL SlopeX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
134 _RL SlopeY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
135 _RL dSigmaDx(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
136 _RL dSigmaDy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
137 _RL tfluxX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
138 _RL tfluxY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
139 _RL cori(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
140 _RL coriU(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
141 _RL coriV(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
142 _RL fCoriU(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
143 _RL fCoriV(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
144 _RL surfkz(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
145
146 C gradqx :: Potential vorticity gradient in x direction
147 C gradqy :: Potential vorticity gradient in y direction
148 C XimX :: Vertical integral of phi_m*K*gradqx
149 C XimY :: Vertical integral of phi_m*K*gradqy
150 C cDopp :: Quasi-Doppler shifted long Rossby wave speed (m/s)
151 C umc :: ubar-c (m/s)
152 C eady :: Eady growth rate (1/s)
153 C urms :: the rms eddy velocity (m/s)
154 C supp :: The suppression factor
155 C ustar :: The eddy induced velocity in the x direction
156 C vstar :: The eddy induced velocity in the y direction
157 C Xix :: Xi in the x direction (m/s**2)
158 C Xiy :: Xi in the y direction (m/s**2)
159 _RL gradqx(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
160 _RL gradqy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
161 _RL XimX(GM_K3D_NModes,1-Olx:sNx+Olx,1-Oly:sNy+Oly)
162 _RL XimY(GM_K3D_NModes,1-Olx:sNx+Olx,1-Oly:sNy+Oly)
163 _RL cDopp(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
164 _RL umc( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
165 _RL eady( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
166 _RL urms( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
167 _RL supp( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
168 _RL ustar(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
169 _RL vstar(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
170 _RL Xix( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
171 _RL Xiy( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
172 #ifdef GM_K3D_PASSIVE
173 C psistar :: eddy induced streamfunction in the y direction
174 _RL psistar(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
175 #endif
176
177
178 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
179
180 C ======================================
181 C Initialise some variables
182 C ======================================
183 small = TINY(zeroRL)
184 update_modes=.FALSE.
185 IF ( DIFFERENT_MULTIPLE(GM_K3D_vecFreq,myTime,deltaTClock) )
186 & update_modes=.TRUE.
187
188 DO j=1-Oly,sNy+Oly
189 DO i=1-Olx,sNx+Olx
190 kLow_C(i,j) = kLowC(i,j,bi,bj)
191 ENDDO
192 ENDDO
193 DO j=1-Oly,sNy+Oly
194 DO i=1-Olx+1,sNx+Olx
195 kLow_U(i,j) = MIN( kLow_C(i,j), kLow_C(i-1,j) )
196 ENDDO
197 ENDDO
198 DO j=1-Oly+1,sNy+Oly
199 DO i=1-Olx,sNx+Olx
200 kLow_V(i,j) = MIN( kLow_C(i,j), kLow_C(i,j-1) )
201 ENDDO
202 ENDDO
203
204 C Dummy values for the edges
205 C This avoids weirdness in gmredi_calc_eigs
206 i=1-Olx
207 DO j=1-Oly,sNy+Oly
208 kLow_U(i,j) = kLow_C(i,j)
209 ENDDO
210 j=1-Oly
211 DO i=1-Olx,sNx+Olx
212 kLow_V(i,j) = kLow_C(i,j)
213 ENDDO
214
215 g_reciprho_sq = (gravity*recip_rhoConst)**2
216 C Gradient of Coriolis
217 DO j=1-Oly+1,sNy+Oly
218 DO i=1-Olx+1,sNx+Olx
219 dfdx(i,j) = ( fCori(i,j,bi,bj)-fCori(i-1,j,bi,bj) )
220 & *recip_dxC(i,j,bi,bj)
221 dfdy(i,j) = ( fCori(i,j,bi,bj)-fCori(i,j-1,bi,bj) )
222 & *recip_dyC(i,j,bi,bj)
223 ENDDO
224 ENDDO
225
226 C Coriolis at C points enforcing a minimum value so
227 C that it is defined at the equator
228 DO j=1-Oly,sNy+Oly
229 DO i=1-Olx,sNx+Olx
230 cori(i,j) = SIGN( MAX( ABS(fCori(i,j,bi,bj)),GM_K3D_minCori ),
231 & fCori(i,j,bi,bj) )
232 ENDDO
233 ENDDO
234 C Coriolis at U and V points
235 DO j=1-Oly+1,sNy+Oly
236 DO i=1-Olx+1,sNx+Olx
237 C Limited so that the inverse is defined at the equator
238 coriU(i,j) = op5*( cori(i,j)+cori(i-1,j) )
239 coriU(i,j) = SIGN( MAX( ABS(coriU(i,j)),GM_K3D_minCori ),
240 & coriU(i,j) )
241
242 coriV(i,j) = op5*( cori(i,j)+cori(i,j-1) )
243 coriV(i,j) = SIGN( MAX( ABS(coriV(i,j)),GM_K3D_minCori ),
244 & coriV(i,j) )
245
246 C Not limited
247 fCoriU(i,j) = op5*( fCori(i,j,bi,bj)+fCori(i-1,j,bi,bj) )
248 fCoriV(i,j) = op5*( fCori(i,j,bi,bj)+fCori(i,j-1,bi,bj) )
249 ENDDO
250 ENDDO
251 DO j=1-Oly,sNy+Oly
252 DO i=1-Olx,sNx+Olx
253 gradf(i,j) = SQRT( dfdx(i,j)**2 + dfdy(i,j)**2 )
254 ENDDO
255 ENDDO
256
257 C Zeroing some cumulative fields
258 DO j=1-Oly,sNy+Oly
259 DO i=1-Olx,sNx+Olx
260 eady(i,j) = zeroRL
261 BVint(i,j) = zeroRL
262 Ubaro(i,j) = zeroRL
263 ENDDO
264 ENDDO
265
266 C Find the zonal velocity at the cell centre
267 C The logicals here are, in order: 1/ go from grid to north/east directions
268 C 2/ go from C to A grid and 3/ apply the mask
269 CALL rotate_uv2en_rl(uVel, vVel, ubar, vbar, .TRUE., .TRUE.,
270 & .TRUE.,Nr,mythid)
271
272 C Square of the buoyancy frequency at the top of a grid cell
273 DO k=2,Nr
274 DO j=1-Oly,sNy+Oly
275 DO i=1-Olx,sNx+Olx
276 N2(i,j,k) = -gravity*recip_rhoConst*sigmaR(i,j,k)
277 ENDDO
278 ENDDO
279 ENDDO
280 C N2(k=1) is always zero
281 k=1
282 DO j=1-Oly,sNy+Oly
283 DO i=1-Olx,sNx+Olx
284 N2(i,j,k) = 0.0
285 N(i,j,k) = 0.0
286 ENDDO
287 ENDDO
288 C Enforce a minimum N2
289 DO k=2,Nr
290 DO j=1-Oly,sNy+Oly
291 DO i=1-Olx,sNx+Olx
292 IF (N2(i,j,k).LT.GM_K3D_minN2) N2(i,j,k)=GM_K3D_minN2
293 N(i,j,k) = SQRT(N2(i,j,k))
294 ENDDO
295 ENDDO
296 ENDDO
297 C Calculate the minimum drho/dz
298 maxDRhoDz = -rhoConst*GM_K3D_minN2/gravity
299
300 C Calculate the barotropic velocity by vertically integrating
301 C and the dividing by the depth of the water column
302 C Note that Ubaro is on the U grid.
303 DO k=1,Nr
304 DO j=1-Oly,sNy+Oly
305 DO i=1-Olx,sNx+Olx
306 Ubaro(i,j) = Ubaro(i,j) +
307 & maskW(i,j,k,bi,bj)*drF(k)*hfacC(i,j,k,bi,bj)
308 & *ubar(i,j,k,bi,bj)
309 ENDDO
310 ENDDO
311 ENDDO
312 DO j=1-Oly,sNy+Oly
313 DO i=1-Olx,sNx+Olx
314 IF (kLow_C(i,j).GT.0) THEN
315 C The minus sign is because r_Low<0
316 Ubaro(i,j) = -Ubaro(i,j)/r_Low(i,j,bi,bj)
317 ENDIF
318 ENDDO
319 ENDDO
320
321 C Integrate the buoyancy frequency vertically using the trapezoidal method.
322 DO k=1,Nr
323 DO j=1-Oly,sNy+Oly
324 DO i=1-Olx,sNx+Olx
325 IF (k.LT.kLow_C(i,j)) THEN
326 BVint(i,j) = BVint(i,j) + hFacC(i,j,k,bi,bj)*drF(k)
327 & *(N(i,j,k)+N(i,j,k+1))
328 ELSEIF (k.EQ.kLow_C(i,j)) THEN
329 C Assume that N(z=-H)=0
330 BVint(i,j) = BVint(i,j) + hFacC(i,j,k,bi,bj)*drF(k)*N(i,j,k)
331 ENDIF
332 ENDDO
333 ENDDO
334 ENDDO
335 DO j=1-Oly,sNy+Oly
336 DO i=1-Olx,sNx+Olx
337 BVint(i,j) = op5*BVint(i,j)
338 ENDDO
339 ENDDO
340
341 C Calculate the eigenvalues and eigenvectors
342 IF (update_modes) THEN
343 CALL GMREDI_CALC_EIGS(
344 I iMin,iMax,jMin,jMax,bi,bj,N2,myThid,
345 I kLow_C, maskC(:,:,:,bi,bj),
346 I hfacC(:,:,:,bi,bj), recip_hfacC(:,:,:,bi,bj),
347 I R_Low(:,:,bi,bj), 1, .TRUE.,
348 O Rmid, modesC(:,:,:,:,bi,bj))
349
350 C Calculate the Rossby Radius
351 DO j=1-Oly+1,sNy+Oly
352 DO i=1-Olx+1,sNx+Olx
353 Req = SQRT(BVint(i,j)/(2*pi*gradf(i,j)))
354 Rdef(i,j,bi,bj) = MIN(Rmid(i,j),Req)
355 ENDDO
356 ENDDO
357 ENDIF
358
359 C Average dsigma/dx and dsigma/dy onto the centre points
360
361 DO k=1,Nr
362 DO j=1-Oly,sNy+Oly-1
363 DO i=1-Olx,sNx+Olx-1
364 dSigmaDx(i,j,k) = op5*(sigmaX(i,j,k)+sigmaX(i+1,j,k))
365 dSigmaDy(i,j,k) = op5*(sigmaY(i,j,k)+sigmaY(i,j+1,k))
366 ENDDO
367 ENDDO
368 ENDDO
369
370 C ===============================
371 C Calculate the Eady growth rate
372 C ===============================
373 DO k=1,Nr
374
375 C The bottom of the grid cell is shallower than the top
376 C integration level, so, advance the depth.
377 IF (-rF(k+1).LE. GM_K3D_EadyMinDepth) CYCLE
378
379 C Don't bother going any deeper since the top of the
380 C cell is deeper than the bottom integration level
381 IF (-rF(k).GE.GM_K3D_EadyMaxDepth) EXIT
382
383 C We are in the integration depth range
384 DO j=1-Oly,sNy+Oly-1
385 DO i=1-Olx,sNx+Olx-1
386 IF (kLow_C(i,j).GE.k) THEN
387 IF (k.NE.kLow_C(i,j)) THEN
388 N2loc = op5*(N2(i,j,k)+N2(i,j,k+1))
389 ELSE
390 N2loc = op5*N2(i,j,k)
391 ENDIF
392 M4loc = g_reciprho_sq*( dSigmaDx(i,j,k)**2
393 & +dSigmaDy(i,j,k)**2 )
394 slope = SQRT(SQRT(M4loc)/N2loc)
395
396 C Limit the slope. Note, this is not all the Eady calculations.
397 IF (slope.LE.GM_K3D_maxSlope) THEN
398 eady(i,j) = eady(i,j)
399 & + hfacC(i,j,k,bi,bj)*drF(k)*M4loc/(N2loc)
400 ELSE
401 eady(i,j) = eady(i,j)
402 & + hfacC(i,j,k,bi,bj)*drF(k)*SQRT(M4loc)
403 & *GM_K3D_maxSlope*GM_K3D_maxSlope
404 ENDIF
405 ENDIF
406 ENDDO
407 ENDDO
408 ENDDO
409
410 DO j=1-Oly,sNy+Oly
411 DO i=1-Olx,sNx+Olx
412 C If the minimum depth for the integration is deeper than ocean
413 C bottom then give the eady growth rate a dummy, non-zero value
414 C to avoid floating point exceptions. These points are taken care
415 C of by setting K3D=GM_K3D_smallK later.
416 IF (kLow_C(i,j).NE.0
417 & .AND. -r_Low(i,j,bi,bj).LT.GM_K3D_EadyMinDepth) THEN
418 eady(i,j) = small
419
420 C Otherwise, multiply eady by the various constants to get the
421 C growth rate.
422 ELSE
423 deltaH = MIN(-r_low(i,j,bi,bj),GM_K3D_EadyMaxDepth)
424 deltaH = deltaH - GM_K3D_EadyMinDepth
425 eady(i,j) = SQRT(eady(i,j)/deltaH)
426
427 ENDIF
428
429 ENDDO
430 ENDDO
431
432 C ======================================
433 C Calculate the diffusivity
434 C ======================================
435 DO j=1-Oly+1,sNy+Oly
436 DO i=1-Olx+1,sNx+Olx-1
437 C Calculate the Visbeck velocity
438 Rurms = MIN(Rdef(i,j,bi,bj),GM_K3D_maxLurms)
439 urms(i,j,1) = GM_K3D_Lambda*eady(i,j)*Rurms
440 C Set the bottom urms to zero
441 k=kLow_C(i,j)
442 IF (k.GT.0) urms(i,j,k) = 0.0
443
444 C Calculate the Rhines scale
445 RRhines(i,j) = SQRT(urms(i,j,1)/gradf(i,j))
446
447 C Calculate the estimated length scale
448 Rmix(i,j) = MIN(Rdef(i,j,bi,bj), RRhines(i,j))
449
450 C Calculate the Doppler shifted long Rossby wave speed
451 C Ubaro is on the U grid so we must average onto the M grid.
452 cDopp(i,j) = op5*( Ubaro(i,j)+Ubaro(i+1,j) )
453 & - gradf(i,j)*Rdef(i,j,bi,bj)*Rdef(i,j,bi,bj)
454 C Limit the wave speed to the namelist variable GM_K3D_maxC
455 IF (ABS(cDopp(i,j)).GT.GM_K3D_maxC) THEN
456 cDopp(i,j) = MAX(GM_K3D_maxC, cDopp(i,j))
457 ENDIF
458
459 ENDDO
460 ENDDO
461
462 C Project the surface urms to the subsurface using the first baroclinic mode
463 CALL GMREDI_CALC_URMS(
464 I iMin,iMax,jMin,jMax,bi,bj,N2,myThid,
465 U urms)
466
467 C Calculate the diffusivity (on the mass grid)
468 DO k=1,Nr
469 DO j=1-Oly,sNy+Oly
470 DO i=1-Olx,sNx+Olx
471 IF (k.LE.kLow_C(i,j)) THEN
472 IF (-r_Low(i,j,bi,bj).LT.GM_K3D_EadyMinDepth) THEN
473 K3D(i,j,k,bi,bj) = GM_K3D_smallK
474 ELSE
475 IF (urms(i,j,k).EQ.0.0) THEN
476 K3D(i,j,k,bi,bj) = GM_K3D_smallK
477 ELSE
478 umc(i,j,k) = ubar(i,j,k,bi,bj) - cDopp(i,j)
479 supp(i,j,k) = 1/( 1 + 4*umc(i,j,k)**2/urms(i,j,1)**2 )
480 K3D(i,j,k,bi,bj) = GM_K3D_gamma*urms(i,j,k)
481 & *Rmix(i,j)*supp(i,j,k)
482 ENDIF
483
484 C Enforce lower and upper bounds on the diffusivity
485 IF (K3D(i,j,k,bi,bj).LT.GM_K3D_smallK)
486 & K3D(i,j,k,bi,bj) = GM_K3D_smallK
487 IF (K3D(i,j,k,bi,bj).GT.GM_maxK3D)
488 & K3D(i,j,k,bi,bj) = GM_maxK3D
489 ENDIF
490 ENDIF
491 ENDDO
492 ENDDO
493 ENDDO
494
495 C ======================================
496 C Find the PV gradient
497 C ======================================
498 C Calculate the surface layer thickness.
499 C Use hMixLayer (calculated in model/src/calc_oce_mxlayer)
500 C for the mixed layer depth.
501
502 C Enforce a minimum surface layer depth
503 DO j=1-Oly,sNy+Oly
504 DO i=1-Olx,sNx+Olx
505 surfkz(i,j) = MIN(-GM_K3D_surfMinDepth,-hMixLayer(i,j,bi,bj))
506 surfkz(i,j) = MAX(surfkz(i,j),R_low(i,j,bi,bj))
507 IF(maskC(i,j,1,bi,bj).EQ.0.0) surfkz(i,j)=0.0
508 surfk(i,j) = 0
509 ENDDO
510 ENDDO
511 DO k=1,Nr
512 DO j=1-Oly,sNy+Oly
513 DO i=1-Olx,sNx+Olx
514 IF (rF(k).GT.surfkz(i,j) .AND. surfkz(i,j).GE.rF(k+1))
515 & surfk(i,j) = k
516 ENDDO
517 ENDDO
518 ENDDO
519
520 C Calculate the isopycnal slope
521 DO j=1-Oly,sNy+Oly-1
522 DO i=1-Olx,sNx+Olx-1
523 SlopeX(i,j,1) = zeroRL
524 SlopeY(i,j,1) = zeroRL
525 ENDDO
526 ENDDO
527 DO k=2,Nr
528 DO j=1-Oly+1,sNy+Oly
529 DO i=1-Olx+1,sNx+Olx
530 IF(surfk(i,j).GE.kLowC(i,j,bi,bj)) THEN
531 C If the surface layer is thinner than the water column
532 C the set the slope to zero to avoid problems.
533 SlopeX(i,j,k) = zeroRL
534 SlopeY(i,j,k) = zeroRL
535
536 ELSE
537 C Calculate the zonal slope at the western cell face (U grid)
538 sigz = MIN( op5*(sigmaR(i,j,k)+sigmaR(i-1,j,k)), maxDRhoDz )
539 sigx = op5*( sigmaX(i,j,k)+sigmaX(i,j,k-1) )
540 slope = sigx/sigz
541 C IF(ABS(slope).GT.GM_K3D_maxSlope)
542 C & slope = SIGN(GM_K3D_maxSlope,slope)
543 IF(ABS(slope).GT.GM_maxSlope)
544 & slope = SIGN(GM_maxSlope,slope)
545 SlopeX(i,j,k)=-maskW(i,j,k-1,bi,bj)*maskW(i,j,k,bi,bj)*slope
546
547 C Calculate the meridional slope at the southern cell face (V grid)
548 sigz = MIN( op5*(sigmaR(i,j,k)+sigmaR(i,j-1,k)), maxDRhoDz )
549 sigy = op5*( sigmaY(i,j,k) + sigmaY(i,j,k-1) )
550 slope = sigy/sigz
551 C IF(ABS(slope).GT.GM_K3D_maxSlope)
552 C & slope = SIGN(GM_K3D_maxSlope,slope)
553 IF(ABS(slope).GT.GM_maxSlope)
554 & slope = SIGN(GM_maxSlope,slope)
555 SlopeY(i,j,k)=-maskS(i,j,k-1,bi,bj)*maskS(i,j,k,bi,bj)*slope
556 ENDIF
557 ENDDO
558 ENDDO
559 ENDDO
560
561 C Calculate the thickness flux
562 C Enforce a zero slope bottom boundary condition for the bottom most cells (k=Nr)
563 k=Nr
564 DO j=1-Oly,sNy+Oly
565 DO i=1-Olx,sNx+Olx
566 C Zonal thickness flux at the western cell face
567 tfluxX(i,j,k) = -fCoriU(i,j)*SlopeX(i,j,k)
568 & *recip_drF(k)*recip_hFacW(i,j,k,bi,bj)
569 C Meridional thickness flux at the southern cell face
570 tfluxY(i,j,k) = -fCoriV(i,j)*SlopeY(i,j,k)
571 & *recip_drF(k)*recip_hFacS(i,j,k,bi,bj)
572 ENDDO
573 ENDDO
574
575 C Calculate the thickness flux for other cells (k<Nr)
576 C SlopeX and SlopeY are zero in dry cells, so, the bottom boundary
577 C condition that the slope is zero is taken care of.
578 C We still need to give special treatment for the surface layer however.
579 DO k=Nr-1,1,-1
580 DO j=1-Oly,sNy+Oly-1
581 DO i=1-Olx,sNx+Olx-1
582 IF(k.LE.surfk(i,j) .AND. .NOT. GM_K3D_likeGM) THEN
583 C We are in the surface layer, so set the thickness flux
584 C based on the average slope over the surface layer
585 C If we are on the edge of a "cliff" the surface layer at the
586 C centre of the grid point could be deeper than the U or V point.
587 C So, we ensure that we always take a sensible slope.
588 IF(kLow_U(i,j).LT.surfk(i,j)) THEN
589 kk=kLow_U(i,j)
590 hsurf = -rLowW(i,j,bi,bj)
591 ELSE
592 kk=surfk(i,j)
593 hsurf = -surfkz(i,j)
594 ENDIF
595 IF(kk.GT.0) THEN
596 IF(kk.EQ.Nr) THEN
597 tfluxX(i,j,k) = -fCoriU(i,j)*maskW(i,j,k,bi,bj)
598 & *SlopeX(i,j,kk)/hsurf
599 ELSE
600 tfluxX(i,j,k) = -fCoriU(i,j)*maskW(i,j,k,bi,bj)
601 & *( SlopeX(i,j,kk)-SlopeX(i,j,kk+1) )/hsurf
602 ENDIF
603 ELSE
604 tfluxX(i,j,k) = zeroRL
605 ENDIF
606
607 IF(kLow_V(i,j).LT.surfk(i,j)) THEN
608 kk=kLow_V(i,j)
609 hsurf = -rLowS(i,j,bi,bj)
610 ELSE
611 kk=surfk(i,j)
612 hsurf = -surfkz(i,j)
613 ENDIF
614 IF(kk.GT.0) THEN
615 IF(kk.EQ.Nr) THEN
616 tfluxY(i,j,k) = -fCoriV(i,j)*maskS(i,j,k,bi,bj)
617 & *SlopeY(i,j,kk)/hsurf
618 ELSE
619 tfluxY(i,j,k) = -fCoriV(i,j)*maskS(i,j,k,bi,bj)
620 & *( SlopeY(i,j,kk)-SlopeY(i,j,kk+1) )/hsurf
621 ENDIF
622 ELSE
623 tfluxY(i,j,k) = zeroRL
624 ENDIF
625
626 ELSE
627 C We are not in the surface layer, so calculate the thickness
628 C flux based on the local isopyncal slope
629
630 C Zonal thickness flux at the western cell face
631 tfluxX(i,j,k)=-fCoriU(i,j)*(SlopeX(i,j,k)-SlopeX(i,j,k+1))
632 & *recip_drF(k)*recip_hFacW(i,j,k,bi,bj)
633 & *maskW(i,j,k,bi,bj)
634
635 C Meridional thickness flux at the southern cell face
636 tfluxY(i,j,k)=-fCoriV(i,j)*(SlopeY(i,j,k)-SlopeY(i,j,k+1))
637 & *recip_drF(k)*recip_hFacS(i,j,k,bi,bj)
638 & *maskS(i,j,k,bi,bj)
639 ENDIF
640 ENDDO
641 ENDDO
642 ENDDO
643
644 C Calculate gradq
645 IF (GM_K3D_likeGM .OR. GM_K3D_beta_eq_0) THEN
646 C Ignore beta in the calculation of grad(q)
647 DO k=1,Nr
648 DO j=1-Oly+1,sNy+Oly
649 DO i=1-Olx+1,sNx+Olx
650 gradqx(i,j,k) = maskW(i,j,k,bi,bj)*tfluxX(i,j,k)
651 gradqy(i,j,k) = maskS(i,j,k,bi,bj)*tfluxY(i,j,k)
652 ENDDO
653 ENDDO
654 ENDDO
655
656 ELSE
657 C Do not ignore beta
658 DO k=1,Nr
659 DO j=1-Oly+1,sNy+Oly
660 DO i=1-Olx+1,sNx+Olx
661 gradqx(i,j,k) = maskW(i,j,k,bi,bj)*(dfdx(i,j)+tfluxX(i,j,k))
662 gradqy(i,j,k) = maskS(i,j,k,bi,bj)*(dfdy(i,j)+tfluxY(i,j,k))
663 ENDDO
664 ENDDO
665 ENDDO
666 ENDIF
667
668 C ======================================
669 C Find Xi and the eddy induced velocities
670 C ======================================
671 C Find the buoyancy frequency at the west and south faces of a cell
672 C This is necessary to find the eigenvectors at those points
673 DO k=1,Nr
674 DO j=1-Oly+1,sNy+Oly
675 DO i=1-Olx+1,sNx+Olx
676 N2W(i,j,k) = maskW(i,j,k,bi,bj)
677 & *( N2(i,j,k)+N2(i-1,j,k) )
678 N2S(i,j,k) = maskS(i,j,k,bi,bj)
679 & *( N2(i,j,k)+N2(i,j-1,k) )
680 ENDDO
681 ENDDO
682 C This fudge is necessary to avoid division by zero in gmredi_calc_eigs.
683 C It does not affect the end result since it's in the overlap region.
684 j=1-Oly
685 DO i=1-Olx,sNx+Olx
686 N2W(i,j,k) = GM_K3D_minN2
687 N2S(i,j,k) = GM_K3D_minN2
688 ENDDO
689 i=1-Olx
690 DO j=1-Oly,sNy+Oly
691 N2W(i,j,k) = GM_K3D_minN2
692 N2S(i,j,k) = GM_K3D_minN2
693 ENDDO
694 ENDDO
695
696 IF(GM_K3D_likeGM) THEN
697 DO k=1,Nr
698 DO j=1-Oly,sNy+Oly
699 DO i=1-Olx,sNx+Olx
700 KPV(i,j,k) = GM_K3D_constK
701 ENDDO
702 ENDDO
703 ENDDO
704 ELSE
705 DO k=1,Nr
706 DO j=1-Oly,sNy+Oly
707 DO i=1-Olx,sNx+Olx
708 KPV(i,j,k) = K3D(i,j,k,bi,bj)
709 ENDDO
710 ENDDO
711 ENDDO
712 ENDIF
713
714 IF (.NOT. GM_K3D_smooth) THEN
715 C Do not expand K grad(q) => no smoothing
716 C May only be done with a constant K, otherwise the
717 C integral constraint is violated.
718 DO k=1,Nr
719 DO j=1-Oly,sNy+Oly
720 DO i=1-Olx,sNx+Olx
721 Xix(i,j,k) = -maskW(i,j,k,bi,bj)*KPV(i,j,k)*gradqx(i,j,k)
722 Xiy(i,j,k) = -maskS(i,j,k,bi,bj)*KPV(i,j,k)*gradqy(i,j,k)
723 ENDDO
724 ENDDO
725 ENDDO
726
727 ELSE
728 C Expand K grad(q) in terms of baroclinic modes to smooth
729 C and satisfy the integral constraint
730
731 C Start with the X direction
732 C ------------------------------
733 C Calculate the eigenvectors at the West face of a cell
734 IF (update_modes) THEN
735 CALL GMREDI_CALC_EIGS(
736 I iMin,iMax,jMin,jMax,bi,bj,N2W,myThid,
737 I kLow_U,maskW(:,:,:,bi,bj),
738 I hfacW(:,:,:,bi,bj),recip_hfacW(:,:,:,bi,bj),
739 I rLowW(:,:,bi,bj),GM_K3D_NModes,.FALSE.,
740 O dummy,modesW(:,:,:,:,bi,bj))
741 ENDIF
742
743 C Calculate Xi_m at the west face of a cell
744 DO j=1-Oly,sNy+Oly
745 DO i=1-Olx,sNx+Olx
746 DO m=1,GM_K3D_NModes
747 XimX(m,i,j) = zeroRL
748 ENDDO
749 ENDDO
750 ENDDO
751 DO k=1,Nr
752 DO j=1-Oly,sNy+Oly
753 DO i=1-Olx,sNx+Olx
754 DO m=1,GM_K3D_NModes
755 XimX(m,i,j) = XimX(m,i,j)
756 & - maskW(i,j,k,bi,bj)*drF(k)*hfacW(i,j,k,bi,bj)
757 & *KPV(i,j,k)*gradqx(i,j,k)*modesW(m,i,j,k,bi,bj)
758 ENDDO
759 ENDDO
760 ENDDO
761 ENDDO
762
763 C Calculate Xi in the X direction at the west face
764 DO k=1,Nr
765 DO j=1-Oly,sNy+Oly
766 DO i=1-Olx,sNx+Olx
767 Xix(i,j,k) = zeroRL
768 ENDDO
769 ENDDO
770 ENDDO
771 DO k=1,Nr
772 DO j=1-Oly,sNy+Oly
773 DO i=1-Olx,sNx+Olx
774 DO m=1,GM_K3D_NModes
775 Xix(i,j,k) = Xix(i,j,k)
776 & + maskW(i,j,k,bi,bj)*XimX(m,i,j)*modesW(m,i,j,k,bi,bj)
777 ENDDO
778 ENDDO
779 ENDDO
780 ENDDO
781
782 C Now the Y direction
783 C ------------------------------
784 C Calculate the eigenvectors at the West face of a cell
785 IF (update_modes) THEN
786 CALL GMREDI_CALC_EIGS(
787 I iMin,iMax,jMin,jMax,bi,bj,N2S,myThid,
788 I kLow_V,maskS(:,:,:,bi,bj),
789 I hfacS(:,:,:,bi,bj),recip_hfacS(:,:,:,bi,bj),
790 I rLowS(:,:,bi,bj), GM_K3D_NModes, .FALSE.,
791 O dummy,modesS(:,:,:,:,bi,bj))
792 ENDIF
793
794 DO j=1-Oly,sNy+Oly
795 DO i=1-Olx,sNx+Olx
796 DO m=1,GM_K3D_NModes
797 XimY(m,i,j) = zeroRL
798 ENDDO
799 ENDDO
800 ENDDO
801 DO k=1,Nr
802 DO j=1-Oly,sNy+Oly
803 DO i=1-Olx,sNx+Olx
804 DO m=1,GM_K3D_NModes
805 XimY(m,i,j) = XimY(m,i,j)
806 & - drF(k)*hfacS(i,j,k,bi,bj)
807 & *KPV(i,j,k)*gradqy(i,j,k)*modesS(m,i,j,k,bi,bj)
808 ENDDO
809 ENDDO
810 ENDDO
811 ENDDO
812
813 C Calculate Xi for Y direction at the south face
814 DO k=1,Nr
815 DO j=1-Oly,sNy+Oly
816 DO i=1-Olx,sNx+Olx
817 Xiy(i,j,k) = zeroRL
818 ENDDO
819 ENDDO
820 ENDDO
821 DO k=1,Nr
822 DO j=1-Oly,sNy+Oly
823 DO i=1-Olx,sNx+Olx
824 DO m=1,GM_K3D_NModes
825 Xiy(i,j,k) = Xiy(i,j,k)
826 & + maskS(i,j,k,bi,bj)*XimY(m,i,j)*modesS(m,i,j,k,bi,bj)
827 ENDDO
828 ENDDO
829 ENDDO
830 ENDDO
831
832 C ENDIF GM_K3D_likeGM
833 ENDIF
834
835
836 C Calculate the eddy induced velocity in the X direction at the west face
837 DO k=1,Nr
838 DO j=1-Oly+1,sNy+Oly
839 DO i=1-Olx+1,sNx+Olx
840 ustar(i,j,k) = -Xix(i,j,k)/coriU(i,j)
841 ENDDO
842 ENDDO
843 ENDDO
844
845 C Calculate the eddy induced velocity in the Y direction at the south face
846 DO k=1,Nr
847 DO j=1-Oly+1,sNy+Oly
848 DO i=1-Olx+1,sNx+Olx
849 vstar(i,j,k) = -Xiy(i,j,k)/coriV(i,j)
850 ENDDO
851 ENDDO
852 ENDDO
853
854 C ======================================
855 C Calculate the eddy induced overturning streamfunction
856 C ======================================
857 #ifdef GM_K3D_PASSIVE
858 k=Nr
859 DO j=1-Oly,sNy+Oly
860 DO i=1-Olx,sNx+Olx
861 psistar(i,j,Nr) = -hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
862 ENDDO
863 ENDDO
864 DO k=Nr-1,1,-1
865 DO j=1-Oly,sNy+Oly
866 DO i=1-Olx,sNx+Olx
867 psistar(i,j,k) = psistar(i,j,k+1)
868 & - hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
869 ENDDO
870 ENDDO
871 ENDDO
872
873 #else
874
875 IF (GM_AdvForm) THEN
876 k=Nr
877 DO j=1-Oly+1,sNy+1
878 DO i=1-Olx+1,sNx+1
879 GM_PsiX(i,j,k,bi,bj) = -hfacW(i,j,k,bi,bj)*drF(k)*ustar(i,j,k)
880 GM_PsiY(i,j,k,bi,bj) = -hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
881 ENDDO
882 ENDDO
883 DO k=Nr-1,1,-1
884 DO j=1-Oly+1,sNy+1
885 DO i=1-Olx+1,sNx+1
886 GM_PsiX(i,j,k,bi,bj) = GM_PsiX(i,j,k+1,bi,bj)
887 & - hfacW(i,j,k,bi,bj)*drF(k)*ustar(i,j,k)
888 GM_PsiY(i,j,k,bi,bj) = GM_PsiY(i,j,k+1,bi,bj)
889 & - hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
890 ENDDO
891 ENDDO
892 ENDDO
893
894 ENDIF
895 #endif
896
897 #ifdef ALLOW_DIAGNOSTICS
898 C Diagnostics
899 IF ( useDiagnostics ) THEN
900 CALL DIAGNOSTICS_FILL(K3D, 'GM_K3D ',0,Nr,0,1,1,myThid)
901 CALL DIAGNOSTICS_FILL(urms, 'GM_URMS ',0,Nr,0,1,1,myThid)
902 CALL DIAGNOSTICS_FILL(Rdef, 'GM_RDEF ',0, 1,0,1,1,myThid)
903 CALL DIAGNOSTICS_FILL(RRhines,'GM_LRHNS',0, 1,0,1,1,myThid)
904 CALL DIAGNOSTICS_FILL(Rmix, 'GM_RMIX ',0, 1,0,1,1,myThid)
905 CALL DIAGNOSTICS_FILL(supp, 'GM_SUPP ',0,Nr,0,1,1,myThid)
906 CALL DIAGNOSTICS_FILL(Xix, 'GM_Xix ',0,Nr,0,1,1,myThid)
907 CALL DIAGNOSTICS_FILL(Xiy, 'GM_Xiy ',0,Nr,0,1,1,myThid)
908 CALL DIAGNOSTICS_FILL(cDopp, 'GM_C ',0, 1,0,1,1,myThid)
909 CALL DIAGNOSTICS_FILL(Ubaro, 'GM_UBARO',0, 1,0,1,1,myThid)
910 CALL DIAGNOSTICS_FILL(eady, 'GM_EADY ',0, 1,0,1,1,myThid)
911 CALL DIAGNOSTICS_FILL(SlopeX, 'GM_Sx ',0,Nr,0,1,1,myThid)
912 CALL DIAGNOSTICS_FILL(SlopeY, 'GM_Sy ',0,Nr,0,1,1,myThid)
913 CALL DIAGNOSTICS_FILL(tfluxX, 'GM_TFLXX',0,Nr,0,1,1,myThid)
914 CALL DIAGNOSTICS_FILL(tfluxY, 'GM_TFLXY',0,Nr,0,1,1,myThid)
915 CALL DIAGNOSTICS_FILL(gradqx, 'GM_dqdx ',0,Nr,0,1,1,myThid)
916 CALL DIAGNOSTICS_FILL(gradqy, 'GM_dqdy ',0,Nr,0,1,1,myThid)
917 CALL DIAGNOSTICS_FILL(surfkz, 'GM_SFLYR',0, 1,0,1,1,myThid)
918 CALL DIAGNOSTICS_FILL(ustar, 'GM_USTAR',0,Nr,0,1,1,myThid)
919 CALL DIAGNOSTICS_FILL(vstar, 'GM_VSTAR',0,Nr,0,1,1,myThid)
920 CALL DIAGNOSTICS_FILL(umc, 'GM_UMC ',0,Nr,0,1,1,myThid)
921 CALL DIAGNOSTICS_FILL(ubar, 'GM_UBAR ',0,Nr,0,1,1,myThid)
922 CALL DIAGNOSTICS_FILL(modesC(1,:,:,:,bi,bj),
923 & 'GM_MODEC',0,Nr,0,1,1,myThid)
924 ENDIF
925 #endif
926
927 #endif /* GM_K3D */
928 RETURN
929 END

  ViewVC Help
Powered by ViewVC 1.1.22