fix: Fixes for linter action and code style (#869)
[FMS.git] / mpp / include / mpp_do_updateV.h
blob1f8c9aa5ba6707eb1193d4939eca95eddcd7be7e
1 ! -*-f90-*-
2 !***********************************************************************
3 !* GNU Lesser General Public License
4 !*
5 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !*
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
15 !* for more details.
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
40 integer :: tMe, dir
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(:)))
48 pointer(ptr,buffer )
49 integer :: buffer_pos
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
55 outunit = stdout()
56 update_flags = XUPDATE+YUPDATE !default
57 if( PRESENT(flags) ) then
58 update_flags = flags
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
65 end if
66 end if
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)
72 recv = .false.
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
79 recv(1) = .true.
80 recv(3) = .true.
81 recv(5) = .true.
82 recv(7) = .true.
83 endif
84 else
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)
89 endif
91 send = recv
93 l_size = size(f_addrsx,1)
94 nlist = size(domain%list(:))
95 ptr = LOC(mpp_domains_stack)
97 !recv
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) )
105 msg1 = 0
106 msg2 = 0
107 msg3 = 0
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)
111 msgsize = 0
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)
116 if(recv(dir)) then
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)
120 end if
121 end do
122 ind_x = ind_x+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
126 else
127 rank_x = -1
128 endif
129 endif
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)
134 if(recv(dir)) then
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)
138 end if
139 end do
140 ind_y = ind_y+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
144 else
145 rank_y = -1
146 endif
147 endif
148 cur_rank = max(rank_x, rank_y)
149 m = from_pe-mpp_root_pe()
150 msg2(m) = msgsize
151 end do
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)
155 msgsize = 0
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)
160 if( send(dir) ) then
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)
164 end if
165 end do
166 ind_x = ind_x+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
170 else
171 rank_x = nlist+1
172 endif
173 endif
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)
178 if( send(dir) ) then
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)
182 end if
183 end do
184 ind_y = ind_y+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
188 else
189 rank_y = nlist+1
190 endif
191 endif
192 m = to_pe-mpp_root_pe()
193 msg3(m) = msgsize
194 cur_rank = min(rank_x, rank_y)
195 enddo
196 call mpp_alltoall(msg3, 1, msg1, 1)
197 ! call mpp_sync_self(check=EVENT_RECV)
198 do m = 0, nlist-1
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")
203 endif
204 enddo
206 ! call mpp_sync_self()
207 write(outunit,*)"NOTE from mpp_do_updateV: message sizes are matched between send and recv for domain " &
208 //trim(domain%name)
209 deallocate(msg1, msg2, msg3)
210 endif
212 !--- recv
213 buffer_pos = 0
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)
217 msgsize = 0
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)
224 if(recv(dir)) then
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)
229 end if
230 end do
231 msgsize = msgsize*2
232 ind_x = ind_x+1
233 ind_y = ind_x
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
237 else
238 rank_x = -1
239 endif
240 rank_y = rank_x
241 endif
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)
247 if(recv(dir)) then
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)
251 end if
252 end do
253 ind_x = ind_x+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
257 else
258 rank_x = -1
259 endif
260 endif
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)
265 if(recv(dir)) then
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)
269 end if
270 end do
271 ind_y = ind_y+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
275 else
276 rank_y = -1
277 endif
278 endif
279 end select
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.' )
289 end if
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
292 end if
293 end do
294 call mpp_clock_end(recv_clock)
295 buffer_recv_size = buffer_pos
296 send_start_pos = buffer_pos
298 !--- send
299 cur_rank = get_rank_send(domain, update_x, update_y, rank_x, rank_y, ind_x, ind_y)
300 nsend = 0
301 call mpp_clock_begin(pack_clock)
302 do while (ind_x .LE. nsend_x .OR. ind_y .LE. nsend_y)
303 pos = buffer_pos
304 !--- make sure the domain stack size is big enough
305 msgsize = 0
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)
310 enddo
311 endif
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)
316 enddo
317 endif
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.')
326 end if
327 end if
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)
334 if( send(dir) ) then
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) )
340 case(ZERO)
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)
344 do k = 1,ke
345 do j = js, je
346 do i = is, ie
347 pos = pos + 2
348 buffer(pos-1) = fieldx(i,j,k)
349 buffer(pos) = fieldy(i,j,k)
350 end do
351 end do
352 end do
353 end do
354 case( MINUS_NINETY )
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)
359 do k = 1,ke
360 do i = is, ie
361 do j = je, js, -1
362 pos = pos + 2
363 buffer(pos-1) = fieldy(i,j,k)
364 buffer(pos) = fieldx(i,j,k)
365 end do
366 end do
367 end do
368 end do
369 else
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)
373 do k = 1,ke
374 do i = is, ie
375 do j = je, js, -1
376 pos = pos + 2
377 buffer(pos-1) = -fieldy(i,j,k)
378 buffer(pos) = fieldx(i,j,k)
379 end do
380 end do
381 end do
382 end do
383 end if
384 case( NINETY )
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)
389 do k = 1,ke
390 do i = ie, is, -1
391 do j = js, je
392 pos = pos + 2
393 buffer(pos-1) = fieldy(i,j,k)
394 buffer(pos) = fieldx(i,j,k)
395 end do
396 end do
397 end do
398 end do
399 else
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)
403 do k = 1,ke
404 do i = ie, is, -1
405 do j = js, je
406 pos = pos + 2
407 buffer(pos-1) = fieldy(i,j,k)
408 buffer(pos) = -fieldx(i,j,k)
409 end do
410 end do
411 end do
412 end do
413 end if
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)
419 do k = 1,ke
420 do j = je, js, -1
421 do i = ie, is, -1
422 pos = pos + 2
423 buffer(pos-1) = fieldx(i,j,k)
424 buffer(pos) = fieldy(i,j,k)
425 end do
426 end do
427 end do
428 end do
429 else
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)
433 do k = 1,ke
434 do j = je, js, -1
435 do i = ie, is, -1
436 pos = pos + 2
437 buffer(pos-1) = -fieldx(i,j,k)
438 buffer(pos) = -fieldy(i,j,k)
439 end do
440 end do
441 end do
442 end do
443 end if
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
447 ind_x = ind_x+1
448 ind_y = ind_x
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
452 else
453 rank_x = nlist+1
454 endif
455 rank_y = rank_x
456 endif
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)
462 if( send(dir) ) then
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) )
467 case(ZERO)
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)
471 do k = 1,ke
472 do j = js, je
473 do i = is, ie
474 pos = pos + 1
475 buffer(pos) = fieldx(i,j,k)
476 end do
477 end do
478 end do
479 end do
480 case(MINUS_NINETY)
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)
485 do k = 1,ke
486 do i = is, ie
487 do j = je, js, -1
488 pos = pos + 1
489 buffer(pos) = fieldy(i,j,k)
490 end do
491 end do
492 end do
493 end do
494 else
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)
498 do k = 1,ke
499 do i = is, ie
500 do j = je, js, -1
501 pos = pos + 1
502 buffer(pos) = -fieldy(i,j,k)
503 end do
504 end do
505 end do
506 end do
507 end if
508 case(NINETY)
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)
512 do k = 1, ke
513 do i = ie, is, -1
514 do j = js, je
515 pos = pos + 1
516 buffer(pos) = fieldy(i,j,k)
517 end do
518 end do
519 end do
520 end do
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)
526 do k = 1,ke
527 do j = je, js, -1
528 do i = ie, is, -1
529 pos = pos + 1
530 buffer(pos) = fieldx(i,j,k)
531 end do
532 end do
533 end do
534 end do
535 else
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)
539 do k = 1,ke
540 do j = je, js, -1
541 do i = ie, is, -1
542 pos = pos + 1
543 buffer(pos) = -fieldx(i,j,k)
544 end do
545 end do
546 end do
547 end do
548 end if
549 end select
550 end if
551 end do
552 ind_x = ind_x+1
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
556 else
557 rank_x = nlist+1
558 endif
559 endif
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)
564 if( send(dir) ) then
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) )
569 case(ZERO)
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)
573 do k = 1,ke
574 do j = js, je
575 do i = is, ie
576 pos = pos + 1
577 buffer(pos) = fieldy(i,j,k)
578 end do
579 end do
580 end do
581 end do
582 case(MINUS_NINETY)
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)
586 do k = 1,ke
587 do i = is, ie
588 do j = je, js, -1
589 pos = pos + 1
590 buffer(pos) = fieldx(i,j,k)
591 end do
592 end do
593 end do
594 end do
595 case(NINETY)
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)
600 do k = 1,ke
601 do i = ie, is, -1
602 do j = js, je
603 pos = pos + 1
604 buffer(pos) = fieldx(i,j,k)
605 end do
606 end do
607 end do
608 end do
609 else
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)
613 do k = 1,ke
614 do i = ie, is, -1
615 do j = js, je
616 pos = pos + 1
617 buffer(pos) = -fieldx(i,j,k)
618 end do
619 end do
620 end do
621 end do
622 end if
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)
628 do k = 1,ke
629 do j = je, js, -1
630 do i = ie, is, -1
631 pos = pos + 1
632 buffer(pos) = fieldy(i,j,k)
633 end do
634 end do
635 end do
636 end do
637 else
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)
641 do k = 1,ke
642 do j = je, js, -1
643 do i = ie, is, -1
644 pos = pos + 1
645 buffer(pos) = -fieldy(i,j,k)
646 end do
647 end do
648 end do
649 end do
650 end if
651 end select
652 endif
653 enddo
654 ind_y = ind_y+1
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
658 else
659 rank_y = nlist+1
660 endif
661 endif
662 end select
663 cur_rank = min(rank_x, rank_y)
664 nsend = nsend + 1
665 send_pe(nsend) = to_pe
666 send_msgsize(nsend) = pos - buffer_pos
667 buffer_pos = pos
668 end do
670 buffer_pos = send_start_pos
671 call mpp_clock_end(pack_clock)
672 call mpp_clock_begin(send_clock)
673 do m = 1, nsend
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
678 end if
679 end do
680 call mpp_clock_end(send_clock)
682 !unpack recv
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)
692 pos = buffer_pos
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)
698 if( recv(dir) ) then
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
704 buffer_pos = pos
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)
708 do k = 1,ke
709 do j = js, je
710 do i = is, ie
711 pos = pos + 2
712 fieldx(i,j,k) = buffer(pos-1)
713 fieldy(i,j,k) = buffer(pos)
714 end do
715 end do
716 end do
717 end do
718 end if ! end if( recv(dir) )
719 end do ! do dir=8,1,-1
720 ind_x = ind_x-1
721 ind_y = ind_x
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
725 else
726 rank_x = nlist+1
727 endif
728 rank_y = rank_x
729 endif
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)
734 if( recv(dir) ) then
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
740 buffer_pos = pos
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)
744 do k = 1,ke
745 do j = js, je
746 do i = is, ie
747 pos = pos + 1
748 fieldy(i,j,k) = buffer(pos)
749 end do
750 end do
751 end do
752 end do
753 end if
754 end do
755 ind_y = ind_y-1
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
759 else
760 rank_y = nlist+1
761 endif
762 endif
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)
766 if( recv(dir) ) then
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
772 buffer_pos = pos
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)
776 do k = 1,ke
777 do j = js, je
778 do i = is, ie
779 pos = pos + 1
780 fieldx(i,j,k) = buffer(pos)
781 end do
782 end do
783 end do
784 end do
785 end if
786 end do
787 ind_x = ind_x-1
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
791 else
792 rank_x = nlist+1
793 endif
794 endif
795 end select
796 cur_rank = min(rank_x, rank_y)
797 end do
798 call mpp_clock_end(unpk_clock)
800 ! ---northern boundary fold
801 shift = 0
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
819 do l=1,l_size
820 ptr_fieldx = f_addrsx(l, 1)
821 ptr_fieldy = f_addrsy(l, 1)
822 do k = 1,ke
823 fieldx(i,j,k) = 0.
824 fieldy(i,j,k) = 0.
825 end do
826 end do
827 end if
828 end do
829 endif
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)
837 case(BGRID_NE)
838 if(domain%symmetry) then
839 is = domain%x(1)%global%begin
840 else
841 is = domain%x(1)%global%begin - 1
842 end if
843 if( is.GT.isd )then
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.' )
846 do l=1,l_size
847 ptr_fieldx = f_addrsx(l, 1)
848 ptr_fieldy = f_addrsy(l, 1)
849 do k = 1,ke
850 do i = isd,is-1
851 fieldx(i,j,k) = fieldx(2*is-i,j,k)
852 fieldy(i,j,k) = fieldy(2*is-i,j,k)
853 end do
854 end do
855 end do
856 end if
857 case(CGRID_NE)
858 is = domain%x(1)%global%begin
859 if( is.GT.isd )then
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.' )
862 do l=1,l_size
863 ptr_fieldy = f_addrsy(l, 1)
864 do k = 1,ke
865 do i = isd,is-1
866 fieldy(i,j,k) = fieldy(2*is-i-1,j,k)
867 end do
868 end do
869 end do
870 end if
871 end select
872 end if
874 !off east edge
875 is = domain%x(1)%global%end
876 if(domain%x(1)%cyclic .AND. is.LT.ied )then
877 ie = ied
878 is = is + 1
879 select case(gridtype)
880 case(BGRID_NE)
881 is = is + shift
882 ie = ie + shift
883 do l=1,l_size
884 ptr_fieldx = f_addrsx(l, 1)
885 ptr_fieldy = f_addrsy(l, 1)
886 do k = 1,ke
887 do i = is,ie
888 fieldx(i,j,k) = -fieldx(i,j,k)
889 fieldy(i,j,k) = -fieldy(i,j,k)
890 end do
891 end do
892 end do
893 case(CGRID_NE)
894 do l=1,l_size
895 ptr_fieldy = f_addrsy(l, 1)
896 do k = 1,ke
897 do i = is, ie
898 fieldy(i,j,k) = -fieldy(i,j,k)
899 end do
900 end do
901 end do
902 end select
903 end if
904 end if
905 else if( BTEST(domain%fold,SOUTH) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then ! ---southern
906 !! boundary fold
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
917 do l=1,l_size
918 ptr_fieldx = f_addrsx(l, 1)
919 ptr_fieldy = f_addrsy(l, 1)
920 do k = 1,ke
921 fieldx(i,j,k) = 0.
922 fieldy(i,j,k) = 0.
923 end do
924 end do
925 end if
926 end do
927 endif
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)
935 case(BGRID_NE)
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.' )
941 do l=1,l_size
942 ptr_fieldx = f_addrsx(l, 1)
943 ptr_fieldy = f_addrsy(l, 1)
944 do k = 1,ke
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)
948 end do
949 end do
950 end do
951 end if
952 case(CGRID_NE)
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.' )
957 do l=1,l_size
958 ptr_fieldy = f_addrsy(l, 1)
959 do k = 1,ke
960 do i = domain%x(1)%data%begin,is-1
961 fieldy(i,j,k) = fieldy(2*is-i-1,j,k)
962 end do
963 end do
964 end do
965 end if
966 end select
967 end if
969 !off east edge
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
973 is = is + 1
974 select case(gridtype)
975 case(BGRID_NE)
976 is = is + shift
977 ie = ie + shift
978 do l=1,l_size
979 ptr_fieldx = f_addrsx(l, 1)
980 ptr_fieldy = f_addrsy(l, 1)
981 do k = 1,ke
982 do i = is,ie
983 fieldx(i,j,k) = -fieldx(i,j,k)
984 fieldy(i,j,k) = -fieldy(i,j,k)
985 end do
986 end do
987 end do
988 case(CGRID_NE)
989 do l=1,l_size
990 ptr_fieldy = f_addrsy(l, 1)
991 do k = 1,ke
992 do i = is, ie
993 fieldy(i,j,k) = -fieldy(i,j,k)
994 end do
995 end do
996 end do
997 end select
998 end if
999 end if
1000 else if( BTEST(domain%fold,WEST) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then ! ---eastern
1001 !! boundary fold
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
1012 do l=1,l_size
1013 ptr_fieldx = f_addrsx(l, 1)
1014 ptr_fieldy = f_addrsy(l, 1)
1015 do k = 1,ke
1016 fieldx(i,j,k) = 0.
1017 fieldy(i,j,k) = 0.
1018 end do
1019 end do
1020 end if
1021 end do
1022 endif
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)
1030 case(BGRID_NE)
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.' )
1036 do l=1,l_size
1037 ptr_fieldx = f_addrsx(l, 1)
1038 ptr_fieldy = f_addrsy(l, 1)
1039 do k = 1,ke
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)
1043 end do
1044 end do
1045 end do
1046 end if
1047 case(CGRID_NE)
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.' )
1052 do l=1,l_size
1053 ptr_fieldx = f_addrsx(l, 1)
1054 do k = 1,ke
1055 do j = domain%y(1)%data%begin,js-1
1056 fieldx(i,j,k) = fieldx(i, 2*js-j-1,k)
1057 end do
1058 end do
1059 end do
1060 end if
1061 end select
1062 end if
1064 !off north edge
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
1068 js = js + 1
1069 select case(gridtype)
1070 case(BGRID_NE)
1071 js = js + shift
1072 je = je + shift
1073 do l=1,l_size
1074 ptr_fieldx = f_addrsx(l, 1)
1075 ptr_fieldy = f_addrsy(l, 1)
1076 do k = 1,ke
1077 do j = js,je
1078 fieldx(i,j,k) = -fieldx(i,j,k)
1079 fieldy(i,j,k) = -fieldy(i,j,k)
1080 end do
1081 end do
1082 end do
1083 case(CGRID_NE)
1084 do l=1,l_size
1085 ptr_fieldx = f_addrsx(l, 1)
1086 do k = 1,ke
1087 do j = js, je
1088 fieldx(i,j,k) = -fieldx(i,j,k)
1089 end do
1090 end do
1091 end do
1092 end select
1093 end if
1094 end if
1095 else if( BTEST(domain%fold,EAST) .AND. (.NOT.BTEST(update_flags,SCALAR_BIT)) )then ! ---eastern
1096 !! boundary fold
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
1107 do l=1,l_size
1108 ptr_fieldx = f_addrsx(l, 1)
1109 ptr_fieldy = f_addrsy(l, 1)
1110 do k = 1,ke
1111 fieldx(i,j,k) = 0.
1112 fieldy(i,j,k) = 0.
1113 end do
1114 end do
1115 end if
1116 end do
1117 endif
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)
1125 case(BGRID_NE)
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.' )
1131 do l=1,l_size
1132 ptr_fieldx = f_addrsx(l, 1)
1133 ptr_fieldy = f_addrsy(l, 1)
1134 do k = 1,ke
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)
1138 end do
1139 end do
1140 end do
1141 end if
1142 case(CGRID_NE)
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.' )
1147 do l=1,l_size
1148 ptr_fieldx = f_addrsx(l, 1)
1149 do k = 1,ke
1150 do j = domain%y(1)%data%begin,js-1
1151 fieldx(i,j,k) = fieldx(i, 2*js-j-1,k)
1152 end do
1153 end do
1154 end do
1155 end if
1156 end select
1157 end if
1159 !off north edge
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
1163 js = js + 1
1164 select case(gridtype)
1165 case(BGRID_NE)
1166 js = js + shift
1167 je = je + shift
1168 do l=1,l_size
1169 ptr_fieldx = f_addrsx(l, 1)
1170 ptr_fieldy = f_addrsy(l, 1)
1171 do k = 1,ke
1172 do j = js,je
1173 fieldx(i,j,k) = -fieldx(i,j,k)
1174 fieldy(i,j,k) = -fieldy(i,j,k)
1175 end do
1176 end do
1177 end do
1178 case(CGRID_NE)
1179 do l=1,l_size
1180 ptr_fieldx = f_addrsx(l, 1)
1181 do k = 1,ke
1182 do j = js, je
1183 fieldx(i,j,k) = -fieldx(i,j,k)
1184 end do
1185 end do
1186 end do
1187 end select
1188 end if
1189 end if
1190 end if
1192 call mpp_clock_begin(wait_clock)
1193 call mpp_sync_self( )
1194 call mpp_clock_end(wait_clock)
1196 return
1198 end subroutine MPP_DO_UPDATE_3D_V_