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