fix: Fixes for linter action and code style (#869)
[FMS.git] / mpp / include / mpp_group_update.h
blob924a2c7911982b41f808dfd37debd894d0ef7ccf
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 ! -*-f90-*-
20 subroutine MPP_CREATE_GROUP_UPDATE_2D_(group, field, domain, flags, position, &
21 whalo, ehalo, shalo, nhalo)
22 type(mpp_group_update_type), intent(inout) :: group
23 MPP_TYPE_, intent(inout) :: field(:,:)
24 type(domain2D), intent(inout) :: domain
25 integer, intent(in), optional :: flags
26 integer, intent(in), optional :: position
27 integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
29 MPP_TYPE_ :: field3D(size(field,1),size(field,2),1)
30 pointer( ptr, field3D )
31 ptr = LOC(field)
33 call mpp_create_group_update(group, field3D, domain, flags, position, whalo, ehalo, shalo, nhalo)
35 return
37 end subroutine MPP_CREATE_GROUP_UPDATE_2D_
39 subroutine MPP_CREATE_GROUP_UPDATE_3D_(group, field, domain, flags, position, whalo, ehalo, shalo, nhalo)
40 type(mpp_group_update_type), intent(inout) :: group
41 MPP_TYPE_, intent(inout) :: field(:,:,:)
42 type(domain2D), intent(inout) :: domain
43 integer, intent(in), optional :: flags
44 integer, intent(in), optional :: position
45 integer, intent(in), optional :: whalo, ehalo, shalo, nhalo ! specify halo region to be updated.
47 integer :: update_position, update_whalo, update_ehalo, update_shalo, update_nhalo
48 integer :: update_flags, isize, jsize, ksize
49 integer :: nscalar
50 character(len=3) :: text
51 logical :: set_mismatch, update_edge_only
52 logical :: recv(8)
54 if(group%initialized) then
55 call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE_3D: group is already initialized")
56 endif
58 if(present(whalo)) then
59 update_whalo = whalo
60 if(abs(update_whalo) > domain%whalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE: "// &
61 "optional argument whalo should not be larger than the whalo when define domain.")
62 else
63 update_whalo = domain%whalo
64 end if
65 if(present(ehalo)) then
66 update_ehalo = ehalo
67 if(abs(update_ehalo) > domain%ehalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE: "// &
68 "optional argument ehalo should not be larger than the ehalo when define domain.")
69 else
70 update_ehalo = domain%ehalo
71 end if
72 if(present(shalo)) then
73 update_shalo = shalo
74 if(abs(update_shalo) > domain%shalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE: "// &
75 "optional argument shalo should not be larger than the shalo when define domain.")
76 else
77 update_shalo = domain%shalo
78 end if
79 if(present(nhalo)) then
80 update_nhalo = nhalo
81 if(abs(update_nhalo) > domain%nhalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE: "// &
82 "optional argument nhalo should not be larger than the nhalo when define domain.")
83 else
84 update_nhalo = domain%nhalo
85 end if
86 update_position = CENTER
87 !--- when there is NINETY or MINUS_NINETY rotation for some contact, the salar data can not be on E or N-cell,
88 if(present(position)) then
89 update_position = position
90 if(domain%rotated_ninety .AND. ( position == EAST .OR. position == NORTH ) ) &
91 call mpp_error(FATAL, 'MPP_CREATE_GROUP_UPDATE_3D: hen there is NINETY or MINUS_NINETY rotation, ' // &
92 'can not use scalar version update_domain for data on E or N-cell' )
93 end if
95 if( domain%max_ntile_pe > 1 ) then
96 call mpp_error(FATAL,'MPP_CREATE_GROUP_UPDATE: do not support multiple tile per processor')
97 endif
99 update_flags = XUPDATE+YUPDATE
100 if(present(flags)) update_flags = flags
102 group%nscalar = group%nscalar + 1
103 nscalar = group%nscalar
104 if( nscalar > MAX_DOMAIN_FIELDS)then
105 write( text,'(i2)' ) MAX_DOMAIN_FIELDS
106 call mpp_error(FATAL,'MPP_CREATE_GROUP_UPDATE: MAX_DOMAIN_FIELDS='//text//' exceeded for group update.' )
107 endif
109 isize = size(field,1); jsize=size(field,2); ksize = size(field,3)
111 group%addrs_s(nscalar) = LOC(field)
112 if( group%nscalar == 1 ) then
113 group%flags_s = update_flags
114 group%whalo_s = update_whalo
115 group%ehalo_s = update_ehalo
116 group%shalo_s = update_shalo
117 group%nhalo_s = update_nhalo
118 group%position = update_position
119 group%isize_s = isize
120 group%jsize_s = jsize
121 group%ksize_s = ksize
122 call mpp_get_memory_domain(domain, group%is_s, group%ie_s, group%js_s, group%je_s, position=position)
124 update_edge_only = BTEST(update_flags, EDGEONLY)
125 recv(1) = BTEST(update_flags,EAST)
126 recv(3) = BTEST(update_flags,SOUTH)
127 recv(5) = BTEST(update_flags,WEST)
128 recv(7) = BTEST(update_flags,NORTH)
129 if( update_edge_only ) then
130 recv(2) = .false.
131 recv(4) = .false.
132 recv(6) = .false.
133 recv(8) = .false.
134 if( .NOT. (recv(1) .OR. recv(3) .OR. recv(5) .OR. recv(7)) ) then
135 recv(1) = .true.
136 recv(3) = .true.
137 recv(5) = .true.
138 recv(7) = .true.
139 endif
140 else
141 recv(2) = recv(1) .AND. recv(3)
142 recv(4) = recv(3) .AND. recv(5)
143 recv(6) = recv(5) .AND. recv(7)
144 recv(8) = recv(7) .AND. recv(1)
145 endif
146 group%recv_s = recv
147 else
148 set_mismatch = .false.
149 set_mismatch = set_mismatch .OR. (group%flags_s .NE. update_flags)
150 set_mismatch = set_mismatch .OR. (group%whalo_s .NE. update_whalo)
151 set_mismatch = set_mismatch .OR. (group%ehalo_s .NE. update_ehalo)
152 set_mismatch = set_mismatch .OR. (group%shalo_s .NE. update_shalo)
153 set_mismatch = set_mismatch .OR. (group%nhalo_s .NE. update_nhalo)
154 set_mismatch = set_mismatch .OR. (group%position .NE. update_position)
155 set_mismatch = set_mismatch .OR. (group%isize_s .NE. isize)
156 set_mismatch = set_mismatch .OR. (group%jsize_s .NE. jsize)
157 set_mismatch = set_mismatch .OR. (group%ksize_s .NE. ksize)
159 if(set_mismatch)then
160 write( text,'(i2)' ) nscalar
161 call mpp_error(FATAL,'MPP_CREATE_GROUP_UPDATE_3D: Incompatible field at count '//text//' for group update.' )
162 endif
163 endif
165 return
167 end subroutine MPP_CREATE_GROUP_UPDATE_3D_
170 subroutine MPP_CREATE_GROUP_UPDATE_4D_(group, field, domain, flags, position, &
171 whalo, ehalo, shalo, nhalo)
172 type(mpp_group_update_type), intent(inout) :: group
173 MPP_TYPE_, intent(inout) :: field(:,:,:,:)
174 type(domain2D), intent(inout) :: domain
175 integer, intent(in), optional :: flags
176 integer, intent(in), optional :: position
177 integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
179 MPP_TYPE_ :: field3D(size(field,1),size(field,2),size(field,3)*size(field,4))
180 pointer( ptr, field3D )
181 ptr = LOC(field)
183 call mpp_create_group_update(group, field3D, domain, flags, position, whalo, ehalo, shalo, nhalo)
185 return
187 end subroutine MPP_CREATE_GROUP_UPDATE_4D_
189 subroutine MPP_CREATE_GROUP_UPDATE_2D_V_( group, fieldx, fieldy, domain, flags, gridtype, &
190 whalo, ehalo, shalo, nhalo)
192 type(mpp_group_update_type), intent(inout) :: group
193 MPP_TYPE_, intent(inout) :: fieldx(:,:), fieldy(:,:)
194 type(domain2D), intent(inout) :: domain
195 integer, intent(in), optional :: flags, gridtype
196 integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
197 MPP_TYPE_ :: field3Dx(size(fieldx,1),size(fieldx,2),1)
198 MPP_TYPE_ :: field3Dy(size(fieldy,1),size(fieldy,2),1)
199 pointer( ptrx, field3Dx )
200 pointer( ptry, field3Dy )
201 ptrx = LOC(fieldx)
202 ptry = LOC(fieldy)
204 call mpp_create_group_update(group, field3Dx, field3Dy, domain, flags, gridtype, &
205 whalo, ehalo, shalo, nhalo)
207 return
209 end subroutine MPP_CREATE_GROUP_UPDATE_2D_V_
213 subroutine MPP_CREATE_GROUP_UPDATE_3D_V_( group, fieldx, fieldy, domain, flags, gridtype, &
214 whalo, ehalo, shalo, nhalo)
215 type(mpp_group_update_type), intent(inout) :: group
216 MPP_TYPE_, intent(inout) :: fieldx(:,:,:), fieldy(:,:,:)
217 type(domain2D), intent(inout) :: domain
218 integer, intent(in), optional :: flags, gridtype
219 integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
221 integer :: update_whalo, update_ehalo, update_shalo, update_nhalo
222 integer :: update_flags, isize_x, jsize_x, ksize_x, isize_y, jsize_y, ksize_y
223 integer :: nvector, update_gridtype, position_x, position_y
224 character(len=3) :: text
225 logical :: set_mismatch, update_edge_only
226 logical :: recv(8)
229 if(group%initialized) then
230 call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE_V: group is already initialized")
231 endif
233 if(present(whalo)) then
234 update_whalo = whalo
235 if(abs(update_whalo) > domain%whalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE_V: "// &
236 "optional argument whalo should not be larger than the whalo when define domain.")
237 else
238 update_whalo = domain%whalo
239 end if
240 if(present(ehalo)) then
241 update_ehalo = ehalo
242 if(abs(update_ehalo) > domain%ehalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE_V: "// &
243 "optional argument ehalo should not be larger than the ehalo when define domain.")
244 else
245 update_ehalo = domain%ehalo
246 end if
247 if(present(shalo)) then
248 update_shalo = shalo
249 if(abs(update_shalo) > domain%shalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE_V: "// &
250 "optional argument shalo should not be larger than the shalo when define domain.")
251 else
252 update_shalo = domain%shalo
253 end if
254 if(present(nhalo)) then
255 update_nhalo = nhalo
256 if(abs(update_nhalo) > domain%nhalo ) call mpp_error(FATAL, "MPP_CREATE_GROUP_UPDATE_V: "// &
257 "optional argument nhalo should not be larger than the nhalo when define domain.")
258 else
259 update_nhalo = domain%nhalo
260 end if
262 update_gridtype = AGRID
263 if(PRESENT(gridtype)) update_gridtype = gridtype
265 if( domain%max_ntile_pe > 1 ) then
266 call mpp_error(FATAL,'MPP_CREATE_GROUP_UPDATE_V: do not support multiple tile per processor')
267 endif
269 update_flags = XUPDATE+YUPDATE !default
270 if( PRESENT(flags) )update_flags = flags
271 ! The following test is so that SCALAR_PAIR can be used alone with the
272 ! same default update pattern as without.
273 if (BTEST(update_flags,SCALAR_BIT)) then
274 if (.NOT.(BTEST(update_flags,WEST) .OR. BTEST(update_flags,EAST) &
275 .OR. BTEST(update_flags,NORTH) .OR. BTEST(update_flags,SOUTH))) &
276 update_flags = update_flags + XUPDATE+YUPDATE !default with SCALAR_PAIR
277 end if
279 group%nvector = group%nvector + 1
280 nvector = group%nvector
281 if( nvector > MAX_DOMAIN_FIELDS)then
282 write( text,'(i2)' ) MAX_DOMAIN_FIELDS
283 call mpp_error(FATAL,'MPP_CREATE_GROUP_UPDATE_V: MAX_DOMAIN_FIELDS='//text//' exceeded for group update.' )
284 endif
286 isize_x = size(fieldx,1); jsize_x = size(fieldx,2); ksize_x = size(fieldx,3)
287 isize_y = size(fieldy,1); jsize_y = size(fieldy,2); ksize_y = size(fieldy,3)
289 if(ksize_x .NE. ksize_y) call mpp_error(FATAL, &
290 'MPP_CREATE_GROUP_UPDATE_V: mismatch of ksize between fieldx and fieldy')
292 group%addrs_x(nvector) = LOC(fieldx)
293 group%addrs_y(nvector) = LOC(fieldy)
295 if( group%nvector == 1 ) then
296 group%flags_v = update_flags
297 group%whalo_v = update_whalo
298 group%ehalo_v = update_ehalo
299 group%shalo_v = update_shalo
300 group%nhalo_v = update_nhalo
301 group%gridtype = update_gridtype
302 group%isize_x = isize_x
303 group%jsize_x = jsize_x
304 group%isize_y = isize_y
305 group%jsize_y = jsize_y
306 group%ksize_v = ksize_x
307 update_edge_only = BTEST(update_flags, EDGEONLY)
308 group%nonsym_edge = .false.
310 recv(1) = BTEST(update_flags,EAST)
311 recv(3) = BTEST(update_flags,SOUTH)
312 recv(5) = BTEST(update_flags,WEST)
313 recv(7) = BTEST(update_flags,NORTH)
314 if( update_edge_only ) then
315 recv(2) = .false.
316 recv(4) = .false.
317 recv(6) = .false.
318 recv(8) = .false.
319 if( .NOT. (recv(1) .OR. recv(3) .OR. recv(5) .OR. recv(7)) ) then
320 recv(1) = .true.
321 recv(3) = .true.
322 recv(5) = .true.
323 recv(7) = .true.
324 endif
325 else
326 recv(2) = recv(1) .AND. recv(3)
327 recv(4) = recv(3) .AND. recv(5)
328 recv(6) = recv(5) .AND. recv(7)
329 recv(8) = recv(7) .AND. recv(1)
330 endif
331 group%recv_x = recv
332 group%recv_y = recv
334 !--- NONSYMEDGE is only for non-symmetric domain and CGRID/DGRID
335 if( .not. domain%symmetry .and. (update_gridtype==CGRID_NE .OR. update_gridtype==DGRID_NE)) then
336 group%nonsym_edge = BTEST(update_flags, NONSYMEDGE)
337 endif
338 if( group%nonsym_edge ) then
339 group%recv_x(2:8:2) = .false.
340 group%recv_y(2:8:2) = .false.
341 if(update_gridtype==CGRID_NE) then
342 group%recv_x(3) = .false.
343 group%recv_x(7) = .false.
344 group%recv_y(1) = .false.
345 group%recv_y(5) = .false.
346 else if(update_gridtype==DGRID_NE) then
347 group%recv_x(1) = .false.
348 group%recv_x(5) = .false.
349 group%recv_y(3) = .false.
350 group%recv_y(7) = .false.
351 endif
352 endif
354 select case(group%gridtype)
355 case (AGRID)
356 position_x = CENTER
357 position_y = CENTER
358 case (BGRID_NE, BGRID_SW)
359 position_x = CORNER
360 position_y = CORNER
361 case (CGRID_NE, CGRID_SW)
362 position_x = EAST
363 position_y = NORTH
364 case (DGRID_NE, DGRID_SW)
365 position_x = NORTH
366 position_y = EAST
367 case default
368 call mpp_error(FATAL, "mpp_CREATE_GROUP_UPDATE_V: invalid value of gridtype")
369 end select
371 call mpp_get_memory_domain(domain, group%is_x, group%ie_x, group%js_x, group%je_x, position=position_x)
372 call mpp_get_memory_domain(domain, group%is_y, group%ie_y, group%js_y, group%je_y, position=position_y)
373 else
374 set_mismatch = .false.
375 set_mismatch = set_mismatch .OR. (group%flags_v .NE. update_flags)
376 set_mismatch = set_mismatch .OR. (group%whalo_v .NE. update_whalo)
377 set_mismatch = set_mismatch .OR. (group%ehalo_v .NE. update_ehalo)
378 set_mismatch = set_mismatch .OR. (group%shalo_v .NE. update_shalo)
379 set_mismatch = set_mismatch .OR. (group%nhalo_v .NE. update_nhalo)
380 set_mismatch = set_mismatch .OR. (group%gridtype .NE. update_gridtype)
381 set_mismatch = set_mismatch .OR. (group%isize_x .NE. isize_x)
382 set_mismatch = set_mismatch .OR. (group%jsize_x .NE. jsize_x)
383 set_mismatch = set_mismatch .OR. (group%isize_y .NE. isize_y)
384 set_mismatch = set_mismatch .OR. (group%jsize_y .NE. jsize_y)
385 set_mismatch = set_mismatch .OR. (group%ksize_v .NE. ksize_x)
387 if(set_mismatch)then
388 write( text,'(i2)' ) nvector
389 call mpp_error(FATAL,'MPP_CREATE_GROUP_UPDATE_V: Incompatible field at count '//text//' for group update.' )
390 endif
391 endif
393 return
395 end subroutine MPP_CREATE_GROUP_UPDATE_3D_V_
397 subroutine MPP_CREATE_GROUP_UPDATE_4D_V_( group, fieldx, fieldy, domain, flags, gridtype, &
398 whalo, ehalo, shalo, nhalo)
400 type(mpp_group_update_type), intent(inout) :: group
401 MPP_TYPE_, intent(inout) :: fieldx(:,:,:,:), fieldy(:,:,:,:)
402 type(domain2D), intent(inout) :: domain
403 integer, intent(in), optional :: flags, gridtype
404 integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
405 MPP_TYPE_ :: field3Dx(size(fieldx,1),size(fieldx,2),size(fieldx,3)*size(fieldx,4))
406 MPP_TYPE_ :: field3Dy(size(fieldy,1),size(fieldy,2),size(fieldy,3)*size(fieldy,4))
407 pointer( ptrx, field3Dx )
408 pointer( ptry, field3Dy )
409 ptrx = LOC(fieldx)
410 ptry = LOC(fieldy)
412 call mpp_create_group_update(group, field3Dx, field3Dy, domain, flags, gridtype, &
413 whalo, ehalo, shalo, nhalo)
415 return
417 end subroutine MPP_CREATE_GROUP_UPDATE_4D_V_
420 subroutine MPP_DO_GROUP_UPDATE_(group, domain, d_type)
421 type(mpp_group_update_type), intent(inout) :: group
422 type(domain2D), intent(inout) :: domain
423 MPP_TYPE_, intent(in) :: d_type
425 integer :: nscalar, nvector, nlist
426 logical :: recv_y(8)
427 integer :: nsend, nrecv, flags_v
428 integer :: msgsize
429 integer :: from_pe, to_pe, buffer_pos, pos
430 integer :: ksize, is, ie, js, je
431 integer :: n, l, m, i, j, k, buffer_start_pos, nk
432 integer :: shift, gridtype, midpoint
433 integer :: npack, nunpack, rotation, isd
435 MPP_TYPE_ :: buffer(mpp_domains_stack_size)
436 MPP_TYPE_ :: field (group%is_s:group%ie_s,group%js_s:group%je_s, group%ksize_s)
437 MPP_TYPE_ :: fieldx(group%is_x:group%ie_x,group%js_x:group%je_x, group%ksize_v)
438 MPP_TYPE_ :: fieldy(group%is_y:group%ie_y,group%js_y:group%je_y, group%ksize_v)
439 pointer(ptr, buffer )
440 pointer(ptr_field, field)
441 pointer(ptr_fieldx, fieldx)
442 pointer(ptr_fieldy, fieldy)
444 nscalar = group%nscalar
445 nvector = group%nvector
446 nlist = size(domain%list(:))
447 gridtype = group%gridtype
449 !--- ksize_s must equal ksize_v
450 if(nvector > 0 .AND. nscalar > 0) then
451 if(group%ksize_s .NE. group%ksize_v) then
452 call mpp_error(FATAL, "MPP_DO_GROUP_UPDATE: ksize_s and ksize_v are not equal")
453 endif
454 ksize = group%ksize_s
455 else if (nscalar > 0) then
456 ksize = group%ksize_s
457 else if (nvector > 0) then
458 ksize = group%ksize_v
459 else
460 call mpp_error(FATAL, "MPP_DO_GROUP_UPDATE: nscalar and nvector are all 0")
461 endif
462 if(nvector > 0) recv_y = group%recv_y
464 ptr = LOC(mpp_domains_stack)
466 !--- set reset_index_s and reset_index_v to 0
467 group%reset_index_s = 0
468 group%reset_index_v = 0
470 if(.not. group%initialized) call set_group_update(group,domain)
472 nrecv = group%nrecv
473 nsend = group%nsend
475 !---pre-post receive.
476 call mpp_clock_begin(group_recv_clock)
477 do m = 1, nrecv
478 msgsize = group%recv_size(m)
479 from_pe = group%from_pe(m)
480 if( msgsize .GT. 0 )then
481 buffer_pos = group%buffer_pos_recv(m)
482 call mpp_recv( buffer(buffer_pos+1), glen=msgsize, from_pe=from_pe, block=.false., &
483 tag=COMM_TAG_1)
484 end if
485 end do
487 !pack the data
488 call mpp_clock_end(group_recv_clock)
490 flags_v = group%flags_v
491 npack = group%npack
493 call mpp_clock_begin(group_pack_clock)
494 !pack the data
495 buffer_start_pos = 0
496 #include <group_update_pack.inc>
497 call mpp_clock_end(group_pack_clock)
499 call mpp_clock_begin(group_send_clock)
500 do n = 1, nsend
501 msgsize = group%send_size(n)
502 if( msgsize .GT. 0 )then
503 buffer_pos = group%buffer_pos_send(n)
504 to_pe = group%to_pe(n)
505 call mpp_send( buffer(buffer_pos+1), plen=msgsize, to_pe=to_pe, tag=COMM_TAG_1)
506 endif
507 enddo
508 call mpp_clock_end(group_send_clock)
510 if(nrecv>0) then
511 call mpp_clock_begin(group_wait_clock)
512 call mpp_sync_self(check=EVENT_RECV)
513 call mpp_clock_end(group_wait_clock)
514 endif
516 !---unpack the buffer
517 nunpack = group%nunpack
518 call mpp_clock_begin(group_unpk_clock)
519 #include <group_update_unpack.inc>
520 call mpp_clock_end(group_unpk_clock)
522 ! ---northern boundary fold
523 shift = 0
524 if(domain%symmetry) shift = 1
525 if( nvector >0 .AND. BTEST(domain%fold,NORTH) .AND. (.NOT.BTEST(flags_v,SCALAR_BIT)) )then
526 j = domain%y(1)%global%end+shift
527 if( domain%y(1)%data%begin.LE.j .AND. j.LE.domain%y(1)%data%end+shift )then !fold is within domain
528 !poles set to 0: BGRID only
529 if( gridtype.EQ.BGRID_NE )then
530 midpoint = (domain%x(1)%global%begin+domain%x(1)%global%end-1+shift)/2
531 j = domain%y(1)%global%end+shift
532 is = domain%x(1)%global%begin; ie = domain%x(1)%global%end+shift
533 if( .NOT. domain%symmetry ) is = is - 1
534 do i = is ,ie, midpoint
535 if( domain%x(1)%data%begin.LE.i .AND. i.LE. domain%x(1)%data%end+shift )then
536 do l=1,nvector
537 ptr_fieldx = group%addrs_x(l)
538 ptr_fieldy = group%addrs_y(l)
539 do k = 1,ksize
540 fieldx(i,j,k) = 0.
541 fieldy(i,j,k) = 0.
542 end do
543 end do
544 end if
545 end do
546 endif
547 ! the following code code block correct an error where the data in your halo coming from
548 ! other half may have the wrong sign
549 !off west edge, when update north or west direction
550 j = domain%y(1)%global%end+shift
551 if ( recv_y(7) .OR. recv_y(5) ) then
552 select case(gridtype)
553 case(BGRID_NE)
554 if(domain%symmetry) then
555 is = domain%x(1)%global%begin
556 else
557 is = domain%x(1)%global%begin - 1
558 end if
559 if( is.GT.domain%x(1)%data%begin )then
561 if( 2*is-domain%x(1)%data%begin.GT.domain%x(1)%data%end+shift ) &
562 call mpp_error( FATAL, 'MPP_DO_UPDATE_V: folded-north BGRID_NE west edge ubound error.' )
563 do l=1,nvector
564 ptr_fieldx = group%addrs_x(l)
565 ptr_fieldy = group%addrs_y(l)
566 do k = 1,ksize
567 do i = domain%x(1)%data%begin,is-1
568 fieldx(i,j,k) = fieldx(2*is-i,j,k)
569 fieldy(i,j,k) = fieldy(2*is-i,j,k)
570 end do
571 end do
572 end do
573 end if
574 case(CGRID_NE)
575 is = domain%x(1)%global%begin
576 isd = domain%x(1)%compute%begin - group%whalo_v
577 if( is.GT.isd )then
578 if( 2*is-domain%x(1)%data%begin-1.GT.domain%x(1)%data%end ) &
579 call mpp_error( FATAL, 'MPP_DO_UPDATE_V: folded-north CGRID_NE west edge ubound error.' )
580 do l=1,nvector
581 ptr_fieldy = group%addrs_y(l)
582 do k = 1,ksize
583 do i = isd,is-1
584 fieldy(i,j,k) = fieldy(2*is-i-1,j,k)
585 end do
586 end do
587 end do
588 end if
589 end select
590 end if
591 !off east edge
592 is = domain%x(1)%global%end
593 if(domain%x(1)%cyclic .AND. is.LT.domain%x(1)%data%end )then
594 ie = domain%x(1)%compute%end+group%ehalo_v
595 is = is + 1
596 select case(gridtype)
597 case(BGRID_NE)
598 is = is + shift
599 ie = ie + shift
600 do l=1,nvector
601 ptr_fieldx = group%addrs_x(l)
602 ptr_fieldy = group%addrs_y(l)
603 do k = 1,ksize
604 do i = is,ie
605 fieldx(i,j,k) = -fieldx(i,j,k)
606 fieldy(i,j,k) = -fieldy(i,j,k)
607 end do
608 end do
609 end do
610 case(CGRID_NE)
611 do l=1,nvector
612 ptr_fieldy = group%addrs_y(l)
613 do k = 1,ksize
614 do i = is, ie
615 fieldy(i,j,k) = -fieldy(i,j,k)
616 end do
617 end do
618 end do
619 end select
620 end if
621 end if
622 else if( BTEST(domain%fold,SOUTH) .OR. BTEST(domain%fold,WEST) .OR. BTEST(domain%fold,EAST) ) then
623 call mpp_error(FATAL, "MPP_DO_GROUP_UPDATE: this interface does not support folded_south, " // &
624 "folded_west of folded_east, contact developer")
625 endif
627 if(nsend>0) then
628 call mpp_clock_begin(group_wait_clock)
629 call mpp_sync_self( )
630 call mpp_clock_end(group_wait_clock)
631 endif
633 end subroutine MPP_DO_GROUP_UPDATE_
636 subroutine MPP_START_GROUP_UPDATE_(group, domain, d_type, reuse_buffer)
637 type(mpp_group_update_type), intent(inout) :: group
638 type(domain2D), intent(inout) :: domain
639 MPP_TYPE_, intent(in) :: d_type
640 logical, optional, intent(in) :: reuse_buffer
642 integer :: nscalar, nvector
643 integer :: nsend, nrecv, flags_v
644 integer :: msgsize, npack, rotation
645 integer :: from_pe, to_pe, buffer_pos, pos
646 integer :: ksize, is, ie, js, je
647 integer :: n, l, m, i, j, k, buffer_start_pos, nk
648 logical :: reuse_buf_pos
649 character(len=8) :: text
651 MPP_TYPE_ :: buffer(size(mpp_domains_stack_nonblock(:)))
652 MPP_TYPE_ :: field (group%is_s:group%ie_s,group%js_s:group%je_s, group%ksize_s)
653 MPP_TYPE_ :: fieldx(group%is_x:group%ie_x,group%js_x:group%je_x, group%ksize_v)
654 MPP_TYPE_ :: fieldy(group%is_y:group%ie_y,group%js_y:group%je_y, group%ksize_v)
655 pointer( ptr, buffer )
656 pointer(ptr_field, field)
657 pointer(ptr_fieldx, fieldx)
658 pointer(ptr_fieldy, fieldy)
660 nscalar = group%nscalar
661 nvector = group%nvector
663 if(nscalar>0) then
664 ksize = group%ksize_s
665 else
666 ksize = group%ksize_v
667 endif
669 !--- set reset_index_s and reset_index_v to 0
670 group%reset_index_s = 0
671 group%reset_index_v = 0
673 reuse_buf_pos = .FALSE.
674 if (PRESENT(reuse_buffer)) reuse_buf_pos = reuse_buffer
676 if (.not. group%initialized) then
677 call set_group_update(group,domain)
678 endif
680 if (.not. reuse_buf_pos) then
681 group%buffer_start_pos = nonblock_group_buffer_pos
682 nonblock_group_buffer_pos = nonblock_group_buffer_pos + group%tot_msgsize
683 mpp_domains_stack_hwm = nonblock_group_buffer_pos + 1
684 if( mpp_domains_stack_hwm .GT. mpp_domains_stack_size )then
685 write( text,'(i8)' )mpp_domains_stack_hwm
686 call mpp_error( FATAL, 'set_group_update: mpp_domains_stack overflow, '// &
687 'call mpp_domains_set_stack_size('//trim(text)//') from all PEs.' )
688 end if
690 else if( group%buffer_start_pos < 0 ) then
691 call mpp_error(FATAL, "MPP_START_GROUP_UPDATE: group%buffer_start_pos is not set")
692 endif
694 nrecv = group%nrecv
695 nsend = group%nsend
697 ptr = LOC(mpp_domains_stack_nonblock)
699 ! Make sure it is not in the middle of the old version of non-blocking halo update.
700 if(num_update>0) call mpp_error(FATAL, "MPP_START_GROUP_UPDATE: can not be called in the middle of "// &
701 "mpp_start_update_domains/mpp_complete_update_domains call")
703 num_nonblock_group_update = num_nonblock_group_update + 1
705 !---pre-post receive.
706 call mpp_clock_begin(nonblock_group_recv_clock)
707 do m = 1, nrecv
708 msgsize = group%recv_size(m)
709 from_pe = group%from_pe(m)
710 if( msgsize .GT. 0 )then
711 buffer_pos = group%buffer_pos_recv(m) + group%buffer_start_pos
712 call mpp_recv( buffer(buffer_pos+1), glen=msgsize, from_pe=from_pe, block=.false., &
713 tag=COMM_TAG_1, request=group%request_recv(m))
714 #ifdef use_libMPI
715 group%type_recv(m) = MPI_TYPE_
716 #endif
717 end if
718 end do
719 call mpp_clock_end(nonblock_group_recv_clock)
721 flags_v = group%flags_v
723 !pack the data
724 call mpp_clock_begin(nonblock_group_pack_clock)
725 npack = group%npack
726 buffer_start_pos = group%buffer_start_pos
727 #include <group_update_pack.inc>
728 call mpp_clock_end(nonblock_group_pack_clock)
730 call mpp_clock_begin(nonblock_group_send_clock)
731 do n = 1, nsend
732 msgsize = group%send_size(n)
733 if( msgsize .GT. 0 )then
734 buffer_pos = group%buffer_pos_send(n) + group%buffer_start_pos
735 to_pe = group%to_pe(n)
736 call mpp_send( buffer(buffer_pos+1), plen=msgsize, to_pe=to_pe, tag=COMM_TAG_1, &
737 request=group%request_send(n))
738 endif
739 enddo
740 call mpp_clock_end(nonblock_group_send_clock)
742 end subroutine MPP_START_GROUP_UPDATE_
744 subroutine MPP_COMPLETE_GROUP_UPDATE_(group, domain, d_type)
745 type(mpp_group_update_type), intent(inout) :: group
746 type(domain2D), intent(inout) :: domain
747 MPP_TYPE_, intent(in) :: d_type
749 integer :: nsend, nrecv, nscalar, nvector
750 integer :: k, buffer_pos, pos, m, n, l
751 integer :: is, ie, js, je, ksize, i, j
752 integer :: shift, gridtype, midpoint, flags_v
753 integer :: nunpack, rotation, buffer_start_pos, nk, isd
754 logical :: recv_y(8)
755 MPP_TYPE_ :: buffer(size(mpp_domains_stack_nonblock(:)))
756 MPP_TYPE_ :: field (group%is_s:group%ie_s,group%js_s:group%je_s, group%ksize_s)
757 MPP_TYPE_ :: fieldx(group%is_x:group%ie_x,group%js_x:group%je_x, group%ksize_v)
758 MPP_TYPE_ :: fieldy(group%is_y:group%ie_y,group%js_y:group%je_y, group%ksize_v)
759 pointer(ptr, buffer )
760 pointer(ptr_field, field)
761 pointer(ptr_fieldx, fieldx)
762 pointer(ptr_fieldy, fieldy)
764 gridtype = group%gridtype
765 flags_v = group%flags_v
766 nscalar = group%nscalar
767 nvector = group%nvector
768 nrecv = group%nrecv
769 nsend = group%nsend
770 if(nscalar>0) then
771 ksize = group%ksize_s
772 else
773 ksize = group%ksize_v
774 endif
775 if(nvector > 0) recv_y = group%recv_y
776 ptr = LOC(mpp_domains_stack_nonblock)
778 if(num_nonblock_group_update < 1) call mpp_error(FATAL, &
779 'mpp_start_group_update must be called before calling mpp_end_group_update')
780 num_nonblock_group_update = num_nonblock_group_update - 1
781 complete_group_update_on = .true.
783 if(nrecv>0) then
784 call mpp_clock_begin(nonblock_group_wait_clock)
785 call mpp_sync_self(check=EVENT_RECV, request=group%request_recv(1:nrecv), &
786 msg_size=group%recv_size(1:nrecv), msg_type=group%type_recv(1:nrecv))
787 call mpp_clock_end(nonblock_group_wait_clock)
788 endif
790 !---unpack the buffer
791 nunpack = group%nunpack
793 call mpp_clock_begin(nonblock_group_unpk_clock)
794 buffer_start_pos = group%buffer_start_pos
795 #include <group_update_unpack.inc>
796 call mpp_clock_end(nonblock_group_unpk_clock)
798 ! ---northern boundary fold
799 shift = 0
800 if(domain%symmetry) shift = 1
801 if( nvector >0 .AND. BTEST(domain%fold,NORTH) .AND. (.NOT.BTEST(flags_v,SCALAR_BIT)) )then
802 j = domain%y(1)%global%end+shift
803 if( domain%y(1)%data%begin.LE.j .AND. j.LE.domain%y(1)%data%end+shift )then !fold is within domain
804 !poles set to 0: BGRID only
805 if( gridtype.EQ.BGRID_NE )then
806 midpoint = (domain%x(1)%global%begin+domain%x(1)%global%end-1+shift)/2
807 j = domain%y(1)%global%end+shift
808 is = domain%x(1)%global%begin; ie = domain%x(1)%global%end+shift
809 if( .NOT. domain%symmetry ) is = is - 1
810 do i = is ,ie, midpoint
811 if( domain%x(1)%data%begin.LE.i .AND. i.LE. domain%x(1)%data%end+shift )then
812 do l=1,nvector
813 ptr_fieldx = group%addrs_x(l)
814 ptr_fieldy = group%addrs_y(l)
815 do k = 1,ksize
816 fieldx(i,j,k) = 0.
817 fieldy(i,j,k) = 0.
818 end do
819 end do
820 end if
821 end do
822 endif
823 ! the following code code block correct an error where the data in your halo coming from
824 ! other half may have the wrong sign
825 !off west edge, when update north or west direction
826 j = domain%y(1)%global%end+shift
827 if ( recv_y(7) .OR. recv_y(5) ) then
828 select case(gridtype)
829 case(BGRID_NE)
830 if(domain%symmetry) then
831 is = domain%x(1)%global%begin
832 else
833 is = domain%x(1)%global%begin - 1
834 end if
835 if( is.GT.domain%x(1)%data%begin )then
837 if( 2*is-domain%x(1)%data%begin.GT.domain%x(1)%data%end+shift ) &
838 call mpp_error( FATAL, 'MPP_DO_UPDATE_V: folded-north BGRID_NE west edge ubound error.' )
839 do l=1,nvector
840 ptr_fieldx = group%addrs_x(l)
841 ptr_fieldy = group%addrs_y(l)
842 do k = 1,ksize
843 do i = domain%x(1)%data%begin,is-1
844 fieldx(i,j,k) = fieldx(2*is-i,j,k)
845 fieldy(i,j,k) = fieldy(2*is-i,j,k)
846 end do
847 end do
848 end do
849 end if
850 case(CGRID_NE)
851 is = domain%x(1)%global%begin
852 isd = domain%x(1)%compute%begin - group%whalo_v
853 if( is.GT.isd)then
854 if( 2*is-domain%x(1)%data%begin-1.GT.domain%x(1)%data%end ) &
855 call mpp_error( FATAL, 'MPP_DO_UPDATE_V: folded-north CGRID_NE west edge ubound error.' )
856 do l=1,nvector
857 ptr_fieldy = group%addrs_y(l)
858 do k = 1,ksize
859 do i = isd,is-1
860 fieldy(i,j,k) = fieldy(2*is-i-1,j,k)
861 end do
862 end do
863 end do
864 end if
865 end select
866 end if
867 !off east edge
868 is = domain%x(1)%global%end
869 if(domain%x(1)%cyclic .AND. is.LT.domain%x(1)%data%end )then
870 ie = domain%x(1)%compute%end+group%ehalo_v
871 is = is + 1
872 select case(gridtype)
873 case(BGRID_NE)
874 is = is + shift
875 ie = ie + shift
876 do l=1,nvector
877 ptr_fieldx = group%addrs_x(l)
878 ptr_fieldy = group%addrs_y(l)
879 do k = 1,ksize
880 do i = is,ie
881 fieldx(i,j,k) = -fieldx(i,j,k)
882 fieldy(i,j,k) = -fieldy(i,j,k)
883 end do
884 end do
885 end do
886 case(CGRID_NE)
887 do l=1,nvector
888 ptr_fieldy = group%addrs_y(l)
889 do k = 1,ksize
890 do i = is, ie
891 fieldy(i,j,k) = -fieldy(i,j,k)
892 end do
893 end do
894 end do
895 end select
896 end if
897 end if
898 else if( BTEST(domain%fold,SOUTH) .OR. BTEST(domain%fold,WEST) .OR. BTEST(domain%fold,EAST) ) then
899 call mpp_error(FATAL, "MPP_COMPLETE_GROUP_UPDATE: this interface does not support folded_south, " // &
900 "folded_west of folded_east, contact developer")
901 endif
903 if(nsend>0) then
904 call mpp_clock_begin(nonblock_group_wait_clock)
905 call mpp_sync_self(check=EVENT_SEND, request=group%request_send(1:nsend) )
906 call mpp_clock_end(nonblock_group_wait_clock)
907 endif
909 if( num_nonblock_group_update == 0) then
910 nonblock_group_buffer_pos = 0
911 endif
913 end subroutine MPP_COMPLETE_GROUP_UPDATE_
915 subroutine MPP_RESET_GROUP_UPDATE_FIELD_2D_(group, field)
916 type(mpp_group_update_type), intent(inout) :: group
917 MPP_TYPE_, intent(in) :: field(:,:)
919 group%reset_index_s = group%reset_index_s + 1
921 if(group%reset_index_s > group%nscalar) &
922 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_2D_: group%reset_index_s > group%nscalar")
923 if(size(field,1) .NE. group%isize_s .OR. size(field,2) .NE. group%jsize_s .OR. group%ksize_s .NE. 1) &
924 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_2D_: size of field does not match the size stored in group")
926 group%addrs_s(group%reset_index_s) = LOC(field)
928 end subroutine MPP_RESET_GROUP_UPDATE_FIELD_2D_
930 subroutine MPP_RESET_GROUP_UPDATE_FIELD_3D_(group, field)
931 type(mpp_group_update_type), intent(inout) :: group
932 MPP_TYPE_, intent(in) :: field(:,:,:)
934 group%reset_index_s = group%reset_index_s + 1
936 if(group%reset_index_s > group%nscalar) &
937 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_3D_: group%reset_index_s > group%nscalar")
938 if(size(field,1) .NE. group%isize_s .OR. size(field,2) .NE. group%jsize_s .OR. size(field,3) .NE. group%ksize_s) &
939 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_3D_: size of field does not match the size stored in group")
941 group%addrs_s(group%reset_index_s) = LOC(field)
943 end subroutine MPP_RESET_GROUP_UPDATE_FIELD_3D_
945 subroutine MPP_RESET_GROUP_UPDATE_FIELD_4D_(group, field)
946 type(mpp_group_update_type), intent(inout) :: group
947 MPP_TYPE_, intent(in) :: field(:,:,:,:)
949 group%reset_index_s = group%reset_index_s + 1
951 if(group%reset_index_s > group%nscalar) &
952 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_4D_: group%reset_index_s > group%nscalar")
953 if(size(field,1) .NE. group%isize_s .OR. size(field,2) .NE. group%jsize_s .OR. &
954 size(field,3)*size(field,4) .NE. group%ksize_s) &
955 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_4D_: size of field does not match the size stored in group")
957 group%addrs_s(group%reset_index_s) = LOC(field)
959 end subroutine MPP_RESET_GROUP_UPDATE_FIELD_4D_
962 subroutine MPP_RESET_GROUP_UPDATE_FIELD_2D_V_(group, fieldx, fieldy)
963 type(mpp_group_update_type), intent(inout) :: group
964 MPP_TYPE_, intent(in) :: fieldx(:,:), fieldy(:,:)
966 group%reset_index_v = group%reset_index_v + 1
968 if(group%reset_index_v > group%nvector) &
969 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_2D_V_: group%reset_index_v > group%nvector")
970 if(size(fieldx,1) .NE. group%isize_x .OR. size(fieldx,2) .NE. group%jsize_x .OR. group%ksize_v .NE. 1) &
971 call mpp_error(FATAL, &
972 & "MPP_RESET_GROUP_UPDATE_FIELD_2D_V_: size of fieldx does not match the size stored in group")
973 if(size(fieldy,1) .NE. group%isize_y .OR. size(fieldy,2) .NE. group%jsize_y ) &
974 call mpp_error(FATAL, &
975 & "MPP_RESET_GROUP_UPDATE_FIELD_2D_V_: size of fieldy does not match the size stored in group")
977 group%addrs_x(group%reset_index_v) = LOC(fieldx)
978 group%addrs_y(group%reset_index_v) = LOC(fieldy)
980 end subroutine MPP_RESET_GROUP_UPDATE_FIELD_2D_V_
983 subroutine MPP_RESET_GROUP_UPDATE_FIELD_3D_V_(group, fieldx, fieldy)
984 type(mpp_group_update_type), intent(inout) :: group
985 MPP_TYPE_, intent(in) :: fieldx(:,:,:), fieldy(:,:,:)
987 group%reset_index_v = group%reset_index_v + 1
989 if(group%reset_index_v > group%nvector) &
990 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_3D_V_: group%reset_index_v > group%nvector")
991 if(size(fieldx,1) .NE. group%isize_x .OR. size(fieldx,2) .NE. group%jsize_x .OR. size(fieldx,3) .NE. group%ksize_v) &
992 call mpp_error(FATAL, &
993 & "MPP_RESET_GROUP_UPDATE_FIELD_3D_V_: size of fieldx does not match the size stored in group")
994 if(size(fieldy,1) .NE. group%isize_y .OR. size(fieldy,2) .NE. group%jsize_y .OR. size(fieldy,3) .NE. group%ksize_v) &
995 call mpp_error(FATAL, &
996 & "MPP_RESET_GROUP_UPDATE_FIELD_3D_V_: size of fieldy does not match the size stored in group")
998 group%addrs_x(group%reset_index_v) = LOC(fieldx)
999 group%addrs_y(group%reset_index_v) = LOC(fieldy)
1001 end subroutine MPP_RESET_GROUP_UPDATE_FIELD_3D_V_
1004 subroutine MPP_RESET_GROUP_UPDATE_FIELD_4D_V_(group, fieldx, fieldy)
1005 type(mpp_group_update_type), intent(inout) :: group
1006 MPP_TYPE_, intent(in) :: fieldx(:,:,:,:), fieldy(:,:,:,:)
1008 group%reset_index_v = group%reset_index_v + 1
1010 if(group%reset_index_v > group%nvector) &
1011 call mpp_error(FATAL, "MPP_RESET_GROUP_UPDATE_FIELD_4D_V_: group%reset_index_v > group%nvector")
1012 if(size(fieldx,1) .NE. group%isize_x .OR. size(fieldx,2) .NE. group%jsize_x .OR. &
1013 size(fieldx,3)*size(fieldx,4) .NE. group%ksize_v) &
1014 call mpp_error(FATAL, &
1015 & "MPP_RESET_GROUP_UPDATE_FIELD_4D_V_: size of fieldx does not match the size stored in group")
1016 if(size(fieldy,1) .NE. group%isize_y .OR. size(fieldy,2) .NE. group%jsize_y .OR. &
1017 size(fieldy,3)*size(fieldy,4) .NE. group%ksize_v) &
1018 call mpp_error(FATAL, &
1019 & "MPP_RESET_GROUP_UPDATE_FIELD_4D_V_: size of fieldy does not match the size stored in group")
1021 group%addrs_x(group%reset_index_v) = LOC(fieldx)
1022 group%addrs_y(group%reset_index_v) = LOC(fieldy)
1024 end subroutine MPP_RESET_GROUP_UPDATE_FIELD_4D_V_