30 #include <fms_platform.h> 33 mpp_declare_pelist, mpp_get_current_pelist, mpp_set_current_pelist, &
34 mpp_init, comm_tag_1, comm_tag_2, comm_tag_3, mpp_sync_self
45 #include <mpp/shmem.fh> 46 integer :: psync(shmem_barrier_sync_size)
47 pointer( p_psync, psync )
75 integer :: next_in_pset, prev_in_pset
76 integer :: root_in_pset
81 real, _allocatable :: stack(:) _null
82 integer, _allocatable :: pelist(:) _null
83 integer, _allocatable :: root_pelist(:) _null
84 integer, _allocatable :: pset(:) _null
85 integer(POINTER_KIND) :: p_stack
86 integer :: lstack, maxstack, hiwm
88 character(len=32) :: name
89 logical :: initialized=.false.
107 character(len=8) ::
value 108 integer :: lenname, lenval
112 call shmem_barrier_all()
113 call shpalloc( p_psync, shmem_barrier_sync_size, err, -1 )
114 call shmem_barrier_all()
116 call pxfgetenv(
'SMA_GLOBAL_ALLOC', 0,
value, lenval, err )
118 'The environment variable SMA_GLOBAL_ALLOC must be set on Irix.' )
128 integer,
intent(in) :: npset
130 integer,
intent(in),
optional :: stacksize
131 integer,
intent(in),
optional :: pelist(:)
132 integer,
intent(in),
optional ::
commid 134 integer :: npes, my_commid
135 integer :: i, j, k, out_unit, errunit
136 integer,
allocatable :: my_pelist(:), root_pelist(:)
147 if(
present(pelist))
then 148 npes =
size(pelist(:))
152 if( mod(npes,npset).NE.0 )
then 153 write(
text,
'(a,2i6)' ) &
154 'MPP_PSET_CREATE: PSET size (npset) must divide npes exactly:'// &
155 ' npset, npes=', npset, npes
160 allocate(my_pelist(0:npes-1) )
161 allocate(root_pelist(0:npes/npset-1) )
162 if(
present(pelist))
then 164 'MPP_PSET_CREATE: when pelist is present, commID should also be present')
168 call mpp_get_current_pelist(my_pelist,
commid = my_commid)
170 do i = 0,npes/npset-1
171 root_pelist(i) = my_pelist(npset*i)
173 write( out_unit,
'(a,i6)' )
'MPP_PSET_CREATE creating PSETs... npset=', npset
174 if(any(my_pelist == pe) )
then 175 if( pset%initialized )
call mpp_error( fatal, &
176 'MPP_PSET_CREATE: PSET already initialized!' )
178 allocate( pset%pelist(0:npes-1) )
179 allocate( pset%root_pelist(0:npes/npset-1) )
180 pset%commID = my_commid
181 pset%pelist = my_pelist
183 pset%root_pelist = root_pelist
184 allocate( pset%pset(0:npset-1) )
185 do i = 0,npes/npset-1
189 if( pe.EQ.pset%pelist(k+j) )
then 190 pset%pset(:) = pset%pelist(k:k+npset-1)
192 pset%root_in_pset = pset%root_pelist(i)
194 pset%prev_in_pset = pset%pelist(k+npset-1)
196 pset%prev_in_pset = pset%pelist(k+j-1)
198 if( j.EQ.npset-1 )
then 199 pset%next_in_pset = pset%pelist(k)
201 pset%next_in_pset = pset%pelist(k+j+1)
207 pset%root = pe.EQ.pset%root_in_pset
211 pset%maxstack = 1000000
212 if(
PRESENT(stacksize) )pset%maxstack = stacksize
213 write( out_unit,
'(a,i8)' ) &
214 'MPP_PSET_CREATE: setting stacksize=', pset%maxstack
216 allocate( pset%stack(pset%maxstack) )
217 #ifdef use_CRI_pointers 218 pset%p_stack = loc(pset%stack)
221 pset%initialized = .true.
225 call mpp_declare_pelist(root_pelist)
228 write( errunit,
'(a,4i6)' )
'MPP_PSET_CREATE: pe, root, next, prev=', &
229 pe, pset%root_in_pset, pset%next_in_pset, pset%prev_in_pset
230 write( errunit,* )
'PE ', pe,
' pset=', pset%pset(:)
231 write( out_unit,* )
'root pelist=', pset%root_pelist(:)
240 if( .NOT.pset%initialized )
call mpp_error( fatal, &
241 'MPP_PSET_DELETE: called with uninitialized PSET.' )
243 deallocate( pset%pelist )
244 deallocate( pset%root_pelist )
245 deallocate( pset%pset )
246 if( pset%root )
deallocate( pset%stack )
247 write( out_unit,
'(a,i10)' ) &
248 'Deleting PSETs... stack high-water-mark=', pset%hiWM
250 pset%initialized = .false.
254 integer(POINTER_KIND),
intent(in) :: ptr
255 integer,
intent(in) :: pe
259 call mpp_send( ptr, pe, tag=comm_tag_1 )
263 integer(POINTER_KIND),
intent(in) :: ptr(:)
264 integer,
intent(in) :: pe
268 call mpp_send( ptr,
size(ptr), pe, tag=comm_tag_2 )
272 integer(POINTER_KIND),
intent(inout) :: ptr
273 integer,
intent(in) :: pe
275 call mpp_recv( ptr, pe, tag=comm_tag_1 )
281 integer(POINTER_KIND),
intent(inout) :: ptr(:)
282 integer,
intent(in) :: pe
285 call mpp_recv( ptr,
size(ptr), pe, tag=comm_tag_2 )
294 integer(POINTER_KIND),
intent(inout) :: ptr
295 integer,
intent(in) :: pe
311 integer(POINTER_KIND) :: length=1
320 call mpi_sgi_globalptr( dummy, length, pe,
commid, ptr, ierror )
321 if( ierror.EQ.-1 )
call mpp_error( fatal, &
322 'MPP_TRANSLATE_REMOTE_PTR: unknown MPI_SGI_GLOBALPTR error.' )
325 'MPP_TRANSLATE_REMOTE_PTR now only works under -Duse_libMPI' )
336 if( .NOT.pset%initialized )
call mpp_error( fatal, &
337 'MPP_PSET_SYNC: called with uninitialized PSET.' )
340 call shmem_barrier( pset%root_in_pset, 0, pset%npset, psync )
351 real,
intent(inout) :: a
354 if( .NOT.pset%initialized )
call mpp_error( fatal, &
355 'MPP_PSET_BROADCAST: called with uninitialized PSET.' )
357 do i = 1,pset%npset-1
358 call mpp_send( a, pset%pset(i), tag=comm_tag_3 )
361 call mpp_recv( a, pset%root_in_pset, tag=comm_tag_3 )
370 type(mpp_pset_type),
intent(in) :: pset
371 integer(POINTER_KIND),
intent(inout) :: ptr
374 if( .NOT.pset%initialized )
call mpp_error( fatal, &
375 'MPP_PSET_BROADCAST_PTR: called with uninitialized PSET.' )
378 do i = 1,pset%npset-1
391 type(mpp_pset_type),
intent(in) :: pset
392 integer(POINTER_KIND),
intent(inout) :: ptr(:)
395 if( .NOT.pset%initialized )
call mpp_error( fatal, &
396 'MPP_PSET_BROADCAST_PTR: called with uninitialized PSET.' )
399 do i = 1,pset%npset-1
412 #ifdef use_CRI_pointers 414 pointer( ptr, dummy )
416 integer(POINTER_KIND),
intent(in) :: ptr
419 integer(POINTER_KIND) :: p
421 if( .NOT.pset%initialized )
call mpp_error( fatal, &
422 'MPP_PSET_CHECK_PTR: called with uninitialized PSET.' )
427 do i = 1,pset%npset-1
435 'MPP_PSET_CHECK_PTR: pointers do not match!' )
447 integer,
intent(in) :: ls, le
448 integer,
intent(out) :: lsp, lep
451 if( .NOT.pset%initialized )
call mpp_error( fatal, &
452 'MPP_PSET_SEGMENT_ARRAY: called with uninitialized PSET.' )
454 if( le-ls+1.LT.pset%npset )
then 455 write(
text,
'(3(a,i6))' ) &
456 'MPP_PSET_ARRAY_SEGMENT: parallel range (', ls,
',', le, &
457 ') is smaller than the number of threads:', pset%npset
464 lep = lsp + ceiling(
REAL(le-lsp+1)/(pset%npset-i) ) - 1
474 integer,
intent(in) :: len
475 #ifdef use_CRI_pointers 477 pointer( ptr, dummy )
478 real :: stack(pset%maxstack)
481 if( .NOT.pset%initialized )
call mpp_error( fatal, &
482 'MPP_PSET_STACK_PUSH: called with uninitialized PSET.' )
483 if( pset%lstack+len.GT.pset%maxstack )
then 484 write(
text,
'(a,3i12)' ) &
485 'MPP_PSET_STACK_PUSH: mpp_pset_stack overflow: '// &
486 .GT.
'len+lstackmaxstack. len, lstack, maxstack=', &
487 len, pset%lstack, pset%maxstack
491 ptr = loc( stack(pset%lstack+1) )
493 pset%lstack = pset%lstack + len
494 pset%hiWM =
max( pset%hiWM, pset%lstack )
496 integer(POINTER_KIND),
intent(out) :: ptr
498 'MPP_PSET_STACK_PUSH only works with Cray pointers.' )
513 if( .NOT.pset%initialized )
call mpp_error( fatal, &
514 'MPP_PSET_STACK_RESET: called with uninitialized PSET.' )
522 type(mpp_pset_type),
intent(in) :: pset
523 character(len=*),
intent(in) :: caller
524 real,
intent(in) :: array(:)
529 integer(LONG_KIND) :: chksum
531 if( .NOT.pset%initialized )
call mpp_error( fatal, &
532 'MPP_PSET_PRINT_CHKSUM: called with uninitialized PSET.' )
536 do_print = pe.EQ.mpp_root_pe()
537 call mpp_set_current_pelist(pset%root_pelist)
540 write( errunit,
'(a,z18)' )trim(caller)//
' chksum=', chksum
542 call mpp_set_current_pelist(pset%pelist)
548 type(mpp_pset_type),
intent(in) :: pset
549 character(len=*),
intent(in) :: caller
550 real,
intent(in) :: array(:,:)
551 real :: array1D( size(array) )
552 #ifdef use_CRI_pointers 553 pointer( p, array1d )
556 array1d = transfer( array, array1d )
562 type(mpp_pset_type),
intent(in) :: pset
563 character(len=*),
intent(in) :: caller
564 real,
intent(in) :: array(:,:,:)
565 real :: array1D( size(array) )
566 #ifdef use_CRI_pointers 567 pointer( p, array1d )
570 array1d = transfer( array, array1d )
576 type(mpp_pset_type),
intent(in) :: pset
577 character(len=*),
intent(in) :: caller
578 real,
intent(in) :: array(:,:,:,:)
579 real :: array1D( size(array) )
580 #ifdef use_CRI_pointers 581 pointer( p, array1d )
584 array1d = transfer( array, array1d )
591 character(len=*),
intent(in) :: caller
593 if( .NOT.pset%initialized )
call mpp_error( fatal, &
594 'MPP_PSET_PRINT_STACK_CHKSUM: called with uninitialized PSET.' )
596 pset%stack(1:pset%lstack) )
604 if( .NOT.pset%initialized )
call mpp_error( fatal, &
605 'MPP_PSET_ROOT: called with uninitialized PSET.' )
614 if( .NOT.pset%initialized )
call mpp_error( fatal, &
615 'MPP_PSET_NUMROOTS: called with uninitialized PSET.' )
621 integer,
intent(out) :: pelist(:)
622 integer,
intent(out),
optional ::
commid 624 if( .NOT.pset%initialized )
call mpp_error( fatal, &
625 'MPP_PSET_GET_ROOT_PELIST: called with uninitialized PSET.' )
626 if(
size(pelist).NE.
size(pset%root_pelist) )
then 627 write(
text,
'(a,2i6)' ) &
628 'pelist argument has wrong size: requested, actual=', &
629 size(pelist),
size(pset%root_pelist)
632 pelist(:) = pset%root_pelist(:)
638 'MPP_PSET_GET_ROOT_PELIST: commID is only defined under -Duse_libMPI.' )
subroutine, public mpp_pset_stack_reset(pset)
subroutine mpp_translate_remote_ptr(ptr, pe)
subroutine mpp_send_ptr_array(ptr, pe)
subroutine, public mpp_pset_print_stack_chksum(pset, caller)
logical module_is_initialized
subroutine, public mpp_pset_segment_array(pset, ls, le, lsp, lep)
subroutine, public mpp_pset_stack_push(pset, ptr, len)
subroutine, public mpp_pset_check_ptr(pset, ptr)
subroutine mpp_send_ptr_scalar(ptr, pe)
subroutine mpp_pset_broadcast_ptr_scalar(pset, ptr)
subroutine, public mpp_pset_create(npset, pset, stacksize, pelist, commID)
subroutine, public mpp_pset_init
subroutine mpp_pset_print_chksum_1d(pset, caller, array)
subroutine, public mpp_pset_delete(pset)
subroutine mpp_recv_ptr_array(ptr, pe)
logical function, public mpp_pset_root(pset)
subroutine mpp_pset_print_chksum_2d(pset, caller, array)
integer function, public mpp_pset_numroots(pset)
subroutine, public mpp_pset_get_root_pelist(pset, pelist, commID)
subroutine, public mpp_pset_sync(pset)
subroutine, public mpp_pset_broadcast(pset, a)
subroutine mpp_recv_ptr_scalar(ptr, pe)
subroutine mpp_pset_print_chksum_3d(pset, caller, array)
subroutine mpp_pset_print_chksum_4d(pset, caller, array)
subroutine mpp_pset_broadcast_ptr_array(pset, ptr)