2 !***********************************************************************
3 !* GNU Lesser General Public License
5 !* This file is part of the GFDL Flexible Modeling
System (FMS
).
7 !* FMS is free software
: you can redistribute it
and/or modify it under
8 !* the terms of the GNU Lesser General Public License as published by
9 !* the Free Software Foundation
, either version
3 of the License
, or (at
10 !* your option
) any later version
.
12 !* FMS is distributed in the hope that it will be useful
, but WITHOUT
13 !* ANY WARRANTY
; without even the implied warranty of MERCHANTABILITY
or
14 !* FITNESS FOR A PARTICULAR PURPOSE
. See the GNU General Public License
17 !* You should have received a copy of the GNU Lesser General Public
18 !* License along with FMS
. If
not, see
<http
://www.gnu.org/licenses/>.
19 !***********************************************************************
20 subroutine
MPP_DO_UPDATE_3D_V_(f_addrsx
,f_addrsy
, domain
, update_x
, update_y
, &
21 d_type
, ke
, gridtype
, flags
)
22 !updates data domain of
3D field whose computational domains have been computed
23 integer(i8_kind
), intent(in
) :: f_addrsx(:,:), f_addrsy(:,:)
24 type(domain2d
), intent(in
) :: domain
25 type(overlapSpec
), intent(in
) :: update_x
, update_y
26 integer
, intent(in
) :: ke
27 MPP_TYPE_
, intent(in
) :: d_type
! creates unique interface
28 integer
, intent(in
) :: gridtype
29 integer
, intent(in
), optional :: flags
31 MPP_TYPE_ :: fieldx(update_x
%xbegin
:update_x
%xend
, update_x
%ybegin
:update_x
%yend
,ke
)
32 MPP_TYPE_ :: fieldy(update_y
%xbegin
:update_y
%xend
, update_y
%ybegin
:update_y
%yend
,ke
)
33 pointer(ptr_fieldx
, fieldx
)
34 pointer(ptr_fieldy
, fieldy
)
36 integer :: update_flags
37 integer :: l_size
, l
, i
, j
, k
, is
, ie
, js
, je
, n
, m
38 integer :: pos
, nlist
, msgsize
, isd
, ied
, jsd
, jed
39 integer :: to_pe
, from_pe
, midpoint
42 integer :: send_start_pos
, nsend
43 integer :: send_msgsize(2*MAXLIST
)
44 integer :: send_pe(2*MAXLIST
)
45 integer
, allocatable :: msg1(:), msg2(:), msg3(:)
46 logical :: send(8), recv(8), update_edge_only
47 MPP_TYPE_ :: buffer(size(mpp_domains_stack(:)))
50 character(len
=8) :: text
51 integer :: buffer_recv_size
, shift
52 integer :: rank_x
, rank_y
, ind_x
, ind_y
, cur_rank
53 integer :: nsend_x
, nsend_y
, nrecv_x
, nrecv_y
, outunit
56 update_flags
= XUPDATE
+YUPDATE
!default
57 if( PRESENT(flags
) ) then
59 ! The following test is so that SCALAR_PAIR can be used alone with the
60 ! same
default update pattern as without
.
61 if (BTEST(update_flags
,SCALAR_BIT
)) then
62 if (.NOT
.(BTEST(update_flags
,WEST
) .OR
. BTEST(update_flags
,EAST
) &
63 .OR
. BTEST(update_flags
,NORTH
) .OR
. BTEST(update_flags
,SOUTH
))) &
64 update_flags
= update_flags
+ XUPDATE
+YUPDATE
!default with SCALAR_PAIR
68 if( BTEST(update_flags
,NORTH
) .AND
. BTEST(domain
%fold
,NORTH
) .AND
. BTEST(gridtype
,SOUTH
) ) &
69 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: Incompatible grid offset and fold.' )
71 update_edge_only
= BTEST(update_flags
, EDGEONLY
)
73 recv(1) = BTEST(update_flags
,EAST
)
74 recv(3) = BTEST(update_flags
,SOUTH
)
75 recv(5) = BTEST(update_flags
,WEST
)
76 recv(7) = BTEST(update_flags
,NORTH
)
77 if( update_edge_only
) then
78 if( .NOT
. (recv(1) .OR
. recv(3) .OR
. recv(5) .OR
. recv(7)) ) then
85 recv(2) = recv(1) .AND
. recv(3)
86 recv(4) = recv(3) .AND
. recv(5)
87 recv(6) = recv(5) .AND
. recv(7)
88 recv(8) = recv(7) .AND
. recv(1)
93 l_size
= size(f_addrsx
,1)
94 nlist
= size(domain
%list(:))
95 ptr
= LOC(mpp_domains_stack
)
98 nsend_x
= update_x
%nsend
99 nsend_y
= update_y
%nsend
100 nrecv_x
= update_x
%nrecv
101 nrecv_y
= update_y
%nrecv
103 if(debug_message_passing
) then
104 allocate(msg1(0:nlist
-1), msg2(0:nlist
-1), msg3(0:nlist
-1) )
108 cur_rank
= get_rank_recv(domain
, update_x
, update_y
, rank_x
, rank_y
, ind_x
, ind_y
)
110 do while (ind_x
.LE
. nrecv_x
.OR
. ind_y
.LE
. nrecv_y
)
112 if(cur_rank
== rank_x
) then
113 from_pe
= update_x
%recv(ind_x
)%pe
114 do n
= 1, update_x
%recv(ind_x
)%count
115 dir
= update_x
%recv(ind_x
)%dir(n
)
117 is
= update_x
%recv(ind_x
)%is(n
); ie
= update_x
%recv(ind_x
)%ie(n
)
118 js
= update_x
%recv(ind_x
)%js(n
); je
= update_x
%recv(ind_x
)%je(n
)
119 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
123 if(ind_x
.LE
. nrecv_x
) then
124 rank_x
= update_x
%recv(ind_x
)%pe
- domain
%pe
125 if(rank_x
.LE
.0) rank_x
= rank_x
+ nlist
130 if(cur_rank
== rank_y
) then
131 from_pe
= update_y
%recv(ind_y
)%pe
132 do n
= 1, update_y
%recv(ind_y
)%count
133 dir
= update_y
%recv(ind_y
)%dir(n
)
135 is
= update_y
%recv(ind_y
)%is(n
); ie
= update_y
%recv(ind_y
)%ie(n
)
136 js
= update_y
%recv(ind_y
)%js(n
); je
= update_y
%recv(ind_y
)%je(n
)
137 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
141 if(ind_y
.LE
. nrecv_y
) then
142 rank_y
= update_y
%recv(ind_y
)%pe
- domain
%pe
143 if(rank_y
.LE
.0) rank_y
= rank_y
+ nlist
148 cur_rank
= max(rank_x
, rank_y
)
149 m
= from_pe
-mpp_root_pe()
153 cur_rank
= get_rank_send(domain
, update_x
, update_y
, rank_x
, rank_y
, ind_x
, ind_y
)
154 do while (ind_x
.LE
. nsend_x
.OR
. ind_y
.LE
. nsend_y
)
156 if(cur_rank
== rank_x
) then
157 to_pe
= update_x
%send(ind_x
)%pe
158 do n
= 1, update_x
%send(ind_x
)%count
159 dir
= update_x
%send(ind_x
)%dir(n
)
161 is
= update_x
%send(ind_x
)%is(n
); ie
= update_x
%send(ind_x
)%ie(n
)
162 js
= update_x
%send(ind_x
)%js(n
); je
= update_x
%send(ind_x
)%je(n
)
163 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
167 if(ind_x
.LE
. nsend_x
) then
168 rank_x
= update_x
%send(ind_x
)%pe
- domain
%pe
169 if(rank_x
.LT
.0) rank_x
= rank_x
+ nlist
174 if(cur_rank
== rank_y
) then
175 to_pe
= update_y
%send(ind_y
)%pe
176 do n
= 1, update_y
%send(ind_y
)%count
177 dir
= update_y
%send(ind_y
)%dir(n
)
179 is
= update_y
%send(ind_y
)%is(n
); ie
= update_y
%send(ind_y
)%ie(n
)
180 js
= update_y
%send(ind_y
)%js(n
); je
= update_y
%send(ind_y
)%je(n
)
181 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
185 if(ind_y
.LE
. nsend_y
) then
186 rank_y
= update_y
%send(ind_y
)%pe
- domain
%pe
187 if(rank_y
.LT
.0) rank_y
= rank_y
+ nlist
192 m
= to_pe
-mpp_root_pe()
194 cur_rank
= min(rank_x
, rank_y
)
196 call
mpp_alltoall(msg3
, 1, msg1
, 1)
197 ! call
mpp_sync_self(check
=EVENT_RECV
)
199 if(msg1(m
) .NE
. msg2(m
)) then
200 print
*, "My pe = ", mpp_pe(), ",domain name =", trim(domain
%name
), ",from pe=", &
201 domain
%list(m
)%pe
, ":send size = ", msg1(m
), ", recv size = ", msg2(m
)
202 call
mpp_error(FATAL
, "mpp_do_updateV: mismatch on send and recv size")
206 ! call
mpp_sync_self()
207 write(outunit
,*)"NOTE from mpp_do_updateV: message sizes are matched between send and recv for domain " &
209 deallocate(msg1
, msg2
, msg3
)
214 cur_rank
= get_rank_recv(domain
, update_x
, update_y
, rank_x
, rank_y
, ind_x
, ind_y
)
215 call
mpp_clock_begin(recv_clock
)
216 do while (ind_x
.LE
. nrecv_x
.OR
. ind_y
.LE
. nrecv_y
)
218 select
case(gridtype
)
219 case(BGRID_NE
, BGRID_SW
, AGRID
)
220 if(cur_rank
== rank_x
) then
221 from_pe
= update_x
%recv(ind_x
)%pe
222 do n
= 1, update_x
%recv(ind_x
)%count
223 dir
= update_x
%recv(ind_x
)%dir(n
)
225 tMe
= update_x
%recv(ind_x
)%tileMe(n
)
226 is
= update_x
%recv(ind_x
)%is(n
); ie
= update_x
%recv(ind_x
)%ie(n
)
227 js
= update_x
%recv(ind_x
)%js(n
); je
= update_x
%recv(ind_x
)%je(n
)
228 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
234 if(ind_x
.LE
. nrecv_x
) then
235 rank_x
= update_x
%recv(ind_x
)%pe
- domain
%pe
236 if(rank_x
.LE
.0) rank_x
= rank_x
+ nlist
242 case(CGRID_NE
, CGRID_SW
)
243 if(cur_rank
== rank_x
) then
244 from_pe
= update_x
%recv(ind_x
)%pe
245 do n
= 1, update_x
%recv(ind_x
)%count
246 dir
= update_x
%recv(ind_x
)%dir(n
)
248 is
= update_x
%recv(ind_x
)%is(n
); ie
= update_x
%recv(ind_x
)%ie(n
)
249 js
= update_x
%recv(ind_x
)%js(n
); je
= update_x
%recv(ind_x
)%je(n
)
250 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
254 if(ind_x
.LE
. nrecv_x
) then
255 rank_x
= update_x
%recv(ind_x
)%pe
- domain
%pe
256 if(rank_x
.LE
.0) rank_x
= rank_x
+ nlist
261 if(cur_rank
== rank_y
) then
262 from_pe
= update_y
%recv(ind_y
)%pe
263 do n
= 1, update_y
%recv(ind_y
)%count
264 dir
= update_y
%recv(ind_y
)%dir(n
)
266 is
= update_y
%recv(ind_y
)%is(n
); ie
= update_y
%recv(ind_y
)%ie(n
)
267 js
= update_y
%recv(ind_y
)%js(n
); je
= update_y
%recv(ind_y
)%je(n
)
268 msgsize
= msgsize
+ (ie
-is
+1)*(je
-js
+1)
272 if(ind_y
.LE
. nrecv_y
) then
273 rank_y
= update_y
%recv(ind_y
)%pe
- domain
%pe
274 if(rank_y
.LE
.0) rank_y
= rank_y
+ nlist
280 cur_rank
= max(rank_x
, rank_y
)
281 msgsize
= msgsize
*ke
*l_size
283 if( msgsize
.GT
.0 )then
284 mpp_domains_stack_hwm
= max( mpp_domains_stack_hwm
, buffer_pos
+msgsize
)
285 if( mpp_domains_stack_hwm
.GT
.mpp_domains_stack_size
)then
286 write( text
,'(i8)' )mpp_domains_stack_hwm
287 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: mpp_domains_stack overflow, '// &
288 'call mpp_domains_set_stack_size('//trim(text)//') from all PEs.' )
290 call
mpp_recv( buffer(buffer_pos
+1), glen
=msgsize
, from_pe
=from_pe
, block
=.false., tag
=COMM_TAG_2
)
291 buffer_pos
= buffer_pos
+ msgsize
294 call
mpp_clock_end(recv_clock
)
295 buffer_recv_size
= buffer_pos
296 send_start_pos
= buffer_pos
299 cur_rank
= get_rank_send(domain
, update_x
, update_y
, rank_x
, rank_y
, ind_x
, ind_y
)
301 call
mpp_clock_begin(pack_clock
)
302 do while (ind_x
.LE
. nsend_x
.OR
. ind_y
.LE
. nsend_y
)
304 !--- make sure the domain stack size is big enough
306 if(cur_rank
== rank_x
) then
307 do n
= 1, update_x
%send(ind_x
)%count
308 dir
= update_x
%send(ind_x
)%dir(n
)
309 if( send(dir
) ) msgsize
= msgsize
+ update_x
%send(ind_x
)%msgsize(n
)
312 if(cur_rank
== rank_y
) then
313 do n
= 1, update_y
%send(ind_y
)%count
314 dir
= update_y
%send(ind_y
)%dir(n
)
315 if( send(dir
) ) msgsize
= msgsize
+ update_y
%send(ind_y
)%msgsize(n
)
319 if( msgsize
.GT
.0 )then
320 msgsize
= msgsize
*ke
*l_size
321 mpp_domains_stack_hwm
= max( mpp_domains_stack_hwm
, pos
+msgsize
)
322 if( mpp_domains_stack_hwm
.GT
.mpp_domains_stack_size
)then
323 write( text
,'(i8)' )mpp_domains_stack_hwm
324 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: mpp_domains_stack overflow, ' // &
325 'call mpp_domains_set_stack_size('//trim(text)//') from all PEs.')
328 select
case( gridtype
)
329 case(BGRID_NE
, BGRID_SW
, AGRID
)
330 if(cur_rank
== rank_x
) then
331 to_pe
= update_x
%send(ind_x
)%pe
332 do n
= 1, update_x
%send(ind_x
)%count
333 dir
= update_x
%send(ind_x
)%dir(n
)
335 tMe
= update_x
%send(ind_x
)%tileMe(n
)
337 is
= update_x
%send(ind_x
)%is(n
); ie
= update_x
%send(ind_x
)%ie(n
)
338 js
= update_x
%send(ind_x
)%js(n
); je
= update_x
%send(ind_x
)%je(n
)
339 select
case( update_x
%send(ind_x
)%rotation(n
) )
341 do l
=1,l_size
! loop over number of fields
342 ptr_fieldx
= f_addrsx(l
,tMe
)
343 ptr_fieldy
= f_addrsy(l
,tMe
)
348 buffer(pos
-1) = fieldx(i
,j
,k
)
349 buffer(pos
) = fieldy(i
,j
,k
)
355 if( BTEST(update_flags
,SCALAR_BIT
) ) then
356 do l
=1,l_size
! loop over number of fields
357 ptr_fieldx
= f_addrsx(l
,tMe
)
358 ptr_fieldy
= f_addrsy(l
,tMe
)
363 buffer(pos
-1) = fieldy(i
,j
,k
)
364 buffer(pos
) = fieldx(i
,j
,k
)
370 do l
=1,l_size
! loop over number of fields
371 ptr_fieldx
= f_addrsx(l
,tMe
)
372 ptr_fieldy
= f_addrsy(l
,tMe
)
377 buffer(pos
-1) = -fieldy(i
,j
,k
)
378 buffer(pos
) = fieldx(i
,j
,k
)
385 if( BTEST(update_flags
,SCALAR_BIT
) ) then
386 do l
=1,l_size
! loop over number of fields
387 ptr_fieldx
= f_addrsx(l
,tMe
)
388 ptr_fieldy
= f_addrsy(l
,tMe
)
393 buffer(pos
-1) = fieldy(i
,j
,k
)
394 buffer(pos
) = fieldx(i
,j
,k
)
400 do l
=1,l_size
! loop over number of fields
401 ptr_fieldx
= f_addrsx(l
,tMe
)
402 ptr_fieldy
= f_addrsy(l
,tMe
)
407 buffer(pos
-1) = fieldy(i
,j
,k
)
408 buffer(pos
) = -fieldx(i
,j
,k
)
414 case( ONE_HUNDRED_EIGHTY
)
415 if( BTEST(update_flags
,SCALAR_BIT
) ) then
416 do l
=1,l_size
! loop over number of fields
417 ptr_fieldx
= f_addrsx(l
,tMe
)
418 ptr_fieldy
= f_addrsy(l
,tMe
)
423 buffer(pos
-1) = fieldx(i
,j
,k
)
424 buffer(pos
) = fieldy(i
,j
,k
)
430 do l
=1,l_size
! loop over number of fields
431 ptr_fieldx
= f_addrsx(l
,tMe
)
432 ptr_fieldy
= f_addrsy(l
,tMe
)
437 buffer(pos
-1) = -fieldx(i
,j
,k
)
438 buffer(pos
) = -fieldy(i
,j
,k
)
444 end select
! select
case( rotation(n
) )
445 end
if ! if( send(dir
) )
446 end
do ! do n
= 1, update_x
%send(ind_x
)%count
449 if(ind_x
.LE
. nsend_x
) then
450 rank_x
= update_x
%send(ind_x
)%pe
- domain
%pe
451 if(rank_x
.LT
.0) rank_x
= rank_x
+ nlist
457 case(CGRID_NE
, CGRID_SW
)
458 if(cur_rank
== rank_x
) then
459 to_pe
= update_x
%send(ind_x
)%pe
460 do n
= 1, update_x
%send(ind_x
)%count
461 dir
= update_x
%send(ind_x
)%dir(n
)
463 tMe
= update_x
%send(ind_x
)%tileMe(n
)
464 is
= update_x
%send(ind_x
)%is(n
); ie
= update_x
%send(ind_x
)%ie(n
)
465 js
= update_x
%send(ind_x
)%js(n
); je
= update_x
%send(ind_x
)%je(n
)
466 select
case( update_x
%send(ind_x
)%rotation(n
) )
468 do l
=1,l_size
! loop over number of fields
469 ptr_fieldx
= f_addrsx(l
, tMe
)
470 ptr_fieldy
= f_addrsy(l
, tMe
)
475 buffer(pos
) = fieldx(i
,j
,k
)
481 if( BTEST(update_flags
,SCALAR_BIT
) ) then
482 do l
=1,l_size
! loop over number of fields
483 ptr_fieldx
= f_addrsx(l
, tMe
)
484 ptr_fieldy
= f_addrsy(l
, tMe
)
489 buffer(pos
) = fieldy(i
,j
,k
)
495 do l
=1,l_size
! loop over number of fields
496 ptr_fieldx
= f_addrsx(l
, tMe
)
497 ptr_fieldy
= f_addrsy(l
, tMe
)
502 buffer(pos
) = -fieldy(i
,j
,k
)
509 do l
=1,l_size
! loop over number of fields
510 ptr_fieldx
= f_addrsx(l
, tMe
)
511 ptr_fieldy
= f_addrsy(l
, tMe
)
516 buffer(pos
) = fieldy(i
,j
,k
)
521 case(ONE_HUNDRED_EIGHTY
)
522 if( BTEST(update_flags
,SCALAR_BIT
) ) then
523 do l
=1,l_size
! loop over number of fields
524 ptr_fieldx
= f_addrsx(l
, tMe
)
525 ptr_fieldy
= f_addrsy(l
, tMe
)
530 buffer(pos
) = fieldx(i
,j
,k
)
536 do l
=1,l_size
! loop over number of fields
537 ptr_fieldx
= f_addrsx(l
, tMe
)
538 ptr_fieldy
= f_addrsy(l
, tMe
)
543 buffer(pos
) = -fieldx(i
,j
,k
)
553 if(ind_x
.LE
. nsend_x
) then
554 rank_x
= update_x
%send(ind_x
)%pe
- domain
%pe
555 if(rank_x
.LT
.0) rank_x
= rank_x
+ nlist
560 if(cur_rank
== rank_y
) then
561 to_pe
= update_y
%send(ind_y
)%pe
562 do n
= 1, update_y
%send(ind_y
)%count
563 dir
= update_y
%send(ind_y
)%dir(n
)
565 tMe
= update_y
%send(ind_y
)%tileMe(n
)
566 is
= update_y
%send(ind_y
)%is(n
); ie
= update_y
%send(ind_y
)%ie(n
)
567 js
= update_y
%send(ind_y
)%js(n
); je
= update_y
%send(ind_y
)%je(n
)
568 select
case( update_y
%send(ind_y
)%rotation(n
) )
570 do l
=1,l_size
! loop over number of fields
571 ptr_fieldx
= f_addrsx(l
, tMe
)
572 ptr_fieldy
= f_addrsy(l
, tMe
)
577 buffer(pos
) = fieldy(i
,j
,k
)
583 do l
=1,l_size
! loop over number of fields
584 ptr_fieldx
= f_addrsx(l
, tMe
)
585 ptr_fieldy
= f_addrsy(l
, tMe
)
590 buffer(pos
) = fieldx(i
,j
,k
)
596 if( BTEST(update_flags
,SCALAR_BIT
) ) then
597 do l
=1,l_size
! loop over number of fields
598 ptr_fieldx
= f_addrsx(l
, tMe
)
599 ptr_fieldy
= f_addrsy(l
, tMe
)
604 buffer(pos
) = fieldx(i
,j
,k
)
610 do l
=1,l_size
! loop over number of fields
611 ptr_fieldx
= f_addrsx(l
, tMe
)
612 ptr_fieldy
= f_addrsy(l
, tMe
)
617 buffer(pos
) = -fieldx(i
,j
,k
)
623 case(ONE_HUNDRED_EIGHTY
)
624 if( BTEST(update_flags
,SCALAR_BIT
) ) then
625 do l
=1,l_size
! loop over number of fields
626 ptr_fieldx
= f_addrsx(l
, tMe
)
627 ptr_fieldy
= f_addrsy(l
, tMe
)
632 buffer(pos
) = fieldy(i
,j
,k
)
638 do l
=1,l_size
! loop over number of fields
639 ptr_fieldx
= f_addrsx(l
, tMe
)
640 ptr_fieldy
= f_addrsy(l
, tMe
)
645 buffer(pos
) = -fieldy(i
,j
,k
)
655 if(ind_y
.LE
. nsend_y
) then
656 rank_y
= update_y
%send(ind_y
)%pe
- domain
%pe
657 if(rank_y
.LT
.0) rank_y
= rank_y
+ nlist
663 cur_rank
= min(rank_x
, rank_y
)
665 send_pe(nsend
) = to_pe
666 send_msgsize(nsend
) = pos
- buffer_pos
670 buffer_pos
= send_start_pos
671 call
mpp_clock_end(pack_clock
)
672 call
mpp_clock_begin(send_clock
)
674 msgsize
= send_msgsize(m
)
675 if( msgsize
.GT
.0 )then
676 call
mpp_send( buffer(buffer_pos
+1), plen
=msgsize
, to_pe
=send_pe(m
), tag
=COMM_TAG_2
)
677 buffer_pos
= buffer_pos
+ msgsize
680 call
mpp_clock_end(send_clock
)
683 !unpack halos in reverse order
684 call
mpp_clock_begin(wait_clock
)
685 call
mpp_sync_self(check
=EVENT_RECV
)
686 call
mpp_clock_end(wait_clock
)
687 buffer_pos
= buffer_recv_size
688 cur_rank
= get_rank_unpack(domain
, update_x
, update_y
, rank_x
, rank_y
, ind_x
, ind_y
)
690 call
mpp_clock_begin(unpk_clock
)
691 do while (ind_x
> 0 .OR
. ind_y
> 0)
693 select
case ( gridtype
)
694 case(BGRID_NE
, BGRID_SW
, AGRID
)
695 if(cur_rank
== rank_x
) then
696 do n
= update_x
%recv(ind_x
)%count
, 1, -1
697 dir
= update_x
%recv(ind_x
)%dir(n
)
699 tMe
= update_x
%recv(ind_x
)%tileMe(n
)
700 is
= update_x
%recv(ind_x
)%is(n
); ie
= update_x
%recv(ind_x
)%ie(n
)
701 js
= update_x
%recv(ind_x
)%js(n
); je
= update_x
%recv(ind_x
)%je(n
)
702 msgsize
= (ie
-is
+1)*(je
-js
+1)*ke
*2*l_size
703 pos
= buffer_pos
- msgsize
705 do l
=1, l_size
! loop over number of fields
706 ptr_fieldx
= f_addrsx(l
, tMe
)
707 ptr_fieldy
= f_addrsy(l
, tMe
)
712 fieldx(i
,j
,k
) = buffer(pos
-1)
713 fieldy(i
,j
,k
) = buffer(pos
)
718 end
if ! end
if( recv(dir
) )
719 end
do ! do dir
=8,1,-1
722 if(ind_x
.GT
. 0) then
723 rank_x
= update_x
%recv(ind_x
)%pe
- domain
%pe
724 if(rank_x
.LE
.0) rank_x
= rank_x
+ nlist
730 case(CGRID_NE
, CGRID_SW
)
731 if(cur_rank
== rank_y
) then
732 do n
= update_y
%recv(ind_y
)%count
, 1, -1
733 dir
= update_y
%recv(ind_y
)%dir(n
)
735 tMe
= update_y
%recv(ind_y
)%tileMe(n
)
736 is
= update_y
%recv(ind_y
)%is(n
); ie
= update_y
%recv(ind_y
)%ie(n
)
737 js
= update_y
%recv(ind_y
)%js(n
); je
= update_y
%recv(ind_y
)%je(n
)
738 msgsize
= (ie
-is
+1)*(je
-js
+1)*ke
*l_size
739 pos
= buffer_pos
- msgsize
741 do l
=1,l_size
! loop over number of fields
742 ptr_fieldx
= f_addrsx(l
, tMe
)
743 ptr_fieldy
= f_addrsy(l
, tMe
)
748 fieldy(i
,j
,k
) = buffer(pos
)
756 if(ind_y
.GT
. 0) then
757 rank_y
= update_y
%recv(ind_y
)%pe
- domain
%pe
758 if(rank_y
.LE
.0) rank_y
= rank_y
+ nlist
763 if(cur_rank
== rank_x
) then
764 do n
= update_x
%recv(ind_x
)%count
, 1, -1
765 dir
= update_x
%recv(ind_x
)%dir(n
)
767 tMe
= update_x
%recv(ind_x
)%tileMe(n
)
768 is
= update_x
%recv(ind_x
)%is(n
); ie
= update_x
%recv(ind_x
)%ie(n
)
769 js
= update_x
%recv(ind_x
)%js(n
); je
= update_x
%recv(ind_x
)%je(n
)
770 msgsize
= (ie
-is
+1)*(je
-js
+1)*ke
*l_size
771 pos
= buffer_pos
- msgsize
773 do l
=1,l_size
! loop over number of fields
774 ptr_fieldx
= f_addrsx(l
, tMe
)
775 ptr_fieldy
= f_addrsy(l
, tMe
)
780 fieldx(i
,j
,k
) = buffer(pos
)
788 if(ind_x
.GT
. 0) then
789 rank_x
= update_x
%recv(ind_x
)%pe
- domain
%pe
790 if(rank_x
.LE
.0) rank_x
= rank_x
+ nlist
796 cur_rank
= min(rank_x
, rank_y
)
798 call
mpp_clock_end(unpk_clock
)
800 ! ---northern boundary fold
802 if(domain
%symmetry
) shift
= 1
803 if( BTEST(domain
%fold
,NORTH
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
804 isd
= domain
%x(1)%compute
%begin
- update_x
%whalo
;
805 ied
= domain
%x(1)%compute
%end
+ update_x
%ehalo
;
806 jsd
= domain
%y(1)%compute
%begin
- update_y
%shalo
;
807 jed
= domain
%y(1)%compute
%end
+ update_y
%nhalo
;
809 j
= domain
%y(1)%global
%end
+shift
810 if( jsd
.LE
. j
.AND
. j
.LE
.jed
+shift
)then
!fold is within domain
811 !poles set to
0: BGRID only
812 if( gridtype
.EQ
.BGRID_NE
)then
813 midpoint
= (domain
%x(1)%global
%begin
+domain
%x(1)%global
%end
-1+shift
)/2
814 j
= domain
%y(1)%global
%end
+shift
815 is
= domain
%x(1)%global
%begin
; ie
= domain
%x(1)%global
%end
+shift
816 if( .NOT
. domain
%symmetry
) is
= is
- 1
817 do i
= is
,ie
, midpoint
818 if( isd
.LE
.i
.AND
. i
.LE
. ied
+shift
)then
820 ptr_fieldx
= f_addrsx(l
, 1)
821 ptr_fieldy
= f_addrsy(l
, 1)
831 ! the following code code block correct an error where the data in your halo coming from
832 ! other half may have the wrong sign
833 !off west edge
, when update north
or west direction
834 j
= domain
%y(1)%global
%end
+shift
835 if ( recv(7) .OR
. recv(5) ) then
836 select
case(gridtype
)
838 if(domain
%symmetry
) then
839 is
= domain
%x(1)%global
%begin
841 is
= domain
%x(1)%global
%begin
- 1
844 if( 2*is
-domain
%x(1)%data
%begin
.GT
.domain
%x(1)%data
%end
+shift
) &
845 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-north BGRID_NE west edge ubound error.' )
847 ptr_fieldx
= f_addrsx(l
, 1)
848 ptr_fieldy
= f_addrsy(l
, 1)
851 fieldx(i
,j
,k
) = fieldx(2*is
-i
,j
,k
)
852 fieldy(i
,j
,k
) = fieldy(2*is
-i
,j
,k
)
858 is
= domain
%x(1)%global
%begin
860 if( 2*is
-domain
%x(1)%data
%begin
-1.GT
.domain
%x(1)%data
%end
) &
861 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-north CGRID_NE west edge ubound error.' )
863 ptr_fieldy
= f_addrsy(l
, 1)
866 fieldy(i
,j
,k
) = fieldy(2*is
-i
-1,j
,k
)
875 is
= domain
%x(1)%global
%end
876 if(domain
%x(1)%cyclic
.AND
. is
.LT
.ied
)then
879 select
case(gridtype
)
884 ptr_fieldx
= f_addrsx(l
, 1)
885 ptr_fieldy
= f_addrsy(l
, 1)
888 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
889 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
895 ptr_fieldy
= f_addrsy(l
, 1)
898 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
905 else if( BTEST(domain
%fold
,SOUTH
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
! ---southern
907 ! NOTE
: symmetry is assumed
for fold
-south boundary
908 j
= domain
%y(1)%global
%begin
909 if( domain
%y(1)%data
%begin
.LE
.j
.AND
. j
.LE
.domain
%y(1)%data
%end
+shift
)then
!fold is within domain
910 midpoint
= (domain
%x(1)%global
%begin
+domain
%x(1)%global
%end
-1+shift
)/2
911 !poles set to
0: BGRID only
912 if( gridtype
.EQ
.BGRID_NE
)then
913 j
= domain
%y(1)%global
%begin
914 is
= domain
%x(1)%global
%begin
; ie
= domain
%x(1)%global
%end
+shift
915 do i
= is
,ie
, midpoint
916 if( domain
%x(1)%data
%begin
.LE
.i
.AND
. i
.LE
. domain
%x(1)%data
%end
+shift
)then
918 ptr_fieldx
= f_addrsx(l
, 1)
919 ptr_fieldy
= f_addrsy(l
, 1)
929 ! the following code code block correct an error where the data in your halo coming from
930 ! other half may have the wrong sign
931 !off west edge
, when update north
or west direction
932 j
= domain
%y(1)%global
%begin
933 if ( recv(3) .OR
. recv(5) ) then
934 select
case(gridtype
)
936 is
= domain
%x(1)%global
%begin
937 if( is
.GT
.domain
%x(1)%data
%begin
)then
939 if( 2*is
-domain
%x(1)%data
%begin
.GT
.domain
%x(1)%data
%end
+shift
) &
940 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-south BGRID_NE west edge ubound error.' )
942 ptr_fieldx
= f_addrsx(l
, 1)
943 ptr_fieldy
= f_addrsy(l
, 1)
945 do i
= domain
%x(1)%data
%begin
,is
-1
946 fieldx(i
,j
,k
) = fieldx(2*is
-i
,j
,k
)
947 fieldy(i
,j
,k
) = fieldy(2*is
-i
,j
,k
)
953 is
= domain
%x(1)%global
%begin
954 if( is
.GT
.domain
%x(1)%data
%begin
)then
955 if( 2*is
-domain
%x(1)%data
%begin
-1.GT
.domain
%x(1)%data
%end
) &
956 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-south CGRID_NE west edge ubound error.' )
958 ptr_fieldy
= f_addrsy(l
, 1)
960 do i
= domain
%x(1)%data
%begin
,is
-1
961 fieldy(i
,j
,k
) = fieldy(2*is
-i
-1,j
,k
)
970 is
= domain
%x(1)%global
%end
971 if(domain
%x(1)%cyclic
.AND
. is
.LT
.domain
%x(1)%data
%end
)then
972 ie
= domain
%x(1)%data
%end
974 select
case(gridtype
)
979 ptr_fieldx
= f_addrsx(l
, 1)
980 ptr_fieldy
= f_addrsy(l
, 1)
983 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
984 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
990 ptr_fieldy
= f_addrsy(l
, 1)
993 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
1000 else if( BTEST(domain
%fold
,WEST
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
! ---eastern
1002 ! NOTE
: symmetry is assumed
for fold
-west boundary
1003 i
= domain
%x(1)%global
%begin
1004 if( domain
%x(1)%data
%begin
.LE
.i
.AND
. i
.LE
.domain
%x(1)%data
%end
+shift
)then
!fold is within domain
1005 midpoint
= (domain
%y(1)%global
%begin
+domain
%y(1)%global
%end
-1+shift
)/2
1006 !poles set to
0: BGRID only
1007 if( gridtype
.EQ
.BGRID_NE
)then
1008 i
= domain
%x(1)%global
%begin
1009 js
= domain
%y(1)%global
%begin
; je
= domain
%y(1)%global
%end
+shift
1010 do j
= js
,je
, midpoint
1011 if( domain
%y(1)%data
%begin
.LE
.j
.AND
. j
.LE
. domain
%y(1)%data
%end
+shift
)then
1013 ptr_fieldx
= f_addrsx(l
, 1)
1014 ptr_fieldy
= f_addrsy(l
, 1)
1024 ! the following code code block correct an error where the data in your halo coming from
1025 ! other half may have the wrong sign
1026 !off south edge
, when update south
or west direction
1027 i
= domain
%x(1)%global
%begin
1028 if ( recv(3) .OR
. recv(5) ) then
1029 select
case(gridtype
)
1031 js
= domain
%y(1)%global
%begin
1032 if( js
.GT
.domain
%y(1)%data
%begin
)then
1034 if( 2*js
-domain
%y(1)%data
%begin
.GT
.domain
%y(1)%data
%end
+shift
) &
1035 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-west BGRID_NE west edge ubound error.' )
1037 ptr_fieldx
= f_addrsx(l
, 1)
1038 ptr_fieldy
= f_addrsy(l
, 1)
1040 do j
= domain
%y(1)%data
%begin
,js
-1
1041 fieldx(i
,j
,k
) = fieldx(i
,2*js
-j
,k
)
1042 fieldy(i
,j
,k
) = fieldy(i
,2*js
-j
,k
)
1048 js
= domain
%y(1)%global
%begin
1049 if( js
.GT
.domain
%y(1)%data
%begin
)then
1050 if( 2*js
-domain
%y(1)%data
%begin
-1.GT
.domain
%y(1)%data
%end
) &
1051 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-west CGRID_NE west edge ubound error.' )
1053 ptr_fieldx
= f_addrsx(l
, 1)
1055 do j
= domain
%y(1)%data
%begin
,js
-1
1056 fieldx(i
,j
,k
) = fieldx(i
, 2*js
-j
-1,k
)
1065 js
= domain
%y(1)%global
%end
1066 if(domain
%y(1)%cyclic
.AND
. js
.LT
.domain
%y(1)%data
%end
)then
1067 je
= domain
%y(1)%data
%end
1069 select
case(gridtype
)
1074 ptr_fieldx
= f_addrsx(l
, 1)
1075 ptr_fieldy
= f_addrsy(l
, 1)
1078 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
1079 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
1085 ptr_fieldx
= f_addrsx(l
, 1)
1088 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
1095 else if( BTEST(domain
%fold
,EAST
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
! ---eastern
1097 ! NOTE
: symmetry is assumed
for fold
-west boundary
1098 i
= domain
%x(1)%global
%end
+shift
1099 if( domain
%x(1)%data
%begin
.LE
.i
.AND
. i
.LE
.domain
%x(1)%data
%end
+shift
)then
!fold is within domain
1100 midpoint
= (domain
%y(1)%global
%begin
+domain
%y(1)%global
%end
-1+shift
)/2
1101 !poles set to
0: BGRID only
1102 if( gridtype
.EQ
.BGRID_NE
)then
1103 i
= domain
%x(1)%global
%end
+shift
1104 js
= domain
%y(1)%global
%begin
; je
= domain
%y(1)%global
%end
+shift
1105 do j
= js
,je
, midpoint
1106 if( domain
%y(1)%data
%begin
.LE
.j
.AND
. j
.LE
. domain
%y(1)%data
%end
+shift
)then
1108 ptr_fieldx
= f_addrsx(l
, 1)
1109 ptr_fieldy
= f_addrsy(l
, 1)
1119 ! the following code code block correct an error where the data in your halo coming from
1120 ! other half may have the wrong sign
1121 !off south edge
, when update south
or west direction
1122 i
= domain
%x(1)%global
%end
+shift
1123 if ( recv(3) .OR
. recv(1) ) then
1124 select
case(gridtype
)
1126 js
= domain
%y(1)%global
%begin
1127 if( js
.GT
.domain
%y(1)%data
%begin
)then
1129 if( 2*js
-domain
%y(1)%data
%begin
.GT
.domain
%y(1)%data
%end
+shift
) &
1130 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-east BGRID_NE west edge ubound error.' )
1132 ptr_fieldx
= f_addrsx(l
, 1)
1133 ptr_fieldy
= f_addrsy(l
, 1)
1135 do j
= domain
%y(1)%data
%begin
,js
-1
1136 fieldx(i
,j
,k
) = fieldx(i
,2*js
-j
,k
)
1137 fieldy(i
,j
,k
) = fieldy(i
,2*js
-j
,k
)
1143 js
= domain
%y(1)%global
%begin
1144 if( js
.GT
.domain
%y(1)%data
%begin
)then
1145 if( 2*js
-domain
%y(1)%data
%begin
-1.GT
.domain
%y(1)%data
%end
) &
1146 call
mpp_error( FATAL
, 'MPP_DO_UPDATE_V: folded-east CGRID_NE west edge ubound error.' )
1148 ptr_fieldx
= f_addrsx(l
, 1)
1150 do j
= domain
%y(1)%data
%begin
,js
-1
1151 fieldx(i
,j
,k
) = fieldx(i
, 2*js
-j
-1,k
)
1160 js
= domain
%y(1)%global
%end
1161 if(domain
%y(1)%cyclic
.AND
. js
.LT
.domain
%y(1)%data
%end
)then
1162 je
= domain
%y(1)%data
%end
1164 select
case(gridtype
)
1169 ptr_fieldx
= f_addrsx(l
, 1)
1170 ptr_fieldy
= f_addrsy(l
, 1)
1173 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
1174 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
1180 ptr_fieldx
= f_addrsx(l
, 1)
1183 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
1192 call
mpp_clock_begin(wait_clock
)
1193 call
mpp_sync_self( )
1194 call
mpp_clock_end(wait_clock
)
1198 end subroutine MPP_DO_UPDATE_3D_V_