! -*-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 .
!***********************************************************************
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! ROUTINES TO INITIALIZE/FINALIZE MPP MODULE: mpp_init, mpp_exit !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! subroutine mpp_init( flags, in, out, err, log )
! integer, optional, intent(in) :: flags, in, out, err, log
subroutine mpp_init( flags, localcomm )
integer, optional, intent(in) :: flags
integer, optional, intent(in) :: localcomm
integer :: my_pe, num_pes, len, i, iunit
logical :: opened, existed
integer :: unit_begin, unit_end, unit_nml, io_status
character(len=5) :: this_pe
type(mpp_type), pointer :: dtype
if( module_is_initialized )return
call MPI_INITIALIZED( opened, error ) !in case called from another MPI package
if(opened .and. .NOT. PRESENT(localcomm)) call mpp_error( FATAL, 'MPP_INIT: communicator is required' )
if( .NOT.opened ) then
call MPI_INIT(error)
mpp_comm_private = MPI_COMM_WORLD
else
mpp_comm_private = localcomm
endif
call MPI_COMM_RANK( mpp_comm_private, pe, error )
call MPI_COMM_SIZE( mpp_comm_private, npes, error )
#ifdef use_MPI_SMA
call SHMEM_BARRIER_ALL()
call SHPALLOC( p_pSync, SHMEM_BARRIER_SYNC_SIZE, error, -1 )
call SHMEM_BARRIER_ALL()
#endif
module_is_initialized = .TRUE.
!PEsets: make defaults illegal
allocate(peset(0:current_peset_max))
peset(:)%count = -1
peset(:)%id = -1
peset(:)%group = -1
peset(:)%start = -1
peset(:)%log2stride = -1
peset(:)%name = " "
!0=single-PE, initialized so that count returns 1
peset(0)%count = 1
allocate( peset(0)%list(1) )
peset(0)%list = pe
current_peset_num = 0
peset(0)%id = mpp_comm_private
call MPI_COMM_GROUP( mpp_comm_private, peset(0)%group, error )
world_peset_num = get_peset( (/(i,i=0,npes-1)/) )
current_peset_num = world_peset_num !initialize current PEset to world
!initialize clocks
call SYSTEM_CLOCK( count=tick0, count_rate=ticks_per_sec, count_max=max_ticks )
tick_rate = 1./ticks_per_sec
clock0 = mpp_clock_id( 'Total runtime', flags=MPP_CLOCK_SYNC )
! Create the bytestream (default) mpp_datatype
mpp_byte%counter = 1
mpp_byte%ndims = 0
allocate(mpp_byte%sizes(0))
allocate(mpp_byte%subsizes(0))
allocate(mpp_byte%starts(0))
mpp_byte%etype = MPI_BYTE
mpp_byte%id = MPI_BYTE
mpp_byte%prev => null()
mpp_byte%next => null()
! Initialize datatype list with mpp_byte
datatypes%head => mpp_byte
datatypes%tail => mpp_byte
datatypes%length = 0
if( PRESENT(flags) )then
debug = flags.EQ.MPP_DEBUG
verbose = flags.EQ.MPP_VERBOSE .OR. debug
end if
call mpp_init_logfile()
call read_input_nml
!--- read namelist
#ifdef INTERNAL_FILE_NML
read (input_nml_file, mpp_nml, iostat=io_status)
#else
unit_begin = 103
unit_end = 512
do unit_nml = unit_begin, unit_end
inquire( unit_nml,OPENED=opened )
if( .NOT.opened )exit
end do
open(unit_nml,file='input.nml', iostat=io_status)
read(unit_nml,mpp_nml,iostat=io_status)
close(unit_nml)
#endif
if (io_status > 0) then
call mpp_error(FATAL,'=>mpp_init: Error reading input.nml')
endif
if(sync_all_clocks .AND. mpp_pe()==mpp_root_pe() ) call mpp_error(NOTE, &
"mpp_mod: mpp_nml variable sync_all_clocks is set to .true., all clocks are synchronized in mpp_clock_begin.")
! non-root pe messages written to other location than stdout()
if(etc_unit_is_stderr) then
etc_unit = stderr()
else
! 9 is reserved for etc_unit
etc_unit=9
inquire(unit=etc_unit,opened=opened)
if(opened) call mpp_error(FATAL,'Unit 9 is already in use (etc_unit) in mpp_comm_mpi')
if (trim(etcfile) /= '/dev/null') then
write( etcfile,'(a,i6.6)' )trim(etcfile)//'.', pe
endif
inquire(file=etcfile, exist=existed)
if(existed) then
open( unit=etc_unit, file=trim(etcfile), status='REPLACE' )
else
open( unit=etc_unit, file=trim(etcfile) )
endif
endif
! max_request is set to maximum of npes * REQUEST_MULTIPLY ( default is 20) and MAX_REQUEST_MIN ( default 10000)
max_request = max(MAX_REQUEST_MIN, mpp_npes()*REQUEST_MULTIPLY)
allocate( request_send(max_request) )
allocate( request_recv(max_request) )
allocate( size_recv(max_request) )
allocate( type_recv(max_request) )
request_send(:) = MPI_REQUEST_NULL
request_recv(:) = MPI_REQUEST_NULL
size_recv(:) = 0
type_recv(:) = 0
!if optional argument logunit=stdout, write messages to stdout instead.
!if specifying non-defaults, you must specify units not yet in use.
! if( PRESENT(in) )then
! inquire( unit=in, opened=opened )
! if( opened )call mpp_error( FATAL, 'MPP_INIT: unable to open stdin.' )
! in_unit=in
! end if
! if( PRESENT(out) )then
! inquire( unit=out, opened=opened )
! if( opened )call mpp_error( FATAL, 'MPP_INIT: unable to open stdout.' )
! out_unit=out
! end if
! if( PRESENT(err) )then
! inquire( unit=err, opened=opened )
! if( opened )call mpp_error( FATAL, 'MPP_INIT: unable to open stderr.' )
! err_unit=err
! end if
! log_unit=get_unit()
! if( PRESENT(log) )then
! inquire( unit=log, opened=opened )
! if( opened .AND. log.NE.out_unit )call mpp_error( FATAL, 'MPP_INIT: unable to open stdlog.' )
! log_unit=log
! end if
!!log_unit can be written to only from root_pe, all others write to stdout
! if( log_unit.NE.out_unit )then
! inquire( unit=log_unit, opened=opened )
! if( opened )call mpp_error( FATAL, 'MPP_INIT: specified unit for stdlog already in use.' )
! if( pe.EQ.root_pe )open( unit=log_unit, file=trim(configfile), status='REPLACE' )
! call mpp_sync()
! if( pe.NE.root_pe )open( unit=log_unit, file=trim(configfile), status='OLD' )
! end if
!messages
iunit = stdlog() ! workaround for lf95.
if( verbose )call mpp_error( NOTE, 'MPP_INIT: initializing MPP module...' )
if( pe.EQ.root_pe )then
write( iunit,'(/a)' )'MPP module '//trim(version)
write( iunit,'(a,i6)' )'MPP started with NPES=', npes
write( iunit,'(a)' )'Using MPI library for message passing...'
write( iunit, '(a,es12.4,a,i10,a)' ) &
'Realtime clock resolution=', tick_rate, ' sec (', ticks_per_sec, ' ticks/sec)'
write( iunit, '(a,es12.4,a,i20,a)' ) &
'Clock rolls over after ', max_ticks*tick_rate, ' sec (', max_ticks, ' ticks)'
write( iunit,'(/a)' )'MPP Parameter module '//trim(mpp_parameter_version)
write( iunit,'(/a)' )'MPP Data module '//trim(mpp_data_version)
end if
stdout_unit = stdout()
call mpp_clock_begin(clock0)
return
end subroutine mpp_init
!#######################################################################
!to be called at the end of a run
subroutine mpp_exit()
integer :: i, j, k, n, nmax, istat, out_unit, log_unit
real :: t, tmin, tmax, tavg, tstd
real :: m, mmin, mmax, mavg, mstd, t_total
logical :: opened
type(mpp_type), pointer :: dtype
if( .NOT.module_is_initialized )return
call mpp_set_current_pelist()
call mpp_clock_end(clock0)
t_total = clocks(clock0)%total_ticks*tick_rate
out_unit = stdout()
log_unit = stdlog()
if( clock_num.GT.0 )then
if( ANY(clocks(1:clock_num)%detailed) )then
call sum_clock_data; call dump_clock_summary
end if
call mpp_sync()
call FLUSH( out_unit )
if( pe.EQ.root_pe )then
write( out_unit,'(/a,i6,a)' ) 'Tabulating mpp_clock statistics across ', npes, ' PEs...'
if( ANY(clocks(1:clock_num)%detailed) ) &
write( out_unit,'(a)' )' ... see mpp_clock.out.#### for details on individual PEs.'
write( out_unit,'(/32x,a)' ) ' tmin tmax tavg tstd tfrac grain pemin pemax'
end if
write( log_unit,'(/37x,a)' ) 'time'
call FLUSH( out_unit )
call mpp_sync()
do i = 1,clock_num
if( .NOT.ANY(peset(clocks(i)%peset_num)%list(:).EQ.pe) )cycle
call mpp_set_current_pelist( peset(clocks(i)%peset_num)%list )
out_unit = stdout()
log_unit = stdlog()
!times between mpp_clock ticks
t = clocks(i)%total_ticks*tick_rate
tmin = t; call mpp_min(tmin)
tmax = t; call mpp_max(tmax)
tavg = t; call mpp_sum(tavg); tavg = tavg/mpp_npes()
tstd = (t-tavg)**2; call mpp_sum(tstd); tstd = sqrt( tstd/mpp_npes() )
if( pe.EQ.root_pe )write( out_unit,'(a32,4f14.6,f7.3,3i6)' ) &
clocks(i)%name, tmin, tmax, tavg, tstd, tavg/t_total, &
clocks(i)%grain, minval(peset(clocks(i)%peset_num)%list), &
maxval(peset(clocks(i)%peset_num)%list)
write(log_unit,'(a32,f14.6)') clocks(i)%name, clocks(i)%total_ticks*tick_rate
end do
if( ANY(clocks(1:clock_num)%detailed) .AND. pe.EQ.root_pe )write( out_unit,'(/32x,a)' ) &
' tmin tmax tavg tstd mmin mmax mavg mstd mavg/tavg'
do i = 1,clock_num
!messages: bytelengths and times
if( .NOT.clocks(i)%detailed )cycle
if( .NOT.ANY(peset(clocks(i)%peset_num)%list(:).EQ.pe) )cycle
call mpp_set_current_pelist( peset(clocks(i)%peset_num)%list )
out_unit = stdout()
do j = 1,MAX_EVENT_TYPES
n = clocks(i)%events(j)%calls; nmax = n
call mpp_max(nmax)
if( nmax.NE.0 )then
!don't divide by n because n might be 0
m = 0
if( n.GT.0 )m = sum(clocks(i)%events(j)%bytes(1:n))
mmin = m; call mpp_min(mmin)
mmax = m; call mpp_max(mmax)
mavg = m; call mpp_sum(mavg); mavg = mavg/mpp_npes()
mstd = (m-mavg)**2; call mpp_sum(mstd); mstd = sqrt( mstd/mpp_npes() )
t = 0
if( n.GT.0 )t = sum(clocks(i)%events(j)%ticks(1:n))*tick_rate
tmin = t; call mpp_min(tmin)
tmax = t; call mpp_max(tmax)
tavg = t; call mpp_sum(tavg); tavg = tavg/mpp_npes()
tstd = (t-tavg)**2; call mpp_sum(tstd); tstd = sqrt( tstd/mpp_npes() )
if( pe.EQ.root_pe )write( out_unit,'(a32,4f11.3,5es11.3)' ) &
trim(clocks(i)%name)//' '//trim(clocks(i)%events(j)%name), &
tmin, tmax, tavg, tstd, mmin, mmax, mavg, mstd, mavg/tavg
end if
end do
end do
end if
call FLUSH( out_unit )
! close down etc_unit: 9
inquire(unit=etc_unit, opened=opened)
if (opened) then
call FLUSH (etc_unit)
close(etc_unit)
endif
! Clear derived data types (skipping list head, mpp_byte)
dtype => datatypes%head
do while (.not. associated(dtype))
dtype => dtype%next
dtype%counter = 1 ! Force deallocation
call mpp_type_free(dtype)
end do
call mpp_set_current_pelist()
call mpp_sync()
call mpp_max(mpp_stack_hwm)
if( pe.EQ.root_pe )write( out_unit,* )'MPP_STACK high water mark=', mpp_stack_hwm
if(mpp_comm_private == MPI_COMM_WORLD ) call MPI_FINALIZE(error)
return
end subroutine mpp_exit
!#######################################################################
subroutine mpp_malloc( ptr, newlen, len )
integer, intent(in) :: newlen
integer, intent(inout) :: len
#ifdef use_MPI_SMA
real :: dummy
pointer( ptr, dummy )
if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_MALLOC: You must first call mpp_init.' )
!use existing allocation if it is enough
if( newlen.LE.len )return
call SHMEM_BARRIER_ALL()
!if the pointer is already allocated, deallocate
! if( len.NE.0 )call SHPDEALLC( ptr, error_8, -1 ) !BWA: error_8 instead of error, see PV 682618 (fixed in mpt.1.3.0.1)
if( len.NE.0 )call SHPDEALLC( ptr, error, -1 )
!allocate new length: assume that the array is KIND=8
call SHPALLOC( ptr, newlen, error, -1 )
len = newlen
call SHMEM_BARRIER_ALL()
#else
integer(POINTER_KIND), intent(in) :: ptr
call mpp_error( FATAL, 'mpp_malloc: requires use_MPI_SMA' )
#endif
return
end subroutine mpp_malloc
#ifdef use_MPI_GSM
!#######################################################################
!--- routine to perform GSM allocations
!
subroutine mpp_gsm_malloc( ptr, len )
use mpp_data_mod, only : peset, current_peset_num
integer(KIND=MPI_ADDRESS_KIND), intent(in) :: len
real :: dummy
integer :: ierror
!argument ptr is a cray pointer, points to a dummy argument in this routine
pointer( ptr, dummy )
include "mpi_gsmf.h"
if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_GSM_MALLOC: You must first call mpp_init.' )
call gsm_alloc(len, GSM_SINGLERANK, 0, peset(current_peset_num)%id, ptr, ierror)
if (ierror .eq. -1) call mpp_error( FATAL, 'MPP_GSM_MALLOC: gsm_alloc failed.' )
return
end subroutine mpp_gsm_malloc
!#######################################################################
!--- routine to free GSM allocations
!
subroutine mpp_gsm_free( ptr )
use mpp_data_mod, only : peset, current_peset_num
real :: dummy
integer :: ierror
!argument ptr is a cray pointer, points to a dummy argument in this routine
pointer( ptr, dummy )
include "mpi_gsmf.h"
call gsm_free(dummy, peset(current_peset_num)%id, ierror)
if (ierror .eq. -1) call mpp_error( FATAL, 'MPP_GSM_FREE: gsm_free failed.' )
return
end subroutine mpp_gsm_free
#endif
!#######################################################################
!set the mpp_stack variable to be at least n LONG words long
subroutine mpp_set_stack_size(n)
integer, intent(in) :: n
character(len=8) :: text
if( n.GT.mpp_stack_size .AND. allocated(mpp_stack) )deallocate(mpp_stack)
if( .NOT.allocated(mpp_stack) )then
allocate( mpp_stack(n) )
mpp_stack_size = n
end if
write( text,'(i8)' )n
if( pe.EQ.root_pe )call mpp_error( NOTE, 'MPP_SET_STACK_SIZE: stack size set to '//text//'.' )
return
end subroutine mpp_set_stack_size
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! BASIC MESSAGE PASSING ROUTINE: mpp_transmit !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine mpp_broadcast_char(data, length, from_pe, pelist )
character(len=*), intent(inout) :: data(:)
integer, intent(in) :: length, from_pe
integer, intent(in), optional :: pelist(:)
integer :: n, i, from_rank, out_unit
character :: str1D(length*size(data(:)))
pointer(lptr, str1D)
if( .NOT.module_is_initialized )call mpp_error( FATAL, 'mpp_broadcast_text: You must first call mpp_init.' )
n = get_peset(pelist); if( peset(n)%count.EQ.1 )return
if( debug )then
call SYSTEM_CLOCK(tick)
if(mpp_pe() == mpp_root_pe()) then
write( stdout_unit,'(a,i18,a,i5,a,2i5,2i8)' )&
'T=',tick, ' PE=',pe, 'mpp_broadcast_text begin: from_pe, length=', from_pe, length
endif
end if
if( .NOT.ANY(from_pe.EQ.peset(current_peset_num)%list) ) &
call mpp_error( FATAL, 'mpp_broadcast_text: broadcasting from invalid PE.' )
if( debug .and. (current_clock.NE.0) )call SYSTEM_CLOCK(start_tick)
! find the rank of from_pe in the pelist.
do i = 1, mpp_npes()
if(peset(n)%list(i) == from_pe) then
from_rank = i - 1
exit
endif
enddo
lptr = LOC (data)
if( mpp_npes().GT.1 ) call MPI_BCAST( data, length*size(data(:)), MPI_CHARACTER, from_rank, peset(n)%id, error )
if( debug .and. (current_clock.NE.0) )call increment_current_clock( EVENT_BROADCAST, length )
return
end subroutine mpp_broadcast_char
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_real8
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_real8_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_real8_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_real8_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_real8_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_real8_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_real8
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_real8_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_real8_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_real8_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_real8_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_real8_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_real8
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_real8_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_real8_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_real8_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_real8_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_real8_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_real8
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_real8_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_real8_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_real8_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_real8_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_real8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ real(DOUBLE_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL8
#include
#ifdef OVERLOAD_C8
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_cmplx8
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_cmplx8_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_cmplx8_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_cmplx8_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_cmplx8_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_cmplx8_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_cmplx8
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_cmplx8_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_cmplx8_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_cmplx8_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_cmplx8_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_cmplx8_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_cmplx8
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_cmplx8_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_cmplx8_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_cmplx8_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_cmplx8_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_cmplx8_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_cmplx8
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_cmplx8_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_cmplx8_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_cmplx8_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_cmplx8_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_cmplx8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ complex(DOUBLE_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 16
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_DOUBLE_COMPLEX
#include
#endif
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_real4
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_real4_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_real4_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_real4_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_real4_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_real4_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_real4
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_real4_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_real4_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_real4_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_real4_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_real4_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_real4
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_real4_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_real4_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_real4_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_real4_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_real4_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_real4
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_real4_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_real4_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_real4_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_real4_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_real4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ real(FLOAT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL4
#include
#ifdef OVERLOAD_C4
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_cmplx4
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_cmplx4_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_cmplx4_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_cmplx4_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_cmplx4_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_cmplx4_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_cmplx4
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_cmplx4_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_cmplx4_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_cmplx4_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_cmplx4_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_cmplx4_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_cmplx4
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_cmplx4_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_cmplx4_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_cmplx4_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_cmplx4_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_cmplx4_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_cmplx4
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_cmplx4_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_cmplx4_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_cmplx4_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_cmplx4_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_cmplx4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ complex(FLOAT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_COMPLEX
#include
#endif
#ifndef no_8byte_integers
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_int8
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_int8_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_int8_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_int8_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_int8_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_int8_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_int8
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_int8_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_int8_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_int8_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_int8_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_int8_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_int8
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_int8_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_int8_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_int8_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_int8_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_int8_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_int8
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_int8_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_int8_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_int8_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_int8_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_int8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(LONG_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER8
#include
#endif
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_int4
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_int4_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_int4_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_int4_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_int4_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_int4_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_int4
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_int4_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_int4_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_int4_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_int4_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_int4_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_int4
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_int4_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_int4_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_int4_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_int4_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_int4_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_int4
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_int4_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_int4_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_int4_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_int4_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_int4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(INT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER4
#include
#ifndef no_8byte_integers
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_logical8
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_logical8_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_logical8_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_logical8_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_logical8_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_logical8_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_logical8
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_logical8_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_logical8_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_logical8_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_logical8_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_logical8_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_logical8
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_logical8_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_logical8_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_logical8_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_logical8_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_logical8_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_logical8
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_logical8_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_logical8_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_logical8_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_logical8_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_logical8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ logical(LONG_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER8
#include
#endif
#undef MPP_TRANSMIT_
#define MPP_TRANSMIT_ mpp_transmit_logical4
#undef MPP_TRANSMIT_SCALAR_
#define MPP_TRANSMIT_SCALAR_ mpp_transmit_logical4_scalar
#undef MPP_TRANSMIT_2D_
#define MPP_TRANSMIT_2D_ mpp_transmit_logical4_2d
#undef MPP_TRANSMIT_3D_
#define MPP_TRANSMIT_3D_ mpp_transmit_logical4_3d
#undef MPP_TRANSMIT_4D_
#define MPP_TRANSMIT_4D_ mpp_transmit_logical4_4d
#undef MPP_TRANSMIT_5D_
#define MPP_TRANSMIT_5D_ mpp_transmit_logical4_5d
#undef MPP_RECV_
#define MPP_RECV_ mpp_recv_logical4
#undef MPP_RECV_SCALAR_
#define MPP_RECV_SCALAR_ mpp_recv_logical4_scalar
#undef MPP_RECV_2D_
#define MPP_RECV_2D_ mpp_recv_logical4_2d
#undef MPP_RECV_3D_
#define MPP_RECV_3D_ mpp_recv_logical4_3d
#undef MPP_RECV_4D_
#define MPP_RECV_4D_ mpp_recv_logical4_4d
#undef MPP_RECV_5D_
#define MPP_RECV_5D_ mpp_recv_logical4_5d
#undef MPP_SEND_
#define MPP_SEND_ mpp_send_logical4
#undef MPP_SEND_SCALAR_
#define MPP_SEND_SCALAR_ mpp_send_logical4_scalar
#undef MPP_SEND_2D_
#define MPP_SEND_2D_ mpp_send_logical4_2d
#undef MPP_SEND_3D_
#define MPP_SEND_3D_ mpp_send_logical4_3d
#undef MPP_SEND_4D_
#define MPP_SEND_4D_ mpp_send_logical4_4d
#undef MPP_SEND_5D_
#define MPP_SEND_5D_ mpp_send_logical4_5d
#undef MPP_BROADCAST_
#define MPP_BROADCAST_ mpp_broadcast_logical4
#undef MPP_BROADCAST_SCALAR_
#define MPP_BROADCAST_SCALAR_ mpp_broadcast_logical4_scalar
#undef MPP_BROADCAST_2D_
#define MPP_BROADCAST_2D_ mpp_broadcast_logical4_2d
#undef MPP_BROADCAST_3D_
#define MPP_BROADCAST_3D_ mpp_broadcast_logical4_3d
#undef MPP_BROADCAST_4D_
#define MPP_BROADCAST_4D_ mpp_broadcast_logical4_4d
#undef MPP_BROADCAST_5D_
#define MPP_BROADCAST_5D_ mpp_broadcast_logical4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ logical(INT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER4
#include
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! GLOBAL REDUCTION ROUTINES: mpp_max, mpp_sum, mpp_min !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_max_real8_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_max_real8_1d
#undef MPP_TYPE_
#define MPP_TYPE_ real(DOUBLE_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL8
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MAX
#include
#ifdef OVERLOAD_R4
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_max_real4_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_max_real4_1d
#undef MPP_TYPE_
#define MPP_TYPE_ real(FLOAT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL4
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MAX
#include
#endif
#ifndef no_8byte_integers
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_max_int8_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_max_int8_1d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(LONG_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER8
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MAX
#include
#endif
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_max_int4_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_max_int4_1d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(INT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER4
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MAX
#include
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_min_real8_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_min_real8_1d
#undef MPP_TYPE_
#define MPP_TYPE_ real(DOUBLE_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL8
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MIN
#include
#ifdef OVERLOAD_R4
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_min_real4_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_min_real4_1d
#undef MPP_TYPE_
#define MPP_TYPE_ real(FLOAT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL4
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MIN
#include
#endif
#ifndef no_8byte_integers
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_min_int8_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_min_int8_1d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(LONG_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER8
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MIN
#include
#endif
#undef MPP_REDUCE_0D_
#define MPP_REDUCE_0D_ mpp_min_int4_0d
#undef MPP_REDUCE_1D_
#define MPP_REDUCE_1D_ mpp_min_int4_1d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(INT_KIND)
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER4
#undef MPI_REDUCE_
#define MPI_REDUCE_ MPI_MIN
#include
#undef MPP_SUM_
#define MPP_SUM_ mpp_sum_real8
#undef MPP_SUM_SCALAR_
#define MPP_SUM_SCALAR_ mpp_sum_real8_scalar
#undef MPP_SUM_2D_
#define MPP_SUM_2D_ mpp_sum_real8_2d
#undef MPP_SUM_3D_
#define MPP_SUM_3D_ mpp_sum_real8_3d
#undef MPP_SUM_4D_
#define MPP_SUM_4D_ mpp_sum_real8_4d
#undef MPP_SUM_5D_
#define MPP_SUM_5D_ mpp_sum_real8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ real(DOUBLE_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL8
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#include
#ifdef OVERLOAD_C8
#undef MPP_SUM_
#define MPP_SUM_ mpp_sum_cmplx8
#undef MPP_SUM_SCALAR_
#define MPP_SUM_SCALAR_ mpp_sum_cmplx8_scalar
#undef MPP_SUM_2D_
#define MPP_SUM_2D_ mpp_sum_cmplx8_2d
#undef MPP_SUM_3D_
#define MPP_SUM_3D_ mpp_sum_cmplx8_3d
#undef MPP_SUM_4D_
#define MPP_SUM_4D_ mpp_sum_cmplx8_4d
#undef MPP_SUM_5D_
#define MPP_SUM_5D_ mpp_sum_cmplx8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ complex(DOUBLE_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_DOUBLE_COMPLEX
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 16
#include
#endif
#ifdef OVERLOAD_R4
#undef MPP_SUM_
#define MPP_SUM_ mpp_sum_real4
#undef MPP_SUM_SCALAR_
#define MPP_SUM_SCALAR_ mpp_sum_real4_scalar
#undef MPP_SUM_2D_
#define MPP_SUM_2D_ mpp_sum_real4_2d
#undef MPP_SUM_3D_
#define MPP_SUM_3D_ mpp_sum_real4_3d
#undef MPP_SUM_4D_
#define MPP_SUM_4D_ mpp_sum_real4_4d
#undef MPP_SUM_5D_
#define MPP_SUM_5D_ mpp_sum_real4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ real(FLOAT_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL4
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#include
#endif
#ifdef OVERLOAD_C4
#undef MPP_SUM_
#define MPP_SUM_ mpp_sum_cmplx4
#undef MPP_SUM_SCALAR_
#define MPP_SUM_SCALAR_ mpp_sum_cmplx4_scalar
#undef MPP_SUM_2D_
#define MPP_SUM_2D_ mpp_sum_cmplx4_2d
#undef MPP_SUM_3D_
#define MPP_SUM_3D_ mpp_sum_cmplx4_3d
#undef MPP_SUM_4D_
#define MPP_SUM_4D_ mpp_sum_cmplx4_4d
#undef MPP_SUM_5D_
#define MPP_SUM_5D_ mpp_sum_cmplx4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ complex(FLOAT_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_COMPLEX
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#include
#endif
#ifndef no_8byte_integers
#undef MPP_SUM_
#define MPP_SUM_ mpp_sum_int8
#undef MPP_SUM_SCALAR_
#define MPP_SUM_SCALAR_ mpp_sum_int8_scalar
#undef MPP_SUM_2D_
#define MPP_SUM_2D_ mpp_sum_int8_2d
#undef MPP_SUM_3D_
#define MPP_SUM_3D_ mpp_sum_int8_3d
#undef MPP_SUM_4D_
#define MPP_SUM_4D_ mpp_sum_int8_4d
#undef MPP_SUM_5D_
#define MPP_SUM_5D_ mpp_sum_int8_5d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(LONG_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER8
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#include
#endif
#undef MPP_SUM_
#define MPP_SUM_ mpp_sum_int4
#undef MPP_SUM_SCALAR_
#define MPP_SUM_SCALAR_ mpp_sum_int4_scalar
#undef MPP_SUM_2D_
#define MPP_SUM_2D_ mpp_sum_int4_2d
#undef MPP_SUM_3D_
#define MPP_SUM_3D_ mpp_sum_int4_3d
#undef MPP_SUM_4D_
#define MPP_SUM_4D_ mpp_sum_int4_4d
#undef MPP_SUM_5D_
#define MPP_SUM_5D_ mpp_sum_int4_5d
#undef MPP_TYPE_
#define MPP_TYPE_ integer(INT_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER4
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#include
!--------------------------------
#undef MPP_SUM_AD_
#define MPP_SUM_AD_ mpp_sum_real8_ad
#undef MPP_SUM_SCALAR_AD_
#define MPP_SUM_SCALAR_AD_ mpp_sum_real8_scalar_ad
#undef MPP_SUM_2D_AD_
#define MPP_SUM_2D_AD_ mpp_sum_real8_2d_ad
#undef MPP_SUM_3D_AD_
#define MPP_SUM_3D_AD_ mpp_sum_real8_3d_ad
#undef MPP_SUM_4D_AD_
#define MPP_SUM_4D_AD_ mpp_sum_real8_4d_ad
#undef MPP_SUM_5D_AD_
#define MPP_SUM_5D_AD_ mpp_sum_real8_5d_ad
#undef MPP_TYPE_
#define MPP_TYPE_ real(DOUBLE_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL8
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#include
#ifdef OVERLOAD_C8
#undef MPP_SUM_AD_
#define MPP_SUM_AD_ mpp_sum_cmplx8_ad
#undef MPP_SUM_SCALAR_AD_
#define MPP_SUM_SCALAR_AD_ mpp_sum_cmplx8_scalar_ad
#undef MPP_SUM_2D_AD_
#define MPP_SUM_2D_AD_ mpp_sum_cmplx8_2d_ad
#undef MPP_SUM_3D_AD_
#define MPP_SUM_3D_AD_ mpp_sum_cmplx8_3d_ad
#undef MPP_SUM_4D_AD_
#define MPP_SUM_4D_AD_ mpp_sum_cmplx8_4d_ad
#undef MPP_SUM_5D_AD_
#define MPP_SUM_5D_AD_ mpp_sum_cmplx8_5d_ad
#undef MPP_TYPE_
#define MPP_TYPE_ complex(DOUBLE_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_DOUBLE_COMPLEX
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 16
#include
#endif
#ifdef OVERLOAD_R4
#undef MPP_SUM_AD_
#define MPP_SUM_AD_ mpp_sum_real4_ad
#undef MPP_SUM_SCALAR_AD_
#define MPP_SUM_SCALAR_AD_ mpp_sum_real4_scalar_ad
#undef MPP_SUM_2D_AD_
#define MPP_SUM_2D_AD_ mpp_sum_real4_2d_ad
#undef MPP_SUM_3D_AD_
#define MPP_SUM_3D_AD_ mpp_sum_real4_3d_ad
#undef MPP_SUM_4D_AD_
#define MPP_SUM_4D_AD_ mpp_sum_real4_4d_ad
#undef MPP_SUM_5D_AD_
#define MPP_SUM_5D_AD_ mpp_sum_real4_5d_ad
#undef MPP_TYPE_
#define MPP_TYPE_ real(FLOAT_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_REAL4
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#include
#endif
#ifdef OVERLOAD_C4
#undef MPP_SUM_AD_
#define MPP_SUM_AD_ mpp_sum_cmplx4_ad
#undef MPP_SUM_SCALAR_AD_
#define MPP_SUM_SCALAR_AD_ mpp_sum_cmplx4_scalar_ad
#undef MPP_SUM_2D_AD_
#define MPP_SUM_2D_AD_ mpp_sum_cmplx4_2d_ad
#undef MPP_SUM_3D_AD_
#define MPP_SUM_3D_AD_ mpp_sum_cmplx4_3d_ad
#undef MPP_SUM_4D_AD_
#define MPP_SUM_4D_AD_ mpp_sum_cmplx4_4d_ad
#undef MPP_SUM_5D_AD_
#define MPP_SUM_5D_AD_ mpp_sum_cmplx4_5d_ad
#undef MPP_TYPE_
#define MPP_TYPE_ complex(FLOAT_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_COMPLEX
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#include
#endif
#ifndef no_8byte_integers
#undef MPP_SUM_AD_
#define MPP_SUM_AD_ mpp_sum_int8_ad
#undef MPP_SUM_SCALAR_AD_
#define MPP_SUM_SCALAR_AD_ mpp_sum_int8_scalar_ad
#undef MPP_SUM_2D_AD_
#define MPP_SUM_2D_AD_ mpp_sum_int8_2d_ad
#undef MPP_SUM_3D_AD_
#define MPP_SUM_3D_AD_ mpp_sum_int8_3d_ad
#undef MPP_SUM_4D_AD_
#define MPP_SUM_4D_AD_ mpp_sum_int8_4d_ad
#undef MPP_SUM_5D_AD_
#define MPP_SUM_5D_AD_ mpp_sum_int8_5d_ad
#undef MPP_TYPE_
#define MPP_TYPE_ integer(LONG_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER8
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 8
#include
#endif
#undef MPP_SUM_AD_
#define MPP_SUM_AD_ mpp_sum_int4_ad
#undef MPP_SUM_SCALAR_AD_
#define MPP_SUM_SCALAR_AD_ mpp_sum_int4_scalar_ad
#undef MPP_SUM_2D_AD_
#define MPP_SUM_2D_AD_ mpp_sum_int4_2d_ad
#undef MPP_SUM_3D_AD_
#define MPP_SUM_3D_AD_ mpp_sum_int4_3d_ad
#undef MPP_SUM_4D_AD_
#define MPP_SUM_4D_AD_ mpp_sum_int4_4d_ad
#undef MPP_SUM_5D_AD_
#define MPP_SUM_5D_AD_ mpp_sum_int4_5d_ad
#undef MPP_TYPE_
#define MPP_TYPE_ integer(INT_KIND)
#undef MPI_TYPE_
#define MPI_TYPE_ MPI_INTEGER4
#undef MPP_TYPE_BYTELEN_
#define MPP_TYPE_BYTELEN_ 4
#include
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! SCATTER AND GATHER ROUTINES: mpp_alltoall !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
#undef MPP_ALLTOALL_
#undef MPP_ALLTOALLV_
#undef MPP_ALLTOALLW_
#undef MPP_TYPE_
#undef MPP_TYPE_BYTELEN_
#undef MPI_TYPE_
#define MPP_ALLTOALL_ mpp_alltoall_int4
#define MPP_ALLTOALLV_ mpp_alltoall_int4_v
#define MPP_ALLTOALLW_ mpp_alltoall_int4_w
#define MPP_TYPE_ integer(INT_KIND)
#define MPP_TYPE_BYTELEN_ 4
#define MPI_TYPE_ MPI_INTEGER4
#include
#undef MPP_ALLTOALL_
#undef MPP_ALLTOALLV_
#undef MPP_ALLTOALLW_
#undef MPP_TYPE_
#undef MPP_TYPE_BYTELEN_
#undef MPI_TYPE_
#define MPP_ALLTOALL_ mpp_alltoall_int8
#define MPP_ALLTOALLV_ mpp_alltoall_int8_v
#define MPP_ALLTOALLW_ mpp_alltoall_int8_w
#define MPP_TYPE_ integer(LONG_KIND)
#define MPP_TYPE_BYTELEN_ 8
#define MPI_TYPE_ MPI_INTEGER8
#include
#undef MPP_ALLTOALL_
#undef MPP_ALLTOALLV_
#undef MPP_ALLTOALLW_
#undef MPP_TYPE_
#undef MPP_TYPE_BYTELEN_
#undef MPI_TYPE_
#define MPP_ALLTOALL_ mpp_alltoall_real4
#define MPP_ALLTOALLV_ mpp_alltoall_real4_v
#define MPP_ALLTOALLW_ mpp_alltoall_real4_w
#define MPP_TYPE_ real(FLOAT_KIND)
#define MPP_TYPE_BYTELEN_ 4
#define MPI_TYPE_ MPI_REAL4
#include
#undef MPP_ALLTOALL_
#undef MPP_ALLTOALLV_
#undef MPP_ALLTOALLW_
#undef MPP_TYPE_
#undef MPP_TYPE_BYTELEN_
#undef MPI_TYPE_
#define MPP_ALLTOALL_ mpp_alltoall_real8
#define MPP_ALLTOALLV_ mpp_alltoall_real8_v
#define MPP_ALLTOALLW_ mpp_alltoall_real8_w
#define MPP_TYPE_ real(DOUBLE_KIND)
#define MPP_TYPE_BYTELEN_ 8
#define MPI_TYPE_ MPI_REAL8
#include
#undef MPP_ALLTOALL_
#undef MPP_ALLTOALLV_
#undef MPP_ALLTOALLW_
#undef MPP_TYPE_
#undef MPP_TYPE_BYTELEN_
#undef MPI_TYPE_
#define MPP_ALLTOALL_ mpp_alltoall_logical4
#define MPP_ALLTOALLV_ mpp_alltoall_logical4_v
#define MPP_ALLTOALLW_ mpp_alltoall_logical4_w
#define MPP_TYPE_ logical(INT_KIND)
#define MPP_TYPE_BYTELEN_ 4
#define MPI_TYPE_ MPI_INTEGER4
#include
#undef MPP_ALLTOALL_
#undef MPP_ALLTOALLV_
#undef MPP_ALLTOALLW_
#undef MPP_TYPE_
#undef MPP_TYPE_BYTELEN_
#undef MPI_TYPE_
#define MPP_ALLTOALL_ mpp_alltoall_logical8
#define MPP_ALLTOALLV_ mpp_alltoall_logical8_v
#define MPP_ALLTOALLW_ mpp_alltoall_logical8_w
#define MPP_TYPE_ logical(LONG_KIND)
#define MPP_TYPE_BYTELEN_ 8
#define MPI_TYPE_ MPI_INTEGER8
#include
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! DATA TRANSFER TYPES: mpp_type_create, mpp_type_free !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
#define MPP_TYPE_CREATE_ mpp_type_create_int4
#define MPP_TYPE_ integer(INT_KIND)
#define MPI_TYPE_ MPI_INTEGER4
#include
#define MPP_TYPE_CREATE_ mpp_type_create_int8
#define MPP_TYPE_ integer(LONG_KIND)
#define MPI_TYPE_ MPI_INTEGER8
#include
#define MPP_TYPE_CREATE_ mpp_type_create_real4
#define MPP_TYPE_ real(FLOAT_KIND)
#define MPI_TYPE_ MPI_REAL4
#include
#define MPP_TYPE_CREATE_ mpp_type_create_real8
#define MPP_TYPE_ real(DOUBLE_KIND)
#define MPI_TYPE_ MPI_REAL8
#include
#define MPP_TYPE_CREATE_ mpp_type_create_logical4
#define MPP_TYPE_ logical(INT_KIND)
#define MPI_TYPE_ MPI_INTEGER4
#include
#define MPP_TYPE_CREATE_ mpp_type_create_logical8
#define MPP_TYPE_ logical(LONG_KIND)
#define MPI_TYPE_ MPI_INTEGER8
#include
! Clear preprocessor flags
#undef MPI_TYPE_
#undef MPP_TYPE_
#undef MPP_TYPE_CREATE_
! NOTE: This should probably not take a pointer, but for now we do this.
subroutine mpp_type_free(dtype)
type(mpp_type), pointer, intent(inout) :: dtype
if (.NOT. module_is_initialized) &
call mpp_error(FATAL, 'MPP_TYPE_FREE: You must first call mpp_init.')
if (current_clock .NE. 0) &
call SYSTEM_CLOCK(start_tick)
if (verbose) &
call mpp_error(NOTE, 'MPP_TYPE_FREE: using MPI_Type_free...')
! Decrement the reference counter
dtype%counter = dtype%counter - 1
if (dtype%counter < 1) then
! De-register the datatype in MPI runtime
call MPI_Type_free(dtype%id, error)
! Remove from list
dtype%prev => dtype%next
! Remove from memory
if (allocated(dtype%sizes)) deallocate(dtype%sizes)
if (allocated(dtype%subsizes)) deallocate(dtype%subsizes)
if (allocated(dtype%starts)) deallocate(dtype%starts)
deallocate(dtype)
datatypes%length = datatypes%length - 1
end if
if (current_clock .NE. 0) &
call increment_current_clock(EVENT_TYPE_FREE, MPP_TYPE_BYTELEN_)
end subroutine mpp_type_free