4 !***********************************************************************
5 !* GNU Lesser General Public License
7 !* This file is part of the GFDL Flexible Modeling
System (FMS
).
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
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
45 integer
, allocatable :: msg1(:), msg2(:)
46 logical :: send(8), recv(8), update_edge_only
47 MPP_TYPE_ :: buffer(size(mpp_domains_stack(:)))
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
55 update_flags
= XUPDATE
+YUPDATE
!default
56 if( PRESENT(flags
) ) then
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
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
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)
91 l_size
= size(f_addrsx
,1)
92 nlist
= size(domain
%list(:))
93 ptr
= LOC(mpp_domains_stack
)
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) )
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
)
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
)
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)
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
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
)
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)
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
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
)
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
)
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
)
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)
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
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
)
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)
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
190 cur_rank
= min(rank_x
, rank_y
)
191 call
mpp_send( msgsize
, plen
=1, to_pe
=to_pe
, tag
=COMM_TAG_1
)
193 call
mpp_sync_self(check
=EVENT_RECV
)
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")
203 write(outunit
,*)"NOTE from mpp_do_updateV: message sizes are matched between send and recv for domain " &
205 deallocate(msg1
, msg2
)
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
)
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
)
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)
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
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
)
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)
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
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
)
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)
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
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
282 buffer_recv_size
= buffer_pos
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)
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
)
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
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
)
309 buffer(pos
-1) = fieldx(i
,j
,k
)
310 buffer(pos
) = fieldy(i
,j
,k
)
317 end
if ! end
if( recv(dir
) )
318 end
do ! do dir
=8,1,-1
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
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
)
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
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
)
347 buffer(pos
) = fieldy(i
,j
,k
)
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
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
)
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
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
)
380 buffer(pos
) = fieldx(i
,j
,k
)
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
397 cur_rank
= min(rank_x
, rank_y
)
400 ! ---northern boundary fold
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
415 ptr_fieldx
= f_addrsx(l
, 1)
416 ptr_fieldy
= f_addrsy(l
, 1)
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
)
433 if(domain
%symmetry
) then
434 is
= domain
%x(1)%global
%begin
436 is
= domain
%x(1)%global
%begin
- 1
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.' )
443 ptr_fieldx
= f_addrsx(l
, 1)
444 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
459 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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
475 select
case(gridtype
)
480 ptr_fieldx
= f_addrsx(l
, 1)
481 ptr_fieldy
= f_addrsy(l
, 1)
484 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
485 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
491 ptr_fieldy
= f_addrsy(l
, 1)
494 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
501 else if( BTEST(domain
%fold
,SOUTH
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
! ---southern
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
514 ptr_fieldx
= f_addrsx(l
, 1)
515 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
538 ptr_fieldx
= f_addrsx(l
, 1)
539 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
554 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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
570 select
case(gridtype
)
575 ptr_fieldx
= f_addrsx(l
, 1)
576 ptr_fieldy
= f_addrsy(l
, 1)
579 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
580 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
586 ptr_fieldy
= f_addrsy(l
, 1)
589 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
596 else if( BTEST(domain
%fold
,WEST
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
! ---eastern
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
609 ptr_fieldx
= f_addrsx(l
, 1)
610 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
633 ptr_fieldx
= f_addrsx(l
, 1)
634 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
649 ptr_fieldx
= f_addrsx(l
, 1)
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
)
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
665 select
case(gridtype
)
670 ptr_fieldx
= f_addrsx(l
, 1)
671 ptr_fieldy
= f_addrsy(l
, 1)
674 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
675 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
681 ptr_fieldx
= f_addrsx(l
, 1)
684 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
691 else if( BTEST(domain
%fold
,EAST
) .AND
. (.NOT
.BTEST(update_flags
,SCALAR_BIT
)) )then
! ---eastern
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
704 ptr_fieldx
= f_addrsx(l
, 1)
705 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
728 ptr_fieldx
= f_addrsx(l
, 1)
729 ptr_fieldy
= f_addrsy(l
, 1)
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
)
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.' )
744 ptr_fieldx
= f_addrsx(l
, 1)
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
)
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
760 select
case(gridtype
)
765 ptr_fieldx
= f_addrsx(l
, 1)
766 ptr_fieldy
= f_addrsy(l
, 1)
769 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
770 fieldy(i
,j
,k
) = -fieldy(i
,j
,k
)
776 ptr_fieldx
= f_addrsx(l
, 1)
779 fieldx(i
,j
,k
) = -fieldx(i
,j
,k
)
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
)
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
)
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)
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
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
)
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)
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
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
)
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)
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
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
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
)
870 !--- make sure the domain stack size is big enough
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
)
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
)
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
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
902 case(CGRID_NE
, CGRID_SW
)
903 if(cur_rank
== rank_x
) then
904 to_pe
= update_x
%send(ind_x
)%pe
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
913 if(cur_rank
== rank_y
) then
914 to_pe
= update_y
%send(ind_y
)%pe
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
924 cur_rank
= min(rank_x
, rank_y
)
926 if( msgsize
.GT
.0 )then
927 msgsize
= msgsize
*ke
*l_size
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
937 call
mpp_sync_self(check
=EVENT_RECV
)
940 buffer_pos
= buffer_recv_size
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
)
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
)
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
) )
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
)
968 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
-1)
969 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
)
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
)
983 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
-1)
984 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
997 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)-buffer(pos
-1)
998 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
1013 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
-1)
1014 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
1027 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
-1)
1028 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)-buffer(pos
)
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
)
1043 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
-1)
1044 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
)
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
)
1057 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)-buffer(pos
-1)
1058 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)-buffer(pos
)
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
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
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
) )
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
)
1095 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
1109 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
)
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
)
1122 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)-buffer(pos
)
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
)
1136 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
)
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
)
1150 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
1163 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)-buffer(pos
)
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
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
) )
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
)
1198 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
)
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
)
1211 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
1225 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)+buffer(pos
)
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
)
1238 fieldx(i
,j
,k
)=fieldx(i
,j
,k
)-buffer(pos
)
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
)
1253 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)+buffer(pos
)
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
)
1266 fieldy(i
,j
,k
)=fieldy(i
,j
,k
)-buffer(pos
)
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
1284 cur_rank
= min(rank_x
, rank_y
)
1285 msgsize
= pos
- buffer_pos
1286 if( msgsize
.GT
.0 )then
1291 call
mpp_sync_self( )
1295 end subroutine MPP_DO_UPDATE_AD_3D_V_