/[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.22 - (show annotations) (download)
Thu Mar 12 15:56:02 2015 UTC (9 years, 3 months ago) by dfer
Branch: MAIN
CVS Tags: checkpoint65k, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m
Changes since 1.21: +51 -53 lines
Fix a few comments + small changes in the code (does not change results of global_ocean.gm_k3d)

1 C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_k3d.F,v 1.21 2015/02/22 01:52:18 m_bates Exp $
2 C $Name: $
3
4 #include "GMREDI_OPTIONS.h"
5
6 C !ROUTINE: GMREDI_K3D
7 C !INTERFACE:
8 SUBROUTINE GMREDI_K3D(
9 I iMin, iMax, jMin, jMax,
10 I sigmaX, sigmaY, sigmaR,
11 I bi, bj, myTime, myThid )
12
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | SUBROUTINE GMREDI_K3D
16 C | o Calculates the 3D diffusivity as per Bates et al. (2013)
17 C *==========================================================*
18 C \ev
19
20 IMPLICIT NONE
21
22 C == Global variables ==
23 #include "SIZE.h"
24 #include "EEPARAMS.h"
25 #include "PARAMS.h"
26 #include "GRID.h"
27 #include "DYNVARS.h"
28 #include "FFIELDS.h"
29 #include "GMREDI.h"
30
31 C !INPUT/OUTPUT PARAMETERS:
32 C == Routine arguments ==
33 C bi, bj :: tile indices
34 C myThid :: My Thread Id. number
35
36 INTEGER iMin,iMax,jMin,jMax
37 _RL sigmaX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
38 _RL sigmaY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
39 _RL sigmaR(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
40 INTEGER bi, bj
41 _RL myTime
42 INTEGER myThid
43
44 #ifdef GM_K3D
45
46 C === Functions ====
47 LOGICAL DIFFERENT_MULTIPLE
48 EXTERNAL DIFFERENT_MULTIPLE
49
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 INTEGER i,j,k,kk,m,kp1
53
54 C update_modes :: Whether to update the eigenmodes
55 LOGICAL update_modes
56
57 C surfk :: index of the depth of the surface layer
58 C kLow_C :: Local version of the index of deepest wet grid cell on tracer grid
59 C kLow_U :: Local version of the index of deepest wet grid cell on U grid
60 C kLow_V :: Local version of the index of deepest wet grid cell on V grid
61 INTEGER surfk(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
62 INTEGER kLow_C(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
63 INTEGER kLow_U(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
64 INTEGER kLow_V(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
65
66 C N2loc :: local N**2
67 C slope :: local slope
68 C Req :: local equatorial deformation radius (m)
69 C deltaH :: local thickness of Eady integration (m)
70 C g_reciprho_sq :: (gravity*recip_rhoConst)**2
71 C M4loc :: local M**4
72 C maxDRhoDz :: maximum value of d(rho)/dz (derived from GM_K3D_minN2)
73 C sigx :: local d(rho)/dx
74 C sigy :: local d(rho)/dy
75 C sigz :: local d(rho)/dz
76 C hsurf :: local surface layer depth
77 C small :: a small number (to avoid floating point exceptions)
78 _RL N2loc(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
79 _RL slope
80 _RL slopeC(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
81 _RL Req
82 _RL deltaH(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
83 _RL g_reciprho_sq
84 _RL M4loc(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
85 _RL M4onN2(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
86 _RL maxDRhoDz
87 _RL sigx, sigy, sigz
88 _RL hsurf, mskp1
89 _RL small
90
91 C dfdy :: gradient of the Coriolis paramater, df/dy, 1/(m*s)
92 C dfdx :: gradient of the Coriolis paramater, df/dx, 1/(m*s)
93 C gradf :: gradient of the Coriolis paramater at a cell centre, 1/(m*s)
94 C Rurms :: a local mixing length used in calculation of urms (m)
95 C RRhines :: The Rhines scale (m)
96 C Rmix :: Mixing length
97 C N2 :: Square of the buoyancy frequency (1/s**2)
98 C N2W :: Square of the buoyancy frequency (1/s**2) averaged to west of grid cell
99 C N2S :: Square of the buoyancy frequency (1/s**2) averaged to south of grid cell
100 C N :: Buoyancy frequency, SQRT(N2)
101 C BVint :: The vertical integral of N (m/s)
102 C ubar :: Zonal velocity on a tracer point (m/s)
103 C vbar :: Meridional velocity on a tracer point (m/s)
104 C Ubaro :: Barotropic velocity (m/s)
105 _RL dfdy( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
106 _RL dfdx( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
107 _RL gradf( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
108 _RL dummy( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
109 _RL Rurms( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
110 _RL RRhines(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
111 _RL Rmix( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
112 _RL N2( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
113 _RL N2W( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
114 _RL N2S( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
115 _RL N( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
116 _RL BVint( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
117 _RL Ubaro( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
118 _RL ubar( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
119
120 _RL tmpU( 1-Olx:sNx+Olx,1-Oly:sNy+Oly )
121 _RL tmpV( 1-Olx:sNx+Olx,1-Oly:sNy+Oly )
122 _RL uFldX( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr )
123 _RL vFldY( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr )
124
125 C Rmid :: Rossby radius (m)
126 C KPV :: Diffusivity (m**2/s)
127 C Kdqdx :: diffusivity multiplied by zonal PV gradient
128 C Kdqdy :: diffusivity multiplied by meridional PV gradient
129 C SlopeX :: isopycnal slope in x direction
130 C SlopeY :: isopycnal slope in y direction
131 C dSigmaDx :: sigmaX averaged onto tracer grid
132 C dSigmaDy :: sigmaY averaged onto tracer grid
133 C tfluxX :: thickness flux in x direction
134 C tfluxY :: thickness flux in y direction
135 C fCoriU :: Coriolis parameter averaged to U points
136 C fCoriV :: Coriolis parameter averaged to V points
137 C cori :: Coriolis parameter forced to be finite near the equator
138 C coriU :: As for cori, but, at U point
139 C coriV :: As for cori, but, at V point
140 C surfkz :: Depth of surface layer (in r units)
141 _RL Rmid(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
142 _RL KPV(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
143 _RL Kdqdy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
144 _RL Kdqdx(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
145 _RL SlopeX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
146 _RL SlopeY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
147 _RL dSigmaDx(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
148 _RL dSigmaDy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
149 _RL tfluxX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
150 _RL tfluxY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
151 _RL cori(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
152 _RL coriU(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
153 _RL coriV(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
154 _RL fCoriU(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
155 _RL fCoriV(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
156 _RL surfkz(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
157
158 C centreX,centreY :: used for calculating averages at centre of cell
159 C numerator,denominator :: of the renormalisation factor
160 C uInt :: column integral of u velocity (sum u*dz)
161 C vInt :: column integral of v velocity (sum v*dz)
162 C KdqdxInt :: column integral of K*dqdx (sum K*dqdx*dz)
163 C KdqdyInt :: column integral of K*dqdy (sum K*dqdy*dz)
164 C uKdqdyInt :: column integral of u*K*dqdy (sum u*K*dqdy*dz)
165 C vKdqdxInt :: column integral of v*K*dqdx (sum v*K*dqdx*dz)
166 C uXiyInt :: column integral of u*Xiy (sum u*Xiy*dz)
167 C vXixInt :: column integral of v*Xix (sum v*Xix*dz)
168 C Renorm :: renormalisation factor at the centre of a cell
169 C RenormU :: renormalisation factor at the western face of a cell
170 C RenormV :: renormalisation factor at the southern face of a cell
171 _RL centreX, centreY
172 _RL numerator, denominator
173 _RL uInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
174 _RL vInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
175 _RL KdqdxInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
176 _RL KdqdyInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
177 _RL uKdqdyInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
178 _RL vKdqdxInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
179 _RL uXiyInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
180 _RL vXixInt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
181 _RL Renorm(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
182 _RL RenormU(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
183 _RL RenormV(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
184
185 C gradqx :: Potential vorticity gradient in x direction
186 C gradqy :: Potential vorticity gradient in y direction
187 C XimX :: Vertical integral of phi_m*K*gradqx
188 C XimY :: Vertical integral of phi_m*K*gradqy
189 C cDopp :: Quasi-Doppler shifted long Rossby wave speed (m/s)
190 C umc :: ubar-c (m/s)
191 C eady :: Eady growth rate (1/s)
192 C urms :: the rms eddy velocity (m/s)
193 C supp :: The suppression factor
194 C ustar :: The eddy induced velocity in the x direction
195 C vstar :: The eddy induced velocity in the y direction
196 C Xix :: Xi in the x direction (m/s**2)
197 C Xiy :: Xi in the y direction (m/s**2)
198 _RL gradqx(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
199 _RL gradqy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
200 _RL XimX(GM_K3D_NModes,1-Olx:sNx+Olx,1-Oly:sNy+Oly)
201 _RL XimY(GM_K3D_NModes,1-Olx:sNx+Olx,1-Oly:sNy+Oly)
202 _RL cDopp(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
203 _RL umc( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
204 _RL eady( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
205 _RL urms( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
206 _RL supp( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
207 _RL ustar(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
208 _RL vstar(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
209 _RL Xix( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
210 _RL Xiy( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
211 #ifdef GM_K3D_PASSIVE
212 C psistar :: eddy induced streamfunction in the y direction
213 _RL psistar(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nr)
214 #endif
215
216 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
217
218 C ======================================
219 C Initialise some variables
220 C ======================================
221 small = TINY(zeroRL)
222 update_modes=.FALSE.
223 IF ( DIFFERENT_MULTIPLE(GM_K3D_vecFreq,myTime,deltaTClock) )
224 & update_modes=.TRUE.
225
226 DO j=1-Oly,sNy+Oly
227 DO i=1-Olx,sNx+Olx
228 kLow_C(i,j) = kLowC(i,j,bi,bj)
229 ENDDO
230 ENDDO
231 DO j=1-Oly,sNy+Oly
232 DO i=1-Olx+1,sNx+Olx
233 kLow_U(i,j) = MIN( kLow_C(i,j), kLow_C(i-1,j) )
234 ENDDO
235 ENDDO
236 DO j=1-Oly+1,sNy+Oly
237 DO i=1-Olx,sNx+Olx
238 kLow_V(i,j) = MIN( kLow_C(i,j), kLow_C(i,j-1) )
239 ENDDO
240 ENDDO
241
242 C Dummy values for the edges. This does not affect the results
243 C but avoids problems when solving for the eigenvalues.
244 i=1-Olx
245 DO j=1-Oly,sNy+Oly
246 kLow_U(i,j) = 0
247 ENDDO
248 j=1-Oly
249 DO i=1-Olx,sNx+Olx
250 kLow_V(i,j) = 0
251 ENDDO
252
253 g_reciprho_sq = (gravity*recip_rhoConst)**2
254 C Gradient of Coriolis
255 DO j=1-Oly+1,sNy+Oly
256 DO i=1-Olx+1,sNx+Olx
257 dfdx(i,j) = ( fCori(i,j,bi,bj)-fCori(i-1,j,bi,bj) )
258 & *recip_dxC(i,j,bi,bj)
259 dfdy(i,j) = ( fCori(i,j,bi,bj)-fCori(i,j-1,bi,bj) )
260 & *recip_dyC(i,j,bi,bj)
261 ENDDO
262 ENDDO
263
264 C Coriolis at C points enforcing a minimum value so
265 C that it is defined at the equator
266 DO j=1-Oly,sNy+Oly
267 DO i=1-Olx,sNx+Olx
268 cori(i,j) = SIGN( MAX( ABS(fCori(i,j,bi,bj)),GM_K3D_minCori ),
269 & fCori(i,j,bi,bj) )
270 ENDDO
271 ENDDO
272 C Coriolis at U and V points
273 DO j=1-Oly,sNy+Oly
274 DO i=1-Olx+1,sNx+Olx
275 C Limited so that the inverse is defined at the equator
276 coriU(i,j) = op5*( cori(i,j)+cori(i-1,j) )
277 coriU(i,j) = SIGN( MAX( ABS(coriU(i,j)),GM_K3D_minCori ),
278 & coriU(i,j) )
279
280 C Not limited
281 fCoriU(i,j) = op5*( fCori(i,j,bi,bj)+fCori(i-1,j,bi,bj) )
282 ENDDO
283 ENDDO
284 DO j=1-Oly+1,sNy+Oly
285 DO i=1-Olx,sNx+Olx
286 C Limited so that the inverse is defined at the equator
287 coriV(i,j) = op5*( cori(i,j)+cori(i,j-1) )
288 coriV(i,j) = SIGN( MAX( ABS(coriV(i,j)),GM_K3D_minCori ),
289 & coriV(i,j) )
290
291 C Not limited
292 fCoriV(i,j) = op5*( fCori(i,j,bi,bj)+fCori(i,j-1,bi,bj) )
293 ENDDO
294 ENDDO
295 C Computing beta
296 IF ( selectCoriMap.EQ.1 ) THEN
297 DO j=1-Oly,sNy+Oly
298 DO i=1-Olx,sNx+Olx
299 gradf(i,j) = beta
300 ENDDO
301 ENDDO
302 ELSEIF ( selectCoriMap.EQ.2 ) THEN
303 DO j=1-Oly,sNy+Oly
304 DO i=1-Olx,sNx+Olx
305 gradf(i,j) = recip_rSphere*fCoriCos(i,j,bi,bj)
306 ENDDO
307 ENDDO
308 ENDIF
309 C Some dummy values at the edges
310 i=1-Olx
311 DO j=1-Oly,sNy+Oly
312 coriU(i,j)=cori(i,j)
313 fCoriU(i,j)=fCori(i,j,bi,bj)
314 ENDDO
315 j=1-Oly
316 DO i=1-Olx,sNx+Olx
317 coriV(i,j)=cori(i,j)
318 fCoriV(i,j)=fCori(i,j,bi,bj)
319 ENDDO
320
321 C Zeroing some cumulative fields
322 DO j=1-Oly,sNy+Oly
323 DO i=1-Olx,sNx+Olx
324 eady(i,j) = zeroRL
325 BVint(i,j) = zeroRL
326 Ubaro(i,j) = zeroRL
327 deltaH(i,j) = zeroRL
328 ENDDO
329 ENDDO
330 DO k=1,Nr
331 DO j=1-Oly,sNy+Oly
332 DO i=1-Olx,sNx+Olx
333 slopeC(i,j,k)=zeroRL
334 ENDDO
335 ENDDO
336 ENDDO
337
338 C initialise remaining 2d variables
339 DO j=1-Oly,sNy+Oly
340 DO i=1-Olx,sNx+Olx
341 dfdy(i,j)=zeroRL
342 dfdy(i,j)=zeroRL
343 Rurms(i,j)=zeroRL
344 RRhines(i,j)=zeroRL
345 Rmix(i,j)=zeroRL
346 ENDDO
347 ENDDO
348 C initialise remaining 3d variables
349 DO k=1,Nr
350 DO j=1-Oly,sNy+Oly
351 DO i=1-Olx,sNx+Olx
352 N2loc(i,j,k)=GM_K3D_minN2
353 N2W(i,j,k) = GM_K3D_minN2
354 N2S(i,j,k) = GM_K3D_minN2
355 M4loc(i,j,k)=zeroRL
356 M4onN2(i,j,k)=zeroRL
357 urms(i,j,k)=zeroRL
358 SlopeX(i,j,k)=zeroRL
359 SlopeY(i,j,k)=zeroRL
360 dSigmaDx(i,j,k)=zeroRL
361 dSigmaDy(i,j,k)=zeroRL
362 gradqx(i,j,k)=zeroRL
363 gradqy(i,j,k)=zeroRL
364 ENDDO
365 ENDDO
366 ENDDO
367
368 C Find the zonal velocity at the cell centre
369 #ifdef ALLOW_EDDYPSI
370 IF (GM_InMomAsStress) THEN
371 DO k=1,Nr
372 DO i = 1-olx,snx+olx
373 DO j = 1-oly,sny+oly
374 uFldX(i,j,k) = uEulerMean(i,j,k,bi,bj)
375 vFldY(i,j,k) = vEulerMean(i,j,k,bi,bj)
376 ENDDO
377 ENDDO
378 ENDDO
379 ELSE
380 #endif
381 DO k=1,Nr
382 DO i = 1-olx,snx+olx
383 DO j = 1-oly,sny+oly
384 uFldX(i,j,k) = uVel(i,j,k,bi,bj)
385 vFldY(i,j,k) = vVel(i,j,k,bi,bj)
386 ENDDO
387 ENDDO
388 ENDDO
389 #ifdef ALLOW_EDDYPSI
390 ENDIF
391 #endif
392
393 C The following comes from rotate_uv2en_rl
394 C This code does two things:
395 C 1) go from C grid velocity points to A grid velocity points
396 C 2) go from model grid directions to east/west directions
397 DO k = 1,Nr
398
399 DO i = 1-Olx,sNx+Olx
400 j=sNy+Oly
401 tmpU(i,j)=zeroRL
402 tmpV(i,j)=zeroRL
403 ENDDO
404 DO j = 1-Oly,sNy+Oly-1
405 i=sNx+Olx
406 tmpU(i,j)=zeroRL
407 tmpV(i,j)=zeroRL
408 DO i = 1-Olx,sNx+Olx-1
409 tmpU(i,j) = 0.5 _d 0
410 & *( uFldX(i+1,j,k) + uFldX(i,j,k) )
411 tmpV(i,j) = 0.5 _d 0
412 & *( vFldY(i,j+1,k) + vFldY(i,j,k) )
413
414 tmpU(i,j) = tmpU(i,j) * maskC(i,j,k,bi,bj)
415 tmpV(i,j) = tmpV(i,j) * maskC(i,j,k,bi,bj)
416 ENDDO
417 ENDDO
418
419 DO j = 1-oly,sny+oly
420 DO i = 1-olx,snx+olx
421 ENDDO
422 ENDDO
423
424 C rotation
425 DO j = 1-oly,sny+oly
426 DO i = 1-olx,snx+olx
427 ubar(i,j,k) =
428 & angleCosC(i,j,bi,bj)*tmpU(i,j)
429 & -angleSinC(i,j,bi,bj)*tmpV(i,j)
430 ENDDO
431 ENDDO
432 ENDDO
433
434 C Square of the buoyancy frequency at the top of a grid cell
435 C Enforce a minimum N2
436 C Mask N2, so it is zero at bottom
437 DO k=2,Nr
438 DO j=1-Oly,sNy+Oly
439 DO i=1-Olx,sNx+Olx
440 N2(i,j,k) = -gravity*recip_rhoConst*sigmaR(i,j,k)
441 N2(i,j,k) = MAX(N2(i,j,k),GM_K3D_minN2)*maskC(i,j,k,bi,bj)
442 N(i,j,k) = SQRT(N2(i,j,k))
443 ENDDO
444 ENDDO
445 ENDDO
446 C N2(k=1) is always zero
447 DO j=1-Oly,sNy+Oly
448 DO i=1-Olx,sNx+Olx
449 N2(i,j,1) = zeroRL
450 N(i,j,1) = zeroRL
451 ENDDO
452 ENDDO
453 C Calculate the minimum drho/dz
454 maxDRhoDz = -rhoConst*GM_K3D_minN2/gravity
455
456 C Calculate the barotropic velocity by vertically integrating
457 C and the dividing by the depth of the water column
458 C Note that Ubaro is at the C-point.
459 DO k=1,Nr
460 DO j=1-Oly,sNy+Oly
461 DO i=1-Olx,sNx+Olx
462 Ubaro(i,j) = Ubaro(i,j) +
463 & drF(k)*hfacC(i,j,k,bi,bj)*ubar(i,j,k)
464 ENDDO
465 ENDDO
466 ENDDO
467 DO j=1-Oly,sNy+Oly
468 DO i=1-Olx,sNx+Olx
469 IF (kLow_C(i,j).GT.0) THEN
470 C The minus sign is because r_Low<0
471 Ubaro(i,j) = -Ubaro(i,j)/r_Low(i,j,bi,bj)
472 ENDIF
473 ENDDO
474 ENDDO
475
476 C Integrate the buoyancy frequency vertically using the trapezoidal method.
477 C Assume that N(z=-H)=0
478 DO k=1,Nr
479 kp1 = min(k+1,Nr)
480 mskp1 = oneRL
481 IF ( k.EQ.Nr ) mskp1 = zeroRL
482 DO j=1-Oly,sNy+Oly
483 DO i=1-Olx,sNx+Olx
484 BVint(i,j) = BVint(i,j) + hFacC(i,j,k,bi,bj)*drF(k)
485 & *op5*(N(i,j,k)+mskp1*N(i,j,kp1))
486 ENDDO
487 ENDDO
488 ENDDO
489
490 C Calculate the eigenvalues and eigenvectors
491 IF (update_modes) THEN
492 CALL GMREDI_CALC_EIGS(
493 I iMin,iMax,jMin,jMax,bi,bj,N2,myThid,
494 I kLow_C, maskC(:,:,:,bi,bj),
495 I hfacC(:,:,:,bi,bj), recip_hfacC(:,:,:,bi,bj),
496 I R_Low(:,:,bi,bj), 1, .TRUE.,
497 O Rmid, modesC(:,:,:,:,bi,bj))
498
499 C Calculate the Rossby Radius
500 DO j=1-Oly+1,sNy+Oly
501 DO i=1-Olx+1,sNx+Olx
502 Req = SQRT(BVint(i,j)/(2*pi*gradf(i,j)))
503 Rdef(i,j,bi,bj) = MIN(Rmid(i,j),Req)
504 ENDDO
505 ENDDO
506 ENDIF
507
508 C Average dsigma/dx and dsigma/dy onto the centre points
509
510 DO k=1,Nr
511 DO j=1-Oly,sNy+Oly-1
512 DO i=1-Olx,sNx+Olx-1
513 dSigmaDx(i,j,k) = op5*(sigmaX(i,j,k)+sigmaX(i+1,j,k))
514 dSigmaDy(i,j,k) = op5*(sigmaY(i,j,k)+sigmaY(i,j+1,k))
515 ENDDO
516 ENDDO
517 ENDDO
518
519 C ===============================
520 C Calculate the Eady growth rate
521 C ===============================
522 DO k=1,Nr
523
524 kp1 = min(k+1,Nr)
525 mskp1 = oneRL
526 IF ( k.EQ.Nr ) mskp1 = zeroRL
527
528 DO j=1-Oly,sNy+Oly-1
529 DO i=1-Olx,sNx+Olx-1
530 M4loc(i,j,k) = g_reciprho_sq*( dSigmaDx(i,j,k)**2
531 & +dSigmaDy(i,j,k)**2 )
532 N2loc(i,j,k) = op5*(N2(i,j,k)+mskp1*N2(i,j,kp1))
533 ENDDO
534 ENDDO
535 C The bottom of the grid cell is shallower than the top
536 C integration level, so, advance the depth.
537 IF (-rF(k+1) .LE. GM_K3D_EadyMinDepth) CYCLE
538
539 C Do not bother going any deeper since the top of the
540 C cell is deeper than the bottom integration level
541 IF (-rF(k).GE.GM_K3D_EadyMaxDepth) EXIT
542
543 C We are in the integration depth range
544 DO j=1-Oly,sNy+Oly-1
545 DO i=1-Olx,sNx+Olx-1
546 IF ( (kLow_C(i,j).GE.k) .AND.
547 & (-hMixLayer(i,j,bi,bj).LE.-rC(k)) ) THEN
548
549 slopeC(i,j,k) = SQRT(M4loc(i,j,k))/N2loc(i,j,k)
550 C Limit the slope. Note, this is not all the Eady calculations.
551 IF (slopeC(i,j,k).LE.GM_maxSlope) THEN
552 M4onN2(i,j,k) = M4loc(i,j,k)/N2loc(i,j,k)
553 ELSE
554 slopeC(i,j,k) = GM_maxslope
555 M4onN2(i,j,k) = SQRT(M4loc(i,j,k))*GM_maxslope
556 ENDIF
557 eady(i,j) = eady(i,j)
558 & + hfacC(i,j,k,bi,bj)*drF(k)*M4onN2(i,j,k)
559 deltaH(i,j) = deltaH(i,j) + drF(k)
560 ENDIF
561 ENDDO
562 ENDDO
563 ENDDO
564
565 DO j=1-Oly,sNy+Oly
566 DO i=1-Olx,sNx+Olx
567 C If the minimum depth for the integration is deeper than the ocean
568 C bottom OR the mixed layer is deeper than the maximum depth of
569 C integration, we set the Eady growth rate to something small
570 C to avoid floating point exceptions.
571 C Later, these areas will be given a small diffusivity.
572 IF (deltaH(i,j).EQ.zeroRL) THEN
573 eady(i,j) = small
574
575 C Otherwise, divide over the integration and take the square root
576 C to actually find the Eady growth rate.
577 ELSE
578 eady(i,j) = SQRT(eady(i,j)/deltaH(i,j))
579
580 ENDIF
581
582 ENDDO
583 ENDDO
584
585 C ======================================
586 C Calculate the diffusivity
587 C ======================================
588 DO j=1-Oly+1,sNy+Oly
589 DO i=1-Olx+1,sNx+Olx-1
590 C Calculate the Visbeck velocity
591 Rurms(i,j) = MIN(Rdef(i,j,bi,bj),GM_K3D_Rmax)
592 urms(i,j,1) = GM_K3D_Lambda*eady(i,j)*Rurms(i,j)
593 C Set the bottom urms to zero
594 k=kLow_C(i,j)
595 IF (k.GT.0) urms(i,j,k) = 0.0
596
597 C Calculate the Rhines scale
598 RRhines(i,j) = SQRT(urms(i,j,1)/gradf(i,j))
599
600 C Calculate the estimated length scale
601 Rmix(i,j) = MIN(Rdef(i,j,bi,bj), RRhines(i,j))
602 Rmix(i,j) = MAX(Rmix(i,j),GM_K3D_Rmin)
603
604 C Calculate the Doppler shifted long Rossby wave speed
605 C Ubaro is at the C-point.
606 cDopp(i,j) = Ubaro(i,j)
607 & - gradf(i,j)*Rdef(i,j,bi,bj)*Rdef(i,j,bi,bj)
608 C Limit the wave speed to the namelist variable GM_K3D_maxC
609 IF (ABS(cDopp(i,j)).GT.GM_K3D_maxC) THEN
610 cDopp(i,j) = MAX(GM_K3D_maxC, cDopp(i,j))
611 ENDIF
612
613 ENDDO
614 ENDDO
615
616 C Project the surface urms to the subsurface using the first baroclinic mode
617 CALL GMREDI_CALC_URMS(
618 I iMin,iMax,jMin,jMax,bi,bj,N2,myThid,
619 U urms)
620
621 C Calculate the diffusivity (on the mass grid)
622 DO k=1,Nr
623 DO j=1-Oly,sNy+Oly
624 DO i=1-Olx,sNx+Olx
625 IF (k.LE.kLow_C(i,j)) THEN
626 IF (deltaH(i,j).EQ.zeroRL) THEN
627 K3D(i,j,k,bi,bj) = GM_K3D_smallK
628 ELSE
629 IF (urms(i,j,k).EQ.0.0) THEN
630 K3D(i,j,k,bi,bj) = GM_K3D_smallK
631 ELSE
632 umc(i,j,k) =ubar(i,j,k) - cDopp(i,j)
633 supp(i,j,k)=1./(1.+GM_K3D_b1*umc(i,j,k)**2/urms(i,j,1)**2)
634 C 2*Rmix gives the diameter
635 K3D(i,j,k,bi,bj) = GM_K3D_gamma*urms(i,j,k)
636 & *2.*Rmix(i,j)*supp(i,j,k)
637 ENDIF
638
639 C Enforce lower and upper bounds on the diffusivity
640 K3D(i,j,k,bi,bj) = MIN(K3D(i,j,k,bi,bj),GM_maxK3D)
641 K3D(i,j,k,bi,bj) = MAX(K3D(i,j,k,bi,bj),GM_K3D_smallK)
642 ENDIF
643 ENDIF
644 ENDDO
645 ENDDO
646 ENDDO
647
648 C ======================================
649 C Find the PV gradient
650 C ======================================
651 C Calculate the surface layer thickness.
652 C Use hMixLayer (calculated in model/src/calc_oce_mxlayer)
653 C for the mixed layer depth.
654
655 C Enforce a minimum surface layer depth
656 DO j=1-Oly,sNy+Oly
657 DO i=1-Olx,sNx+Olx
658 surfkz(i,j) = MIN(-GM_K3D_surfMinDepth,-hMixLayer(i,j,bi,bj))
659 surfkz(i,j) = MAX(surfkz(i,j),R_low(i,j,bi,bj))
660 IF(maskC(i,j,1,bi,bj).EQ.0.0) surfkz(i,j)=0.0
661 surfk(i,j) = 0
662 ENDDO
663 ENDDO
664 DO k=1,Nr
665 DO j=1-Oly,sNy+Oly
666 DO i=1-Olx,sNx+Olx
667 IF (rF(k).GT.surfkz(i,j) .AND. surfkz(i,j).GE.rF(k+1))
668 & surfk(i,j) = k
669 ENDDO
670 ENDDO
671 ENDDO
672
673 C Calculate the isopycnal slope
674 DO j=1-Oly,sNy+Oly-1
675 DO i=1-Olx,sNx+Olx-1
676 SlopeX(i,j,1) = zeroRL
677 SlopeY(i,j,1) = zeroRL
678 ENDDO
679 ENDDO
680 DO k=2,Nr
681 DO j=1-Oly+1,sNy+Oly
682 DO i=1-Olx+1,sNx+Olx
683 IF(surfk(i,j).GE.kLowC(i,j,bi,bj)) THEN
684 C If the surface layer is thinner than the water column
685 C the set the slope to zero to avoid problems.
686 SlopeX(i,j,k) = zeroRL
687 SlopeY(i,j,k) = zeroRL
688
689 ELSE
690 C Calculate the zonal slope at the western cell face (U grid)
691 sigz = MIN( op5*(sigmaR(i,j,k)+sigmaR(i-1,j,k)), maxDRhoDz )
692 sigx = op5*( sigmaX(i,j,k)+sigmaX(i,j,k-1) )
693 slope = sigx/sigz
694 IF(ABS(slope).GT.GM_maxSlope)
695 & slope = SIGN(GM_maxSlope,slope)
696 SlopeX(i,j,k)=-maskW(i,j,k-1,bi,bj)*maskW(i,j,k,bi,bj)*slope
697
698 C Calculate the meridional slope at the southern cell face (V grid)
699 sigz = MIN( op5*(sigmaR(i,j,k)+sigmaR(i,j-1,k)), maxDRhoDz )
700 sigy = op5*( sigmaY(i,j,k) + sigmaY(i,j,k-1) )
701 slope = sigy/sigz
702 IF(ABS(slope).GT.GM_maxSlope)
703 & slope = SIGN(GM_maxSlope,slope)
704 SlopeY(i,j,k)=-maskS(i,j,k-1,bi,bj)*maskS(i,j,k,bi,bj)*slope
705 ENDIF
706 ENDDO
707 ENDDO
708 ENDDO
709
710 C Calculate the thickness flux and diffusivity. These may be altered later
711 C depending on namelist options.
712 C Enforce a zero slope bottom boundary condition for the bottom most cells (k=Nr)
713 k=Nr
714 DO j=1-Oly,sNy+Oly
715 DO i=1-Olx,sNx+Olx
716 C Zonal thickness flux at the western cell face
717 tfluxX(i,j,k) = -fCoriU(i,j)*SlopeX(i,j,k)
718 & *recip_drF(k)*recip_hFacW(i,j,k,bi,bj)
719 C Meridional thickness flux at the southern cell face
720 tfluxY(i,j,k) = -fCoriV(i,j)*SlopeY(i,j,k)
721 & *recip_drF(k)*recip_hFacS(i,j,k,bi,bj)
722
723 C Use the interior diffusivity. Note: if we are using a
724 C constant diffusivity KPV is overwritten later
725 KPV(i,j,k) = K3D(i,j,k,bi,bj)
726
727 ENDDO
728 ENDDO
729
730 C Calculate the thickness flux and diffusivity and for other cells (k<Nr)
731 DO k=Nr-1,1,-1
732 DO j=1-Oly,sNy+Oly
733 DO i=1-Olx,sNx+Olx
734 C Zonal thickness flux at the western cell face
735 tfluxX(i,j,k)=-fCoriU(i,j)*(SlopeX(i,j,k)-SlopeX(i,j,k+1))
736 & *recip_drF(k)*recip_hFacW(i,j,k,bi,bj)
737 & *maskW(i,j,k,bi,bj)
738
739 C Meridional thickness flux at the southern cell face
740 tfluxY(i,j,k)=-fCoriV(i,j)*(SlopeY(i,j,k)-SlopeY(i,j,k+1))
741 & *recip_drF(k)*recip_hFacS(i,j,k,bi,bj)
742 & *maskS(i,j,k,bi,bj)
743
744 C Use the interior diffusivity. Note: if we are using a
745 C constant diffusivity KPV is overwritten later
746 KPV(i,j,k) = K3D(i,j,k,bi,bj)
747 ENDDO
748 ENDDO
749 ENDDO
750
751 C Special treatment for the surface layer (if necessary), which overwrites
752 C values in the previous loops.
753 IF (GM_K3D_ThickSheet .OR. GM_K3D_surfK) THEN
754 DO k=Nr-1,1,-1
755 DO j=1-Oly,sNy+Oly
756 DO i=1-Olx,sNx+Olx
757 IF(k.LE.surfk(i,j)) THEN
758 C We are in the surface layer. Change the thickness flux
759 C and diffusivity as necessary.
760
761 IF (GM_K3D_ThickSheet) THEN
762 C We are in the surface layer, so set the thickness flux
763 C based on the average slope over the surface layer
764 C If we are on the edge of a "cliff" the surface layer at the
765 C centre of the grid point could be deeper than the U or V point.
766 C So, we ensure that we always take a sensible slope.
767 IF(kLow_U(i,j).LT.surfk(i,j)) THEN
768 kk=kLow_U(i,j)
769 hsurf = -rLowW(i,j,bi,bj)
770 ELSE
771 kk=surfk(i,j)
772 hsurf = -surfkz(i,j)
773 ENDIF
774 IF(kk.GT.0) THEN
775 IF(kk.EQ.Nr) THEN
776 tfluxX(i,j,k) = -fCoriU(i,j)*maskW(i,j,k,bi,bj)
777 & *SlopeX(i,j,kk)/hsurf
778 ELSE
779 tfluxX(i,j,k) = -fCoriU(i,j)*maskW(i,j,k,bi,bj)
780 & *( SlopeX(i,j,kk)-SlopeX(i,j,kk+1) )/hsurf
781 ENDIF
782 ELSE
783 tfluxX(i,j,k) = zeroRL
784 ENDIF
785
786 IF(kLow_V(i,j).LT.surfk(i,j)) THEN
787 kk=kLow_V(i,j)
788 hsurf = -rLowS(i,j,bi,bj)
789 ELSE
790 kk=surfk(i,j)
791 hsurf = -surfkz(i,j)
792 ENDIF
793 IF(kk.GT.0) THEN
794 IF(kk.EQ.Nr) THEN
795 tfluxY(i,j,k) = -fCoriV(i,j)*maskS(i,j,k,bi,bj)
796 & *SlopeY(i,j,kk)/hsurf
797 ELSE
798 tfluxY(i,j,k) = -fCoriV(i,j)*maskS(i,j,k,bi,bj)
799 & *( SlopeY(i,j,kk)-SlopeY(i,j,kk+1) )/hsurf
800 ENDIF
801 ELSE
802 tfluxY(i,j,k) = zeroRL
803 ENDIF
804 ENDIF
805
806 IF (GM_K3D_surfK) THEN
807 C Use a constant K in the surface layer.
808 KPV(i,j,k) = GM_K3D_constK
809 ENDIF
810 ENDIF
811 ENDDO
812 ENDDO
813 ENDDO
814 ENDIF
815
816 C Calculate gradq
817 IF (GM_K3D_beta_eq_0) THEN
818 C Ignore beta in the calculation of grad(q)
819 DO k=1,Nr
820 DO j=1-Oly+1,sNy+Oly
821 DO i=1-Olx+1,sNx+Olx
822 gradqx(i,j,k) = maskW(i,j,k,bi,bj)*tfluxX(i,j,k)
823 gradqy(i,j,k) = maskS(i,j,k,bi,bj)*tfluxY(i,j,k)
824 ENDDO
825 ENDDO
826 ENDDO
827
828 ELSE
829 C Do not ignore beta
830 DO k=1,Nr
831 DO j=1-Oly+1,sNy+Oly
832 DO i=1-Olx+1,sNx+Olx
833 gradqx(i,j,k) = maskW(i,j,k,bi,bj)*(dfdx(i,j)+tfluxX(i,j,k))
834 gradqy(i,j,k) = maskS(i,j,k,bi,bj)*(dfdy(i,j)+tfluxY(i,j,k))
835 ENDDO
836 ENDDO
837 ENDDO
838 ENDIF
839
840 C ======================================
841 C Find Xi and the eddy induced velocities
842 C ======================================
843 C Find the buoyancy frequency at the west and south faces of a cell
844 C This is necessary to find the eigenvectors at those points
845 DO k=1,Nr
846 DO j=1-Oly+1,sNy+Oly
847 DO i=1-Olx+1,sNx+Olx
848 N2W(i,j,k) = maskW(i,j,k,bi,bj)
849 & *( N2(i,j,k)+N2(i-1,j,k) )
850 N2S(i,j,k) = maskS(i,j,k,bi,bj)
851 & *( N2(i,j,k)+N2(i,j-1,k) )
852 ENDDO
853 ENDDO
854 ENDDO
855
856 C If GM_K3D_use_constK=.TRUE., the diffusivity for the eddy transport is
857 C set to a constant equal to GM_K3D_constK.
858 C If the diffusivity is constant the method here is the same as GM.
859 C If GM_K3D_constRedi=.TRUE. K3D will be set equal to GM_K3D_constK later.
860 IF(GM_K3D_use_constK) THEN
861 DO k=1,Nr
862 DO j=1-Oly,sNy+Oly
863 DO i=1-Olx,sNx+Olx
864 KPV(i,j,k) = GM_K3D_constK
865 ENDDO
866 ENDDO
867 ENDDO
868 ENDIF
869
870 IF (.NOT. GM_K3D_smooth) THEN
871 C Do not expand K grad(q) => no smoothing
872 C May only be done with a constant K, otherwise the
873 C integral constraint is violated.
874 DO k=1,Nr
875 DO j=1-Oly,sNy+Oly
876 DO i=1-Olx,sNx+Olx
877 Xix(i,j,k) = -maskW(i,j,k,bi,bj)*KPV(i,j,k)*gradqx(i,j,k)
878 Xiy(i,j,k) = -maskS(i,j,k,bi,bj)*KPV(i,j,k)*gradqy(i,j,k)
879 ENDDO
880 ENDDO
881 ENDDO
882
883 ELSE
884 C Expand K grad(q) in terms of baroclinic modes to smooth
885 C and satisfy the integral constraint
886
887 C Start with the X direction
888 C ------------------------------
889 C Calculate the eigenvectors at the West face of a cell
890 IF (update_modes) THEN
891 CALL GMREDI_CALC_EIGS(
892 I iMin,iMax,jMin,jMax,bi,bj,N2W,myThid,
893 I kLow_U,maskW(:,:,:,bi,bj),
894 I hfacW(:,:,:,bi,bj),recip_hfacW(:,:,:,bi,bj),
895 I rLowW(:,:,bi,bj),GM_K3D_NModes,.FALSE.,
896 O dummy,modesW(:,:,:,:,bi,bj))
897 ENDIF
898
899 C Calculate Xi_m at the west face of a cell
900 DO j=1-Oly,sNy+Oly
901 DO i=1-Olx,sNx+Olx
902 DO m=1,GM_K3D_NModes
903 XimX(m,i,j) = zeroRL
904 ENDDO
905 ENDDO
906 ENDDO
907 DO k=1,Nr
908 DO j=1-Oly,sNy+Oly
909 DO i=1-Olx,sNx+Olx
910 DO m=1,GM_K3D_NModes
911 Kdqdx(i,j,k) = KPV(i,j,k)*gradqx(i,j,k)
912 XimX(m,i,j) = XimX(m,i,j)
913 & - maskW(i,j,k,bi,bj)*drF(k)*hfacW(i,j,k,bi,bj)
914 & *Kdqdx(i,j,k)*modesW(m,i,j,k,bi,bj)
915 ENDDO
916 ENDDO
917 ENDDO
918 ENDDO
919
920 C Calculate Xi in the X direction at the west face
921 DO k=1,Nr
922 DO j=1-Oly,sNy+Oly
923 DO i=1-Olx,sNx+Olx
924 Xix(i,j,k) = zeroRL
925 ENDDO
926 ENDDO
927 ENDDO
928 DO k=1,Nr
929 DO j=1-Oly,sNy+Oly
930 DO i=1-Olx,sNx+Olx
931 DO m=1,GM_K3D_NModes
932 Xix(i,j,k) = Xix(i,j,k)
933 & + maskW(i,j,k,bi,bj)*XimX(m,i,j)*modesW(m,i,j,k,bi,bj)
934 ENDDO
935 ENDDO
936 ENDDO
937 ENDDO
938
939 C Now the Y direction
940 C ------------------------------
941 C Calculate the eigenvectors at the West face of a cell
942 IF (update_modes) THEN
943 CALL GMREDI_CALC_EIGS(
944 I iMin,iMax,jMin,jMax,bi,bj,N2S,myThid,
945 I kLow_V,maskS(:,:,:,bi,bj),
946 I hfacS(:,:,:,bi,bj),recip_hfacS(:,:,:,bi,bj),
947 I rLowS(:,:,bi,bj), GM_K3D_NModes, .FALSE.,
948 O dummy,modesS(:,:,:,:,bi,bj))
949 ENDIF
950
951 DO j=1-Oly,sNy+Oly
952 DO i=1-Olx,sNx+Olx
953 DO m=1,GM_K3D_NModes
954 XimY(m,i,j) = zeroRL
955 ENDDO
956 ENDDO
957 ENDDO
958 DO k=1,Nr
959 DO j=1-Oly,sNy+Oly
960 DO i=1-Olx,sNx+Olx
961 DO m=1,GM_K3D_NModes
962 Kdqdy(i,j,k) = KPV(i,j,k)*gradqy(i,j,k)
963 XimY(m,i,j) = XimY(m,i,j)
964 & - drF(k)*hfacS(i,j,k,bi,bj)
965 & *Kdqdy(i,j,k)*modesS(m,i,j,k,bi,bj)
966 ENDDO
967 ENDDO
968 ENDDO
969 ENDDO
970
971 C Calculate Xi for Y direction at the south face
972 DO k=1,Nr
973 DO j=1-Oly,sNy+Oly
974 DO i=1-Olx,sNx+Olx
975 Xiy(i,j,k) = zeroRL
976 ENDDO
977 ENDDO
978 ENDDO
979 DO k=1,Nr
980 DO j=1-Oly,sNy+Oly
981 DO i=1-Olx,sNx+Olx
982 DO m=1,GM_K3D_NModes
983 Xiy(i,j,k) = Xiy(i,j,k)
984 & + maskS(i,j,k,bi,bj)*XimY(m,i,j)*modesS(m,i,j,k,bi,bj)
985 ENDDO
986 ENDDO
987 ENDDO
988 ENDDO
989
990 C ENDIF (.NOT. GM_K3D_smooth)
991 ENDIF
992
993 C Calculate the renormalisation factor
994 DO j=1-Oly,sNy+Oly
995 DO i=1-Olx,sNx+Olx
996 uInt(i,j)=zeroRL
997 vInt(i,j)=zeroRL
998 KdqdyInt(i,j)=zeroRL
999 KdqdxInt(i,j)=zeroRL
1000 uKdqdyInt(i,j)=zeroRL
1001 vKdqdxInt(i,j)=zeroRL
1002 uXiyInt(i,j)=zeroRL
1003 vXixInt(i,j)=zeroRL
1004 Renorm(i,j)=oneRL
1005 RenormU(i,j)=oneRL
1006 RenormV(i,j)=oneRL
1007 ENDDO
1008 ENDDO
1009 DO k=1,Nr
1010 DO j=1-Oly,sNy+Oly-1
1011 DO i=1-Olx,sNx+Olx-1
1012 centreX = op5*(uVel(i,j,k,bi,bj)+uVel(i+1,j,k,bi,bj))
1013 centreY = op5*(Kdqdy(i,j,k) +Kdqdy(i,j+1,k) )
1014 C For the numerator
1015 uInt(i,j) = uInt(i,j)
1016 & + centreX*hfacC(i,j,k,bi,bj)*drF(k)
1017 KdqdyInt(i,j) = KdqdyInt(i,j)
1018 & + centreY*hfacC(i,j,k,bi,bj)*drF(k)
1019 uKdqdyInt(i,j) = uKdqdyInt(i,j)
1020 & + centreX*centreY*hfacC(i,j,k,bi,bj)*drF(k)
1021 C For the denominator
1022 centreY = op5*(Xiy(i,j,k) + Xiy(i,j+1,k))
1023 uXiyInt(i,j) = uXiyInt(i,j)
1024 & + centreX*centreY*hfacC(i,j,k,bi,bj)*drF(k)
1025
1026 centreX = op5*(Kdqdx(i,j,k) +Kdqdx(i+1,j,k))
1027 centreY = op5*(vVel(i,j,k,bi,bj)+vVel(i,j+1,k,bi,bj) )
1028 C For the numerator
1029 vInt(i,j) = vInt(i,j)
1030 & + centreY*hfacC(i,j,k,bi,bj)*drF(k)
1031 KdqdxInt(i,j) = KdqdxInt(i,j)
1032 & + CentreX*hfacC(i,j,k,bi,bj)*drF(k)
1033 vKdqdxInt(i,j) = vKdqdxInt(i,j)
1034 & + centreY*centreX*hfacC(i,j,k,bi,bj)*drF(k)
1035 C For the denominator
1036 centreX = op5*(Xix(i,j,k) + Xix(i+1,j,k))
1037 vXixInt(i,j) = vXixInt(i,j)
1038 & + centreY*centreX*hfacC(i,j,k,bi,bj)*drF(k)
1039
1040 ENDDO
1041 ENDDO
1042 ENDDO
1043
1044 DO j=1-Oly,sNy+Oly-1
1045 DO i=1-Olx,sNx+Olx-1
1046 IF (kLowC(i,j,bi,bj).GT.0) THEN
1047 numerator =
1048 & (uKdqdyInt(i,j)-uInt(i,j)*KdqdyInt(i,j)/R_low(i,j,bi,bj))
1049 & -(vKdqdxInt(i,j)-vInt(i,j)*KdqdxInt(i,j)/R_low(i,j,bi,bj))
1050 denominator = uXiyInt(i,j) - vXixInt(i,j)
1051 C We can have troubles with floating point exceptions if the denominator
1052 C of the renormalisation if the ocean is resting (e.g. intial conditions).
1053 C So we make the renormalisation factor one if the denominator is very small
1054 C The renormalisation factor is supposed to correct the error in the extraction of
1055 C potential energy associated with the truncation of the expansion. Thus, we
1056 C enforce a minimum value for the renormalisation factor.
1057 C We also enforce a maximum renormalisation factor.
1058 IF (denominator.GT.small) THEN
1059 Renorm(i,j) = ABS(numerator/denominator)
1060 Renorm(i,j) = MAX(Renorm(i,j),GM_K3D_minRenorm)
1061 Renorm(i,j) = MIN(Renorm(i,j),GM_K3D_maxRenorm)
1062 ENDIF
1063 ENDIF
1064 ENDDO
1065 ENDDO
1066 C Now put it back on to the velocity grids
1067 DO j=1-Oly+1,sNy+Oly-1
1068 DO i=1-Olx+1,sNx+Olx-1
1069 RenormU(i,j) = op5*(Renorm(i-1,j)+Renorm(i,j))
1070 RenormV(i,j) = op5*(Renorm(i,j-1)+Renorm(i,j))
1071 ENDDO
1072 ENDDO
1073
1074 C Calculate the eddy induced velocity in the X direction at the west face
1075 DO k=1,Nr
1076 DO j=1-Oly+1,sNy+Oly
1077 DO i=1-Olx+1,sNx+Olx
1078 ustar(i,j,k) = -RenormU(i,j)*Xix(i,j,k)/coriU(i,j)
1079 ENDDO
1080 ENDDO
1081 ENDDO
1082
1083 C Calculate the eddy induced velocity in the Y direction at the south face
1084 DO k=1,Nr
1085 DO j=1-Oly+1,sNy+Oly
1086 DO i=1-Olx+1,sNx+Olx
1087 vstar(i,j,k) = -RenormV(i,j)*Xiy(i,j,k)/coriV(i,j)
1088 ENDDO
1089 ENDDO
1090 ENDDO
1091
1092 C ======================================
1093 C Calculate the eddy induced overturning streamfunction
1094 C ======================================
1095 #ifdef GM_K3D_PASSIVE
1096 k=Nr
1097 DO j=1-Oly,sNy+Oly
1098 DO i=1-Olx,sNx+Olx
1099 psistar(i,j,Nr) = -hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
1100 ENDDO
1101 ENDDO
1102 DO k=Nr-1,1,-1
1103 DO j=1-Oly,sNy+Oly
1104 DO i=1-Olx,sNx+Olx
1105 psistar(i,j,k) = psistar(i,j,k+1)
1106 & - hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
1107 ENDDO
1108 ENDDO
1109 ENDDO
1110
1111 #else
1112
1113 IF (GM_AdvForm) THEN
1114 k=Nr
1115 DO j=1-Oly+1,sNy+1
1116 DO i=1-Olx+1,sNx+1
1117 GM_PsiX(i,j,k,bi,bj) = -hfacW(i,j,k,bi,bj)*drF(k)*ustar(i,j,k)
1118 GM_PsiY(i,j,k,bi,bj) = -hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
1119 ENDDO
1120 ENDDO
1121 DO k=Nr-1,1,-1
1122 DO j=1-Oly+1,sNy+1
1123 DO i=1-Olx+1,sNx+1
1124 GM_PsiX(i,j,k,bi,bj) = GM_PsiX(i,j,k+1,bi,bj)
1125 & - hfacW(i,j,k,bi,bj)*drF(k)*ustar(i,j,k)
1126 GM_PsiY(i,j,k,bi,bj) = GM_PsiY(i,j,k+1,bi,bj)
1127 & - hfacS(i,j,k,bi,bj)*drF(k)*vstar(i,j,k)
1128 ENDDO
1129 ENDDO
1130 ENDDO
1131
1132 ENDIF
1133 #endif
1134
1135 #ifdef ALLOW_DIAGNOSTICS
1136 C Diagnostics
1137 IF ( useDiagnostics ) THEN
1138 CALL DIAGNOSTICS_FILL(K3D, 'GM_K3D ',0,Nr,1,bi,bj,myThid)
1139 CALL DIAGNOSTICS_FILL(KPV, 'GM_KPV ',0,Nr,2,bi,bj,myThid)
1140 CALL DIAGNOSTICS_FILL(urms, 'GM_URMS ',0,Nr,2,bi,bj,myThid)
1141 CALL DIAGNOSTICS_FILL(Rdef, 'GM_RDEF ',0, 1,1,bi,bj,myThid)
1142 CALL DIAGNOSTICS_FILL(Rurms, 'GM_RURMS',0, 1,2,bi,bj,myThid)
1143 CALL DIAGNOSTICS_FILL(RRhines,'GM_RRHNS',0, 1,2,bi,bj,myThid)
1144 CALL DIAGNOSTICS_FILL(Rmix, 'GM_RMIX ',0, 1,2,bi,bj,myThid)
1145 CALL DIAGNOSTICS_FILL(supp, 'GM_SUPP ',0,Nr,2,bi,bj,myThid)
1146 CALL DIAGNOSTICS_FILL(Xix, 'GM_Xix ',0,Nr,2,bi,bj,myThid)
1147 CALL DIAGNOSTICS_FILL(Xiy, 'GM_Xiy ',0,Nr,2,bi,bj,myThid)
1148 CALL DIAGNOSTICS_FILL(cDopp, 'GM_C ',0, 1,2,bi,bj,myThid)
1149 CALL DIAGNOSTICS_FILL(Ubaro, 'GM_UBARO',0, 1,2,bi,bj,myThid)
1150 CALL DIAGNOSTICS_FILL(eady, 'GM_EADY ',0, 1,2,bi,bj,myThid)
1151 CALL DIAGNOSTICS_FILL(SlopeX, 'GM_Sx ',0,Nr,2,bi,bj,myThid)
1152 CALL DIAGNOSTICS_FILL(SlopeY, 'GM_Sy ',0,Nr,2,bi,bj,myThid)
1153 CALL DIAGNOSTICS_FILL(tfluxX, 'GM_TFLXX',0,Nr,2,bi,bj,myThid)
1154 CALL DIAGNOSTICS_FILL(tfluxY, 'GM_TFLXY',0,Nr,2,bi,bj,myThid)
1155 CALL DIAGNOSTICS_FILL(gradqx, 'GM_dqdx ',0,Nr,2,bi,bj,myThid)
1156 CALL DIAGNOSTICS_FILL(gradqy, 'GM_dqdy ',0,Nr,2,bi,bj,myThid)
1157 CALL DIAGNOSTICS_FILL(Kdqdy, 'GM_Kdqdy',0,Nr,2,bi,bj,myThid)
1158 CALL DIAGNOSTICS_FILL(Kdqdx, 'GM_Kdqdx',0,Nr,2,bi,bj,myThid)
1159 CALL DIAGNOSTICS_FILL(surfkz, 'GM_SFLYR',0, 1,2,bi,bj,myThid)
1160 CALL DIAGNOSTICS_FILL(ustar, 'GM_USTAR',0,Nr,2,bi,bj,myThid)
1161 CALL DIAGNOSTICS_FILL(vstar, 'GM_VSTAR',0,Nr,2,bi,bj,myThid)
1162 CALL DIAGNOSTICS_FILL(umc, 'GM_UMC ',0,Nr,2,bi,bj,myThid)
1163 CALL DIAGNOSTICS_FILL(ubar, 'GM_UBAR ',0,Nr,2,bi,bj,myThid)
1164 CALL DIAGNOSTICS_FILL(modesC, 'GM_MODEC',0,Nr,1,bi,bj,myThid)
1165 CALL DIAGNOSTICS_FILL(M4loc, 'GM_M4 ',0,Nr,2,bi,bj,myThid)
1166 CALL DIAGNOSTICS_FILL(N2loc, 'GM_N2 ',0,Nr,2,bi,bj,myThid)
1167 CALL DIAGNOSTICS_FILL(M4onN2, 'GM_M4_N2',0,Nr,2,bi,bj,myThid)
1168 CALL DIAGNOSTICS_FILL(slopeC, 'GM_SLOPE',0,Nr,2,bi,bj,myThid)
1169 CALL DIAGNOSTICS_FILL(Renorm, 'GM_RENRM',0, 1,2,bi,bj,myThid)
1170 #ifdef GM_K3D_PASSIVE
1171 CALL DIAGNOSTICS_FILL(psistar,'GM_PSTAR',0,Nr,2,bi,bj,myThid)
1172 #endif
1173 ENDIF
1174 #endif
1175
1176 C For the Redi diffusivity, we set K3D to a constant if
1177 C GM_K3D_constRedi=.TRUE.
1178 IF (GM_K3D_constRedi) THEN
1179 DO k=1,Nr
1180 DO j=1-Oly,sNy+Oly
1181 DO i=1-Olx,sNx+Olx
1182 K3D(i,j,k,bi,bj) = GM_K3D_constK
1183 ENDDO
1184 ENDDO
1185 ENDDO
1186 ENDIF
1187
1188 #ifdef ALLOW_DIAGNOSTICS
1189 IF ( useDiagnostics )
1190 & CALL DIAGNOSTICS_FILL(K3D, 'GM_K3D_T',0,Nr,1,bi,bj,myThid)
1191 #endif
1192
1193 #endif /* GM_K3D */
1194 RETURN
1195 END

  ViewVC Help
Powered by ViewVC 1.1.22