1 |
gforget |
1.1 |
C $Header: /u/gcmpack/MITgcm/model/src/rotate_uv2en.F,v 1.1 2010/03/04 03:42:30 gforget Exp $ |
2 |
|
|
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "CPP_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
C-- File rotate_uv2en.F: Routines to handle a vector coordinate system rotation. |
7 |
|
|
C-- Contents |
8 |
|
|
C-- o ROTATE_UV2EN_RL |
9 |
|
|
C-- o ROTATE_UV2EN_RS |
10 |
|
|
|
11 |
|
|
subroutine rotate_uv2en_rl( |
12 |
|
|
U uFldX, vFldY, |
13 |
|
|
U uFldE, vFldN, |
14 |
|
|
I xy2en, switchGrid, applyMask, kSize, mythid |
15 |
|
|
& ) |
16 |
|
|
|
17 |
|
|
c ================================================================== |
18 |
|
|
c SUBROUTINE rotate_uv2en_rl |
19 |
|
|
c ================================================================== |
20 |
|
|
c |
21 |
|
|
c o uFldX/vFldY are in the model grid directions. |
22 |
|
|
c o uFldE/vFldN are eastward/northward. |
23 |
|
|
c o This routine goes from uFldX/vFldY to uFldE/vFldN (for xy2en=.TRUE.) |
24 |
|
|
c or vice versa (for xy2en=.FALSE.). |
25 |
|
|
c o uFldX/vFldY may be at the C grid velocity points, or at the A grid |
26 |
|
|
c velocity points (i.e. the C grid cell center). uFldE/vFldN are always |
27 |
|
|
c at the cell center. If switchGrid=.TRUE. we go from C grid uFldX/vFldY |
28 |
|
|
c to A grid uFldE/vFldN, or vice versa. If switchGrid=.FALSE. we go |
29 |
|
|
c from A grid uFldX/vFldY to A grid uFldE/vFldN, or vice versa. |
30 |
|
|
c o If applyMask=.TRUE. then masks are applied to the output. |
31 |
|
|
c If kSize=1 (resp. nr) we then use the surface (resp. 3D) masks. |
32 |
|
|
c o In any case it is assumed that exchanges are done on the outside. |
33 |
|
|
c |
34 |
|
|
c ================================================================== |
35 |
|
|
c SUBROUTINE rotate_uv2en_rl |
36 |
|
|
c ================================================================== |
37 |
|
|
|
38 |
|
|
implicit none |
39 |
|
|
|
40 |
|
|
c == global variables == |
41 |
|
|
|
42 |
|
|
#include "EEPARAMS.h" |
43 |
|
|
#include "SIZE.h" |
44 |
|
|
#include "PARAMS.h" |
45 |
|
|
#include "GRID.h" |
46 |
|
|
|
47 |
|
|
c == routine arguments == |
48 |
|
|
|
49 |
|
|
integer kSize |
50 |
|
|
logical xy2en, switchGrid, applyMask |
51 |
|
|
_RL uFldX(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
52 |
|
|
_RL vFldY(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
53 |
|
|
_RL uFldE(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
54 |
|
|
_RL vFldN(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
55 |
|
|
|
56 |
|
|
integer mythid |
57 |
|
|
|
58 |
|
|
c == local variables == |
59 |
|
|
|
60 |
|
|
integer bi,bj |
61 |
|
|
integer i,j,k,kk |
62 |
|
|
_RL tmpU(1-olx:snx+olx,1-oly:sny+oly) |
63 |
|
|
_RL tmpV(1-olx:snx+olx,1-oly:sny+oly) |
64 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
65 |
|
|
|
66 |
|
|
c == end of interface == |
67 |
|
|
|
68 |
|
|
if ( (kSize.NE.1).AND.(kSize.NE.nr) |
69 |
|
|
& .AND.(applyMask) ) then |
70 |
|
|
WRITE(msgBuf,'(2A,I4,A)') ' ROTATE_UV2EN: ', |
71 |
|
|
& 'no mask has ',kSize,' levels' |
72 |
|
|
CALL PRINT_ERROR(msgBuf, myThid) |
73 |
|
|
STOP 'ABNROMAL END: S/R ROTATE_UV2EN' |
74 |
|
|
endif |
75 |
|
|
|
76 |
|
|
do bj = mybylo(mythid),mybyhi(mythid) |
77 |
|
|
do bi = mybxlo(mythid),mybxhi(mythid) |
78 |
|
|
do k = 1,kSize |
79 |
|
|
|
80 |
|
|
if ( (kSize.EQ.1).AND.(usingPCoords) ) then |
81 |
|
|
kk=nr |
82 |
|
|
else |
83 |
|
|
kk=k |
84 |
|
|
endif |
85 |
|
|
|
86 |
|
|
if ( xy2en ) then |
87 |
|
|
c go from uFldX/vFldY to uFldE/vFldN |
88 |
|
|
if ( switchGrid ) then |
89 |
|
|
C 1a) go from C grid velocity points to A grid velocity points |
90 |
|
|
do i = 1-olx,snx+olx |
91 |
|
|
tmpU(i,sny+Oly)=0. |
92 |
|
|
tmpV(i,sny+Oly)=0. |
93 |
|
|
enddo |
94 |
|
|
do j = 1-oly,sny+oly-1 |
95 |
|
|
tmpU(snx+Olx,j)=0. |
96 |
|
|
tmpV(snx+Olx,j)=0. |
97 |
|
|
do i = 1-olx,snx+olx-1 |
98 |
|
|
tmpU(i,j) = 0.5 _d 0 |
99 |
|
|
& *( uFldX(i+1,j,k,bi,bj) + uFldX(i,j,k,bi,bj) ) |
100 |
|
|
tmpV(i,j) = 0.5 _d 0 |
101 |
|
|
& *( vFldY(i,j+1,k,bi,bj) + vFldY(i,j,k,bi,bj) ) |
102 |
|
|
if (applyMask) then |
103 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
104 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
105 |
|
|
endif |
106 |
|
|
enddo |
107 |
|
|
enddo |
108 |
|
|
else |
109 |
|
|
C 1b) stay at A grid velocity points (i.e. at the C grid cell center) |
110 |
|
|
do j = 1-oly,sny+oly |
111 |
|
|
do i = 1-olx,snx+olx |
112 |
|
|
tmpU(i,j) = uFldX(i,j,k,bi,bj) |
113 |
|
|
tmpV(i,j) = vFldY(i,j,k,bi,bj) |
114 |
|
|
if (applyMask) then |
115 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
116 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
117 |
|
|
endif |
118 |
|
|
enddo |
119 |
|
|
enddo |
120 |
|
|
endif!if ( switchGrid ) then |
121 |
|
|
|
122 |
|
|
C 2) rotation |
123 |
|
|
do j = 1-oly,sny+oly |
124 |
|
|
do i = 1-olx,snx+olx |
125 |
|
|
uFldE(i,j,k,bi,bj) = |
126 |
|
|
& angleCosC(i,j,bi,bj)*tmpU(i,j) |
127 |
|
|
& -angleSinC(i,j,bi,bj)*tmpV(i,j) |
128 |
|
|
vFldN(i,j,k,bi,bj) = |
129 |
|
|
& angleSinC(i,j,bi,bj)*tmpU(i,j) |
130 |
|
|
& +angleCosC(i,j,bi,bj)*tmpV(i,j) |
131 |
|
|
enddo |
132 |
|
|
enddo |
133 |
|
|
|
134 |
|
|
else!if (xy2en) then |
135 |
|
|
c go from uFldE/vFldN to uFldX/vFldY |
136 |
|
|
|
137 |
|
|
C 1) rotation |
138 |
|
|
do j = 1-oly,sny+oly |
139 |
|
|
do i = 1-olx,snx+olx |
140 |
|
|
tmpU(i,j) = |
141 |
|
|
& angleCosC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
142 |
|
|
& +angleSinC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
143 |
|
|
tmpV(i,j) = |
144 |
|
|
& -angleSinC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
145 |
|
|
& +angleCosC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
146 |
|
|
enddo |
147 |
|
|
enddo |
148 |
|
|
|
149 |
|
|
if ( switchGrid ) then |
150 |
|
|
C 2a) go from A grid velocity points to C grid velocity points |
151 |
|
|
do i = 1-olx,snx+olx |
152 |
|
|
uFldX(i,1,k,bi,bj)=0. |
153 |
|
|
vFldY(i,1,k,bi,bj)=0. |
154 |
|
|
enddo |
155 |
|
|
do j = 1-oly+1,sny+oly |
156 |
|
|
uFldX(1,j,k,bi,bj)=0. |
157 |
|
|
vFldY(1,j,k,bi,bj)=0. |
158 |
|
|
do i = 1-olx+1,snx+olx |
159 |
|
|
uFldX(i,j,k,bi,bj) = 0.5 _d 0 |
160 |
|
|
& *( tmpU(i-1,j) + tmpU(i,j) ) |
161 |
|
|
vFldY(i,j,k,bi,bj) = 0.5 _d 0 |
162 |
|
|
& *( tmpV(i,j-1) + tmpV(i,j) ) |
163 |
|
|
if (applyMask) then |
164 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskW(i,j,kk,bi,bj) |
165 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskS(i,j,kk,bi,bj) |
166 |
|
|
endif |
167 |
|
|
enddo |
168 |
|
|
enddo |
169 |
|
|
else |
170 |
|
|
C 2b) stay at A grid velocity points (i.e. at the C grid cell center) |
171 |
|
|
do j = 1-oly,sny+oly |
172 |
|
|
do i = 1-olx,snx+olx |
173 |
|
|
uFldX(i,j,k,bi,bj) = tmpU(i,j) |
174 |
|
|
vFldY(i,j,k,bi,bj) = tmpV(i,j) |
175 |
|
|
if (applyMask) then |
176 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
177 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
178 |
|
|
endif |
179 |
|
|
enddo |
180 |
|
|
enddo |
181 |
|
|
endif!if ( switchGrid ) then |
182 |
|
|
|
183 |
|
|
endif!if (xy2en) then |
184 |
|
|
|
185 |
|
|
enddo |
186 |
|
|
enddo |
187 |
|
|
enddo |
188 |
|
|
|
189 |
|
|
return |
190 |
|
|
end |
191 |
|
|
|
192 |
|
|
subroutine rotate_uv2en_rs( |
193 |
|
|
U uFldX, vFldY, |
194 |
|
|
U uFldE, vFldN, |
195 |
|
|
I xy2en, switchGrid, applyMask, kSize, mythid |
196 |
|
|
& ) |
197 |
|
|
|
198 |
|
|
c ================================================================== |
199 |
|
|
c SUBROUTINE rotate_uv2en_rs |
200 |
|
|
c ================================================================== |
201 |
|
|
c |
202 |
|
|
c o uFldX/vFldY are in the model grid directions. |
203 |
|
|
c o uFldE/vFldN are eastward/northward. |
204 |
|
|
c o This routine goes from uFldX/vFldY to uFldE/vFldN (for xy2en=.TRUE.) |
205 |
|
|
c or vice versa (for xy2en=.FALSE.). |
206 |
|
|
c o uFldX/vFldY may be at the C grid velocity points, or at the A grid |
207 |
|
|
c velocity points (i.e. the C grid cell center). uFldE/vFldN are always |
208 |
|
|
c at the cell center. If switchGrid=.TRUE. we go from C grid uFldX/vFldY |
209 |
|
|
c to A grid uFldE/vFldN, or vice versa. If switchGrid=.FALSE. we go |
210 |
|
|
c from A grid uFldX/vFldY to A grid uFldE/vFldN, or vice versa. |
211 |
|
|
c o If applyMask=.TRUE. then masks are applied to the output. |
212 |
|
|
c If kSize=1 (resp. nr) we then use the surface (resp. 3D) masks. |
213 |
|
|
c o In any case it is assumed that exchanges are done on the outside. |
214 |
|
|
c |
215 |
|
|
c ================================================================== |
216 |
|
|
c SUBROUTINE rotate_uv2en_rs |
217 |
|
|
c ================================================================== |
218 |
|
|
|
219 |
|
|
implicit none |
220 |
|
|
|
221 |
|
|
c == global variables == |
222 |
|
|
|
223 |
|
|
#include "EEPARAMS.h" |
224 |
|
|
#include "SIZE.h" |
225 |
|
|
#include "PARAMS.h" |
226 |
|
|
#include "GRID.h" |
227 |
|
|
|
228 |
|
|
c == routine arguments == |
229 |
|
|
|
230 |
|
|
integer kSize |
231 |
|
|
logical xy2en, switchGrid, applyMask |
232 |
|
|
_RS uFldX(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
233 |
|
|
_RS vFldY(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
234 |
|
|
_RS uFldE(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
235 |
|
|
_RS vFldN(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
236 |
|
|
|
237 |
|
|
integer mythid |
238 |
|
|
|
239 |
|
|
c == local variables == |
240 |
|
|
|
241 |
|
|
integer bi,bj |
242 |
|
|
integer i,j,k,kk |
243 |
|
|
_RS tmpU(1-olx:snx+olx,1-oly:sny+oly) |
244 |
|
|
_RS tmpV(1-olx:snx+olx,1-oly:sny+oly) |
245 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
246 |
|
|
|
247 |
|
|
c == end of interface == |
248 |
|
|
|
249 |
|
|
if ( (kSize.NE.1).AND.(kSize.NE.nr) |
250 |
|
|
& .AND.(applyMask) ) then |
251 |
|
|
WRITE(msgBuf,'(2A,I4,A)') ' ROTATE_UV2EN: ', |
252 |
|
|
& 'no mask has ',kSize,' levels' |
253 |
|
|
CALL PRINT_ERROR(msgBuf, myThid) |
254 |
|
|
STOP 'ABNROMAL END: S/R ROTATE_UV2EN' |
255 |
|
|
endif |
256 |
|
|
|
257 |
|
|
do bj = mybylo(mythid),mybyhi(mythid) |
258 |
|
|
do bi = mybxlo(mythid),mybxhi(mythid) |
259 |
|
|
do k = 1,kSize |
260 |
|
|
|
261 |
|
|
if ( (kSize.EQ.1).AND.(usingPCoords) ) then |
262 |
|
|
kk=nr |
263 |
|
|
else |
264 |
|
|
kk=k |
265 |
|
|
endif |
266 |
|
|
|
267 |
|
|
if ( xy2en ) then |
268 |
|
|
c go from uFldX/vFldY to uFldE/vFldN |
269 |
|
|
if ( switchGrid ) then |
270 |
|
|
C 1a) go from C grid velocity points to A grid velocity points |
271 |
|
|
do i = 1-olx,snx+olx |
272 |
|
|
tmpU(i,sny+Oly)=0. |
273 |
|
|
tmpV(i,sny+Oly)=0. |
274 |
|
|
enddo |
275 |
|
|
do j = 1-oly,sny+oly-1 |
276 |
|
|
tmpU(snx+Olx,j)=0. |
277 |
|
|
tmpV(snx+Olx,j)=0. |
278 |
|
|
do i = 1-olx,snx+olx-1 |
279 |
|
|
tmpU(i,j) = 0.5 _d 0 |
280 |
|
|
& *( uFldX(i+1,j,k,bi,bj) + uFldX(i,j,k,bi,bj) ) |
281 |
|
|
tmpV(i,j) = 0.5 _d 0 |
282 |
|
|
& *( vFldY(i,j+1,k,bi,bj) + vFldY(i,j,k,bi,bj) ) |
283 |
|
|
if (applyMask) then |
284 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
285 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
286 |
|
|
endif |
287 |
|
|
enddo |
288 |
|
|
enddo |
289 |
|
|
else |
290 |
|
|
C 1b) stay at A grid velocity points (i.e. at the C grid cell center) |
291 |
|
|
do j = 1-oly,sny+oly |
292 |
|
|
do i = 1-olx,snx+olx |
293 |
|
|
tmpU(i,j) = uFldX(i,j,k,bi,bj) |
294 |
|
|
tmpV(i,j) = vFldY(i,j,k,bi,bj) |
295 |
|
|
if (applyMask) then |
296 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
297 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
298 |
|
|
endif |
299 |
|
|
enddo |
300 |
|
|
enddo |
301 |
|
|
endif!if ( switchGrid ) then |
302 |
|
|
|
303 |
|
|
C 2) rotation |
304 |
|
|
do j = 1-oly,sny+oly |
305 |
|
|
do i = 1-olx,snx+olx |
306 |
|
|
uFldE(i,j,k,bi,bj) = |
307 |
|
|
& angleCosC(i,j,bi,bj)*tmpU(i,j) |
308 |
|
|
& -angleSinC(i,j,bi,bj)*tmpV(i,j) |
309 |
|
|
vFldN(i,j,k,bi,bj) = |
310 |
|
|
& angleSinC(i,j,bi,bj)*tmpU(i,j) |
311 |
|
|
& +angleCosC(i,j,bi,bj)*tmpV(i,j) |
312 |
|
|
enddo |
313 |
|
|
enddo |
314 |
|
|
|
315 |
|
|
else!if (xy2en) then |
316 |
|
|
c go from uFldE/vFldN to uFldX/vFldY |
317 |
|
|
|
318 |
|
|
C 1) rotation |
319 |
|
|
do j = 1-oly,sny+oly |
320 |
|
|
do i = 1-olx,snx+olx |
321 |
|
|
tmpU(i,j) = |
322 |
|
|
& angleCosC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
323 |
|
|
& +angleSinC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
324 |
|
|
tmpV(i,j) = |
325 |
|
|
& -angleSinC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
326 |
|
|
& +angleCosC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
327 |
|
|
enddo |
328 |
|
|
enddo |
329 |
|
|
|
330 |
|
|
if ( switchGrid ) then |
331 |
|
|
C 2a) go from A grid velocity points to C grid velocity points |
332 |
|
|
do i = 1-olx,snx+olx |
333 |
|
|
uFldX(i,1,k,bi,bj)=0. |
334 |
|
|
vFldY(i,1,k,bi,bj)=0. |
335 |
|
|
enddo |
336 |
|
|
do j = 1-oly+1,sny+oly |
337 |
|
|
uFldX(1,j,k,bi,bj)=0. |
338 |
|
|
vFldY(1,j,k,bi,bj)=0. |
339 |
|
|
do i = 1-olx+1,snx+olx |
340 |
|
|
uFldX(i,j,k,bi,bj) = 0.5 _d 0 |
341 |
|
|
& *( tmpU(i-1,j) + tmpU(i,j) ) |
342 |
|
|
vFldY(i,j,k,bi,bj) = 0.5 _d 0 |
343 |
|
|
& *( tmpV(i,j-1) + tmpV(i,j) ) |
344 |
|
|
if (applyMask) then |
345 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskW(i,j,kk,bi,bj) |
346 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskS(i,j,kk,bi,bj) |
347 |
|
|
endif |
348 |
|
|
enddo |
349 |
|
|
enddo |
350 |
|
|
else |
351 |
|
|
C 2b) stay at A grid velocity points (i.e. at the C grid cell center) |
352 |
|
|
do j = 1-oly,sny+oly |
353 |
|
|
do i = 1-olx,snx+olx |
354 |
|
|
uFldX(i,j,k,bi,bj) = tmpU(i,j) |
355 |
|
|
vFldY(i,j,k,bi,bj) = tmpV(i,j) |
356 |
|
|
if (applyMask) then |
357 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
358 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
359 |
|
|
endif |
360 |
|
|
enddo |
361 |
|
|
enddo |
362 |
|
|
endif!if ( switchGrid ) then |
363 |
|
|
|
364 |
|
|
endif!if (xy2en) then |
365 |
|
|
|
366 |
|
|
enddo |
367 |
|
|
enddo |
368 |
|
|
enddo |
369 |
|
|
|
370 |
|
|
return |
371 |
|
|
end |
372 |
|
|
|