1 |
ifenty |
1.1 |
C $Header: /u/gcmpack/MITgcm_contrib/gael/verification/global_oce_llc90/code/rotate_uv2en.F,v 1.1 2013/04/09 17:23:18 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 |
|
|
do j = 1-oly,sny+oly |
89 |
|
|
do i = 1-olx,snx+olx |
90 |
|
|
uFldE(i,j,k,bi,bj) = 0. _d 0 |
91 |
|
|
vFldN(i,j,k,bi,bj) = 0. _d 0 |
92 |
|
|
tmpU(i,j) = 0. _d 0 |
93 |
|
|
tmpV(i,j) = 0. _d 0 |
94 |
|
|
enddo |
95 |
|
|
enddo |
96 |
|
|
if ( switchGrid ) then |
97 |
|
|
C 1a) go from C grid velocity points to A grid velocity points |
98 |
|
|
do j = 1,sny |
99 |
|
|
do i = 1,snx |
100 |
|
|
tmpU(i,j) = 0.5 _d 0 |
101 |
|
|
& *( uFldX(i+1,j,k,bi,bj) + uFldX(i,j,k,bi,bj) ) |
102 |
|
|
tmpV(i,j) = 0.5 _d 0 |
103 |
|
|
& *( vFldY(i,j+1,k,bi,bj) + vFldY(i,j,k,bi,bj) ) |
104 |
|
|
if (applyMask) then |
105 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
106 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
107 |
|
|
endif |
108 |
|
|
enddo |
109 |
|
|
enddo |
110 |
|
|
else |
111 |
|
|
C 1b) stay at A grid velocity points (i.e. at the C grid cell center) |
112 |
|
|
do j = 1,sny |
113 |
|
|
do i = 1,snx |
114 |
|
|
tmpU(i,j) = uFldX(i,j,k,bi,bj) |
115 |
|
|
tmpV(i,j) = vFldY(i,j,k,bi,bj) |
116 |
|
|
if (applyMask) then |
117 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
118 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
119 |
|
|
endif |
120 |
|
|
enddo |
121 |
|
|
enddo |
122 |
|
|
endif!if ( switchGrid ) then |
123 |
|
|
|
124 |
|
|
C 2) rotation |
125 |
|
|
do j = 1,sny |
126 |
|
|
do i = 1,snx |
127 |
|
|
uFldE(i,j,k,bi,bj) = |
128 |
|
|
& angleCosC(i,j,bi,bj)*tmpU(i,j) |
129 |
|
|
& -angleSinC(i,j,bi,bj)*tmpV(i,j) |
130 |
|
|
vFldN(i,j,k,bi,bj) = |
131 |
|
|
& angleSinC(i,j,bi,bj)*tmpU(i,j) |
132 |
|
|
& +angleCosC(i,j,bi,bj)*tmpV(i,j) |
133 |
|
|
enddo |
134 |
|
|
enddo |
135 |
|
|
|
136 |
|
|
else!if (xy2en) then |
137 |
|
|
c go from uFldE/vFldN to uFldX/vFldY |
138 |
|
|
do j = 1-oly,sny+oly |
139 |
|
|
do i = 1-olx,snx+olx |
140 |
|
|
uFldX(i,j,k,bi,bj) = 0. _d 0 |
141 |
|
|
vFldY(i,j,k,bi,bj) = 0. _d 0 |
142 |
|
|
tmpU(i,j) = 0. _d 0 |
143 |
|
|
tmpV(i,j) = 0. _d 0 |
144 |
|
|
enddo |
145 |
|
|
enddo |
146 |
|
|
C 1) rotation |
147 |
|
|
do j = 1,sny |
148 |
|
|
do i = 1-olx,snx+olx |
149 |
|
|
tmpU(i,j) = |
150 |
|
|
& angleCosC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
151 |
|
|
& +angleSinC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
152 |
|
|
tmpV(i,j) = |
153 |
|
|
& -angleSinC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
154 |
|
|
& +angleCosC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
155 |
|
|
enddo |
156 |
|
|
enddo |
157 |
|
|
do j = 1-oly,sny+oly |
158 |
|
|
do i = 1,snx |
159 |
|
|
tmpU(i,j) = |
160 |
|
|
& angleCosC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
161 |
|
|
& +angleSinC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
162 |
|
|
tmpV(i,j) = |
163 |
|
|
& -angleSinC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
164 |
|
|
& +angleCosC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
165 |
|
|
enddo |
166 |
|
|
enddo |
167 |
|
|
|
168 |
|
|
if ( switchGrid ) then |
169 |
|
|
C 2a) go from A grid velocity points to C grid velocity points |
170 |
|
|
do j = 1,sny |
171 |
|
|
do i = 1,snx |
172 |
|
|
uFldX(i,j,k,bi,bj) = 0.5 _d 0 |
173 |
|
|
& *( tmpU(i-1,j) + tmpU(i,j) ) |
174 |
|
|
vFldY(i,j,k,bi,bj) = 0.5 _d 0 |
175 |
|
|
& *( tmpV(i,j-1) + tmpV(i,j) ) |
176 |
|
|
if (applyMask) then |
177 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskW(i,j,kk,bi,bj) |
178 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskS(i,j,kk,bi,bj) |
179 |
|
|
endif |
180 |
|
|
enddo |
181 |
|
|
enddo |
182 |
|
|
else |
183 |
|
|
C 2b) stay at A grid velocity points (i.e. at the C grid cell center) |
184 |
|
|
do j = 1,sny |
185 |
|
|
do i = 1,snx |
186 |
|
|
uFldX(i,j,k,bi,bj) = tmpU(i,j) |
187 |
|
|
vFldY(i,j,k,bi,bj) = tmpV(i,j) |
188 |
|
|
if (applyMask) then |
189 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
190 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
191 |
|
|
endif |
192 |
|
|
enddo |
193 |
|
|
enddo |
194 |
|
|
endif!if ( switchGrid ) then |
195 |
|
|
|
196 |
|
|
endif!if (xy2en) then |
197 |
|
|
|
198 |
|
|
enddo |
199 |
|
|
enddo |
200 |
|
|
enddo |
201 |
|
|
|
202 |
|
|
return |
203 |
|
|
end |
204 |
|
|
|
205 |
|
|
subroutine rotate_uv2en_rs( |
206 |
|
|
U uFldX, vFldY, |
207 |
|
|
U uFldE, vFldN, |
208 |
|
|
I xy2en, switchGrid, applyMask, kSize, mythid |
209 |
|
|
& ) |
210 |
|
|
|
211 |
|
|
c ================================================================== |
212 |
|
|
c SUBROUTINE rotate_uv2en_rs |
213 |
|
|
c ================================================================== |
214 |
|
|
c |
215 |
|
|
c o uFldX/vFldY are in the model grid directions. |
216 |
|
|
c o uFldE/vFldN are eastward/northward. |
217 |
|
|
c o This routine goes from uFldX/vFldY to uFldE/vFldN (for xy2en=.TRUE.) |
218 |
|
|
c or vice versa (for xy2en=.FALSE.). |
219 |
|
|
c o uFldX/vFldY may be at the C grid velocity points, or at the A grid |
220 |
|
|
c velocity points (i.e. the C grid cell center). uFldE/vFldN are always |
221 |
|
|
c at the cell center. If switchGrid=.TRUE. we go from C grid uFldX/vFldY |
222 |
|
|
c to A grid uFldE/vFldN, or vice versa. If switchGrid=.FALSE. we go |
223 |
|
|
c from A grid uFldX/vFldY to A grid uFldE/vFldN, or vice versa. |
224 |
|
|
c o If applyMask=.TRUE. then masks are applied to the output. |
225 |
|
|
c If kSize=1 (resp. nr) we then use the surface (resp. 3D) masks. |
226 |
|
|
c o In any case it is assumed that exchanges are done on the outside. |
227 |
|
|
c |
228 |
|
|
c ================================================================== |
229 |
|
|
c SUBROUTINE rotate_uv2en_rs |
230 |
|
|
c ================================================================== |
231 |
|
|
|
232 |
|
|
implicit none |
233 |
|
|
|
234 |
|
|
c == global variables == |
235 |
|
|
|
236 |
|
|
#include "EEPARAMS.h" |
237 |
|
|
#include "SIZE.h" |
238 |
|
|
#include "PARAMS.h" |
239 |
|
|
#include "GRID.h" |
240 |
|
|
|
241 |
|
|
c == routine arguments == |
242 |
|
|
|
243 |
|
|
integer kSize |
244 |
|
|
logical xy2en, switchGrid, applyMask |
245 |
|
|
_RS uFldX(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
246 |
|
|
_RS vFldY(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
247 |
|
|
_RS uFldE(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
248 |
|
|
_RS vFldN(1-olx:snx+olx,1-oly:sny+oly,kSize,nsx,nsy) |
249 |
|
|
|
250 |
|
|
integer mythid |
251 |
|
|
|
252 |
|
|
c == local variables == |
253 |
|
|
|
254 |
|
|
integer bi,bj |
255 |
|
|
integer i,j,k,kk |
256 |
|
|
_RS tmpU(1-olx:snx+olx,1-oly:sny+oly) |
257 |
|
|
_RS tmpV(1-olx:snx+olx,1-oly:sny+oly) |
258 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
259 |
|
|
|
260 |
|
|
c == end of interface == |
261 |
|
|
|
262 |
|
|
if ( (kSize.NE.1).AND.(kSize.NE.nr) |
263 |
|
|
& .AND.(applyMask) ) then |
264 |
|
|
WRITE(msgBuf,'(2A,I4,A)') ' ROTATE_UV2EN: ', |
265 |
|
|
& 'no mask has ',kSize,' levels' |
266 |
|
|
CALL PRINT_ERROR(msgBuf, myThid) |
267 |
|
|
STOP 'ABNROMAL END: S/R ROTATE_UV2EN' |
268 |
|
|
endif |
269 |
|
|
|
270 |
|
|
do bj = mybylo(mythid),mybyhi(mythid) |
271 |
|
|
do bi = mybxlo(mythid),mybxhi(mythid) |
272 |
|
|
do k = 1,kSize |
273 |
|
|
|
274 |
|
|
if ( (kSize.EQ.1).AND.(usingPCoords) ) then |
275 |
|
|
kk=nr |
276 |
|
|
else |
277 |
|
|
kk=k |
278 |
|
|
endif |
279 |
|
|
|
280 |
|
|
if ( xy2en ) then |
281 |
|
|
c go from uFldX/vFldY to uFldE/vFldN |
282 |
|
|
if ( switchGrid ) then |
283 |
|
|
C 1a) go from C grid velocity points to A grid velocity points |
284 |
|
|
do i = 1-olx,snx+olx |
285 |
|
|
tmpU(i,sny+Oly)=0. |
286 |
|
|
tmpV(i,sny+Oly)=0. |
287 |
|
|
enddo |
288 |
|
|
do j = 1-oly,sny+oly-1 |
289 |
|
|
tmpU(snx+Olx,j)=0. |
290 |
|
|
tmpV(snx+Olx,j)=0. |
291 |
|
|
do i = 1-olx,snx+olx-1 |
292 |
|
|
tmpU(i,j) = 0.5 _d 0 |
293 |
|
|
& *( uFldX(i+1,j,k,bi,bj) + uFldX(i,j,k,bi,bj) ) |
294 |
|
|
tmpV(i,j) = 0.5 _d 0 |
295 |
|
|
& *( vFldY(i,j+1,k,bi,bj) + vFldY(i,j,k,bi,bj) ) |
296 |
|
|
if (applyMask) then |
297 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
298 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
299 |
|
|
endif |
300 |
|
|
enddo |
301 |
|
|
enddo |
302 |
|
|
else |
303 |
|
|
C 1b) stay at A grid velocity points (i.e. at the C grid cell center) |
304 |
|
|
do j = 1-oly,sny+oly |
305 |
|
|
do i = 1-olx,snx+olx |
306 |
|
|
tmpU(i,j) = uFldX(i,j,k,bi,bj) |
307 |
|
|
tmpV(i,j) = vFldY(i,j,k,bi,bj) |
308 |
|
|
if (applyMask) then |
309 |
|
|
tmpU(i,j) = tmpU(i,j) * maskC(i,j,kk,bi,bj) |
310 |
|
|
tmpV(i,j) = tmpV(i,j) * maskC(i,j,kk,bi,bj) |
311 |
|
|
endif |
312 |
|
|
enddo |
313 |
|
|
enddo |
314 |
|
|
endif!if ( switchGrid ) then |
315 |
|
|
|
316 |
|
|
C 2) rotation |
317 |
|
|
do j = 1-oly,sny+oly |
318 |
|
|
do i = 1-olx,snx+olx |
319 |
|
|
uFldE(i,j,k,bi,bj) = |
320 |
|
|
& angleCosC(i,j,bi,bj)*tmpU(i,j) |
321 |
|
|
& -angleSinC(i,j,bi,bj)*tmpV(i,j) |
322 |
|
|
vFldN(i,j,k,bi,bj) = |
323 |
|
|
& angleSinC(i,j,bi,bj)*tmpU(i,j) |
324 |
|
|
& +angleCosC(i,j,bi,bj)*tmpV(i,j) |
325 |
|
|
enddo |
326 |
|
|
enddo |
327 |
|
|
|
328 |
|
|
else!if (xy2en) then |
329 |
|
|
c go from uFldE/vFldN to uFldX/vFldY |
330 |
|
|
|
331 |
|
|
C 1) rotation |
332 |
|
|
do j = 1-oly,sny+oly |
333 |
|
|
do i = 1-olx,snx+olx |
334 |
|
|
tmpU(i,j) = |
335 |
|
|
& angleCosC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
336 |
|
|
& +angleSinC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
337 |
|
|
tmpV(i,j) = |
338 |
|
|
& -angleSinC(i,j,bi,bj)*uFldE(i,j,k,bi,bj) |
339 |
|
|
& +angleCosC(i,j,bi,bj)*vFldN(i,j,k,bi,bj) |
340 |
|
|
enddo |
341 |
|
|
enddo |
342 |
|
|
|
343 |
|
|
if ( switchGrid ) then |
344 |
|
|
C 2a) go from A grid velocity points to C grid velocity points |
345 |
|
|
do i = 1-olx,snx+olx |
346 |
|
|
uFldX(i,1,k,bi,bj)=0. |
347 |
|
|
vFldY(i,1,k,bi,bj)=0. |
348 |
|
|
enddo |
349 |
|
|
do j = 1-oly+1,sny+oly |
350 |
|
|
uFldX(1,j,k,bi,bj)=0. |
351 |
|
|
vFldY(1,j,k,bi,bj)=0. |
352 |
|
|
do i = 1-olx+1,snx+olx |
353 |
|
|
uFldX(i,j,k,bi,bj) = 0.5 _d 0 |
354 |
|
|
& *( tmpU(i-1,j) + tmpU(i,j) ) |
355 |
|
|
vFldY(i,j,k,bi,bj) = 0.5 _d 0 |
356 |
|
|
& *( tmpV(i,j-1) + tmpV(i,j) ) |
357 |
|
|
if (applyMask) then |
358 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskW(i,j,kk,bi,bj) |
359 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskS(i,j,kk,bi,bj) |
360 |
|
|
endif |
361 |
|
|
enddo |
362 |
|
|
enddo |
363 |
|
|
else |
364 |
|
|
C 2b) stay at A grid velocity points (i.e. at the C grid cell center) |
365 |
|
|
do j = 1-oly,sny+oly |
366 |
|
|
do i = 1-olx,snx+olx |
367 |
|
|
uFldX(i,j,k,bi,bj) = tmpU(i,j) |
368 |
|
|
vFldY(i,j,k,bi,bj) = tmpV(i,j) |
369 |
|
|
if (applyMask) then |
370 |
|
|
uFldX(i,j,k,bi,bj)=uFldX(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
371 |
|
|
vFldY(i,j,k,bi,bj)=vFldY(i,j,k,bi,bj)*maskC(i,j,kk,bi,bj) |
372 |
|
|
endif |
373 |
|
|
enddo |
374 |
|
|
enddo |
375 |
|
|
endif!if ( switchGrid ) then |
376 |
|
|
|
377 |
|
|
endif!if (xy2en) then |
378 |
|
|
|
379 |
|
|
enddo |
380 |
|
|
enddo |
381 |
|
|
enddo |
382 |
|
|
|
383 |
|
|
return |
384 |
|
|
end |
385 |
|
|
|