fix: Fixes for linter action and code style (#869)
[FMS.git] / mpp / include / mpp_do_updateV_ad.h
blob65a803b6b8b57cfa3cb06fd9b3914ffc12723b7f
1 ! -*-f90-*-
4 !***********************************************************************
5 !* GNU Lesser General Public License
6 !*
7 !* This file is part of the GFDL Flexible Modeling System (FMS).
8 !*
9 !* FMS is free software: you can redistribute it and/or modify it under
10 !* the terms of the GNU Lesser General Public License as published by
11 !* the Free Software Foundation, either version 3 of the License, or (at
12 !* your option) any later version.
14 !* FMS is distributed in the hope that it will be useful, but WITHOUT
15 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
16 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17 !* for more details.
19 !* You should have received a copy of the GNU Lesser General Public
20 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
21 !***********************************************************************
23 subroutine MPP_DO_UPDATE_AD_3D_V_(f_addrsx,f_addrsy, domain, update_x, update_y, &
24 d_type, ke, gridtype, flags)
25 !updates data domain of 3D field whose computational domains have been computed
26 integer(i8_kind), intent(in) :: f_addrsx(:,:), f_addrsy(:,:)
27 type(domain2d), intent(in) :: domain
28 type(overlapSpec), intent(in) :: update_x, update_y
29 integer, intent(in) :: ke
30 MPP_TYPE_, intent(in) :: d_type ! creates unique interface
31 integer, intent(in) :: gridtype
32 integer, intent(in), optional :: flags
34 MPP_TYPE_ :: fieldx(update_x%xbegin:update_x%xend, update_x%ybegin:update_x%yend,ke)
35 MPP_TYPE_ :: fieldy(update_y%xbegin:update_y%xend, update_y%ybegin:update_y%yend,ke)
36 pointer(ptr_fieldx, fieldx)
37 pointer(ptr_fieldy, fieldy)
39 integer :: update_flags
40 integer :: l_size, l, i, j, k, is, ie, js, je, n, m
41 integer :: pos, nlist, msgsize
42 integer :: to_pe, from_pe, midpoint
43 integer :: tMe, dir
45 integer, allocatable :: msg1(:), msg2(:)
46 logical :: send(8), recv(8), update_edge_only
47 MPP_TYPE_ :: buffer(size(mpp_domains_stack(:)))
48 pointer(ptr,buffer )
49 integer :: buffer_pos
50 integer :: buffer_recv_size, shift
51 integer :: rank_x, rank_y, ind_x, ind_y, cur_rank
52 integer :: nsend_x, nsend_y, nrecv_x, nrecv_y, outunit
54 outunit = stdout()
55 update_flags = XUPDATE+YUPDATE !default
56 if( PRESENT(flags) ) then
57 update_flags = flags
58 ! The following test is so that SCALAR_PAIR can be used alone with the
59 ! same default update pattern as without.
60 if (BTEST(update_flags,SCALAR_BIT)) then
61 if (.NOT.(BTEST(update_flags,WEST) .OR. BTEST(update_flags,EAST) &
62 .OR. BTEST(update_flags,NORTH) .OR. BTEST(update_flags,SOUTH))) &
63 update_flags = update_flags + XUPDATE+YUPDATE !default with SCALAR_PAIR
64 end if
65 end if
67 if( BTEST(update_flags,NORTH) .AND. BTEST(domain%fold,NORTH) .AND. BTEST(gridtype,SOUTH) ) &
68 call mpp_error( FATAL, 'MPP_DO_UPDATE_V: Incompatible grid offset and fold.' )
70 update_edge_only = BTEST(update_flags, EDGEONLY)
71 recv(1) = BTEST(update_flags,EAST)
72 recv(3) = BTEST(update_flags,SOUTH)
73 recv(5) = BTEST(update_flags,WEST)
74 recv(7) = BTEST(update_flags,NORTH)
75 if( update_edge_only ) then
76 if( .NOT. (recv(1) .OR. recv(3) .OR. recv(5) .OR. recv(7)) ) then
77 recv(1) = .true.
78 recv(3) = .true.
79 recv(5) = .true.
80 recv(7) = .true.
81 endif
82 else
83 recv(2) = recv(1) .AND. recv(3)
84 recv(4) = recv(3) .AND. recv(5)
85 recv(6) = recv(5) .AND. recv(7)
86 recv(8) = recv(7) .AND. recv(1)
87 endif
89 send = recv
91 l_size = size(f_addrsx,1)
92 nlist = size(domain%list(:))
93 ptr = LOC(mpp_domains_stack)
95 !recv
96 nsend_x = update_x%nsend
97 nsend_y = update_y%nsend
98 nrecv_x = update_x%nrecv
99 nrecv_y = update_y%nrecv
101 if(debug_message_passing) then
102 allocate(msg1(0:nlist-1), msg2(0:nlist-1) )
103 msg1 = 0
104 msg2 = 0
105 cur_rank = get_rank_recv(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
107 do while (ind_x .LE. nrecv_x .OR. ind_y .LE. nrecv_y)
108 msgsize = 0
109 if(cur_rank == rank_x) then
110 from_pe = update_x%recv(ind_x)%pe
111 do n = 1, update_x%recv(ind_x)%count
112 dir = update_x%recv(ind_x)%dir(n)
113 if(recv(dir)) then
114 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
115 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
116 msgsize = msgsize + (ie-is+1)*(je-js+1)
117 end if
118 end do
119 ind_x = ind_x+1
120 if(ind_x .LE. nrecv_x) then
121 rank_x = update_x%recv(ind_x)%pe - domain%pe
122 if(rank_x .LE.0) rank_x = rank_x + nlist
123 else
124 rank_x = -1
125 endif
126 endif
127 if(cur_rank == rank_y) then
128 from_pe = update_y%recv(ind_y)%pe
129 do n = 1, update_y%recv(ind_y)%count
130 dir = update_y%recv(ind_y)%dir(n)
131 if(recv(dir)) then
132 is = update_y%recv(ind_y)%is(n); ie = update_y%recv(ind_y)%ie(n)
133 js = update_y%recv(ind_y)%js(n); je = update_y%recv(ind_y)%je(n)
134 msgsize = msgsize + (ie-is+1)*(je-js+1)
135 end if
136 end do
137 ind_y = ind_y+1
138 if(ind_y .LE. nrecv_y) then
139 rank_y = update_y%recv(ind_y)%pe - domain%pe
140 if(rank_y .LE.0) rank_y = rank_y + nlist
141 else
142 rank_y = -1
143 endif
144 endif
145 cur_rank = max(rank_x, rank_y)
146 m = from_pe-mpp_root_pe()
147 call mpp_recv( msg1(m), glen=1, from_pe=from_pe, block=.FALSE., tag=COMM_TAG_1)
148 msg2(m) = msgsize
149 end do
151 cur_rank = get_rank_send(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
152 do while (ind_x .LE. nsend_x .OR. ind_y .LE. nsend_y)
153 msgsize = 0
154 if(cur_rank == rank_x) then
155 to_pe = update_x%send(ind_x)%pe
156 do n = 1, update_x%send(ind_x)%count
157 dir = update_x%send(ind_x)%dir(n)
158 if( send(dir) ) then
159 is = update_x%send(ind_x)%is(n); ie = update_x%send(ind_x)%ie(n)
160 js = update_x%send(ind_x)%js(n); je = update_x%send(ind_x)%je(n)
161 msgsize = msgsize + (ie-is+1)*(je-js+1)
162 end if
163 end do
164 ind_x = ind_x+1
165 if(ind_x .LE. nsend_x) then
166 rank_x = update_x%send(ind_x)%pe - domain%pe
167 if(rank_x .LT.0) rank_x = rank_x + nlist
168 else
169 rank_x = nlist+1
170 endif
171 endif
172 if(cur_rank == rank_y) then
173 to_pe = update_y%send(ind_y)%pe
174 do n = 1, update_y%send(ind_y)%count
175 dir = update_y%send(ind_y)%dir(n)
176 if( send(dir) ) then
177 is = update_y%send(ind_y)%is(n); ie = update_y%send(ind_y)%ie(n)
178 js = update_y%send(ind_y)%js(n); je = update_y%send(ind_y)%je(n)
179 msgsize = msgsize + (ie-is+1)*(je-js+1)
180 end if
181 end do
182 ind_y = ind_y+1
183 if(ind_y .LE. nsend_y) then
184 rank_y = update_y%send(ind_y)%pe - domain%pe
185 if(rank_y .LT.0) rank_y = rank_y + nlist
186 else
187 rank_y = nlist+1
188 endif
189 endif
190 cur_rank = min(rank_x, rank_y)
191 call mpp_send( msgsize, plen=1, to_pe=to_pe, tag=COMM_TAG_1)
192 enddo
193 call mpp_sync_self(check=EVENT_RECV)
194 do m = 0, nlist-1
195 if(msg1(m) .NE. msg2(m)) then
196 print*, "My pe = ", mpp_pe(), ",domain name =", trim(domain%name), ",from pe=", &
197 domain%list(m)%pe, ":send size = ", msg1(m), ", recv size = ", msg2(m)
198 call mpp_error(FATAL, "mpp_do_updateV: mismatch on send and recv size")
199 endif
200 enddo
202 call mpp_sync_self()
203 write(outunit,*)"NOTE from mpp_do_updateV: message sizes are matched between send and recv for domain " &
204 //trim(domain%name)
205 deallocate(msg1, msg2)
206 endif
208 buffer_pos = 0
209 cur_rank = get_rank_recv(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
211 do while (ind_x .LE. nrecv_x .OR. ind_y .LE. nrecv_y)
212 msgsize = 0
213 select case(gridtype)
214 case(BGRID_NE, BGRID_SW, AGRID)
215 if(cur_rank == rank_x) then
216 from_pe = update_x%recv(ind_x)%pe
217 do n = 1, update_x%recv(ind_x)%count
218 dir = update_x%recv(ind_x)%dir(n)
219 if(recv(dir)) then
220 tMe = update_x%recv(ind_x)%tileMe(n)
221 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
222 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
223 msgsize = msgsize + (ie-is+1)*(je-js+1)
224 end if
225 end do
226 msgsize = msgsize*2
227 ind_x = ind_x+1
228 ind_y = ind_x
229 if(ind_x .LE. nrecv_x) then
230 rank_x = update_x%recv(ind_x)%pe - domain%pe
231 if(rank_x .LE.0) rank_x = rank_x + nlist
232 else
233 rank_x = -1
234 endif
235 rank_y = rank_x
236 endif
237 case(CGRID_NE, CGRID_SW)
238 if(cur_rank == rank_x) then
239 from_pe = update_x%recv(ind_x)%pe
240 do n = 1, update_x%recv(ind_x)%count
241 dir = update_x%recv(ind_x)%dir(n)
242 if(recv(dir)) then
243 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
244 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
245 msgsize = msgsize + (ie-is+1)*(je-js+1)
246 end if
247 end do
248 ind_x = ind_x+1
249 if(ind_x .LE. nrecv_x) then
250 rank_x = update_x%recv(ind_x)%pe - domain%pe
251 if(rank_x .LE.0) rank_x = rank_x + nlist
252 else
253 rank_x = -1
254 endif
255 endif
256 if(cur_rank == rank_y) then
257 from_pe = update_y%recv(ind_y)%pe
258 do n = 1, update_y%recv(ind_y)%count
259 dir = update_y%recv(ind_y)%dir(n)
260 if(recv(dir)) then
261 is = update_y%recv(ind_y)%is(n); ie = update_y%recv(ind_y)%ie(n)
262 js = update_y%recv(ind_y)%js(n); je = update_y%recv(ind_y)%je(n)
263 msgsize = msgsize + (ie-is+1)*(je-js+1)
264 end if
265 end do
266 ind_y = ind_y+1
267 if(ind_y .LE. nrecv_y) then
268 rank_y = update_y%recv(ind_y)%pe - domain%pe
269 if(rank_y .LE.0) rank_y = rank_y + nlist
270 else
271 rank_y = -1
272 endif
273 endif
274 end select
275 cur_rank = max(rank_x, rank_y)
276 msgsize = msgsize*ke*l_size
278 if( msgsize.GT.0 )then
279 buffer_pos = buffer_pos + msgsize
280 end if
281 end do
282 buffer_recv_size = buffer_pos
284 !unpacking
285 buffer_pos = buffer_recv_size
286 cur_rank = get_rank_unpack(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
288 do while (ind_x > 0 .OR. ind_y > 0)
289 pos = buffer_pos
290 select case ( gridtype )
291 case(BGRID_NE, BGRID_SW, AGRID)
292 if(cur_rank == rank_x) then
293 do n = update_x%recv(ind_x)%count, 1, -1
294 dir = update_x%recv(ind_x)%dir(n)
295 if( recv(dir) ) then
296 tMe = update_x%recv(ind_x)%tileMe(n)
297 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
298 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
299 msgsize = (ie-is+1)*(je-js+1)*ke*2*l_size
300 pos = buffer_pos - msgsize
301 buffer_pos = pos
302 do l=1, l_size ! loop over number of fields
303 ptr_fieldx = f_addrsx(l, tMe)
304 ptr_fieldy = f_addrsy(l, tMe)
305 do k = 1,ke
306 do j = js, je
307 do i = is, ie
308 pos = pos + 2
309 buffer(pos-1) = fieldx(i,j,k)
310 buffer(pos) = fieldy(i,j,k)
311 fieldx(i,j,k) = 0.
312 fieldy(i,j,k) = 0.
313 end do
314 end do
315 end do
316 end do
317 end if ! end if( recv(dir) )
318 end do ! do dir=8,1,-1
319 ind_x = ind_x-1
320 ind_y = ind_x
321 if(ind_x .GT. 0) then
322 rank_x = update_x%recv(ind_x)%pe - domain%pe
323 if(rank_x .LE.0) rank_x = rank_x + nlist
324 else
325 rank_x = nlist+1
326 endif
327 rank_y = rank_x
328 endif
329 case(CGRID_NE, CGRID_SW)
330 if(cur_rank == rank_y) then
331 do n = update_y%recv(ind_y)%count, 1, -1
332 dir = update_y%recv(ind_y)%dir(n)
333 if( recv(dir) ) then
334 tMe = update_y%recv(ind_y)%tileMe(n)
335 is = update_y%recv(ind_y)%is(n); ie = update_y%recv(ind_y)%ie(n)
336 js = update_y%recv(ind_y)%js(n); je = update_y%recv(ind_y)%je(n)
337 msgsize = (ie-is+1)*(je-js+1)*ke*l_size
338 pos = buffer_pos - msgsize
339 buffer_pos = pos
340 do l=1,l_size ! loop over number of fields
341 ptr_fieldx = f_addrsx(l, tMe)
342 ptr_fieldy = f_addrsy(l, tMe)
343 do k = 1,ke
344 do j = js, je
345 do i = is, ie
346 pos = pos + 1
347 buffer(pos) = fieldy(i,j,k)
348 fieldy(i,j,k) = 0.
349 end do
350 end do
351 end do
352 end do
353 end if
354 end do
355 ind_y = ind_y-1
356 if(ind_y .GT. 0) then
357 rank_y = update_y%recv(ind_y)%pe - domain%pe
358 if(rank_y .LE.0) rank_y = rank_y + nlist
359 else
360 rank_y = nlist+1
361 endif
362 endif
363 if(cur_rank == rank_x) then
364 do n = update_x%recv(ind_x)%count, 1, -1
365 dir = update_x%recv(ind_x)%dir(n)
366 if( recv(dir) ) then
367 tMe = update_x%recv(ind_x)%tileMe(n)
368 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
369 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
370 msgsize = (ie-is+1)*(je-js+1)*ke*l_size
371 pos = buffer_pos - msgsize
372 buffer_pos = pos
373 do l=1,l_size ! loop over number of fields
374 ptr_fieldx = f_addrsx(l, tMe)
375 ptr_fieldy = f_addrsy(l, tMe)
376 do k = 1,ke
377 do j = js, je
378 do i = is, ie
379 pos = pos + 1
380 buffer(pos) = fieldx(i,j,k)
381 fieldx(i,j,k) = 0.
382 end do
383 end do
384 end do
385 end do
386 end if
387 end do
388 ind_x = ind_x-1
389 if(ind_x .GT. 0) then
390 rank_x = update_x%recv(ind_x)%pe - domain%pe
391 if(rank_x .LE.0) rank_x = rank_x + nlist
392 else
393 rank_x = nlist+1
394 endif
395 endif
396 end select
397 cur_rank = min(rank_x, rank_y)
398 end do
400 ! ---northern boundary fold
401 shift = 0
402 if(domain%symmetry) shift = 1
403 if( BTEST(domain%fold,NORTH) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then
404 j = domain%y(1)%global%end+shift
405 if( domain%y(1)%data%begin.LE.j .AND. j.LE.domain%y(1)%data%end+shift )then !fold is within domain
406 !poles set to 0: BGRID only
407 if( gridtype.EQ.BGRID_NE )then
408 midpoint = (domain%x(1)%global%begin+domain%x(1)%global%end-1+shift)/2
409 j = domain%y(1)%global%end+shift
410 is = domain%x(1)%global%begin; ie = domain%x(1)%global%end+shift
411 if( .NOT. domain%symmetry ) is = is - 1
412 do i = is ,ie, midpoint
413 if( domain%x(1)%data%begin.LE.i .AND. i.LE. domain%x(1)%data%end+shift )then
414 do l=1,l_size
415 ptr_fieldx = f_addrsx(l, 1)
416 ptr_fieldy = f_addrsy(l, 1)
417 do k = 1,ke
418 fieldx(i,j,k) = 0.
419 fieldy(i,j,k) = 0.
420 end do
421 end do
422 end if
423 end do
424 endif
426 ! the following code code block correct an error where the data in your halo coming from
427 ! other half may have the wrong sign
428 !off west edge, when update north or west direction
429 j = domain%y(1)%global%end+shift
430 if ( recv(7) .OR. recv(5) ) then
431 select case(gridtype)
432 case(BGRID_NE)
433 if(domain%symmetry) then
434 is = domain%x(1)%global%begin
435 else
436 is = domain%x(1)%global%begin - 1
437 end if
438 if( is.GT.domain%x(1)%data%begin )then
440 if( 2*is-domain%x(1)%data%begin.GT.domain%x(1)%data%end+shift ) &
441 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-north BGRID_NE west edge ubound error.' )
442 do l=1,l_size
443 ptr_fieldx = f_addrsx(l, 1)
444 ptr_fieldy = f_addrsy(l, 1)
445 do k = 1,ke
446 do i = domain%x(1)%data%begin,is-1
447 fieldx(2*is-i,j,k) = fieldx(2*is-i,j,k) + fieldx(i,j,k)
448 fieldy(2*is-i,j,k) = fieldy(2*is-i,j,k) + fieldy(i,j,k)
449 end do
450 end do
451 end do
452 end if
453 case(CGRID_NE)
454 is = domain%x(1)%global%begin
455 if( is.GT.domain%x(1)%data%begin )then
456 if( 2*is-domain%x(1)%data%begin-1.GT.domain%x(1)%data%end ) &
457 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-north CGRID_NE west edge ubound error.' )
458 do l=1,l_size
459 ptr_fieldy = f_addrsy(l, 1)
460 do k = 1,ke
461 do i = domain%x(1)%data%begin,is-1
462 fieldy(2*is-i-1,j,k) = fieldy(2*is-i-1,j,k) + fieldy(i,j,k)
463 end do
464 end do
465 end do
466 end if
467 end select
468 end if
470 !off east edge
471 is = domain%x(1)%global%end
472 if(domain%x(1)%cyclic .AND. is.LT.domain%x(1)%data%end )then
473 ie = domain%x(1)%data%end
474 is = is + 1
475 select case(gridtype)
476 case(BGRID_NE)
477 is = is + shift
478 ie = ie + shift
479 do l=1,l_size
480 ptr_fieldx = f_addrsx(l, 1)
481 ptr_fieldy = f_addrsy(l, 1)
482 do k = 1,ke
483 do i = is,ie
484 fieldx(i,j,k) = -fieldx(i,j,k)
485 fieldy(i,j,k) = -fieldy(i,j,k)
486 end do
487 end do
488 end do
489 case(CGRID_NE)
490 do l=1,l_size
491 ptr_fieldy = f_addrsy(l, 1)
492 do k = 1,ke
493 do i = is, ie
494 fieldy(i,j,k) = -fieldy(i,j,k)
495 end do
496 end do
497 end do
498 end select
499 end if
500 end if
501 else if( BTEST(domain%fold,SOUTH) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then ! ---southern
502 !! boundary fold
503 ! NOTE: symmetry is assumed for fold-south boundary
504 j = domain%y(1)%global%begin
505 if( domain%y(1)%data%begin.LE.j .AND. j.LE.domain%y(1)%data%end+shift )then !fold is within domain
506 midpoint = (domain%x(1)%global%begin+domain%x(1)%global%end-1+shift)/2
507 !poles set to 0: BGRID only
508 if( gridtype.EQ.BGRID_NE )then
509 j = domain%y(1)%global%begin
510 is = domain%x(1)%global%begin; ie = domain%x(1)%global%end+shift
511 do i = is ,ie, midpoint
512 if( domain%x(1)%data%begin.LE.i .AND. i.LE. domain%x(1)%data%end+shift )then
513 do l=1,l_size
514 ptr_fieldx = f_addrsx(l, 1)
515 ptr_fieldy = f_addrsy(l, 1)
516 do k = 1,ke
517 fieldx(i,j,k) = 0.
518 fieldy(i,j,k) = 0.
519 end do
520 end do
521 end if
522 end do
523 endif
525 ! the following code code block correct an error where the data in your halo coming from
526 ! other half may have the wrong sign
527 !off west edge, when update north or west direction
528 j = domain%y(1)%global%begin
529 if ( recv(3) .OR. recv(5) ) then
530 select case(gridtype)
531 case(BGRID_NE)
532 is = domain%x(1)%global%begin
533 if( is.GT.domain%x(1)%data%begin )then
535 if( 2*is-domain%x(1)%data%begin.GT.domain%x(1)%data%end+shift ) &
536 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-south BGRID_NE west edge ubound error.' )
537 do l=1,l_size
538 ptr_fieldx = f_addrsx(l, 1)
539 ptr_fieldy = f_addrsy(l, 1)
540 do k = 1,ke
541 do i = domain%x(1)%data%begin,is-1
542 fieldx(2*is-i,j,k) = fieldx(2*is-i,j,k) + fieldx(i,j,k)
543 fieldy(2*is-i,j,k) = fieldy(2*is-i,j,k) + fieldy(i,j,k)
544 end do
545 end do
546 end do
547 end if
548 case(CGRID_NE)
549 is = domain%x(1)%global%begin
550 if( is.GT.domain%x(1)%data%begin )then
551 if( 2*is-domain%x(1)%data%begin-1.GT.domain%x(1)%data%end ) &
552 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-south CGRID_NE west edge ubound error.' )
553 do l=1,l_size
554 ptr_fieldy = f_addrsy(l, 1)
555 do k = 1,ke
556 do i = domain%x(1)%data%begin,is-1
557 fieldy(2*is-i-1,j,k) = fieldy(2*is-i-1,j,k) + fieldy(i,j,k)
558 end do
559 end do
560 end do
561 end if
562 end select
563 end if
565 !off east edge
566 is = domain%x(1)%global%end
567 if(domain%x(1)%cyclic .AND. is.LT.domain%x(1)%data%end )then
568 ie = domain%x(1)%data%end
569 is = is + 1
570 select case(gridtype)
571 case(BGRID_NE)
572 is = is + shift
573 ie = ie + shift
574 do l=1,l_size
575 ptr_fieldx = f_addrsx(l, 1)
576 ptr_fieldy = f_addrsy(l, 1)
577 do k = 1,ke
578 do i = is,ie
579 fieldx(i,j,k) = -fieldx(i,j,k)
580 fieldy(i,j,k) = -fieldy(i,j,k)
581 end do
582 end do
583 end do
584 case(CGRID_NE)
585 do l=1,l_size
586 ptr_fieldy = f_addrsy(l, 1)
587 do k = 1,ke
588 do i = is, ie
589 fieldy(i,j,k) = -fieldy(i,j,k)
590 end do
591 end do
592 end do
593 end select
594 end if
595 end if
596 else if( BTEST(domain%fold,WEST) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then ! ---eastern
597 !! boundary fold
598 ! NOTE: symmetry is assumed for fold-west boundary
599 i = domain%x(1)%global%begin
600 if( domain%x(1)%data%begin.LE.i .AND. i.LE.domain%x(1)%data%end+shift )then !fold is within domain
601 midpoint = (domain%y(1)%global%begin+domain%y(1)%global%end-1+shift)/2
602 !poles set to 0: BGRID only
603 if( gridtype.EQ.BGRID_NE )then
604 i = domain%x(1)%global%begin
605 js = domain%y(1)%global%begin; je = domain%y(1)%global%end+shift
606 do j = js ,je, midpoint
607 if( domain%y(1)%data%begin.LE.j .AND. j.LE. domain%y(1)%data%end+shift )then
608 do l=1,l_size
609 ptr_fieldx = f_addrsx(l, 1)
610 ptr_fieldy = f_addrsy(l, 1)
611 do k = 1,ke
612 fieldx(i,j,k) = 0.
613 fieldy(i,j,k) = 0.
614 end do
615 end do
616 end if
617 end do
618 endif
620 ! the following code code block correct an error where the data in your halo coming from
621 ! other half may have the wrong sign
622 !off south edge, when update south or west direction
623 i = domain%x(1)%global%begin
624 if ( recv(3) .OR. recv(5) ) then
625 select case(gridtype)
626 case(BGRID_NE)
627 js = domain%y(1)%global%begin
628 if( js.GT.domain%y(1)%data%begin )then
630 if( 2*js-domain%y(1)%data%begin.GT.domain%y(1)%data%end+shift ) &
631 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-west BGRID_NE west edge ubound error.' )
632 do l=1,l_size
633 ptr_fieldx = f_addrsx(l, 1)
634 ptr_fieldy = f_addrsy(l, 1)
635 do k = 1,ke
636 do j = domain%y(1)%data%begin,js-1
637 fieldx(i,2*js-j,k) = fieldx(i,2*js-j,k) + fieldx(i,j,k)
638 fieldy(i,2*js-j,k) = fieldy(i,2*js-j,k) + fieldy(i,j,k)
639 end do
640 end do
641 end do
642 end if
643 case(CGRID_NE)
644 js = domain%y(1)%global%begin
645 if( js.GT.domain%y(1)%data%begin )then
646 if( 2*js-domain%y(1)%data%begin-1.GT.domain%y(1)%data%end ) &
647 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-west CGRID_NE west edge ubound error.' )
648 do l=1,l_size
649 ptr_fieldx = f_addrsx(l, 1)
650 do k = 1,ke
651 do j = domain%y(1)%data%begin,js-1
652 fieldx(i, 2*js-j-1,k) = fieldx(i, 2*js-j-1,k) + fieldx(i,j,k)
653 end do
654 end do
655 end do
656 end if
657 end select
658 end if
660 !off north edge
661 js = domain%y(1)%global%end
662 if(domain%y(1)%cyclic .AND. js.LT.domain%y(1)%data%end )then
663 je = domain%y(1)%data%end
664 js = js + 1
665 select case(gridtype)
666 case(BGRID_NE)
667 js = js + shift
668 je = je + shift
669 do l=1,l_size
670 ptr_fieldx = f_addrsx(l, 1)
671 ptr_fieldy = f_addrsy(l, 1)
672 do k = 1,ke
673 do j = js,je
674 fieldx(i,j,k) = -fieldx(i,j,k)
675 fieldy(i,j,k) = -fieldy(i,j,k)
676 end do
677 end do
678 end do
679 case(CGRID_NE)
680 do l=1,l_size
681 ptr_fieldx = f_addrsx(l, 1)
682 do k = 1,ke
683 do j = js, je
684 fieldx(i,j,k) = -fieldx(i,j,k)
685 end do
686 end do
687 end do
688 end select
689 end if
690 end if
691 else if( BTEST(domain%fold,EAST) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then ! ---eastern
692 !! boundary fold
693 ! NOTE: symmetry is assumed for fold-west boundary
694 i = domain%x(1)%global%end+shift
695 if( domain%x(1)%data%begin.LE.i .AND. i.LE.domain%x(1)%data%end+shift )then !fold is within domain
696 midpoint = (domain%y(1)%global%begin+domain%y(1)%global%end-1+shift)/2
697 !poles set to 0: BGRID only
698 if( gridtype.EQ.BGRID_NE )then
699 i = domain%x(1)%global%end+shift
700 js = domain%y(1)%global%begin; je = domain%y(1)%global%end+shift
701 do j = js ,je, midpoint
702 if( domain%y(1)%data%begin.LE.j .AND. j.LE. domain%y(1)%data%end+shift )then
703 do l=1,l_size
704 ptr_fieldx = f_addrsx(l, 1)
705 ptr_fieldy = f_addrsy(l, 1)
706 do k = 1,ke
707 fieldx(i,j,k) = 0.
708 fieldy(i,j,k) = 0.
709 end do
710 end do
711 end if
712 end do
713 endif
715 ! the following code code block correct an error where the data in your halo coming from
716 ! other half may have the wrong sign
717 !off south edge, when update south or west direction
718 i = domain%x(1)%global%end+shift
719 if ( recv(3) .OR. recv(1) ) then
720 select case(gridtype)
721 case(BGRID_NE)
722 js = domain%y(1)%global%begin
723 if( js.GT.domain%y(1)%data%begin )then
725 if( 2*js-domain%y(1)%data%begin.GT.domain%y(1)%data%end+shift ) &
726 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-east BGRID_NE west edge ubound error.' )
727 do l=1,l_size
728 ptr_fieldx = f_addrsx(l, 1)
729 ptr_fieldy = f_addrsy(l, 1)
730 do k = 1,ke
731 do j = domain%y(1)%data%begin,js-1
732 fieldx(i,2*js-j,k) = fieldx(i,2*js-j,k) + fieldx(i,j,k)
733 fieldy(i,2*js-j,k) = fieldy(i,2*js-j,k) + fieldy(i,j,k)
734 end do
735 end do
736 end do
737 end if
738 case(CGRID_NE)
739 js = domain%y(1)%global%begin
740 if( js.GT.domain%y(1)%data%begin )then
741 if( 2*js-domain%y(1)%data%begin-1.GT.domain%y(1)%data%end ) &
742 call mpp_error( FATAL, 'MPP_DO_UPDATE_AD_V: folded-east CGRID_NE west edge ubound error.' )
743 do l=1,l_size
744 ptr_fieldx = f_addrsx(l, 1)
745 do k = 1,ke
746 do j = domain%y(1)%data%begin,js-1
747 fieldx(i, 2*js-j-1,k) = fieldx(i, 2*js-j-1,k) + fieldx(i,j,k)
748 end do
749 end do
750 end do
751 end if
752 end select
753 end if
755 !off north edge
756 js = domain%y(1)%global%end
757 if(domain%y(1)%cyclic .AND. js.LT.domain%y(1)%data%end )then
758 je = domain%y(1)%data%end
759 js = js + 1
760 select case(gridtype)
761 case(BGRID_NE)
762 js = js + shift
763 je = je + shift
764 do l=1,l_size
765 ptr_fieldx = f_addrsx(l, 1)
766 ptr_fieldy = f_addrsy(l, 1)
767 do k = 1,ke
768 do j = js,je
769 fieldx(i,j,k) = -fieldx(i,j,k)
770 fieldy(i,j,k) = -fieldy(i,j,k)
771 end do
772 end do
773 end do
774 case(CGRID_NE)
775 do l=1,l_size
776 ptr_fieldx = f_addrsx(l, 1)
777 do k = 1,ke
778 do j = js, je
779 fieldx(i,j,k) = -fieldx(i,j,k)
780 end do
781 end do
782 end do
783 end select
784 end if
785 end if
786 end if
787 !unpacking done
789 !--- recv
790 buffer_pos = 0
791 cur_rank = get_rank_recv(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
793 do while (ind_x .LE. nrecv_x .OR. ind_y .LE. nrecv_y)
794 msgsize = 0
795 select case(gridtype)
796 case(BGRID_NE, BGRID_SW, AGRID)
797 if(cur_rank == rank_x) then
798 from_pe = update_x%recv(ind_x)%pe
799 do n = 1, update_x%recv(ind_x)%count
800 dir = update_x%recv(ind_x)%dir(n)
801 if(recv(dir)) then
802 tMe = update_x%recv(ind_x)%tileMe(n)
803 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
804 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
805 msgsize = msgsize + (ie-is+1)*(je-js+1)
806 end if
807 end do
808 msgsize = msgsize*2
809 ind_x = ind_x+1
810 ind_y = ind_x
811 if(ind_x .LE. nrecv_x) then
812 rank_x = update_x%recv(ind_x)%pe - domain%pe
813 if(rank_x .LE.0) rank_x = rank_x + nlist
814 else
815 rank_x = -1
816 endif
817 rank_y = rank_x
818 endif
819 case(CGRID_NE, CGRID_SW)
820 if(cur_rank == rank_x) then
821 from_pe = update_x%recv(ind_x)%pe
822 do n = 1, update_x%recv(ind_x)%count
823 dir = update_x%recv(ind_x)%dir(n)
824 if(recv(dir)) then
825 is = update_x%recv(ind_x)%is(n); ie = update_x%recv(ind_x)%ie(n)
826 js = update_x%recv(ind_x)%js(n); je = update_x%recv(ind_x)%je(n)
827 msgsize = msgsize + (ie-is+1)*(je-js+1)
828 end if
829 end do
830 ind_x = ind_x+1
831 if(ind_x .LE. nrecv_x) then
832 rank_x = update_x%recv(ind_x)%pe - domain%pe
833 if(rank_x .LE.0) rank_x = rank_x + nlist
834 else
835 rank_x = -1
836 endif
837 endif
838 if(cur_rank == rank_y) then
839 from_pe = update_y%recv(ind_y)%pe
840 do n = 1, update_y%recv(ind_y)%count
841 dir = update_y%recv(ind_y)%dir(n)
842 if(recv(dir)) then
843 is = update_y%recv(ind_y)%is(n); ie = update_y%recv(ind_y)%ie(n)
844 js = update_y%recv(ind_y)%js(n); je = update_y%recv(ind_y)%je(n)
845 msgsize = msgsize + (ie-is+1)*(je-js+1)
846 end if
847 end do
848 ind_y = ind_y+1
849 if(ind_y .LE. nrecv_y) then
850 rank_y = update_y%recv(ind_y)%pe - domain%pe
851 if(rank_y .LE.0) rank_y = rank_y + nlist
852 else
853 rank_y = -1
854 endif
855 endif
856 end select
857 cur_rank = max(rank_x, rank_y)
858 msgsize = msgsize*ke*l_size
860 if( msgsize.GT.0 )then
861 call mpp_send( buffer(buffer_pos+1), plen=msgsize, to_pe=from_pe, tag=COMM_TAG_2 )
862 buffer_pos = buffer_pos + msgsize
863 end if
864 end do
865 buffer_recv_size = buffer_pos
866 cur_rank = get_rank_send(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
868 do while (ind_x .LE. nsend_x .OR. ind_y .LE. nsend_y)
869 pos = buffer_pos
870 !--- make sure the domain stack size is big enough
871 msgsize = 0
873 if(cur_rank == rank_x) then
874 to_pe = update_x%send(ind_x)%pe
875 do n = 1, update_x%send(ind_x)%count
876 dir = update_x%send(ind_x)%dir(n)
877 if( send(dir) ) msgsize = msgsize + update_x%send(ind_x)%msgsize(n)
878 enddo
879 endif
880 if(cur_rank == rank_y) then
881 to_pe = update_y%send(ind_y)%pe
882 do n = 1, update_y%send(ind_y)%count
883 dir = update_y%send(ind_y)%dir(n)
884 if( send(dir) ) msgsize = msgsize + update_y%send(ind_y)%msgsize(n)
885 enddo
886 endif
888 select case( gridtype )
889 case(BGRID_NE, BGRID_SW, AGRID)
890 if(cur_rank == rank_x) then
891 to_pe = update_x%send(ind_x)%pe
892 ind_x = ind_x+1
893 ind_y = ind_x
894 if(ind_x .LE. nsend_x) then
895 rank_x = update_x%send(ind_x)%pe - domain%pe
896 if(rank_x .LT.0) rank_x = rank_x + nlist
897 else
898 rank_x = nlist+1
899 endif
900 rank_y = rank_x
901 endif
902 case(CGRID_NE, CGRID_SW)
903 if(cur_rank == rank_x) then
904 to_pe = update_x%send(ind_x)%pe
905 ind_x = ind_x+1
906 if(ind_x .LE. nsend_x) then
907 rank_x = update_x%send(ind_x)%pe - domain%pe
908 if(rank_x .LT.0) rank_x = rank_x + nlist
909 else
910 rank_x = nlist+1
911 endif
912 endif
913 if(cur_rank == rank_y) then
914 to_pe = update_y%send(ind_y)%pe
915 ind_y = ind_y+1
916 if(ind_y .LE. nsend_y) then
917 rank_y = update_y%send(ind_y)%pe - domain%pe
918 if(rank_y .LT.0) rank_y = rank_y + nlist
919 else
920 rank_y = nlist+1
921 endif
922 endif
923 end select
924 cur_rank = min(rank_x, rank_y)
926 if( msgsize.GT.0 )then
927 msgsize = msgsize*ke*l_size
928 end if
930 if( msgsize.GT.0 )then
931 call mpp_recv( buffer(buffer_pos+1), glen=msgsize, from_pe=to_pe, block=.false., tag=COMM_TAG_2 )
932 buffer_pos = buffer_pos + msgsize
933 end if
935 enddo
937 call mpp_sync_self(check=EVENT_RECV)
939 !--- send
940 buffer_pos = buffer_recv_size
941 pos = buffer_pos
943 cur_rank = get_rank_send(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
945 do while (ind_x .LE. nsend_x .OR. ind_y .LE. nsend_y)
946 pos = buffer_pos
947 select case( gridtype )
948 case(BGRID_NE, BGRID_SW, AGRID)
949 if(cur_rank == rank_x) then
950 to_pe = update_x%send(ind_x)%pe
951 do n = 1, update_x%send(ind_x)%count
952 dir = update_x%send(ind_x)%dir(n)
953 if( send(dir) ) then
954 tMe = update_x%send(ind_x)%tileMe(n)
956 is = update_x%send(ind_x)%is(n); ie = update_x%send(ind_x)%ie(n)
957 js = update_x%send(ind_x)%js(n); je = update_x%send(ind_x)%je(n)
959 select case( update_x%send(ind_x)%rotation(n) )
960 case(ZERO)
961 do l=1,l_size ! loop over number of fields
962 ptr_fieldx = f_addrsx(l,tMe)
963 ptr_fieldy = f_addrsy(l,tMe)
964 do k = 1,ke
965 do j = js, je
966 do i = is, ie
967 pos = pos + 2
968 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos-1)
969 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos)
970 end do
971 end do
972 end do
973 end do
974 case( MINUS_NINETY )
975 if( BTEST(update_flags,SCALAR_BIT) ) then
976 do l=1,l_size ! loop over number of fields
977 ptr_fieldx = f_addrsx(l,tMe)
978 ptr_fieldy = f_addrsy(l,tMe)
979 do k = 1,ke
980 do i = is, ie
981 do j = je, js, -1
982 pos = pos + 2
983 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos-1)
984 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
985 end do
986 end do
987 end do
988 end do
989 else
990 do l=1,l_size ! loop over number of fields
991 ptr_fieldx = f_addrsx(l,tMe)
992 ptr_fieldy = f_addrsy(l,tMe)
993 do k = 1,ke
994 do i = is, ie
995 do j = je, js, -1
996 pos = pos + 2
997 fieldy(i,j,k)=fieldy(i,j,k)-buffer(pos-1)
998 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
999 end do
1000 end do
1001 end do
1002 end do
1003 end if
1004 case( NINETY )
1005 if( BTEST(update_flags,SCALAR_BIT) ) then
1006 do l=1,l_size ! loop over number of fields
1007 ptr_fieldx = f_addrsx(l,tMe)
1008 ptr_fieldy = f_addrsy(l,tMe)
1009 do k = 1,ke
1010 do i = ie, is, -1
1011 do j = js, je
1012 pos = pos + 2
1013 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos-1)
1014 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
1015 end do
1016 end do
1017 end do
1018 end do
1019 else
1020 do l=1,l_size ! loop over number of fields
1021 ptr_fieldx = f_addrsx(l,tMe)
1022 ptr_fieldy = f_addrsy(l,tMe)
1023 do k = 1,ke
1024 do i = ie, is, -1
1025 do j = js, je
1026 pos = pos + 2
1027 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos-1)
1028 fieldx(i,j,k)=fieldx(i,j,k)-buffer(pos)
1029 end do
1030 end do
1031 end do
1032 end do
1033 end if
1034 case( ONE_HUNDRED_EIGHTY )
1035 if( BTEST(update_flags,SCALAR_BIT) ) then
1036 do l=1,l_size ! loop over number of fields
1037 ptr_fieldx = f_addrsx(l,tMe)
1038 ptr_fieldy = f_addrsy(l,tMe)
1039 do k = 1,ke
1040 do j = je, js, -1
1041 do i = ie, is, -1
1042 pos = pos + 2
1043 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos-1)
1044 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos)
1045 end do
1046 end do
1047 end do
1048 end do
1049 else
1050 do l=1,l_size ! loop over number of fields
1051 ptr_fieldx = f_addrsx(l,tMe)
1052 ptr_fieldy = f_addrsy(l,tMe)
1053 do k = 1,ke
1054 do j = je, js, -1
1055 do i = ie, is, -1
1056 pos = pos + 2
1057 fieldx(i,j,k)=fieldx(i,j,k)-buffer(pos-1)
1058 fieldy(i,j,k)=fieldy(i,j,k)-buffer(pos)
1059 end do
1060 end do
1061 end do
1062 end do
1063 end if
1064 end select ! select case( rotation(n) )
1065 end if ! if( send(dir) )
1066 end do ! do n = 1, update_x%send(ind_x)%count
1067 ind_x = ind_x+1
1068 ind_y = ind_x
1069 if(ind_x .LE. nsend_x) then
1070 rank_x = update_x%send(ind_x)%pe - domain%pe
1071 if(rank_x .LT.0) rank_x = rank_x + nlist
1072 else
1073 rank_x = nlist+1
1074 endif
1075 rank_y = rank_x
1076 endif
1077 case(CGRID_NE, CGRID_SW)
1078 if(cur_rank == rank_x) then
1079 to_pe = update_x%send(ind_x)%pe
1080 do n = 1, update_x%send(ind_x)%count
1081 dir = update_x%send(ind_x)%dir(n)
1082 if( send(dir) ) then
1083 tMe = update_x%send(ind_x)%tileMe(n)
1084 is = update_x%send(ind_x)%is(n); ie = update_x%send(ind_x)%ie(n)
1085 js = update_x%send(ind_x)%js(n); je = update_x%send(ind_x)%je(n)
1086 select case( update_x%send(ind_x)%rotation(n) )
1087 case(ZERO)
1088 do l=1,l_size ! loop over number of fields
1089 ptr_fieldx = f_addrsx(l, tMe)
1090 ptr_fieldy = f_addrsy(l, tMe)
1091 do k = 1,ke
1092 do j = js, je
1093 do i = is, ie
1094 pos = pos + 1
1095 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
1096 end do
1097 end do
1098 end do
1099 end do
1100 case(MINUS_NINETY)
1101 if( BTEST(update_flags,SCALAR_BIT) ) then
1102 do l=1,l_size ! loop over number of fields
1103 ptr_fieldx = f_addrsx(l, tMe)
1104 ptr_fieldy = f_addrsy(l, tMe)
1105 do k = 1,ke
1106 do i = is, ie
1107 do j = je, js, -1
1108 pos = pos + 1
1109 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos)
1110 end do
1111 end do
1112 end do
1113 end do
1114 else
1115 do l=1,l_size ! loop over number of fields
1116 ptr_fieldx = f_addrsx(l, tMe)
1117 ptr_fieldy = f_addrsy(l, tMe)
1118 do k = 1,ke
1119 do i = is, ie
1120 do j = je, js, -1
1121 pos = pos + 1
1122 fieldy(i,j,k)=fieldy(i,j,k)-buffer(pos)
1123 end do
1124 end do
1125 end do
1126 end do
1127 end if
1128 case(NINETY)
1129 do l=1,l_size ! loop over number of fields
1130 ptr_fieldx = f_addrsx(l, tMe)
1131 ptr_fieldy = f_addrsy(l, tMe)
1132 do k = 1, ke
1133 do i = ie, is, -1
1134 do j = js, je
1135 pos = pos + 1
1136 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos)
1137 end do
1138 end do
1139 end do
1140 end do
1141 case(ONE_HUNDRED_EIGHTY)
1142 if( BTEST(update_flags,SCALAR_BIT) ) then
1143 do l=1,l_size ! loop over number of fields
1144 ptr_fieldx = f_addrsx(l, tMe)
1145 ptr_fieldy = f_addrsy(l, tMe)
1146 do k = 1,ke
1147 do j = je, js, -1
1148 do i = ie, is, -1
1149 pos = pos + 1
1150 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
1151 end do
1152 end do
1153 end do
1154 end do
1155 else
1156 do l=1,l_size ! loop over number of fields
1157 ptr_fieldx = f_addrsx(l, tMe)
1158 ptr_fieldy = f_addrsy(l, tMe)
1159 do k = 1,ke
1160 do j = je, js, -1
1161 do i = ie, is, -1
1162 pos = pos + 1
1163 fieldx(i,j,k)=fieldx(i,j,k)-buffer(pos)
1164 end do
1165 end do
1166 end do
1167 end do
1168 end if
1169 end select
1170 end if
1171 end do
1172 ind_x = ind_x+1
1173 if(ind_x .LE. nsend_x) then
1174 rank_x = update_x%send(ind_x)%pe - domain%pe
1175 if(rank_x .LT.0) rank_x = rank_x + nlist
1176 else
1177 rank_x = nlist+1
1178 endif
1179 endif
1180 if(cur_rank == rank_y) then
1181 to_pe = update_y%send(ind_y)%pe
1182 do n = 1, update_y%send(ind_y)%count
1183 dir = update_y%send(ind_y)%dir(n)
1184 if( send(dir) ) then
1185 tMe = update_y%send(ind_y)%tileMe(n)
1186 is = update_y%send(ind_y)%is(n); ie = update_y%send(ind_y)%ie(n)
1187 js = update_y%send(ind_y)%js(n); je = update_y%send(ind_y)%je(n)
1189 select case( update_y%send(ind_y)%rotation(n) )
1190 case(ZERO)
1191 do l=1,l_size ! loop over number of fields
1192 ptr_fieldx = f_addrsx(l, tMe)
1193 ptr_fieldy = f_addrsy(l, tMe)
1194 do k = 1,ke
1195 do j = js, je
1196 do i = is, ie
1197 pos = pos + 1
1198 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos)
1199 end do
1200 end do
1201 end do
1202 end do
1203 case(MINUS_NINETY)
1204 do l=1,l_size ! loop over number of fields
1205 ptr_fieldx = f_addrsx(l, tMe)
1206 ptr_fieldy = f_addrsy(l, tMe)
1207 do k = 1,ke
1208 do i = is, ie
1209 do j = je, js, -1
1210 pos = pos + 1
1211 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
1212 end do
1213 end do
1214 end do
1215 end do
1216 case(NINETY)
1217 if( BTEST(update_flags,SCALAR_BIT) ) then
1218 do l=1,l_size ! loop over number of fields
1219 ptr_fieldx = f_addrsx(l, tMe)
1220 ptr_fieldy = f_addrsy(l, tMe)
1221 do k = 1,ke
1222 do i = ie, is, -1
1223 do j = js, je
1224 pos = pos + 1
1225 fieldx(i,j,k)=fieldx(i,j,k)+buffer(pos)
1226 end do
1227 end do
1228 end do
1229 end do
1230 else
1231 do l=1,l_size ! loop over number of fields
1232 ptr_fieldx = f_addrsx(l, tMe)
1233 ptr_fieldy = f_addrsy(l, tMe)
1234 do k = 1,ke
1235 do i = ie, is, -1
1236 do j = js, je
1237 pos = pos + 1
1238 fieldx(i,j,k)=fieldx(i,j,k)-buffer(pos)
1239 end do
1240 end do
1241 end do
1242 end do
1243 end if
1244 case(ONE_HUNDRED_EIGHTY)
1245 if( BTEST(update_flags,SCALAR_BIT) ) then
1246 do l=1,l_size ! loop over number of fields
1247 ptr_fieldx = f_addrsx(l, tMe)
1248 ptr_fieldy = f_addrsy(l, tMe)
1249 do k = 1,ke
1250 do j = je, js, -1
1251 do i = ie, is, -1
1252 pos = pos + 1
1253 fieldy(i,j,k)=fieldy(i,j,k)+buffer(pos)
1254 end do
1255 end do
1256 end do
1257 end do
1258 else
1259 do l=1,l_size ! loop over number of fields
1260 ptr_fieldx = f_addrsx(l, tMe)
1261 ptr_fieldy = f_addrsy(l, tMe)
1262 do k = 1,ke
1263 do j = je, js, -1
1264 do i = ie, is, -1
1265 pos = pos + 1
1266 fieldy(i,j,k)=fieldy(i,j,k)-buffer(pos)
1267 end do
1268 end do
1269 end do
1270 end do
1271 end if
1272 end select
1273 endif
1274 enddo
1275 ind_y = ind_y+1
1276 if(ind_y .LE. nsend_y) then
1277 rank_y = update_y%send(ind_y)%pe - domain%pe
1278 if(rank_y .LT.0) rank_y = rank_y + nlist
1279 else
1280 rank_y = nlist+1
1281 endif
1282 endif
1283 end select
1284 cur_rank = min(rank_x, rank_y)
1285 msgsize = pos - buffer_pos
1286 if( msgsize.GT.0 )then
1287 buffer_pos = pos
1288 end if
1289 end do
1291 call mpp_sync_self( )
1293 return
1295 end subroutine MPP_DO_UPDATE_AD_3D_V_