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

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

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


Revision 1.7 - (show annotations) (download)
Mon Nov 24 10:23:50 2014 UTC (9 years, 6 months ago) by dfer
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, HEAD
Changes since 1.6: +1 -3 lines
Remove  unnecessary for loop.

1 C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_calc_eigs.F,v 1.6 2014/05/20 11:42:28 m_bates Exp $
2 C $Name: $
3
4 #include "GMREDI_OPTIONS.h"
5
6 C !ROUTINE: EIGENVAL
7 C !INTERFACE:
8 SUBROUTINE GMREDI_CALC_EIGS(
9 I iMin, iMax, jMin, jMax,
10 I bi, bj, N2, myThid, kLow,
11 I mask, hfac, recip_hfac,
12 I rlow, nmodes, writediag,
13 O Rmid, vec)
14
15 C !DESCRIPTION: \bv
16 C *==========================================================*
17 C | SUBROUTINE GMREDI_CALC_URMS
18 C | o Calculate the vertical structure of the rms eddy
19 C | velocity based on baroclinic modal decomposition
20 C *==========================================================*
21 C \ev
22
23 IMPLICIT NONE
24
25 C == Global variables ==
26 #include "SIZE.h"
27 #include "GRID.h"
28 #include "EEPARAMS.h"
29 #include "PARAMS.h"
30 #include "GMREDI.h"
31
32 C !INPUT/OUTPUT PARAMETERS:
33 C == Routine arguments ==
34 C bi, bj :: tile indices
35 LOGICAL writediag
36 INTEGER iMin,iMax,jMin,jMax
37 INTEGER bi, bj
38 INTEGER myThid
39 INTEGER nmodes
40 INTEGER kLow(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
41 _RL mask(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
42 _RL hfac(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
43 _RL recip_hfac(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
44 _RL rlow(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
45 _RL N2( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
46 _RL Rmid(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
47 _RL vec(nmodes,1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
48
49 #ifdef GM_K3D
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 INTEGER i,j,k,kk,m
53 # ifdef HAVE_LAPACK
54 C info :: error code from LAPACK
55 C idx :: index used for sorting the eigenvalues
56 C a3d :: lower diagonal of eigenvalue problem
57 C b3d :: diagonal of eigenvalue problem
58 C c3d :: upper diagonal of eigenvalue problem
59 C val :: Eigenvalue (wavenumber) of the first baroclinic mode
60 C eigR :: Real component of all the eigenvalues in a water column
61 C eigI :: Imaginary component of all the eigenvalues in a water column
62 C vecs :: All the eigenvectors of a water column
63 C dummy :: Redundant variable for calling lapack
64 C work :: Work array for lapack
65 C array :: Array containing the matrix with a,b,c
66 C eigval :: Vector containing the eigenvalues
67 INTEGER info
68 INTEGER idx
69 _RL a3d( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
70 _RL b3d( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
71 _RL c3d( 1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
72 _RL val( 1-Olx:sNx+Olx,1-Oly:sNy+Oly)
73 _RL eigR(Nr),eigI(Nr),vecs(Nr,Nr),dummy(1,Nr),work(Nr*Nr)
74 _RL array(Nr,Nr)
75 _RL eigval(0:nmodes)
76 # else
77 C drNr :: distance from bottom of cell to cell centre at Nr
78 C BuoyFreq :: buoyancy frequency, SQRT(N2)
79 C intN :: Vertical integral of BuoyFreq to each grid cell centre
80 C intN0 :: Vertical integral of BuoyFreq to z=0
81 C c1 :: intN0/pi
82 C nEigs :: number of eigenvalues/vectors to calculate
83 _RL drNr
84 _RL BuoyFreq(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr+1)
85 _RL intN(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
86 _RL intN0(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
87 _RL c1(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
88 INTEGER nEigs(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
89 # endif
90 C small :: a small number (used to avoid floating point exceptions)
91 C vecint :: vertical integral of eigenvector and/or square of eigenvector
92 C fCori2 :: square of the Coriolis parameter
93 _RL small
94 _RL vecint(nmodes,1-Olx:sNx+Olx,1-Oly:sNy+Oly)
95 _RL fCori2(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
96 CHARACTER*(MAX_LEN_MBUF) msgBuf
97
98 small = TINY(zeroRL)
99
100 C Square of the Coriolis parameter
101 DO i=1-OLx,sNx+OLx
102 DO j=1-OLy,sNy+OLy
103 fCori2(i,j) = fCori(i,j,bi,bj)*fCori(i,j,bi,bj)
104 ENDDO
105 ENDDO
106
107 DO k=1,Nr
108 DO j=1-Oly,sNy+Oly
109 DO i=1-Olx,sNx+Olx
110 DO m=1,nmodes
111 vec(m,i,j,k) = zeroRL
112 ENDDO
113 ENDDO
114 ENDDO
115 ENDDO
116
117 # ifdef HAVE_LAPACK
118 C Calculate the tridiagonal operator matrix for
119 C f^2 d/dz 1/N^2 d/dz
120 C a3d is the lower off-diagonal, b3d is the diagonal
121 C and c3d is the upper off-diagonal
122 DO k=1,Nr
123 DO j=1-OLy,sNy+OLy
124 DO i=1-OLx,sNx+OLx
125 IF (kLow(i,j) .GT. 0) THEN
126 IF (k.EQ.1) THEN
127 a3d(i,j,k) = zeroRL
128 c3d(i,j,k) = fCori2(i,j)*recip_hFac(i,j,k)
129 & *recip_drC(k+1)*recip_drF(k)/N2(i,j,k+1)
130 b3d(i,j,k) = -c3d(i,j,k)
131
132 ELSEIF (k.LT.kLow(i,j)) THEN
133 a3d(i,j,k) = fCori2(i,j)*recip_hFac(i,j,k)
134 & *recip_drF(k)*recip_drC(k)/N2(i,j,k)
135 c3d(i,j,k) = fCori2(i,j)*recip_hFac(i,j,k)
136 & *recip_drF(k)*recip_drC(k+1)/N2(i,j,k+1)
137 b3d(i,j,k) = -a3d(i,j,k)-c3d(i,j,k)
138
139 ELSEIF (k.EQ.kLow(i,j)) THEN
140 a3d(i,j,k) = fCori2(i,j)*recip_hFac(i,j,k)
141 & *recip_drF(k)*recip_drC(k)/N2(i,j,k)
142 c3d(i,j,k) = zeroRL
143 b3d(i,j,k) = -a3d(i,j,k)
144
145 ELSE
146 a3d(i,j,k) = zeroRL
147 b3d(i,j,k) = zeroRL
148 c3d(i,j,k) = zeroRL
149 ENDIF
150
151 ELSE
152 a3d(i,j,k) = zeroRL
153 b3d(i,j,k) = zeroRL
154 c3d(i,j,k) = zeroRL
155 ENDIF
156 ENDDO
157 ENDDO
158 ENDDO
159
160 DO j=1-Oly,sNy+Oly
161 DO i=1-Olx,sNx+Olx
162 IF (kLow(i,j).GT.0) THEN
163 DO kk=1,Nr
164 DO k=1,Nr
165 array(k,kk) = zeroRL
166 ENDDO
167 ENDDO
168
169 k=1
170 array(k,k) = b3d(i,j,k)
171 array(k,k+1) = c3d(i,j,k)
172 DO k=2,Nr-1
173 array(k,k-1) = a3d(i,j,k)
174 array(k,k) = b3d(i,j,k)
175 array(k,k+1) = c3d(i,j,k)
176 ENDDO
177 k=Nr
178 array(k,k-1) = a3d(i,j,k)
179 array(k,k) = b3d(i,j,k)
180
181 CALL DGEEV('N','V',Nr,array,Nr,eigR,eigI,dummy,1,
182 & vecs,Nr,work,Nr*Nr,info)
183 IF( info.LT.0 ) THEN
184 WRITE(msgBuf,'(A,x,2(A1,I2),A1,x,A,I4)')
185 & 'GMREDI_CALC_EIGS problem with arguments for DGEEV at',
186 & '(',i,',',j,')', 'error code =',info
187 CALL PRINT_ERROR( msgBuf , myThid )
188
189 ELSEIF(info.GT.0 ) THEN
190 WRITE(msgBuf,'(A,x,2(A1,I2),A1,x,A,I4)')
191 & 'GMREDI_CALC_EIGS problems with eigensolver DGEEV at',
192 & '(',i,',',j,')', 'error code =',info
193 CALL PRINT_ERROR( msgBuf , myThid )
194
195 ENDIF
196
197 C Find the second largest eigenvalue (the Rossby radius)
198 C and the first M baroclinic modes (eigenvectors)
199 DO m=0,nmodes
200 eigval(m) = -HUGE(zeroRL)
201 ENDDO
202
203 DO k=1,kLow(i,j)
204 eigval(0) = MAX(eigval(0),eigR(k))
205 ENDDO
206 DO m=1,MIN(nmodes,klow(i,j)-1)
207 DO k=1,kLow(i,j)
208 IF (eigR(k).LT.eigval(m-1)) THEN
209 eigval(m) = MAX(eigval(m),eigR(k))
210 IF (eigval(m).EQ.eigR(k)) idx=k
211 ENDIF
212 ENDDO
213 IF(vecs(1,idx).LT.zeroRL) THEN
214 DO k=1,Nr
215 vec(m,i,j,k) = -vecs(k,idx)
216 ENDDO
217 ELSE
218 DO k=1,Nr
219 vec(m,i,j,k) = vecs(k,idx)
220 ENDDO
221 ENDIF
222 ENDDO
223 val(i,j) = eigval(1)
224 ELSE
225 val(i,j)=zeroRL
226 DO k=1,Nr
227 DO m=1,nmodes
228 vec(m,i,j,k)=zeroRL
229 ENDDO
230 ENDDO
231
232 ENDIF
233 ENDDO
234 ENDDO
235
236 DO j=1-Oly,sNy+Oly
237 DO i=1-Olx,sNx+Olx
238 IF (kLow(i,j).GT.2 .AND. val(i,j).NE.zeroRL) THEN
239 Rmid(i,j) = 1.0/(SQRT(ABS(val(i,j)))+small)
240 ELSE
241 Rmid(i,j) = zeroRL
242 ENDIF
243 ENDDO
244 ENDDO
245
246 # else
247 DO k=1,Nr
248 DO j=1-Oly,sNy+Oly
249 DO i=1-Olx,sNx+Olx
250 BuoyFreq(i,j,k) = mask(i,j,k)*SQRT(N2(i,j,k))
251 ENDDO
252 ENDDO
253 ENDDO
254 k=Nr+1
255 DO j=1-Oly,sNy+Oly
256 DO i=1-Olx,sNx+Olx
257 BuoyFreq(i,j,k) = zeroRL
258 ENDDO
259 ENDDO
260
261 C integrate N using something like Simpson s rule (but not quite)
262 C drC*( (N(k+1)+N(k+2))/2 + (N(k)+N(k+1))/2 )/2
263 C when k=Nr, say that N(k+2)=0 and N(k+1)=0
264 k=Nr
265 C drNr is like drC(Nr+1)/2
266 drNr = rC(Nr)-rF(Nr+1)
267 DO j=1-Oly,sNy+Oly
268 DO i=1-Olx,sNx+Olx
269 intN(i,j,k) = op5*BuoyFreq(i,j,k)*drNr
270 ENDDO
271 ENDDO
272 DO k=Nr-1,1,-1
273 DO j=1-Oly,sNy+Oly
274 DO i=1-Olx,sNx+Olx
275 intN(i,j,k) = intN(i,j,k+1)
276 & + drC(k)*( op25*BuoyFreq(i,j,k+2) + op5*BuoyFreq(i,j,k)
277 & + op25*BuoyFreq(i,j,k+1) )
278 ENDDO
279 ENDDO
280 ENDDO
281
282 C intN integrates to z=rC(1). We want to integrate to z=0.
283 C Assume that N(0)=0 and N(1)=0.
284 C drC(1) is like half a grid cell.
285 DO j=1-Oly,sNy+Oly
286 DO i=1-Olx,sNx+Olx
287 intN0(i,j) = intN(i,j,1)
288 & + drC(1)*op5*BuoyFreq(i,j,2)
289 ENDDO
290 ENDDO
291
292 DO j=1-Oly,sNy+Oly
293 DO i=1-Olx,sNx+Olx
294 c1(i,j) = intN0(i,j)/pi
295 Rmid(i,j) = c1(i,j)/ABS(fCori(i,j,bi,bj))
296 ENDDO
297 ENDDO
298
299 DO j=1-Oly,sNy+Oly
300 DO i=1-Olx,sNx+Olx
301 nEigs(i,j) = MIN(klow(i,j),nmodes)
302 ENDDO
303 ENDDO
304 DO k=1,Nr
305 DO j=1-Oly,sNy+Oly
306 DO i=1-Olx,sNx+Olx
307 IF (mask(i,j,k).NE.0.0) THEN
308 DO m=1,nEigs(i,j)
309 vec(m,i,j,k) = -COS(intN(i,j,k)/(m*c1(i,j)))
310 ENDDO
311 ENDIF
312 ENDDO
313 ENDDO
314 ENDDO
315
316 C The WKB approximation for the baroclinic mode does not always
317 C integrate to zero so we adjust it.
318 DO j=1-Oly,sNy+Oly
319 DO i=1-Olx,sNx+Olx
320 DO m=1,nEigs(i,j)
321 vecint(m,i,j) = zeroRL
322 ENDDO
323 ENDDO
324 ENDDO
325 DO k=1,Nr
326 DO j=1-Oly,sNy+Oly
327 DO i=1-Olx,sNx+Olx
328 DO m=1,nEigs(i,j)
329 vecint(m,i,j) = vecint(m,i,j) + hfac(i,j,k)*vec(m,i,j,k)*drF(k)
330 ENDDO
331 ENDDO
332 ENDDO
333 ENDDO
334 DO j=1-Oly,sNy+Oly
335 DO i=1-Olx,sNx+Olx
336 DO m=1,nEigs(i,j)
337 vecint(m,i,j) = vecint(m,i,j)/(-rlow(i,j)+small)
338 ENDDO
339 ENDDO
340 ENDDO
341 DO k=1,Nr
342 DO j=1-Oly,sNy+Oly
343 DO i=1-Olx,sNx+Olx
344 DO m=1,nEigs(i,j)
345 vec(m,i,j,k) = vec(m,i,j,k) - vecint(m,i,j)
346 ENDDO
347 ENDDO
348 ENDDO
349 ENDDO
350 # endif
351
352 C Normalise the eigenvectors
353 DO j=1-Oly,sNy+Oly
354 DO i=1-Olx,sNx+Olx
355 DO m=1,nmodes
356 vecint(m,i,j) = zeroRL
357 ENDDO
358 ENDDO
359 ENDDO
360
361 DO k=1,Nr
362 DO j=1-Oly,sNy+Oly
363 DO i=1-Olx,sNx+Olx
364 DO m=1,nmodes
365 vecint(m,i,j) = vecint(m,i,j) +
366 & mask(i,j,k)*drF(k)*hfac(i,j,k)
367 & *vec(m,i,j,k)*vec(m,i,j,k)
368 ENDDO
369 ENDDO
370 ENDDO
371 ENDDO
372
373 DO j=1-Oly,sNy+Oly
374 DO i=1-Olx,sNx+Olx
375 DO m=1,nmodes
376 vecint(m,i,j) = SQRT(vecint(m,i,j))
377 ENDDO
378 ENDDO
379 ENDDO
380
381 DO k=1,Nr
382 DO j=1-Oly,sNy+Oly
383 DO i=1-Olx,sNx+Olx
384 DO m=1,nmodes
385 vec(m,i,j,k) = vec(m,i,j,k)/(vecint(m,i,j)+small)
386 ENDDO
387 ENDDO
388 ENDDO
389 ENDDO
390
391 # ifdef ALLOW_DIAGNOSTICS
392 C Diagnostics
393 IF ( useDiagnostics.AND.writediag ) THEN
394 # ifdef HAVE_LAPACK
395 CALL DIAGNOSTICS_FILL(a3d, 'GM_A3D ',0,Nr,0,1,1,myThid)
396 CALL DIAGNOSTICS_FILL(b3d, 'GM_B3D ',0,Nr,0,1,1,myThid)
397 CALL DIAGNOSTICS_FILL(c3d, 'GM_C3D ',0,Nr,0,1,1,myThid)
398 # endif
399 ENDIF
400 # endif
401
402 #endif /* GM_K3D */
403
404 RETURN
405 END

  ViewVC Help
Powered by ViewVC 1.1.22