! -*-f90-*- !*********************************************************************** !* GNU Lesser General Public License !* !* This file is part of the GFDL Flexible Modeling System (FMS). !* !* FMS is free software: you can redistribute it and/or modify it under !* the terms of the GNU Lesser General Public License as published by !* the Free Software Foundation, either version 3 of the License, or (at !* your option) any later version. !* !* FMS is distributed in the hope that it will be useful, but WITHOUT !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License !* for more details. !* !* You should have received a copy of the GNU Lesser General Public !* License along with FMS. If not, see . !*********************************************************************** !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! ! MISCELLANEOUS UTILITIES: mpp_error ! ! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine mpp_error_basic( errortype, errormsg ) !a very basic error handler !uses ABORT and FLUSH calls, may need to use cpp to rename integer, intent(in) :: errortype character(len=*), intent(in), optional :: errormsg character(len=512) :: text logical :: opened integer :: istat, errunit, outunit if( .NOT.module_is_initialized )call ABORT() select case( errortype ) case(NOTE) text = 'NOTE' !just FYI case(WARNING) text = 'WARNING' !probable error case(FATAL) text = 'FATAL' !fatal error case default text = 'WARNING: non-existent errortype (must be NOTE|WARNING|FATAL)' end select if( npes.GT.1 )write( text,'(a,i5)' )trim(text)//' from PE', pe !this is the mpp part if( PRESENT(errormsg) )text = trim(text)//': '//trim(errormsg) outunit = stdout() errunit = stderr() select case( errortype ) case(NOTE) write( outunit,'(a)' )trim(text) case default #ifdef __SX write( errunit, * )trim(text) #else write( errunit,'(/a/)' )trim(text) #endif write( outunit,'(/a/)' )trim(text) if( errortype.EQ.FATAL .OR. warnings_are_fatal )then call FLUSH(outunit) #ifdef sgi_mipspro call TRACE_BACK_STACK_AND_PRINT() #endif call ABORT() !automatically calls traceback on Cray systems end if end select error_state = errortype return end subroutine mpp_error_basic !##################################################################### function get_peset(pelist) integer :: get_peset !makes a PE set out of a PE list !a PE list is an ordered list of PEs !a PE set is a triad (start,log2stride,size) for SHMEM, an a communicator for MPI !if stride is non-uniform or not a power of 2, will return error (not required for MPI but enforced for uniformity) integer, intent(in), optional :: pelist(:) integer :: group integer :: i, n, stride, errunit integer, allocatable :: sorted(:) if( .NOT.PRESENT(pelist) )then !set it to current_peset_num get_peset = current_peset_num; return end if if( size(pelist(:)).EQ.1 .AND. npes.GT.1 )then !collective ops on single PEs should return get_peset = 0; return end if errunit = stderr() !--- first make sure pelist is monotonically increasing. do n = 2, size(pelist(:)) if(pelist(n) <= pelist(n-1)) call mpp_error(FATAL, "GET_PESET: pelist is not monotonically increasing") enddo allocate( sorted(size(pelist(:))) ) sorted = pelist if( debug )write( errunit,* )'pelist=', pelist !find if this array matches any existing peset do i = 1,peset_num if( debug )write( errunit,'(a,3i6)' )'pe, i, peset_num=', pe, i, peset_num if( size(sorted(:)).EQ.size(peset(i)%list(:)) )then if( ALL(sorted.EQ.peset(i)%list) )then deallocate(sorted) get_peset = i; return end if end if end do !not found, so create new peset peset_num = peset_num + 1 if( peset_num > current_peset_max ) call expand_peset() i = peset_num !shorthand !create list allocate( peset(i)%list(size(sorted(:))) ) peset(i)%list(:) = sorted(:) peset(i)%count = size(sorted(:)) peset(i)%start = sorted(1) n = size(sorted(:)) if( size(sorted(:)).GT.1 )then stride = sorted(2)-sorted(1) if( ANY(sorted(2:n)-sorted(1:n-1).NE.stride) ) & call mpp_error( WARNING, 'GET_PESET: pelist must have constant stride.' ) peset(i)%log2stride = nint( log(real(stride))/log(2.) ) if( 2**peset(i)%log2stride.NE.stride )call mpp_error( WARNING, 'GET_PESET: pelist must have power-of-2 stride.' ) else peset(i)%log2stride = 0 end if deallocate(sorted) get_peset = i return end function get_peset !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! ! SYNCHRONIZATION ROUTINES: mpp_sync, mpp_sync_self ! ! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine mpp_sync( pelist, check) !synchronize PEs in list integer, intent(in), optional :: pelist(:) integer, intent(in), optional :: check integer :: n call mpp_sync_self(pelist) n = get_peset(pelist); if( peset(n)%count.EQ.1 )return if( current_clock.NE.0 )call SYSTEM_CLOCK(start_tick) if( n.EQ.world_peset_num )then call SHMEM_BARRIER_ALL() !special call is faster else call SHMEM_BARRIER( peset(n)%start, peset(n)%log2stride, peset(n)%count, sync ) end if if( current_clock.NE.0 )call increment_current_clock(EVENT_WAIT) return end subroutine mpp_sync !##################################################################### ! ! ! Local synchronization. ! ! ! mpp_transmit is implemented as asynchronous ! put/send and synchronous ! get/recv. mpp_sync_self guarantees that outstanding ! asynchronous operations from the calling PE are complete. If ! pelist is supplied, mpp_sync_self checks only for ! outstanding puts to the PEs in pelist. ! ! If pelist is omitted, the context is assumed to be the ! current pelist. This call implies synchronization across the PEs in ! pelist, or the current pelist if pelist is absent. ! ! ! subroutine mpp_sync_self( pelist, check, request ) !this is to check if current PE's outstanding puts are complete !but we can't use shmem_fence because we are actually waiting for !a remote PE to complete its get integer, intent(in), optional :: pelist(:) integer, intent(in), optional :: check integer, intent(inout), optional :: request(:) integer :: i, m, n, stride n = get_peset(pelist) if( current_clock.NE.0 )call SYSTEM_CLOCK(start_tick) #ifdef _CRAYT90 call SHMEM_UDCFLUSH !invalidate data cache #endif do m = 1,peset(n)%count i = peset(n)%list(m) call SHMEM_INT8_WAIT( status(i), MPP_WAIT ) !wait for status.NE.MPP_WAIT end do if( current_clock.NE.0 )call increment_current_clock(EVENT_WAIT) return end subroutine mpp_sync_self