Oasis3 4.0.2
mpp_mod_oa.F90
Go to the documentation of this file.
00001 #ifndef key_noIO
00002 !-----------------------------------------------------------------------
00003 !                 Communication for message-passing codes
00004 !
00005 ! AUTHOR: V. Balaji (vbalaji@noaa.gov)
00006 !         Princeton University/GFDL
00007 !
00008 ! MODIFICATIONS: Reiner Vogelsang (reiner@sgi.com)
00009 !                Sophie Valcke: renamed the routine with _oa suffix
00010 !
00011 ! This program is free software; The author agrees that you can
00012 ! redistribute and/or modify this version of the program under the
00013 ! terms of the Lesser GNU General Public License as published
00014 ! by the Free Software Foundation.
00015 !
00016 ! This program is distributed in the hope that it will be useful,
00017 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
00018 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00019 ! Lesser GNU General Public License for more details
00020 ! (http://www.gnu.org/copyleft/lesser.html).
00021 !-----------------------------------------------------------------------
00022 
00023 !these are used to determine hardware/OS/compiler
00024 !#include <os.h>
00025 
00026 !onlysgi_mipspro one of SMA or MPI can be used
00027 !(though mixing calls is allowed, this module will not)
00028 #ifdef use_libSMA
00029 #undef use_libMPI
00030 #endif
00031 !
00032 #if defined(_CRAYT3E) || defined(_CRAYT3D) || defined(sgi_mipspro)
00033 #define SGICRAY_MPP
00034 #endif
00035 !
00036 !shmalloc is used on MPP SGI/Cray systems for shmem
00037 #if defined(use_libSMA) && defined(SGICRAY_MPP)
00038 #define use_shmalloc
00039 #endif
00040 
00041 module mpp_mod_oa
00042 use mod_kinds_mpp
00043 #include <os.h>
00044 !string BWA is used to tag lines that are bug workarounds and will disappear
00045 !when offending compiler bug is fixed
00046 !a generalized communication package for use with shmem and MPI
00047 !will add: co_array_fortran, MPI2
00048 !Balaji (vb@gfdl.gov) 11 May 1998
00049 #ifdef sgi_mipspro
00050 #ifdef use_libSMA
00051 !  use shmem_interface
00052 #endif
00053 !#ifdef use_libMPI
00054 !  use mpi
00055 !#endif
00056 #endif
00057   implicit none
00058   private
00059   character(len=128), private :: version= 
00060        '$Id: mpp_mod_oa.F90 2826 2010-12-10 11:14:21Z valcke $'
00061   character(len=128), private :: tagname= 
00062        '$Name$'
00063 
00064 !various lengths (see shpalloc) are estimated in "words" which are 32bit on SGI, 64bit on Cray
00065 !these are also the expected sizeof of args to MPI/shmem libraries
00066 #ifdef _CRAY
00067   integer(LONG_KIND), private :: word(1)
00068 #endif
00069 #ifdef sgi_mipspro
00070   integer(INT_KIND), private :: word(1)
00071 #endif
00072 
00073 #ifdef SGICRAY
00074 !see intro_io(3F): to see why these values are used rather than 5,6,0
00075   integer, private :: in_unit=100, out_unit=101, err_unit=102
00076 #else
00077   integer, private :: in_unit=5, out_unit=6, err_unit=0
00078 #endif
00079   integer :: log_unit, etc_unit
00080   logical, private :: module_is_initialized=.FALSE.
00081   integer, private :: pe=0, node=0, npes=1, root_pe=0
00082   integer, private :: error
00083   integer, parameter, private :: MAXPES=2048 !used for dimensioning stuff that might be indexed by pe
00084   character(len=32) :: configfile='logfile.out'
00085   character(len=32) :: etcfile='._mpp.nonrootpe.stdout'
00086   logical,save::logfile_defined=.false.,opened
00087   integer::io_num
00088 
00089 !initialization flags
00090   integer, parameter, public :: MPP_VERBOSE=1, MPP_DEBUG=2
00091   logical, private :: verbose=.FALSE., debug=.FALSE.
00092 
00093 !flags to transmit routines
00094   integer, parameter, public :: ALL_PES=-1, ANY_PE=-2, NULL_PE=-3
00095 
00096 !errortype flags
00097   integer, parameter, public :: NOTE=0, WARNING=1, FATAL=2
00098   logical, private :: warnings_are_fatal = .FALSE.
00099   integer, private :: error_state=0
00100 
00101   integer(LONG_KIND), parameter, private :: MPP_WAIT=-1, MPP_READY=-2
00102 #ifdef use_libSMA
00103 #include <mpp/shmem.fh>
00104   integer :: sync(SHMEM_REDUCE_SYNC_SIZE+SHMEM_BCAST_SYNC_SIZE+SHMEM_BARRIER_SYNC_SIZE)
00105 !status and remote_data_loc are used to synchronize communication is MPP_TRANSMIT
00106 #ifdef use_shmalloc
00107   integer(LONG_KIND), private, dimension(0:MAXPES) :: status, remote_data_loc
00108 #else
00109   integer(LONG_KIND), private, allocatable, dimension(:) :: status, remote_data_loc
00110 #endif
00111   integer, private :: mpp_from_pe !used to announce from where data is coming from
00112 #ifdef use_shmalloc
00113 !we call shpalloc in mpp_init() to ensure all these are remotely accessible
00114 !on PVP where shpalloc doesn't exist, module variables are automatically
00115 !guaranteed to be remotely accessible
00116   pointer( ptr_sync, sync )
00117   pointer( ptr_status, status )
00118   pointer( ptr_from, mpp_from_pe )
00119   pointer( ptr_remote, remote_data_loc )
00120 #endif
00121 #endif /* use_libSMA */
00122 #ifdef use_libMPI
00123 !#ifndef sgi_mipspro
00124 !sgi_mipspro gets this from 'use mpi'
00125 #include <mpif.h>
00126 !!include 'mpif.h'
00127 !#endif
00128 !tag is never used, but must be initialized to non-negative integer
00129   integer, private :: tag=1, stat(MPI_STATUS_SIZE)
00130 !  integer, private, allocatable :: request(:)
00131   integer, public, allocatable :: mpp_request(:)
00132 #ifdef _CRAYT3E
00133 !BWA: mpif.h on t3e currently does not contain MPI_INTEGER8 datatype
00134 !(O2k and t90 do)
00135 !(t3e: fixed on 3.3 I believe)
00136   integer, parameter :: MPI_INTEGER8=MPI_INTEGER
00137 #endif
00138 #endif /* use_libMPI */
00139 
00140 !mpp_stack is used by SHMEM collective ops
00141 !must be SHPALLOC'd on SGICRAY_MPP, but is allocatable on PVP
00142 #ifdef use_shmalloc
00143   real(DOUBLE_KIND), private :: mpp_stack(1)
00144   pointer( ptr_stack, mpp_stack )
00145 #else
00146   real(DOUBLE_KIND), private, allocatable :: mpp_stack(:)
00147 #endif
00148   integer, private :: mpp_stack_size=0, mpp_stack_hwm=0
00149 
00150 !peset hold communicators as SHMEM-compatible triads (start, log2(stride), num)
00151   type, private :: communicator
00152      character(len=32) :: name
00153      integer, pointer :: list(:)
00154      integer :: count
00155 #ifdef use_libSMA
00156      integer :: start, log2stride
00157 #elif use_libMPI
00158      integer :: id, group    !MPI communicator and group id for this PE set
00159 #endif
00160   end type
00161   integer, parameter :: PESET_MAX=32 !should be .LE. max num of MPI communicators
00162   type(communicator) :: peset(0:PESET_MAX) !0 is a dummy used to hold single-PE "self" communicator
00163   integer :: peset_num=0, current_peset_num=0
00164   integer :: world_peset_num !the world communicator
00165 
00166 !performance profiling
00167 !  This profiles every type of MPI/SHMEM call within
00168 !    a specified region of high-level code
00169 !  Initialize or retrieve a clock with
00170 !  id = mpp_clock_id( 'Region identifier name' )
00171 !  Then set caliper points around the region using:
00172 !  call mpp_clock_begin(id)
00173 !  ...
00174 !  call mpp_clock_end(id)
00175 !  mpp_exit will print out the results.
00176 #ifdef __sgi
00177 #define SYSTEM_CLOCK system_clock_sgi
00178 #endif
00179 
00180 #ifdef use_libMPI
00181 #define SYSTEM_CLOCK system_clock_mpi
00182 #endif
00183 
00184 #if defined(__sgi) || defined(use_libMPI)
00185   integer(LONG_KIND), private :: tick, ticks_per_sec, max_ticks, start_tick, end_tick, tick0=0
00186 #else
00187   integer, private :: tick, ticks_per_sec, max_ticks, start_tick, end_tick, tick0=0
00188 #endif
00189   real, private :: tick_rate
00190   integer, private, parameter :: MAX_CLOCKS=100, MAX_EVENT_TYPES=5, MAX_EVENTS=40000
00191 !event types
00192   integer, private, parameter :: EVENT_ALLREDUCE=1, EVENT_BROADCAST=2, EVENT_RECV=3, EVENT_SEND=4, EVENT_WAIT=5
00193   integer, private :: clock_num=0, current_clock=0
00194   integer, private :: clock0    !measures total runtime from mpp_init to mpp_exit
00195   integer, private :: clock_grain=HUGE(1)
00196 !the event contains information for each type of event (e.g SHMEM_PUT)
00197   type, private :: event
00198      character(len=16)  :: name
00199      integer(LONG_KIND) :: ticks(MAX_EVENTS), bytes(MAX_EVENTS)
00200      integer            :: calls
00201   end type event
00202 !a clock contains an array of event profiles for a region
00203   integer, parameter, public :: MPP_CLOCK_SYNC=1, MPP_CLOCK_DETAILED=2
00204   type, private :: clock
00205      character(len=32) :: name
00206 #if defined(__sgi) || defined(use_libMPI)
00207      integer(LONG_KIND) :: tick
00208 #else
00209      integer :: tick
00210 #endif
00211      integer(LONG_KIND) :: total_ticks
00212      integer :: peset_num
00213      logical :: sync_on_begin, detailed
00214      type(event), pointer :: events(:) !if needed, allocate to MAX_EVENT_TYPES
00215   end type
00216   type(clock) :: clocks(MAX_CLOCKS)
00217 
00218   integer,parameter :: MAX_BINS=20
00219   TYPE :: Clock_Data_Summary
00220     character(len=16) :: name
00221     real(DOUBLE_KIND) :: msg_size_sums(MAX_BINS)
00222     real(DOUBLE_KIND) :: msg_time_sums(MAX_BINS)
00223     real(DOUBLE_KIND) :: total_data
00224     real(DOUBLE_KIND) :: total_time
00225     integer(LONG_KIND) :: msg_size_cnts(MAX_BINS)
00226     integer(LONG_KIND) :: total_cnts
00227   END TYPE Clock_Data_Summary
00228 
00229   TYPE :: Summary_Struct
00230     character(len=16)         :: name
00231     type (Clock_Data_Summary) :: event(MAX_EVENT_TYPES)
00232   END TYPE Summary_Struct
00233   type(Summary_Struct) :: clock_summary(MAX_CLOCKS)
00234   
00235 !public interfaces
00236   interface mpp_max
00237      module procedure mpp_max_real8
00238 #ifndef no_8byte_integers
00239      module procedure mpp_max_int8
00240 #endif
00241 #ifndef no_4byte_reals
00242      module procedure mpp_max_real4
00243 #endif
00244      module procedure mpp_max_int4
00245   end interface
00246   interface mpp_min
00247      module procedure mpp_min_real8
00248 #ifndef no_8byte_integers
00249      module procedure mpp_min_int8
00250 #endif
00251 #ifndef no_4byte_reals
00252      module procedure mpp_min_real4
00253 #endif
00254      module procedure mpp_min_int4
00255   end interface
00256   interface mpp_sum
00257 #ifndef no_8byte_integers
00258      module procedure mpp_sum_int8
00259      module procedure mpp_sum_int8_scalar
00260      module procedure mpp_sum_int8_2d
00261      module procedure mpp_sum_int8_3d
00262      module procedure mpp_sum_int8_4d
00263      module procedure mpp_sum_int8_5d
00264 #endif
00265      module procedure mpp_sum_real8
00266      module procedure mpp_sum_real8_scalar
00267      module procedure mpp_sum_real8_2d
00268      module procedure mpp_sum_real8_3d
00269      module procedure mpp_sum_real8_4d
00270      module procedure mpp_sum_real8_5d
00271      module procedure mpp_sum_cmplx8
00272      module procedure mpp_sum_cmplx8_scalar
00273      module procedure mpp_sum_cmplx8_2d
00274      module procedure mpp_sum_cmplx8_3d
00275      module procedure mpp_sum_cmplx8_4d
00276      module procedure mpp_sum_cmplx8_5d
00277      module procedure mpp_sum_int4
00278      module procedure mpp_sum_int4_scalar
00279      module procedure mpp_sum_int4_2d
00280      module procedure mpp_sum_int4_3d
00281      module procedure mpp_sum_int4_4d
00282      module procedure mpp_sum_int4_5d
00283 #ifndef no_4byte_reals
00284      module procedure mpp_sum_real4
00285      module procedure mpp_sum_real4_scalar
00286      module procedure mpp_sum_real4_2d
00287      module procedure mpp_sum_real4_3d
00288      module procedure mpp_sum_real4_4d
00289      module procedure mpp_sum_real4_5d
00290 #endif
00291 #ifndef no_4byte_cmplx
00292      module procedure mpp_sum_cmplx4
00293      module procedure mpp_sum_cmplx4_scalar
00294      module procedure mpp_sum_cmplx4_2d
00295      module procedure mpp_sum_cmplx4_3d
00296      module procedure mpp_sum_cmplx4_4d
00297      module procedure mpp_sum_cmplx4_5d
00298 #endif
00299   end interface
00300   interface mpp_transmit
00301      module procedure mpp_transmit_real8
00302      module procedure mpp_transmit_real8_scalar
00303      module procedure mpp_transmit_real8_2d
00304      module procedure mpp_transmit_real8_3d
00305      module procedure mpp_transmit_real8_4d
00306      module procedure mpp_transmit_real8_5d
00307      module procedure mpp_transmit_cmplx8
00308      module procedure mpp_transmit_cmplx8_scalar
00309      module procedure mpp_transmit_cmplx8_2d
00310      module procedure mpp_transmit_cmplx8_3d
00311      module procedure mpp_transmit_cmplx8_4d
00312      module procedure mpp_transmit_cmplx8_5d
00313 #ifndef no_8byte_integers
00314      module procedure mpp_transmit_int8
00315      module procedure mpp_transmit_int8_scalar
00316      module procedure mpp_transmit_int8_2d
00317      module procedure mpp_transmit_int8_3d
00318      module procedure mpp_transmit_int8_4d
00319      module procedure mpp_transmit_int8_5d
00320      module procedure mpp_transmit_logical8
00321      module procedure mpp_transmit_logical8_scalar
00322      module procedure mpp_transmit_logical8_2d
00323      module procedure mpp_transmit_logical8_3d
00324      module procedure mpp_transmit_logical8_4d
00325      module procedure mpp_transmit_logical8_5d
00326 #endif
00327 #ifndef no_4byte_reals
00328      module procedure mpp_transmit_real4
00329      module procedure mpp_transmit_real4_scalar
00330      module procedure mpp_transmit_real4_2d
00331      module procedure mpp_transmit_real4_3d
00332      module procedure mpp_transmit_real4_4d
00333      module procedure mpp_transmit_real4_5d
00334 #endif
00335 #ifndef no_4byte_cmplx
00336      module procedure mpp_transmit_cmplx4
00337      module procedure mpp_transmit_cmplx4_scalar
00338      module procedure mpp_transmit_cmplx4_2d
00339      module procedure mpp_transmit_cmplx4_3d
00340      module procedure mpp_transmit_cmplx4_4d
00341      module procedure mpp_transmit_cmplx4_5d
00342 #endif
00343      module procedure mpp_transmit_int4
00344      module procedure mpp_transmit_int4_scalar
00345      module procedure mpp_transmit_int4_2d
00346      module procedure mpp_transmit_int4_3d
00347      module procedure mpp_transmit_int4_4d
00348      module procedure mpp_transmit_int4_5d
00349      module procedure mpp_transmit_logical4
00350      module procedure mpp_transmit_logical4_scalar
00351      module procedure mpp_transmit_logical4_2d
00352      module procedure mpp_transmit_logical4_3d
00353      module procedure mpp_transmit_logical4_4d
00354      module procedure mpp_transmit_logical4_5d
00355   end interface
00356   interface mpp_recv
00357      module procedure mpp_recv_real8
00358      module procedure mpp_recv_real8_scalar
00359      module procedure mpp_recv_real8_2d
00360      module procedure mpp_recv_real8_3d
00361      module procedure mpp_recv_real8_4d
00362      module procedure mpp_recv_real8_5d
00363      module procedure mpp_recv_cmplx8
00364      module procedure mpp_recv_cmplx8_scalar
00365      module procedure mpp_recv_cmplx8_2d
00366      module procedure mpp_recv_cmplx8_3d
00367      module procedure mpp_recv_cmplx8_4d
00368      module procedure mpp_recv_cmplx8_5d
00369 #ifndef no_8byte_integers
00370      module procedure mpp_recv_int8
00371      module procedure mpp_recv_int8_scalar
00372      module procedure mpp_recv_int8_2d
00373      module procedure mpp_recv_int8_3d
00374      module procedure mpp_recv_int8_4d
00375      module procedure mpp_recv_int8_5d
00376      module procedure mpp_recv_logical8
00377      module procedure mpp_recv_logical8_scalar
00378      module procedure mpp_recv_logical8_2d
00379      module procedure mpp_recv_logical8_3d
00380      module procedure mpp_recv_logical8_4d
00381      module procedure mpp_recv_logical8_5d
00382 #endif
00383 #ifndef no_4byte_reals
00384      module procedure mpp_recv_real4
00385      module procedure mpp_recv_real4_scalar
00386      module procedure mpp_recv_real4_2d
00387      module procedure mpp_recv_real4_3d
00388      module procedure mpp_recv_real4_4d
00389      module procedure mpp_recv_real4_5d
00390 #endif
00391 #ifndef no_4byte_cmplx
00392      module procedure mpp_recv_cmplx4
00393      module procedure mpp_recv_cmplx4_scalar
00394      module procedure mpp_recv_cmplx4_2d
00395      module procedure mpp_recv_cmplx4_3d
00396      module procedure mpp_recv_cmplx4_4d
00397      module procedure mpp_recv_cmplx4_5d
00398 #endif
00399      module procedure mpp_recv_int4
00400      module procedure mpp_recv_int4_scalar
00401      module procedure mpp_recv_int4_2d
00402      module procedure mpp_recv_int4_3d
00403      module procedure mpp_recv_int4_4d
00404      module procedure mpp_recv_int4_5d
00405      module procedure mpp_recv_logical4
00406      module procedure mpp_recv_logical4_scalar
00407      module procedure mpp_recv_logical4_2d
00408      module procedure mpp_recv_logical4_3d
00409      module procedure mpp_recv_logical4_4d
00410      module procedure mpp_recv_logical4_5d
00411   end interface
00412   interface mpp_send
00413      module procedure mpp_send_real8
00414      module procedure mpp_send_real8_scalar
00415      module procedure mpp_send_real8_2d
00416      module procedure mpp_send_real8_3d
00417      module procedure mpp_send_real8_4d
00418      module procedure mpp_send_real8_5d
00419      module procedure mpp_send_cmplx8
00420      module procedure mpp_send_cmplx8_scalar
00421      module procedure mpp_send_cmplx8_2d
00422      module procedure mpp_send_cmplx8_3d
00423      module procedure mpp_send_cmplx8_4d
00424      module procedure mpp_send_cmplx8_5d
00425 #ifndef no_8byte_integers
00426      module procedure mpp_send_int8
00427      module procedure mpp_send_int8_scalar
00428      module procedure mpp_send_int8_2d
00429      module procedure mpp_send_int8_3d
00430      module procedure mpp_send_int8_4d
00431      module procedure mpp_send_int8_5d
00432      module procedure mpp_send_logical8
00433      module procedure mpp_send_logical8_scalar
00434      module procedure mpp_send_logical8_2d
00435      module procedure mpp_send_logical8_3d
00436      module procedure mpp_send_logical8_4d
00437      module procedure mpp_send_logical8_5d
00438 #endif
00439 #ifndef no_4byte_reals
00440      module procedure mpp_send_real4
00441      module procedure mpp_send_real4_scalar
00442      module procedure mpp_send_real4_2d
00443      module procedure mpp_send_real4_3d
00444      module procedure mpp_send_real4_4d
00445      module procedure mpp_send_real4_5d
00446 #endif
00447 #ifndef no_4byte_cmplx
00448      module procedure mpp_send_cmplx4
00449      module procedure mpp_send_cmplx4_scalar
00450      module procedure mpp_send_cmplx4_2d
00451      module procedure mpp_send_cmplx4_3d
00452      module procedure mpp_send_cmplx4_4d
00453      module procedure mpp_send_cmplx4_5d
00454 #endif
00455      module procedure mpp_send_int4
00456      module procedure mpp_send_int4_scalar
00457      module procedure mpp_send_int4_2d
00458      module procedure mpp_send_int4_3d
00459      module procedure mpp_send_int4_4d
00460      module procedure mpp_send_int4_5d
00461      module procedure mpp_send_logical4
00462      module procedure mpp_send_logical4_scalar
00463      module procedure mpp_send_logical4_2d
00464      module procedure mpp_send_logical4_3d
00465      module procedure mpp_send_logical4_4d
00466      module procedure mpp_send_logical4_5d
00467   end interface
00468 
00469   interface mpp_broadcast
00470      module procedure mpp_broadcast_real8
00471      module procedure mpp_broadcast_real8_scalar
00472      module procedure mpp_broadcast_real8_2d
00473      module procedure mpp_broadcast_real8_3d
00474      module procedure mpp_broadcast_real8_4d
00475      module procedure mpp_broadcast_real8_5d
00476      module procedure mpp_broadcast_cmplx8
00477      module procedure mpp_broadcast_cmplx8_scalar
00478      module procedure mpp_broadcast_cmplx8_2d
00479      module procedure mpp_broadcast_cmplx8_3d
00480      module procedure mpp_broadcast_cmplx8_4d
00481      module procedure mpp_broadcast_cmplx8_5d
00482 #ifndef no_8byte_integers
00483      module procedure mpp_broadcast_int8
00484      module procedure mpp_broadcast_int8_scalar
00485      module procedure mpp_broadcast_int8_2d
00486      module procedure mpp_broadcast_int8_3d
00487      module procedure mpp_broadcast_int8_4d
00488      module procedure mpp_broadcast_int8_5d
00489      module procedure mpp_broadcast_logical8
00490      module procedure mpp_broadcast_logical8_scalar
00491      module procedure mpp_broadcast_logical8_2d
00492      module procedure mpp_broadcast_logical8_3d
00493      module procedure mpp_broadcast_logical8_4d
00494      module procedure mpp_broadcast_logical8_5d
00495 #endif
00496 #ifndef no_4byte_reals
00497      module procedure mpp_broadcast_real4
00498      module procedure mpp_broadcast_real4_scalar
00499      module procedure mpp_broadcast_real4_2d
00500      module procedure mpp_broadcast_real4_3d
00501      module procedure mpp_broadcast_real4_4d
00502      module procedure mpp_broadcast_real4_5d
00503 #endif
00504 #ifndef no_4byte_cmplx
00505      module procedure mpp_broadcast_cmplx4
00506      module procedure mpp_broadcast_cmplx4_scalar
00507      module procedure mpp_broadcast_cmplx4_2d
00508      module procedure mpp_broadcast_cmplx4_3d
00509      module procedure mpp_broadcast_cmplx4_4d
00510      module procedure mpp_broadcast_cmplx4_5d
00511 #endif
00512      module procedure mpp_broadcast_int4
00513      module procedure mpp_broadcast_int4_scalar
00514      module procedure mpp_broadcast_int4_2d
00515      module procedure mpp_broadcast_int4_3d
00516      module procedure mpp_broadcast_int4_4d
00517      module procedure mpp_broadcast_int4_5d
00518      module procedure mpp_broadcast_logical4
00519      module procedure mpp_broadcast_logical4_scalar
00520      module procedure mpp_broadcast_logical4_2d
00521      module procedure mpp_broadcast_logical4_3d
00522      module procedure mpp_broadcast_logical4_4d
00523      module procedure mpp_broadcast_logical4_5d
00524   end interface
00525 
00526   interface mpp_chksum
00527 #ifndef no_8byte_integers
00528      module procedure mpp_chksum_i8_1d
00529      module procedure mpp_chksum_i8_2d
00530      module procedure mpp_chksum_i8_3d
00531      module procedure mpp_chksum_i8_4d
00532 #endif
00533      module procedure mpp_chksum_i4_1d
00534      module procedure mpp_chksum_i4_2d
00535      module procedure mpp_chksum_i4_3d
00536      module procedure mpp_chksum_i4_4d
00537      module procedure mpp_chksum_r8_0d
00538      module procedure mpp_chksum_r8_1d
00539      module procedure mpp_chksum_r8_2d
00540      module procedure mpp_chksum_r8_3d
00541      module procedure mpp_chksum_r8_4d
00542      module procedure mpp_chksum_r8_5d
00543      module procedure mpp_chksum_c8_0d
00544      module procedure mpp_chksum_c8_1d
00545      module procedure mpp_chksum_c8_2d
00546      module procedure mpp_chksum_c8_3d
00547      module procedure mpp_chksum_c8_4d
00548      module procedure mpp_chksum_c8_5d
00549 #ifndef no_4byte_reals
00550      module procedure mpp_chksum_r4_0d
00551      module procedure mpp_chksum_r4_1d
00552      module procedure mpp_chksum_r4_2d
00553      module procedure mpp_chksum_r4_3d
00554      module procedure mpp_chksum_r4_4d
00555      module procedure mpp_chksum_r4_5d
00556 #endif
00557 #ifndef no_4byte_cmplx
00558      module procedure mpp_chksum_c4_0d
00559      module procedure mpp_chksum_c4_1d
00560      module procedure mpp_chksum_c4_2d
00561      module procedure mpp_chksum_c4_3d
00562      module procedure mpp_chksum_c4_4d
00563      module procedure mpp_chksum_c4_5d
00564 #endif
00565   end interface
00566 
00567   interface mpp_error
00568      module procedure mpp_error_basic
00569      module procedure mpp_error_mesg
00570      module procedure mpp_error_noargs
00571   end interface
00572 
00573 #ifdef use_libSMA
00574 !currently SMA contains no generic shmem_wait for different integer kinds:
00575 !I have inserted one here
00576   interface shmem_integer_wait
00577      module procedure shmem_int4_wait_local
00578      module procedure shmem_int8_wait_local
00579   end interface
00580 #endif
00581   public :: mpp_chksum, mpp_max, mpp_min, mpp_sum
00582   public :: mpp_exit, mpp_init
00583   public :: mpp_pe, mpp_node, mpp_npes, mpp_root_pe, mpp_set_root_pe, mpp_set_stack_size
00584   public :: mpp_clock_begin, mpp_clock_end, mpp_clock_id, mpp_clock_set_grain
00585   public :: mpp_error, mpp_error_state, mpp_set_warn_level
00586   public :: mpp_sync, mpp_sync_self
00587   public :: mpp_transmit, mpp_send, mpp_recv, mpp_broadcast
00588   public :: stdin, stdout, stderr, stdlog
00589   public :: mpp_declare_pelist, mpp_get_current_pelist, mpp_set_current_pelist
00590 #ifdef use_shmalloc
00591   public :: mpp_malloc
00592 #endif
00593 
00594   contains
00595 
00596 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00597 !                                                                             !
00598 !       ROUTINES TO INITIALIZE/FINALIZE MPP MODULE: mpp_init, mpp_exit        !
00599 !                                                                             !
00600 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00601 
00602     subroutine mpp_init( flags,mpp_comm ,logfile)
00603       integer, optional, intent(in) :: flags,mpp_comm
00604       character(len=*), optional, intent(in) :: logfile
00605 !    subroutine mpp_init( flags, in, out, err, log )
00606 !      integer, optional, intent(in) :: flags, in, out, err, log
00607       integer :: my_pe, num_pes, len
00608       integer :: i
00609       integer :: comm_intern
00610       logical :: opened
00611 #ifdef _CRAYT3E
00612       intrinsic my_pe
00613 #endif
00614 
00615       if( module_is_initialized )return
00616 
00617 #ifdef use_libSMA
00618       call START_PES(0)         !the argument 0 means extract from environment variable NPES on PVP/SGI, from mpprun -n on t3e
00619       pe = my_pe()
00620       node = pe                 !on an SMP this should return node ID rather than PE number.
00621       npes = num_pes()
00622 #elif use_libMPI
00623       call MPI_INITIALIZED( opened, error ) !in case called from another MPI package
00624       if(PRESENT(mpp_comm)) then
00625         comm_intern=mpp_comm
00626       else
00627         comm_intern=MPI_COMM_WORLD
00628       endif
00629 !
00630       if(present(logfile)) then
00631         etcfile=trim(logfile)
00632         configfile=trim(logfile)
00633         logfile_defined=.true.
00634       endif
00635 !
00636       if( .NOT.opened )then
00637       call MPI_INIT(error)
00638         comm_intern=MPI_COMM_WORLD
00639       endif
00640       call MPI_COMM_RANK(  comm_intern, pe,   error )
00641       call MPI_COMM_SIZE(  comm_intern, npes, error )
00642       allocate( mpp_request(0:npes-1) )
00643       mpp_request(:) = MPI_REQUEST_NULL
00644 #endif
00645       module_is_initialized = .TRUE.
00646 
00647 !PEsets: make defaults illegal
00648       peset(:)%count = -1
00649 #ifdef use_libSMA
00650       peset(:)%start = -1
00651       peset(:)%log2stride = -1
00652 #elif use_libMPI
00653       peset(:)%id = -1
00654       peset(:)%group = -1
00655 #endif
00656 !0=single-PE, initialized so that count returns 1
00657       peset(0)%count = 1
00658       allocate( peset(0)%list(1) )
00659       peset(0)%list = pe
00660 #ifdef use_libMPI
00661       current_peset_num = 0
00662       peset(0)%id = comm_intern
00663       call MPI_COMM_GROUP(  peset(0)%id, peset(0)%group, error )
00664 #endif
00665       world_peset_num = get_peset( (/(i,i=0,npes-1)/) )
00666       current_peset_num = world_peset_num !initialize current PEset to world
00667 
00668 !initialize clocks
00669       call SYSTEM_CLOCK( count=tick0, count_rate=ticks_per_sec, count_max=max_ticks )
00670       tick_rate = 1./ticks_per_sec
00671       clock0 = mpp_clock_id( 'Total runtime', flags=MPP_CLOCK_SYNC )
00672 
00673       if( PRESENT(flags) )then
00674           debug   = flags.EQ.MPP_DEBUG
00675           verbose = flags.EQ.MPP_VERBOSE .OR. debug
00676       end if
00677 
00678 #ifdef use_libSMA
00679 #ifdef use_shmalloc
00680 !we use shpalloc to ensure all these are remotely accessible
00681       len=0; ptr_sync = LOC(pe)   !null initialization
00682       call mpp_malloc( ptr_sync,        size(TRANSFER(sync,word)),            len )
00683       len=0; ptr_status = LOC(pe)  !null initialization
00684       call mpp_malloc( ptr_status, npes*size(TRANSFER(status(0),word)),   len )
00685       len=0; ptr_remote = LOC(pe) !null initialization
00686       call mpp_malloc( ptr_remote, npes*size(TRANSFER(remote_data_loc(0),word)), len )
00687       len=0; ptr_from = LOC(pe)   !null initialization
00688       call mpp_malloc( ptr_from,        size(TRANSFER(mpp_from_pe,word)),     len )
00689 #else
00690       allocate( status(0:npes-1) )
00691       allocate( remote_data_loc(0:npes-1) )
00692 #endif
00693       sync(:) = SHMEM_SYNC_VALUE
00694       status(0:npes-1) = MPP_READY
00695       remote_data_loc(0:npes-1) = MPP_WAIT
00696       call mpp_set_stack_size(32768) !default initial value
00697 #endif
00698 
00699 !logunit: log messages are written to configfile.out by default
00700       etc_unit=get_unit()
00701 !      write( etcfile,'(a,i4.4)' )trim(etcfile)//'.', pe
00702 !rv
00703 !rv Status 'REPLACE' eads to an unpredictable behaviour on SX-6
00704 !rv      if( pe.EQ.root_pe )open( unit=etc_unit, file=trim(etcfile), status='REPLACE' )
00705 !rv
00706       if( pe.EQ.root_pe )then
00707           if(logfile_defined) then
00708             inquire( file=trim(etcfile), opened=opened,number=io_num)
00709             if(opened) then
00710               etc_unit=io_num
00711             else
00712               open( unit=etc_unit, file=trim(etcfile), status='UNKNOWN' )
00713               close(etc_unit)
00714             endif
00715           else
00716             open( unit=etc_unit, file=trim(etcfile), status='UNKNOWN' )
00717             close(etc_unit)
00718           endif
00719       endif
00720 
00721       call mpp_sync()                         
00722 !rv
00723       if(present(logfile)) then
00724         inquire( file=trim(etcfile), opened=opened,number=io_num)
00725         if(opened) then
00726           if( pe.NE.root_pe ) then
00727             etc_unit=io_num
00728           endif
00729         else
00730           if( pe.NE.root_pe ) then
00731             open( unit=etc_unit, file=trim(etcfile), status='UNKNOWN' )
00732             close(etc_unit)
00733           endif
00734         endif
00735       else
00736         if( pe.NE.root_pe ) then
00737           open( unit=etc_unit, file=trim(etcfile), status='OLD' )
00738           close(etc_unit)
00739         endif
00740       endif
00741 !            write(0,*)'etc_unit=',etc_unit
00742 !rv
00743 
00744 !if optional argument logunit=stdout, write messages to stdout instead.
00745 !if specifying non-defaults, you must specify units not yet in use.
00746 !      if( PRESENT(in) )then
00747 !          inquire( unit=in, opened=opened )
00748 !          if( opened )call mpp_error( FATAL, 'MPP_INIT: unable to open stdin.' )
00749 !          in_unit=in
00750 !      end if
00751 !      if( PRESENT(out) )then
00752 !          inquire( unit=out, opened=opened )
00753 !          if( opened )call mpp_error( FATAL, 'MPP_INIT: unable to open stdout.' )
00754 !          out_unit=out
00755 !      end if
00756 !      if( PRESENT(err) )then
00757 !          inquire( unit=err, opened=opened )
00758 !          if( opened )call mpp_error( FATAL, 'MPP_INIT: unable to open stderr.' )
00759 !          err_unit=err
00760 !      end if
00761 !      log_unit=get_unit()
00762 !      if( PRESENT(log) )then
00763 !          inquire( unit=log, opened=opened )
00764 !          if( opened .AND. log.NE.out_unit )call mpp_error( FATAL, 'MPP_INIT: unable to open stdlog.' )
00765 !          log_unit=log
00766 !      end if
00767 !!log_unit can be written to only from root_pe, all others write to stdout
00768 !      if( log_unit.NE.out_unit )then
00769 !          inquire( unit=log_unit, opened=opened )
00770 !          if( opened )call mpp_error( FATAL, 'MPP_INIT: specified unit for stdlog already in use.' )
00771 !          if( pe.EQ.root_pe )open( unit=log_unit, file=trim(configfile), status='REPLACE' )
00772 !          call mpp_sync()
00773 !          if( pe.NE.root_pe )open( unit=log_unit, file=trim(configfile), status='OLD' )
00774 !      end if
00775       if( pe.EQ.root_pe )then
00776           log_unit = get_unit()
00777 !rv Status 'REPLACE' leads to an unpredictable  behaviour on the SX-6.
00778 !rv          open( unit=log_unit, file=trim(configfile), status='REPLACE' )
00779         if(logfile_defined) then
00780           inquire( file=trim(configfile), opened=opened,number=io_num)
00781           if(opened)then
00782             log_unit=io_num
00783           else
00784             open( unit=log_unit, file=trim(configfile), status='UNKNOWN' )
00785             close(log_unit)
00786           endif
00787         else
00788           open( unit=log_unit, file=trim(configfile), status='UNKNOWN' )
00789           close(log_unit)
00790         endif
00791       end if
00792 !messages
00793       if( verbose )call mpp_error( NOTE, 'MPP_INIT: initializing MPP module...' )
00794       if( pe.EQ.root_pe )then
00795           WRITE( nulprt,'(/a)' )'  ---- MPP module '//TRIM(version)//TRIM(tagname)
00796           WRITE( nulprt,'(a,i4)' )'  ---- MPP started with NPES=', npes
00797 #ifdef use_libSMA
00798           WRITE( nulprt,'(a)' )'  ---- Using SMA (shmem) library for message passing...'
00799 #endif
00800 #ifdef use_libMPI
00801           WRITE( nulprt,'(a)' )'  ---- Using MPI library for message passing...'
00802 #endif
00803           write( nulprt, '(a,es12.4,a,i20,a)' ) &
00804                '  ---- Realtime clock resolution=', tick_rate, ' sec (', ticks_per_sec, ' ticks/sec)'
00805           write( nulprt, '(a,es12.4,a,i20,a)' ) &
00806                '  ---- Clock rolls over after ', max_ticks*tick_rate, ' sec (', max_ticks, ' ticks)'
00807       end if
00808       call mpp_clock_begin(clock0)
00809 
00810       return
00811     end subroutine mpp_init
00812 
00813     function stdin()
00814       integer :: stdin
00815       stdin = in_unit
00816       return
00817     end function stdin
00818 
00819     function stdout()
00820       integer :: stdout
00821       integer::tmp_unit
00822       logical::opened
00823       stdout = out_unit
00824       if( pe.NE.root_pe ) then
00825         stdout = etc_unit
00826         inquire( file=trim(etcfile), opened=opened ,number=tmp_unit)
00827         if(.not.opened) open(file=etcfile,unit=etc_unit,status='UNKNOWN')
00828       endif
00829       if(logfile_defined)then
00830 
00831         if( pe.EQ.root_pe )then
00832             inquire( file=trim(etcfile), opened=opened ,number=tmp_unit)
00833             if( opened )then
00834 
00835                 out_unit=tmp_unit
00836                 stdout=out_unit
00837 !              write(0,*) 'stdout:',tmp_unit,etcfile
00838                 call FLUSH(out_unit)
00839 
00840             else
00841                 tmp_unit=get_unit()
00842 !              write(0,*) 'stdout(not):',tmp_unit,etcfile
00843 !ac
00844                 open( unit=tmp_unit, status='UNKNOWN', file=trim(etcfile), err=10 )
00845 !ac
00846               out_unit=tmp_unit
00847             end if
00848 
00849             stdout = out_unit
00850 
00851         endif
00852         
00853       endif
00854       return
00855  10 logfile_defined=.false.
00856     call mpp_error( FATAL, 'STDOUT: unable to open '//trim(etcfile))
00857     end function stdout
00858 
00859     function stderr()
00860       integer :: stderr
00861       integer :: tmp_unit
00862       logical :: opened
00863 
00864       stderr = err_unit
00865 
00866       if(logfile_defined)then
00867 
00868         if( pe.EQ.root_pe )then
00869             inquire( file=trim(configfile), opened=opened ,number=tmp_unit)
00870             if( opened )then
00871                 err_unit=tmp_unit
00872                 stderr=err_unit
00873                 call FLUSH(err_unit)
00874 !              write(0,*) 'stderr:',tmp_unit
00875             else
00876 !              write(0,*) 'stderr(not):',tmp_unit
00877                 tmp_unit=get_unit()
00878 !ac
00879                 open( unit=tmp_unit, status='UNKNOWN', file=trim(configfile), err=10 )
00880 !ac
00881               err_unit=tmp_unit
00882             end if
00883 
00884             stderr = err_unit
00885 
00886         endif
00887         
00888       endif
00889       return
00890 
00891  10 logfile_defined=.false.
00892     call mpp_error( FATAL, 'STDERR: unable to open '//trim(configfile))
00893 
00894     end function stderr
00895 
00896     function stdlog()
00897       integer :: stdlog
00898       logical :: opened
00899       if( pe.EQ.root_pe )then
00900           inquire( file=trim(configfile), opened=opened )
00901           if( opened )then
00902           call FLUSH(log_unit)
00903           else
00904               log_unit=get_unit()
00905 !ac
00906               open( unit=log_unit, status='UNKNOWN', file=trim(configfile), err=10 )
00907 !ac
00908           end if
00909           stdlog = log_unit
00910       else
00911           stdlog = etc_unit
00912       end if
00913       return
00914    10 call mpp_error( FATAL, 'STDLOG: unable to open '//trim(configfile)//'.' )
00915     end function stdlog
00916 
00917     subroutine mpp_exit()
00918 !to be called at the end of a run
00919       integer :: i, j, k, n, nmax
00920       real :: t, tmin, tmax, tavg, tstd
00921       real :: m, mmin, mmax, mavg, mstd
00922       real :: t_total
00923 
00924       if( .NOT.module_is_initialized )return
00925       call mpp_clock_end(clock0)
00926       t_total = clocks(clock0)%total_ticks*tick_rate
00927       if( clock_num.GT.0 )then
00928           if( ANY(clocks(1:clock_num)%detailed) )then
00929               call sum_clock_data; call dump_clock_summary
00930           end if
00931           if( pe.EQ.root_pe )then
00932               write( nulprt,'(/a,i4,a)' ) 'Tabulating mpp_clock statistics across ', npes, ' PEs...'
00933               if( ANY(clocks(1:clock_num)%detailed) ) &
00934                    write( nulprt,'(a)' )'   ... see mpp_clock.out.#### for details on individual PEs.'
00935               write( nulprt,'(/32x,a)' ) '          tmin          tmax          tavg          tstd  tfrac'
00936           end if
00937           do i = 1,clock_num
00938              call mpp_set_current_pelist() !implied global barrier
00939              current_peset_num = clocks(i)%peset_num
00940              if( .NOT.ANY(peset(current_peset_num)%list(:).EQ.pe) )cycle
00941 !times between mpp_clock ticks
00942              t = clocks(i)%total_ticks*tick_rate
00943              tmin = t; call mpp_min(tmin)
00944              tmax = t; call mpp_max(tmax)
00945              tavg = t; call mpp_sum(tavg); tavg = tavg/mpp_npes()
00946              tstd = (t-tavg)**2; call mpp_sum(tstd); tstd = sqrt( tstd/mpp_npes() )
00947              if( pe.EQ.root_pe )write( nulprt,'(a32,4f14.6,f7.3)' ) &
00948                   clocks(i)%name, tmin, tmax, tavg, tstd, tavg/t_total
00949           end do
00950           if( ANY(clocks(1:clock_num)%detailed) .AND. pe.EQ.root_pe )write( nulprt,'(/32x,a)' ) &
00951                '       tmin       tmax       tavg       tstd       mmin       mmax       mavg       mstd  mavg/tavg'
00952           do i = 1,clock_num
00953 !messages: bytelengths and times
00954          if( .NOT.clocks(i)%detailed )cycle
00955              do j = 1,MAX_EVENT_TYPES
00956                 n = clocks(i)%events(j)%calls; nmax = n
00957                 call mpp_max(nmax)
00958                 if( nmax.NE.0 )then
00959 !don't divide by n because n might be 0
00960                     m = 0
00961                     if( n.GT.0 )m = sum(clocks(i)%events(j)%bytes(1:n))
00962                     mmin = m; call mpp_min(mmin)
00963                     mmax = m; call mpp_max(mmax)
00964                     mavg = m; call mpp_sum(mavg); mavg = mavg/mpp_npes()
00965                     mstd = (m-mavg)**2; call mpp_sum(mstd); mstd = sqrt( mstd/mpp_npes() )
00966                     t = 0
00967                     if( n.GT.0 )t = sum(clocks(i)%events(j)%ticks(1:n))*tick_rate
00968                     tmin = t; call mpp_min(tmin)
00969                     tmax = t; call mpp_max(tmax)
00970                     tavg = t; call mpp_sum(tavg); tavg = tavg/mpp_npes()
00971                     tstd = (t-tavg)**2; call mpp_sum(tstd); tstd = sqrt( tstd/mpp_npes() )
00972                     if( pe.EQ.root_pe )write( nulprt,'(a32,4f11.3,5es11.3)' ) &
00973                          trim(clocks(i)%name)//' '//trim(clocks(i)%events(j)%name), &
00974                          tmin, tmax, tavg, tstd, mmin, mmax, mavg, mstd, mavg/tavg
00975                 end if
00976              end do
00977           end do
00978       end if
00979       call mpp_set_current_pelist()
00980       call mpp_sync()
00981       call mpp_max(mpp_stack_hwm)
00982       if( pe.EQ.root_pe )write( nulprt,* )'MPP_STACK high water mark=', mpp_stack_hwm
00983 #ifdef use_libMPI
00984 !reiner Let's do the MPI_finalize outside the mpp environment
00985 !reiner      call MPI_FINALIZE(error)
00986 #endif
00987 
00988       return
00989     end subroutine mpp_exit
00990 
00991     function mpp_pe()
00992       integer :: mpp_pe
00993 
00994       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_PE: You must first call mpp_init.' )
00995       mpp_pe = pe
00996       return
00997     end function mpp_pe
00998 
00999     function mpp_node()
01000       integer :: mpp_node
01001 
01002       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_NODE: You must first call mpp_init.' )
01003       mpp_node = node
01004       return
01005     end function mpp_node
01006 
01007     function mpp_npes()
01008       integer :: mpp_npes
01009 
01010       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_NPES: You must first call mpp_init.' )
01011 !      mpp_npes = npes
01012       mpp_npes = size(peset(current_peset_num)%list)
01013       return
01014     end function mpp_npes
01015 
01016     function mpp_root_pe()
01017       integer :: mpp_root_pe
01018 
01019       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_ROOT_PE: You must first call mpp_init.' )
01020       mpp_root_pe = root_pe
01021       return
01022     end function mpp_root_pe
01023 
01024     subroutine mpp_set_root_pe(num)
01025       integer, intent(in) :: num
01026       logical :: opened
01027 
01028       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_SET_ROOT_PE: You must first call mpp_init.' )
01029       if( .NOT.(ANY(num.EQ.peset(current_peset_num)%list)) ) &
01030            call mpp_error( FATAL, 'MPP_SET_ROOT_PE: you cannot set a root PE outside the current pelist.' )
01031 !actions to take if root_pe has changed:
01032 ! open log_unit on new root_pe, close it on old root_pe and point its log_unit to stdout.
01033 !      if( num.NE.root_pe )then  !root_pe has changed
01034 !          if( pe.EQ.num )then
01035 !!on the new root_pe
01036 !              if( log_unit.NE.out_unit )then
01037 !                  inquire( unit=log_unit, opened=opened )
01038 !                  if( .NOT.opened )open( unit=log_unit, status='OLD', file=trim(configfile), position='APPEND' )
01039 !              end if
01040 !          else if( pe.EQ.root_pe )then
01041 !!on the old root_pe
01042 !              if( log_unit.NE.out_unit )then
01043 !                  inquire( unit=log_unit, opened=opened )
01044 !                  if( opened )close(log_unit)
01045 !                  log_unit = out_unit
01046 !              end if
01047 !          end if
01048 !      end if
01049       root_pe = num
01050       return
01051     end subroutine mpp_set_root_pe
01052 
01053     subroutine mpp_declare_pelist( pelist, name )
01054 !this call is written specifically to accommodate a brain-dead MPI restriction
01055 !that requires a parent communicator to create a child communicator:
01056 !in other words: a pelist cannot go off and declare a communicator, but every PE
01057 !in the parent, including those not in pelist(:), must get together for the
01058 !MPI_COMM_CREATE call. The parent is typically  peset(0)%id, though it could also
01059 !be a subset that includes all PEs in pelist.
01060 !This restriction does not apply to SMA but to have uniform code,
01061 !you may as well call it. It must be placed in a context where all PEs call it.
01062 !Subsequent calls that use the pelist should be called PEs in the pelist only.
01063       integer, intent(in) :: pelist(:)
01064       character(len=*), optional :: name
01065       integer :: i
01066 
01067       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_DECLARE_PELIST: You must first call mpp_init.' )
01068       i = get_peset(pelist)
01069       write( peset(i)%name,'(a,i2.2)' ) 'PElist', i !default name
01070       if( PRESENT(name) )peset(i)%name = name
01071       return
01072     end subroutine mpp_declare_pelist
01073 
01074     subroutine mpp_set_current_pelist( pelist )
01075 !Once we branch off into a PE subset, we want subsequent "global" calls to
01076 !sync only across this subset. This is declared as the current pelist (peset(current_peset_num)%list)
01077 !when current_peset all pelist ops with no pelist should apply the current pelist.
01078 !also, we set the start PE in this pelist to be the root_pe.
01079 !unlike mpp_declare_pelist, this is called by the PEs in the pelist only
01080 !so if the PEset has not been previously declared, this will hang in MPI.
01081 !if pelist is omitted, we reset pelist to the world pelist.
01082       integer, intent(in), optional :: pelist(:)
01083       integer :: i
01084 
01085       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_SET_CURRENT_PELIST: You must first call mpp_init.' )
01086       if( PRESENT(pelist) )then
01087           if( .NOT.ANY(pe.EQ.pelist) )call mpp_error( FATAL, 'MPP_SET_CURRENT_PELIST: pe must be in pelist.' )
01088           current_peset_num = get_peset(pelist)
01089       else
01090           current_peset_num = world_peset_num
01091       end if
01092       call mpp_set_root_pe( MINVAL(peset(current_peset_num)%list) )
01093       call mpp_sync()           !this is called to make sure everyone in the current pelist is here.
01094 !      npes = mpp_npes()
01095       return
01096     end subroutine mpp_set_current_pelist
01097 
01098     subroutine mpp_get_current_pelist( pelist, name )
01099 !this is created for use by mpp_define_domains within a pelist
01100 !will be published but not publicized
01101       integer, intent(out) :: pelist(:)
01102       character(len=*), intent(out), optional :: name
01103 
01104       if( size(pelist).NE.size(peset(current_peset_num)%list) ) &
01105            call mpp_error( FATAL, 'MPP_GET_CURRENT_PELIST: size(pelist) is wrong.' )
01106       pelist(:) = peset(current_peset_num)%list(:)
01107       if( PRESENT(name) )name = peset(current_peset_num)%name
01108 
01109       return
01110     end subroutine mpp_get_current_pelist
01111 
01112     function get_peset(pelist)
01113       integer :: get_peset
01114 !makes a PE set out of a PE list
01115 !a PE list is an ordered list of PEs
01116 !a PE set is a triad (start,log2stride,size) for SHMEM, an a communicator for MPI
01117 !if stride is non-uniform or not a power of 2, will return error (not required for MPI but enforced for uniformity)
01118       integer, intent(in), optional :: pelist(:)
01119       integer :: group
01120       integer :: i, n, stride
01121       integer, allocatable :: sorted(:)
01122 
01123       if( .NOT.PRESENT(pelist) )then !set it to current_peset_num
01124           get_peset = current_peset_num; return
01125       end if
01126       if( size(pelist).EQ.1 .AND. npes.GT.1 )then    !collective ops on single PEs should return
01127           get_peset = 0; return
01128       end if
01129 !make a sorted list
01130       n = 1
01131       if( ascend_sort(pelist).NE.1 )call mpp_error( FATAL, 'GET_PESET: sort error.' )   !result is the array sorted(:)
01132       if( debug )write( stderr(),* )'pelist=', pelist, ' sorted=', sorted
01133 !find if this array matches any existing peset
01134       do i = 1,peset_num
01135          if( debug )write( stderr(),'(a,3i4)' )'pe, i, peset_num=', pe, i, peset_num
01136          if( size(sorted).EQ.size(peset(i)%list) )then
01137              if( ALL(sorted.EQ.peset(i)%list) )then
01138                  deallocate(sorted)
01139                  get_peset = i; return
01140              end if
01141          end if
01142       end do
01143 !not found, so create new peset
01144       peset_num = peset_num + 1
01145       if( peset_num.GE.PESET_MAX )call mpp_error( FATAL, 'GET_PESET: number of PE sets exceeds PESET_MAX.' )
01146       i = peset_num             !shorthand
01147 !create list
01148       allocate( peset(i)%list(size(sorted)) )
01149       peset(i)%list(:) = sorted(:)
01150       peset(i)%count = size(sorted)
01151 #ifdef use_libSMA
01152       peset(i)%start = sorted(1)
01153       if( size(sorted).GT.1 )then
01154           stride = sorted(2)-sorted(1)
01155           if( ANY(sorted(2:n)-sorted(1:n-1).NE.stride) ) &
01156                call mpp_error( WARNING, 'GET_PESET: pelist must have constant stride.' )
01157           peset(i)%log2stride = nint( log(real(stride))/log(2.) )
01158           if( 2**peset(i)%log2stride.NE.stride )call mpp_error( WARNING, 'GET_PESET: pelist must have power-of-2 stride.' )
01159       else
01160           peset(i)%log2stride = 0
01161       end if
01162 #elif use_libMPI
01163       call MPI_GROUP_INCL( peset(current_peset_num)%group, size(sorted), sorted, peset(i)%group, error )
01164       call MPI_COMM_CREATE( peset(current_peset_num)%id, peset(i)%group, peset(i)%id, error )
01165 #endif
01166       deallocate(sorted)
01167       get_peset = i
01168 
01169       return
01170 
01171       contains
01172         
01173         recursive function ascend_sort(a) result(a_sort)
01174           integer :: a_sort
01175           integer, intent(in) :: a(:)
01176           integer :: b, i
01177           if( size(a).EQ.1 .OR. ALL(a.EQ.a(1)) )then
01178               allocate( sorted(n) )
01179               sorted(n) = a(1)
01180               a_sort = n
01181               return
01182           end if
01183           b = minval(a)
01184           n = n + 1
01185           i = ascend_sort( pack(a,mask=a.NE.b) )
01186           a_sort = i - 1
01187           sorted(i-1) = b
01188           return
01189         end function ascend_sort
01190 
01191     end function get_peset
01192 
01193 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01194 !                                                                             !
01195 !                        PERFORMANCE PROFILING CALLS                          !
01196 !                                                                             !
01197 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01198     subroutine mpp_clock_set_grain( grain )
01199       integer, intent(in) :: grain
01200 !set the granularity of times: only clocks whose grain is lower than
01201 !clock_grain are triggered, finer-grained clocks are dormant.
01202 !clock_grain is initialized to HUGE(1), so all clocks are triggered if
01203 !this is never called.   
01204       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_CLOCK_SET_GRAIN: You must first call mpp_init.' )
01205 
01206       clock_grain = grain
01207       return
01208     end subroutine mpp_clock_set_grain
01209 
01210     function mpp_clock_id( name, flags, grain )
01211 !return an ID for a new or existing clock
01212       integer :: mpp_clock_id
01213       character(len=*), intent(in) :: name
01214       integer, intent(in), optional :: flags, grain
01215       integer :: i
01216 
01217       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_CLOCK_ID: You must first call mpp_init.' )
01218 !if grain is present, the clock is only triggered if it
01219 !is low ("coarse") enough: compared to clock_grain
01220 !finer-grained clocks are dormant.
01221 !if grain is absent, clock is triggered.
01222       if( PRESENT(grain) )then
01223           if( grain.GT.clock_grain )then
01224               mpp_clock_id = 0
01225               return
01226           end if
01227       end if
01228       mpp_clock_id = 1
01229       if( clock_num.EQ.0 )then  !first
01230 !         allocate( clocks(MAX_CLOCKS) )
01231           clock_num = mpp_clock_id
01232           call clock_init(mpp_clock_id,name,flags)
01233       else
01234           FIND_CLOCK: do while( trim(name).NE.trim(clocks(mpp_clock_id)%name) )
01235              mpp_clock_id = mpp_clock_id + 1
01236              if( mpp_clock_id.GT.clock_num )then
01237                  if( mpp_clock_id.GT.MAX_CLOCKS )then
01238                      call mpp_error( WARNING, 'MPP_CLOCK_ID: too many clock requests, this one is ignored.' )
01239                  else               !new clock: initialize
01240                      clock_num = mpp_clock_id
01241                      call clock_init(mpp_clock_id,name,flags)
01242                      exit FIND_CLOCK
01243                  end if
01244              end if
01245           end do FIND_CLOCK
01246       endif
01247       return
01248     end function mpp_clock_id
01249 
01250     subroutine mpp_clock_begin(id)
01251       integer, intent(in) :: id
01252 
01253       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_CLOCK_BEGIN: You must first call mpp_init.' )
01254       if( id.EQ.0 )return
01255       if( id.LT.0 .OR. id.GT.clock_num )call mpp_error( FATAL, 'MPP_CLOCK_BEGIN: invalid id.' )
01256 
01257       if( clocks(id)%peset_num.EQ.0 )clocks(id)%peset_num = current_peset_num
01258       if( clocks(id)%peset_num.NE.current_peset_num ) &
01259            call mpp_error( FATAL, 'MPP_CLOCK_BEGIN: cannot change pelist context of a clock.' )
01260       if( clocks(id)%sync_on_begin )then
01261 !do an untimed sync at the beginning of the clock
01262 !this puts all PEs in the current pelist on par, so that measurements begin together
01263 !ending time will be different, thus measuring load imbalance for this clock.
01264           current_clock = 0; call mpp_sync()
01265       end if
01266       current_clock = id
01267       call SYSTEM_CLOCK( clocks(id)%tick )
01268       return
01269     end subroutine mpp_clock_begin
01270 
01271     subroutine mpp_clock_end(id)
01272 !the id argument is not used for anything at present
01273       integer, intent(in), optional :: id
01274       integer(LONG_KIND) :: delta
01275 
01276       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_CLOCK_END: You must first call mpp_init.' )
01277       if( id.EQ.0 )return
01278       if( id.LT.0 .OR. id.GT.clock_num )call mpp_error( FATAL, 'MPP_CLOCK_BEGIN: invalid id.' )
01279       call SYSTEM_CLOCK(end_tick)
01280       if( clocks(id)%peset_num.NE.current_peset_num ) &
01281            call mpp_error( FATAL, 'MPP_CLOCK_END: cannot change pelist context of a clock.' )
01282       delta = end_tick - clocks(id)%tick
01283       if( delta.LT.0 )then
01284           write( stderr(),* )'pe, id, start_tick, end_tick, delta, max_ticks=', pe, id, clocks(id)%tick, end_tick, delta, max_ticks
01285           delta = delta + max_ticks + 1
01286           call mpp_error( WARNING, 'MPP_CLOCK_END: Clock rollover, assumed single roll.' )
01287       end if
01288       clocks(id)%total_ticks = clocks(id)%total_ticks + delta
01289       current_clock = 0
01290       return
01291     end subroutine mpp_clock_end
01292 
01293     subroutine increment_current_clock( event_id, bytes )
01294       integer, intent(in) :: event_id
01295       integer, intent(in), optional :: bytes
01296       integer :: n
01297       integer(LONG_KIND) :: delta
01298 
01299       if( current_clock.EQ.0 )return
01300       if( current_clock.LT.0 .OR. current_clock.GT.clock_num )call mpp_error( FATAL, 'MPP_CLOCK_BEGIN: invalid current_clock.' )
01301       if( .NOT.clocks(current_clock)%detailed )return
01302       call SYSTEM_CLOCK(end_tick)
01303       n = clocks(current_clock)%events(event_id)%calls + 1
01304 
01305       if( n.EQ.MAX_EVENTS )call mpp_error( WARNING, &
01306            'MPP_CLOCK: events exceed MAX_EVENTS, ignore detailed profiling data for clock '//trim(clocks(current_clock)%name) )
01307       if( n.GT.MAX_EVENTS )return
01308 
01309       clocks(current_clock)%events(event_id)%calls = n
01310       delta = end_tick - start_tick
01311       if( delta.LT.0 )then
01312           delta = delta + max_ticks + 1
01313           call mpp_error( WARNING, 'MPP_CLOCK_END: Clock rollover, assumed single roll.' )
01314       end if
01315       clocks(current_clock)%events(event_id)%ticks(n) = delta
01316       if( PRESENT(bytes) )clocks(current_clock)%events(event_id)%bytes(n) = bytes
01317       return
01318     end subroutine increment_current_clock
01319 
01320   subroutine dump_clock_summary()
01321     implicit none
01322 
01323     real :: total_time,total_time_all,total_data
01324     real :: msg_size,eff_BW,s
01325     integer :: SD_UNIT
01326     integer :: total_calls
01327     integer :: i,j,k,ct
01328     integer :: msg_cnt
01329     character(len=2)  :: u
01330     character(len=18) :: filename
01331     character(len=20),dimension(MAX_BINS),save :: bin
01332 
01333     data bin( 1)  /'  0   -    8    B:  '/
01334     data bin( 2)  /'  8   -   16    B:  '/
01335     data bin( 3)  /' 16   -   32    B:  '/
01336     data bin( 4)  /' 32   -   64    B:  '/
01337     data bin( 5)  /' 64   -  128    B:  '/
01338     data bin( 6)  /'128   -  256    B:  '/
01339     data bin( 7)  /'256   -  512    B:  '/
01340     data bin( 8)  /'512   - 1024    B:  '/
01341     data bin( 9)  /'  1.0 -    2.1 KB:  '/
01342     data bin(10)  /'  2.1 -    4.1 KB:  '/
01343     data bin(11)  /'  4.1 -    8.2 KB:  '/
01344     data bin(12)  /'  8.2 -   16.4 KB:  '/
01345     data bin(13)  /' 16.4 -   32.8 KB:  '/
01346     data bin(14)  /' 32.8 -   65.5 KB:  '/
01347     data bin(15)  /' 65.5 -  131.1 KB:  '/
01348     data bin(16)  /'131.1 -  262.1 KB:  '/
01349     data bin(17)  /'262.1 -  524.3 KB:  '/
01350     data bin(18)  /'524.3 - 1048.6 KB:  '/
01351     data bin(19)  /'  1.0 -    2.1 MB:  '/
01352     data bin(20)  /' >2.1          MB:  '/
01353 
01354     if( .NOT.ANY(clocks(1:clock_num)%detailed) )return
01355     write( filename,'(a,i4.4)' )'mpp_clock.out.', pe
01356 
01357     SD_UNIT = get_unit()
01358     open(SD_UNIT,file=trim(filename),form='formatted')
01359 
01360     COMM_TYPE: do ct = 1,clock_num
01361 
01362       if( .NOT.clocks(ct)%detailed )cycle
01363       write(SD_UNIT,*) &
01364           clock_summary(ct)%name(1:15),' Communication Data for PE ',pe
01365 
01366       write(SD_UNIT,*) ' '
01367       write(SD_UNIT,*) ' '
01368 
01369       total_time_all = 0.0
01370       EVENT_TYPE: do k = 1,MAX_EVENT_TYPES-1
01371 
01372         if(clock_summary(ct)%event(k)%total_time == 0.0)cycle
01373 
01374         total_time = clock_summary(ct)%event(k)%total_time
01375         total_time_all = total_time_all + total_time
01376         total_data = clock_summary(ct)%event(k)%total_data
01377         total_calls = clock_summary(ct)%event(k)%total_cnts
01378 
01379         write(SD_UNIT,1000) clock_summary(ct)%event(k)%name(1:9) // ':'
01380 
01381         write(SD_UNIT,1001) 'Total Data: ',total_data*1.0e-6, &
01382                             'MB; Total Time: ', total_time, &
01383                             'secs; Total Calls: ',total_calls
01384 
01385         write(SD_UNIT,*) ' '
01386         write(SD_UNIT,1002) '     Bin            Counts      Avg Size        Eff B/W'
01387         write(SD_UNIT,*) ' '
01388 
01389         BIN_LOOP: do j=1,MAX_BINS
01390 
01391           if(clock_summary(ct)%event(k)%msg_size_cnts(j)==0)cycle
01392 
01393           if(j<=8)then
01394             s = 1.0
01395             u = ' B'
01396           elseif(j<=18)then
01397             s = 1.0e-3
01398             u = 'KB'
01399           else
01400             s = 1.0e-6
01401             u = 'MB'
01402           endif
01403 
01404           msg_cnt = clock_summary(ct)%event(k)%msg_size_cnts(j)
01405           msg_size = &
01406             s*(clock_summary(ct)%event(k)%msg_size_sums(j)/real(msg_cnt))
01407           eff_BW = (1.0e-6)*( clock_summary(ct)%event(k)%msg_size_sums(j) / &
01408                                   clock_summary(ct)%event(k)%msg_time_sums(j) )
01409 
01410           write(SD_UNIT,1003) bin(j),msg_cnt,msg_size,u,eff_BW
01411 
01412         end do BIN_LOOP
01413 
01414         write(SD_UNIT,*) ' '
01415         write(SD_UNIT,*) ' '
01416       end do EVENT_TYPE
01417 
01418    ! "Data-less" WAIT
01419 
01420       if(clock_summary(ct)%event(MAX_EVENT_TYPES)%total_time>0.0)then
01421 
01422         total_time = clock_summary(ct)%event(MAX_EVENT_TYPES)%total_time
01423         total_time_all = total_time_all + total_time
01424         total_calls = clock_summary(ct)%event(MAX_EVENT_TYPES)%total_cnts
01425 
01426         write(SD_UNIT,1000) clock_summary(ct)%event(MAX_EVENT_TYPES)%name(1:9) // ':'
01427 
01428         write(SD_UNIT,1004) 'Total Calls: ',total_calls,'; Total Time: ', &
01429                              total_time,'secs'
01430 
01431       endif
01432 
01433       write(SD_UNIT,*) ' '
01434       write(SD_UNIT,1005) 'Total communication time spent for ' // &
01435                       clock_summary(ct)%name(1:9) // ': ',total_time_all,'secs'
01436       write(SD_UNIT,*) ' '
01437       write(SD_UNIT,*) ' '
01438       write(SD_UNIT,*) ' '
01439 
01440     end do COMM_TYPE
01441 
01442     close(SD_UNIT)
01443 
01444 1000  format(a)
01445 1001  format(a,f8.2,a,f8.2,a,i6)
01446 1002  format(a)
01447 1003  format(a,i6,'    ','  ',f6.1,a,'    ',f7.3,'MB/sec')
01448 1004  format(a,i8,a,f9.2,a)
01449 1005  format(a,f9.2,a)
01450     return
01451   end subroutine dump_clock_summary
01452 
01453       integer function get_unit()
01454         implicit none
01455 
01456         integer,save :: i
01457         logical :: l_open
01458 
01459         i = 10
01460         do i=10,99
01461            inquire(unit=i,opened=l_open)
01462            if(.not.l_open)exit
01463         end do
01464 
01465         if(i==100)then
01466             call mpp_error(FATAL,'Unable to get I/O unit')
01467         else
01468             get_unit = i
01469         endif
01470 
01471         return
01472       end function get_unit
01473 
01474   subroutine sum_clock_data()
01475     implicit none
01476 
01477     integer :: i,j,k,ct,event_size,event_cnt
01478     real    :: msg_time
01479 
01480     CLOCK_TYPE: do ct=1,clock_num
01481       if( .NOT.clocks(ct)%detailed )cycle
01482       EVENT_TYPE: do j=1,MAX_EVENT_TYPES-1
01483         event_cnt = clocks(ct)%events(j)%calls
01484         EVENT_SUMMARY: do i=1,event_cnt
01485 
01486         clock_summary(ct)%event(j)%total_cnts = &
01487               clock_summary(ct)%event(j)%total_cnts + 1
01488 
01489         event_size = clocks(ct)%events(j)%bytes(i)
01490 
01491         k = find_bin(event_size)
01492 
01493         clock_summary(ct)%event(j)%msg_size_cnts(k) = &
01494               clock_summary(ct)%event(j)%msg_size_cnts(k) + 1
01495 
01496         clock_summary(ct)%event(j)%msg_size_sums(k) = &
01497               clock_summary(ct)%event(j)%msg_size_sums(k) &
01498             + clocks(ct)%events(j)%bytes(i)
01499 
01500         clock_summary(ct)%event(j)%total_data = &
01501               clock_summary(ct)%event(j)%total_data &
01502             + clocks(ct)%events(j)%bytes(i)
01503 
01504         msg_time = clocks(ct)%events(j)%ticks(i)
01505         msg_time = tick_rate * real( clocks(ct)%events(j)%ticks(i) )
01506 
01507         clock_summary(ct)%event(j)%msg_time_sums(k) = &
01508               clock_summary(ct)%event(j)%msg_time_sums(k) + msg_time
01509 
01510         clock_summary(ct)%event(j)%total_time = &
01511               clock_summary(ct)%event(j)%total_time + msg_time
01512 
01513         end do EVENT_SUMMARY
01514       end do EVENT_TYPE
01515 
01516       j = MAX_EVENT_TYPES ! WAITs
01517            ! "msg_size_cnts" doesn't really mean anything for WAIT
01518            ! but position will be used to store number of counts for now.
01519 
01520       event_cnt = clocks(ct)%events(j)%calls
01521       clock_summary(ct)%event(j)%msg_size_cnts(1) = event_cnt
01522       clock_summary(ct)%event(j)%total_cnts       = event_cnt
01523 
01524       msg_time = tick_rate * real( sum ( clocks(ct)%events(j)%ticks(1:event_cnt) ) )
01525       clock_summary(ct)%event(j)%msg_time_sums(1) = &
01526               clock_summary(ct)%event(j)%msg_time_sums(1) + msg_time
01527 
01528       clock_summary(ct)%event(j)%total_time = clock_summary(ct)%event(j)%msg_time_sums(1)
01529 
01530     end do CLOCK_TYPE
01531 
01532     return
01533     contains
01534       integer function find_bin(event_size)
01535         implicit none
01536 
01537         integer,intent(in) :: event_size
01538         integer :: k,msg_size
01539 
01540         msg_size = 8
01541         k = 1
01542         do while(event_size>msg_size .and. k<MAX_BINS)
01543            k = k+1
01544            msg_size = msg_size*2
01545         end do
01546         find_bin = k
01547         return
01548       end function find_bin
01549 
01550   end subroutine sum_clock_data
01551 
01552   subroutine clock_init(id,name,flags)
01553     integer, intent(in) :: id
01554     character(len=*), intent(in) :: name
01555     integer, intent(in), optional :: flags
01556     integer :: i
01557 
01558     clocks(id)%name = name
01559     clocks(id)%tick = 0
01560     clocks(id)%total_ticks = 0
01561     clocks(id)%sync_on_begin = .FALSE.
01562     clocks(id)%detailed      = .FALSE.
01563     clocks(id)%peset_num = 0
01564     if( PRESENT(flags) )then
01565         if( BTEST(flags,0) )clocks(id)%sync_on_begin = .TRUE.
01566         if( BTEST(flags,1) )clocks(id)%detailed      = .TRUE.
01567     end if
01568     if( clocks(id)%detailed )then
01569         allocate( clocks(id)%events(MAX_EVENT_TYPES) )
01570         clocks(id)%events(EVENT_ALLREDUCE)%name = 'ALLREDUCE'
01571         clocks(id)%events(EVENT_BROADCAST)%name = 'BROADCAST'
01572         clocks(id)%events(EVENT_RECV)%name = 'RECV'
01573         clocks(id)%events(EVENT_SEND)%name = 'SEND'
01574         clocks(id)%events(EVENT_WAIT)%name = 'WAIT'
01575         do i=1,MAX_EVENT_TYPES
01576            clocks(id)%events(i)%ticks(:) = 0
01577            clocks(id)%events(i)%bytes(:) = 0
01578            clocks(id)%events(i)%calls = 0
01579         end do
01580         clock_summary(id)%name = name
01581         clock_summary(id)%event(EVENT_ALLREDUCE)%name = 'ALLREDUCE'
01582         clock_summary(id)%event(EVENT_BROADCAST)%name = 'BROADCAST'
01583         clock_summary(id)%event(EVENT_RECV)%name = 'RECV'
01584         clock_summary(id)%event(EVENT_SEND)%name = 'SEND'
01585         clock_summary(id)%event(EVENT_WAIT)%name = 'WAIT'
01586         do i=1,MAX_EVENT_TYPES
01587            clock_summary(id)%event(i)%msg_size_sums(:) = 0.0
01588            clock_summary(id)%event(i)%msg_time_sums(:) = 0.0
01589            clock_summary(id)%event(i)%total_data = 0.0
01590            clock_summary(id)%event(i)%total_time = 0.0
01591            clock_summary(id)%event(i)%msg_size_cnts(:) = 0
01592            clock_summary(id)%event(i)%total_cnts = 0
01593         end do
01594     end if
01595     return
01596   end subroutine clock_init
01597 
01598 #ifdef __sgi
01599     subroutine system_clock_sgi( count, count_rate, count_max )
01600 !mimics F90 SYSTEM_CLOCK intrinsic
01601       integer(LONG_KIND), intent(out), optional :: count, count_rate, count_max
01602       integer(LONG_KIND) :: sgi_tick, sgi_ticks_per_sec, sgi_max_tick
01603 !sgi_max_tick currently returns 64
01604 !count must return a number between 0 and count_max
01605       integer(LONG_KIND), save :: maxtick=0
01606       if( maxtick.EQ.0 )then
01607           maxtick = sgi_max_tick() !actually reports #bits in maxtick
01608           if( maxtick.LT.BIT_SIZE(maxtick) )then
01609               maxtick = 2**maxtick
01610           else
01611               maxtick = huge(maxtick)
01612           end if
01613       end if
01614       if( PRESENT(count) )then
01615           count = modulo( sgi_tick()-tick0, maxtick )
01616 !          count = sgi_tick()
01617       end if
01618       if( PRESENT(count_rate) )then
01619           count_rate = sgi_ticks_per_sec()
01620       end if
01621       if( PRESENT(count_max) )then
01622           count_max = maxtick-1
01623       end if
01624       return
01625     end subroutine system_clock_sgi
01626 #endif
01627 
01628 
01629 #ifdef use_libMPI
01630     subroutine system_clock_mpi( count, count_rate, count_max )
01631 !mimics F90 SYSTEM_CLOCK intrinsic
01632       integer(LONG_KIND), intent(out), optional :: count, count_rate, count_max
01633 !count must return a number between 0 and count_max
01634       integer(LONG_KIND), parameter :: maxtick=HUGE(count_max)
01635       logical,           save       :: first_call = .true.
01636       real(DOUBLE_KIND), save       :: count0  ! use to prevent integer overflow
01637       if ( first_call ) count0 = MPI_WTime(); first_call = .false.
01638       if( PRESENT(count) )then
01639           count = (MPI_WTime()-count0)/MPI_WTick()
01640       end if
01641       if( PRESENT(count_rate) )then
01642           count_rate = MPI_Wtick()**(-1)
01643       end if
01644       if( PRESENT(count_max) )then
01645           count_max = maxtick-1
01646       end if
01647       return
01648     end subroutine system_clock_mpi
01649 #endif
01650 
01651 
01652 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01653 !                                                                             !
01654 !                BASIC MESSAGE PASSING ROUTINE: mpp_transmit                  !
01655 !                                                                             !
01656 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01657 
01658 #define MPP_TRANSMIT_ mpp_transmit_real8
01659 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_real8_scalar
01660 #define MPP_TRANSMIT_2D_ mpp_transmit_real8_2d
01661 #define MPP_TRANSMIT_3D_ mpp_transmit_real8_3d
01662 #define MPP_TRANSMIT_4D_ mpp_transmit_real8_4d
01663 #define MPP_TRANSMIT_5D_ mpp_transmit_real8_5d
01664 #define MPP_RECV_ mpp_recv_real8
01665 #define MPP_RECV_SCALAR_ mpp_recv_real8_scalar
01666 #define MPP_RECV_2D_ mpp_recv_real8_2d
01667 #define MPP_RECV_3D_ mpp_recv_real8_3d
01668 #define MPP_RECV_4D_ mpp_recv_real8_4d
01669 #define MPP_RECV_5D_ mpp_recv_real8_5d
01670 #define MPP_SEND_ mpp_send_real8
01671 #define MPP_SEND_SCALAR_ mpp_send_real8_scalar
01672 #define MPP_SEND_2D_ mpp_send_real8_2d
01673 #define MPP_SEND_3D_ mpp_send_real8_3d
01674 #define MPP_SEND_4D_ mpp_send_real8_4d
01675 #define MPP_SEND_5D_ mpp_send_real8_5d
01676 #define MPP_BROADCAST_ mpp_broadcast_real8
01677 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_real8_scalar
01678 #define MPP_BROADCAST_2D_ mpp_broadcast_real8_2d
01679 #define MPP_BROADCAST_3D_ mpp_broadcast_real8_3d
01680 #define MPP_BROADCAST_4D_ mpp_broadcast_real8_4d
01681 #define MPP_BROADCAST_5D_ mpp_broadcast_real8_5d
01682 #define MPP_TYPE_ real(DOUBLE_KIND)
01683 #define MPP_TYPE_BYTELEN_ 8
01684 #ifdef use_LAM_MPI
01685 #define MPI_TYPE_ MPI_DOUBLE_PRECISION
01686 #else 
01687 #define MPI_TYPE_ MPI_REAL8
01688 #endif
01689 #define SHMEM_BROADCAST_ SHMEM_BROADCAST8
01690 #define SHMEM_GET_ SHMEM_GET8
01691 #include <mpp_transmit.h>
01692 
01693 #ifndef no_4byte_reals
01694 #define MPP_TRANSMIT_ mpp_transmit_real4
01695 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_real4_scalar
01696 #define MPP_TRANSMIT_2D_ mpp_transmit_real4_2d
01697 #define MPP_TRANSMIT_3D_ mpp_transmit_real4_3d
01698 #define MPP_TRANSMIT_4D_ mpp_transmit_real4_4d
01699 #define MPP_TRANSMIT_5D_ mpp_transmit_real4_5d
01700 #define MPP_RECV_ mpp_recv_real4
01701 #define MPP_RECV_SCALAR_ mpp_recv_real4_scalar
01702 #define MPP_RECV_2D_ mpp_recv_real4_2d
01703 #define MPP_RECV_3D_ mpp_recv_real4_3d
01704 #define MPP_RECV_4D_ mpp_recv_real4_4d
01705 #define MPP_RECV_5D_ mpp_recv_real4_5d
01706 #define MPP_SEND_ mpp_send_real4
01707 #define MPP_SEND_SCALAR_ mpp_send_real4_scalar
01708 #define MPP_SEND_2D_ mpp_send_real4_2d
01709 #define MPP_SEND_3D_ mpp_send_real4_3d
01710 #define MPP_SEND_4D_ mpp_send_real4_4d
01711 #define MPP_SEND_5D_ mpp_send_real4_5d
01712 #define MPP_BROADCAST_ mpp_broadcast_real4
01713 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_real4_scalar
01714 #define MPP_BROADCAST_2D_ mpp_broadcast_real4_2d
01715 #define MPP_BROADCAST_3D_ mpp_broadcast_real4_3d
01716 #define MPP_BROADCAST_4D_ mpp_broadcast_real4_4d
01717 #define MPP_BROADCAST_5D_ mpp_broadcast_real4_5d
01718 #define MPP_TYPE_ real(FLOAT_KIND)
01719 #define MPP_TYPE_BYTELEN_ 4
01720 #ifdef use_LAM_MPI
01721 #define MPI_TYPE_ MPI_REAL
01722 #else
01723 #define MPI_TYPE_ MPI_REAL4
01724 #endif
01725 #define SHMEM_BROADCAST_ SHMEM_BROADCAST4
01726 #define SHMEM_GET_ SHMEM_GET4
01727 #include <mpp_transmit.h>
01728 #endif
01729 
01730 #define MPP_TRANSMIT_ mpp_transmit_cmplx8
01731 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_cmplx8_scalar
01732 #define MPP_TRANSMIT_2D_ mpp_transmit_cmplx8_2d
01733 #define MPP_TRANSMIT_3D_ mpp_transmit_cmplx8_3d
01734 #define MPP_TRANSMIT_4D_ mpp_transmit_cmplx8_4d
01735 #define MPP_TRANSMIT_5D_ mpp_transmit_cmplx8_5d
01736 #define MPP_RECV_ mpp_recv_cmplx8
01737 #define MPP_RECV_SCALAR_ mpp_recv_cmplx8_scalar
01738 #define MPP_RECV_2D_ mpp_recv_cmplx8_2d
01739 #define MPP_RECV_3D_ mpp_recv_cmplx8_3d
01740 #define MPP_RECV_4D_ mpp_recv_cmplx8_4d
01741 #define MPP_RECV_5D_ mpp_recv_cmplx8_5d
01742 #define MPP_SEND_ mpp_send_cmplx8
01743 #define MPP_SEND_SCALAR_ mpp_send_cmplx8_scalar
01744 #define MPP_SEND_2D_ mpp_send_cmplx8_2d
01745 #define MPP_SEND_3D_ mpp_send_cmplx8_3d
01746 #define MPP_SEND_4D_ mpp_send_cmplx8_4d
01747 #define MPP_SEND_5D_ mpp_send_cmplx8_5d
01748 #define MPP_BROADCAST_ mpp_broadcast_cmplx8
01749 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_cmplx8_scalar
01750 #define MPP_BROADCAST_2D_ mpp_broadcast_cmplx8_2d
01751 #define MPP_BROADCAST_3D_ mpp_broadcast_cmplx8_3d
01752 #define MPP_BROADCAST_4D_ mpp_broadcast_cmplx8_4d
01753 #define MPP_BROADCAST_5D_ mpp_broadcast_cmplx8_5d
01754 #define MPP_TYPE_ complex(DOUBLE_KIND)
01755 #define MPP_TYPE_BYTELEN_ 16
01756 #define MPI_TYPE_ MPI_DOUBLE_COMPLEX
01757 #define SHMEM_BROADCAST_ SHMEM_BROADCAST8
01758 #define SHMEM_GET_ SHMEM_GET128
01759 #include <mpp_transmit.h>
01760 
01761 #ifndef no_4byte_cmplx
01762 #define MPP_TRANSMIT_ mpp_transmit_cmplx4
01763 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_cmplx4_scalar
01764 #define MPP_TRANSMIT_2D_ mpp_transmit_cmplx4_2d
01765 #define MPP_TRANSMIT_3D_ mpp_transmit_cmplx4_3d
01766 #define MPP_TRANSMIT_4D_ mpp_transmit_cmplx4_4d
01767 #define MPP_TRANSMIT_5D_ mpp_transmit_cmplx4_5d
01768 #define MPP_RECV_ mpp_recv_cmplx4
01769 #define MPP_RECV_SCALAR_ mpp_recv_cmplx4_scalar
01770 #define MPP_RECV_2D_ mpp_recv_cmplx4_2d
01771 #define MPP_RECV_3D_ mpp_recv_cmplx4_3d
01772 #define MPP_RECV_4D_ mpp_recv_cmplx4_4d
01773 #define MPP_RECV_5D_ mpp_recv_cmplx4_5d
01774 #define MPP_SEND_ mpp_send_cmplx4
01775 #define MPP_SEND_SCALAR_ mpp_send_cmplx4_scalar
01776 #define MPP_SEND_2D_ mpp_send_cmplx4_2d
01777 #define MPP_SEND_3D_ mpp_send_cmplx4_3d
01778 #define MPP_SEND_4D_ mpp_send_cmplx4_4d
01779 #define MPP_SEND_5D_ mpp_send_cmplx4_5d
01780 #define MPP_BROADCAST_ mpp_broadcast_cmplx4
01781 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_cmplx4_scalar
01782 #define MPP_BROADCAST_2D_ mpp_broadcast_cmplx4_2d
01783 #define MPP_BROADCAST_3D_ mpp_broadcast_cmplx4_3d
01784 #define MPP_BROADCAST_4D_ mpp_broadcast_cmplx4_4d
01785 #define MPP_BROADCAST_5D_ mpp_broadcast_cmplx4_5d
01786 #define MPP_TYPE_ complex(FLOAT_KIND)
01787 #define MPP_TYPE_BYTELEN_ 8
01788 #define MPI_TYPE_ MPI_COMPLEX
01789 #define SHMEM_BROADCAST_ SHMEM_BROADCAST4
01790 #define SHMEM_GET_ SHMEM_GET64
01791 #include <mpp_transmit.h>
01792 #endif
01793 
01794 #ifndef no_8byte_integers
01795 #define MPP_TRANSMIT_ mpp_transmit_int8
01796 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_int8_scalar
01797 #define MPP_TRANSMIT_2D_ mpp_transmit_int8_2d
01798 #define MPP_TRANSMIT_3D_ mpp_transmit_int8_3d
01799 #define MPP_TRANSMIT_4D_ mpp_transmit_int8_4d
01800 #define MPP_TRANSMIT_5D_ mpp_transmit_int8_5d
01801 #define MPP_RECV_ mpp_recv_int8
01802 #define MPP_RECV_SCALAR_ mpp_recv_int8_scalar
01803 #define MPP_RECV_2D_ mpp_recv_int8_2d
01804 #define MPP_RECV_3D_ mpp_recv_int8_3d
01805 #define MPP_RECV_4D_ mpp_recv_int8_4d
01806 #define MPP_RECV_5D_ mpp_recv_int8_5d
01807 #define MPP_SEND_ mpp_send_int8
01808 #define MPP_SEND_SCALAR_ mpp_send_int8_scalar
01809 #define MPP_SEND_2D_ mpp_send_int8_2d
01810 #define MPP_SEND_3D_ mpp_send_int8_3d
01811 #define MPP_SEND_4D_ mpp_send_int8_4d
01812 #define MPP_SEND_5D_ mpp_send_int8_5d
01813 #define MPP_BROADCAST_ mpp_broadcast_int8
01814 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_int8_scalar
01815 #define MPP_BROADCAST_2D_ mpp_broadcast_int8_2d
01816 #define MPP_BROADCAST_3D_ mpp_broadcast_int8_3d
01817 #define MPP_BROADCAST_4D_ mpp_broadcast_int8_4d
01818 #define MPP_BROADCAST_5D_ mpp_broadcast_int8_5d
01819 #define MPP_TYPE_ integer(LONG_KIND)
01820 #define MPP_TYPE_BYTELEN_ 8
01821 #ifdef use_LAM_MPI
01822 #define MPI_TYPE_ MPI_INTEGER
01823 #else
01824 #define MPI_TYPE_ MPI_INTEGER8
01825 #endif
01826 #define SHMEM_BROADCAST_ SHMEM_BROADCAST8
01827 #define SHMEM_GET_ SHMEM_GET8
01828 #include <mpp_transmit.h>
01829 #endif
01830 
01831 #define MPP_TRANSMIT_ mpp_transmit_int4
01832 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_int4_scalar
01833 #define MPP_TRANSMIT_2D_ mpp_transmit_int4_2d
01834 #define MPP_TRANSMIT_3D_ mpp_transmit_int4_3d
01835 #define MPP_TRANSMIT_4D_ mpp_transmit_int4_4d
01836 #define MPP_TRANSMIT_5D_ mpp_transmit_int4_5d
01837 #define MPP_RECV_ mpp_recv_int4
01838 #define MPP_RECV_SCALAR_ mpp_recv_int4_scalar
01839 #define MPP_RECV_2D_ mpp_recv_int4_2d
01840 #define MPP_RECV_3D_ mpp_recv_int4_3d
01841 #define MPP_RECV_4D_ mpp_recv_int4_4d
01842 #define MPP_RECV_5D_ mpp_recv_int4_5d
01843 #define MPP_SEND_ mpp_send_int4
01844 #define MPP_SEND_SCALAR_ mpp_send_int4_scalar
01845 #define MPP_SEND_2D_ mpp_send_int4_2d
01846 #define MPP_SEND_3D_ mpp_send_int4_3d
01847 #define MPP_SEND_4D_ mpp_send_int4_4d
01848 #define MPP_SEND_5D_ mpp_send_int4_5d
01849 #define MPP_BROADCAST_ mpp_broadcast_int4
01850 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_int4_scalar
01851 #define MPP_BROADCAST_2D_ mpp_broadcast_int4_2d
01852 #define MPP_BROADCAST_3D_ mpp_broadcast_int4_3d
01853 #define MPP_BROADCAST_4D_ mpp_broadcast_int4_4d
01854 #define MPP_BROADCAST_5D_ mpp_broadcast_int4_5d
01855 #define MPP_TYPE_ integer(INT_KIND)
01856 #define MPP_TYPE_BYTELEN_ 4
01857 #ifdef use_LAM_MPI
01858 #define MPI_TYPE_ MPI_INTEGER
01859 #else
01860 #define MPI_TYPE_ MPI_INTEGER4
01861 #endif
01862 #define SHMEM_BROADCAST_ SHMEM_BROADCAST4
01863 #define SHMEM_GET_ SHMEM_GET4
01864 #include <mpp_transmit.h>
01865 
01866 #ifndef no_8byte_integers
01867 #define MPP_TRANSMIT_ mpp_transmit_logical8
01868 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_logical8_scalar
01869 #define MPP_TRANSMIT_2D_ mpp_transmit_logical8_2d
01870 #define MPP_TRANSMIT_3D_ mpp_transmit_logical8_3d
01871 #define MPP_TRANSMIT_4D_ mpp_transmit_logical8_4d
01872 #define MPP_TRANSMIT_5D_ mpp_transmit_logical8_5d
01873 #define MPP_RECV_ mpp_recv_logical8
01874 #define MPP_RECV_SCALAR_ mpp_recv_logical8_scalar
01875 #define MPP_RECV_2D_ mpp_recv_logical8_2d
01876 #define MPP_RECV_3D_ mpp_recv_logical8_3d
01877 #define MPP_RECV_4D_ mpp_recv_logical8_4d
01878 #define MPP_RECV_5D_ mpp_recv_logical8_5d
01879 #define MPP_SEND_ mpp_send_logical8
01880 #define MPP_SEND_SCALAR_ mpp_send_logical8_scalar
01881 #define MPP_SEND_2D_ mpp_send_logical8_2d
01882 #define MPP_SEND_3D_ mpp_send_logical8_3d
01883 #define MPP_SEND_4D_ mpp_send_logical8_4d
01884 #define MPP_SEND_5D_ mpp_send_logical8_5d
01885 #define MPP_BROADCAST_ mpp_broadcast_logical8
01886 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_logical8_scalar
01887 #define MPP_BROADCAST_2D_ mpp_broadcast_logical8_2d
01888 #define MPP_BROADCAST_3D_ mpp_broadcast_logical8_3d
01889 #define MPP_BROADCAST_4D_ mpp_broadcast_logical8_4d
01890 #define MPP_BROADCAST_5D_ mpp_broadcast_logical8_5d
01891 #define MPP_TYPE_ logical(LONG_KIND)
01892 #define MPP_TYPE_BYTELEN_ 8
01893 #ifdef use_LAM_MPI
01894 #define MPI_TYPE_ MPI_INTEGER
01895 #else
01896 #define MPI_TYPE_ MPI_INTEGER8
01897 #endif
01898 #define SHMEM_BROADCAST_ SHMEM_BROADCAST8
01899 #define SHMEM_GET_ SHMEM_GET8
01900 #include <mpp_transmit.h>
01901 #endif
01902 
01903 #define MPP_TRANSMIT_ mpp_transmit_logical4
01904 #define MPP_TRANSMIT_SCALAR_ mpp_transmit_logical4_scalar
01905 #define MPP_TRANSMIT_2D_ mpp_transmit_logical4_2d
01906 #define MPP_TRANSMIT_3D_ mpp_transmit_logical4_3d
01907 #define MPP_TRANSMIT_4D_ mpp_transmit_logical4_4d
01908 #define MPP_TRANSMIT_5D_ mpp_transmit_logical4_5d
01909 #define MPP_RECV_ mpp_recv_logical4
01910 #define MPP_RECV_SCALAR_ mpp_recv_logical4_scalar
01911 #define MPP_RECV_2D_ mpp_recv_logical4_2d
01912 #define MPP_RECV_3D_ mpp_recv_logical4_3d
01913 #define MPP_RECV_4D_ mpp_recv_logical4_4d
01914 #define MPP_RECV_5D_ mpp_recv_logical4_5d
01915 #define MPP_SEND_ mpp_send_logical4
01916 #define MPP_SEND_SCALAR_ mpp_send_logical4_scalar
01917 #define MPP_SEND_2D_ mpp_send_logical4_2d
01918 #define MPP_SEND_3D_ mpp_send_logical4_3d
01919 #define MPP_SEND_4D_ mpp_send_logical4_4d
01920 #define MPP_SEND_5D_ mpp_send_logical4_5d
01921 #define MPP_BROADCAST_ mpp_broadcast_logical4
01922 #define MPP_BROADCAST_SCALAR_ mpp_broadcast_logical4_scalar
01923 #define MPP_BROADCAST_2D_ mpp_broadcast_logical4_2d
01924 #define MPP_BROADCAST_3D_ mpp_broadcast_logical4_3d
01925 #define MPP_BROADCAST_4D_ mpp_broadcast_logical4_4d
01926 #define MPP_BROADCAST_5D_ mpp_broadcast_logical4_5d
01927 #define MPP_TYPE_ logical(INT_KIND)
01928 #define MPP_TYPE_BYTELEN_ 4
01929 #ifdef use_LAM_MPI
01930 #define MPI_TYPE_ MPI_INTEGER
01931 #else
01932 #define MPI_TYPE_ MPI_INTEGER4
01933 #endif
01934 #define SHMEM_BROADCAST_ SHMEM_BROADCAST4
01935 #define SHMEM_GET_ SHMEM_GET4
01936 #include <mpp_transmit.h>
01937 
01938 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01939 !                                                                             !
01940 !            GLOBAL REDUCTION ROUTINES: mpp_max, mpp_sum, mpp_min             !
01941 !                                                                             !
01942 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
01943 
01944 #define MPP_REDUCE_ mpp_max_real8
01945 #define MPP_TYPE_ real(DOUBLE_KIND)
01946 #define SHMEM_REDUCE_ SHMEM_REAL8_MAX_TO_ALL
01947 #ifdef use_LAM_MPI
01948 #define MPI_TYPE_ MPI_DOUBLE_PRECISION
01949 #else
01950 #define MPI_TYPE_ MPI_REAL8
01951 #endif
01952 #define MPI_REDUCE_ MPI_MAX
01953 #include <mpp_reduce.h>
01954 
01955 #define MPP_REDUCE_ mpp_max_real4
01956 #define MPP_TYPE_ real(FLOAT_KIND)
01957 #define SHMEM_REDUCE_ SHMEM_REAL4_MAX_TO_ALL
01958 #ifdef use_LAM_MPI
01959 #define MPI_TYPE_ MPI_REAL
01960 #else
01961 #define MPI_TYPE_ MPI_REAL4
01962 #endif
01963 #define MPI_REDUCE_ MPI_MAX
01964 #include <mpp_reduce.h>
01965 
01966 #ifndef no_8byte_integers   
01967 #define MPP_REDUCE_ mpp_max_int8
01968 #define MPP_TYPE_ integer(LONG_KIND)
01969 #define SHMEM_REDUCE_ SHMEM_INT8_MAX_TO_ALL
01970 #ifdef use_LAM_MPI
01971 #define MPI_TYPE_ MPI_INTEGER
01972 #else
01973 #define MPI_TYPE_ MPI_INTEGER8
01974 #endif
01975 #define MPI_REDUCE_ MPI_MAX
01976 #include <mpp_reduce.h>
01977 #endif
01978 
01979 #define MPP_REDUCE_ mpp_max_int4
01980 #define MPP_TYPE_ integer(INT_KIND)
01981 #define SHMEM_REDUCE_ SHMEM_INT4_MAX_TO_ALL
01982 #ifdef use_LAM_MPI
01983 #define MPI_TYPE_ MPI_INTEGER
01984 #else
01985 #define MPI_TYPE_ MPI_INTEGER4
01986 #endif
01987 #define MPI_REDUCE_ MPI_MAX
01988 #include <mpp_reduce.h>
01989 
01990 #define MPP_REDUCE_ mpp_min_real8
01991 #define MPP_TYPE_ real(DOUBLE_KIND)
01992 #define SHMEM_REDUCE_ SHMEM_REAL8_MIN_TO_ALL
01993 #ifdef use_LAM_MPI
01994 #define MPI_TYPE_ MPI_DOUBLE_PRECISION
01995 #else
01996 #define MPI_TYPE_ MPI_REAL8
01997 #endif
01998 #define MPI_REDUCE_ MPI_MIN
01999 #include <mpp_reduce.h>
02000 
02001 #ifndef no_4byte_reals
02002 #define MPP_REDUCE_ mpp_min_real4
02003 #define MPP_TYPE_ real(FLOAT_KIND)
02004 #define SHMEM_REDUCE_ SHMEM_REAL4_MIN_TO_ALL
02005 #ifdef use_LAM_MPI
02006 #define MPI_TYPE_ MPI_REAL
02007 #else
02008 #define MPI_TYPE_ MPI_REAL4
02009 #endif
02010 #define MPI_REDUCE_ MPI_MIN
02011 #include <mpp_reduce.h>
02012 #endif
02013 
02014 #ifndef no_8byte_integers   
02015 #define MPP_REDUCE_ mpp_min_int8
02016 #define MPP_TYPE_ integer(LONG_KIND)
02017 #define SHMEM_REDUCE_ SHMEM_INT8_MIN_TO_ALL
02018 #ifdef use_LAM_MPI
02019 #define MPI_TYPE_ MPI_INTEGER
02020 #else
02021 #define MPI_TYPE_ MPI_INTEGER8
02022 #endif
02023 #define MPI_REDUCE_ MPI_MIN
02024 #include <mpp_reduce.h>
02025 #endif
02026 
02027 #define MPP_REDUCE_ mpp_min_int4
02028 #define MPP_TYPE_ integer(INT_KIND)
02029 #define SHMEM_REDUCE_ SHMEM_INT4_MIN_TO_ALL
02030 #ifdef use_LAM_MPI
02031 #define MPI_TYPE_ MPI_INTEGER
02032 #else
02033 #define MPI_TYPE_ MPI_INTEGER4
02034 #endif
02035 #define MPI_REDUCE_ MPI_MIN
02036 #include <mpp_reduce.h>
02037 
02038 #define MPP_SUM_ mpp_sum_real8
02039 #define MPP_SUM_SCALAR_ mpp_sum_real8_scalar
02040 #define MPP_SUM_2D_ mpp_sum_real8_2d
02041 #define MPP_SUM_3D_ mpp_sum_real8_3d
02042 #define MPP_SUM_4D_ mpp_sum_real8_4d
02043 #define MPP_SUM_5D_ mpp_sum_real8_5d
02044 #define MPP_TYPE_ real(DOUBLE_KIND)
02045 #define SHMEM_SUM_ SHMEM_REAL8_SUM_TO_ALL
02046 #ifdef use_LAM_MPI
02047 #define MPI_TYPE_ MPI_DOUBLE_PRECISION
02048 #else
02049 #define MPI_TYPE_ MPI_REAL8
02050 #endif
02051 #define MPP_TYPE_BYTELEN_ 8
02052 #include <mpp_sum.h>
02053 
02054 #ifndef no_4byte_reals
02055 #define MPP_SUM_ mpp_sum_real4
02056 #define MPP_SUM_SCALAR_ mpp_sum_real4_scalar
02057 #define MPP_SUM_2D_ mpp_sum_real4_2d
02058 #define MPP_SUM_3D_ mpp_sum_real4_3d
02059 #define MPP_SUM_4D_ mpp_sum_real4_4d
02060 #define MPP_SUM_5D_ mpp_sum_real4_5d
02061 #define MPP_TYPE_ real(FLOAT_KIND)
02062 #define SHMEM_SUM_ SHMEM_REAL4_SUM_TO_ALL
02063 #ifdef use_LAM_MPI
02064 #define MPI_TYPE_ MPI_REAL
02065 #else
02066 #define MPI_TYPE_ MPI_REAL4
02067 #endif
02068 #define MPP_TYPE_BYTELEN_ 4
02069 #include <mpp_sum.h>
02070 #endif
02071 
02072 #define MPP_SUM_ mpp_sum_cmplx8
02073 #define MPP_SUM_SCALAR_ mpp_sum_cmplx8_scalar
02074 #define MPP_SUM_2D_ mpp_sum_cmplx8_2d
02075 #define MPP_SUM_3D_ mpp_sum_cmplx8_3d
02076 #define MPP_SUM_4D_ mpp_sum_cmplx8_4d
02077 #define MPP_SUM_5D_ mpp_sum_cmplx8_5d
02078 #define MPP_TYPE_ complex(DOUBLE_KIND)
02079 #define SHMEM_SUM_ SHMEM_COMP8_SUM_TO_ALL
02080 #define MPI_TYPE_ MPI_DOUBLE_COMPLEX
02081 #define MPP_TYPE_BYTELEN_ 16
02082 #include <mpp_sum.h>
02083 
02084 #ifndef no_4byte_cmplx
02085 #define MPP_SUM_ mpp_sum_cmplx4
02086 #define MPP_SUM_SCALAR_ mpp_sum_cmplx4_scalar
02087 #define MPP_SUM_2D_ mpp_sum_cmplx4_2d
02088 #define MPP_SUM_3D_ mpp_sum_cmplx4_3d
02089 #define MPP_SUM_4D_ mpp_sum_cmplx4_4d
02090 #define MPP_SUM_5D_ mpp_sum_cmplx4_5d
02091 #define MPP_TYPE_ complex(FLOAT_KIND)
02092 #define SHMEM_SUM_ SHMEM_COMP4_SUM_TO_ALL
02093 #define MPI_TYPE_ MPI_COMPLEX
02094 #define MPP_TYPE_BYTELEN_ 8
02095 #include <mpp_sum.h>
02096 #endif
02097 
02098 #ifndef no_8byte_integers
02099 #define MPP_SUM_ mpp_sum_int8
02100 #define MPP_SUM_SCALAR_ mpp_sum_int8_scalar
02101 #define MPP_SUM_2D_ mpp_sum_int8_2d
02102 #define MPP_SUM_3D_ mpp_sum_int8_3d
02103 #define MPP_SUM_4D_ mpp_sum_int8_4d
02104 #define MPP_SUM_5D_ mpp_sum_int8_5d
02105 #define MPP_TYPE_ integer(LONG_KIND)
02106 #define SHMEM_SUM_ SHMEM_INT8_SUM_TO_ALL
02107 #ifdef use_LAM_MPI
02108 #define MPI_TYPE_ MPI_INTEGER
02109 #else
02110 #define MPI_TYPE_ MPI_INTEGER8
02111 #endif
02112 #define MPP_TYPE_BYTELEN_ 8
02113 #include <mpp_sum.h>
02114 #endif
02115 
02116 #define MPP_SUM_ mpp_sum_int4
02117 #define MPP_SUM_SCALAR_ mpp_sum_int4_scalar
02118 #define MPP_SUM_2D_ mpp_sum_int4_2d
02119 #define MPP_SUM_3D_ mpp_sum_int4_3d
02120 #define MPP_SUM_4D_ mpp_sum_int4_4d
02121 #define MPP_SUM_5D_ mpp_sum_int4_5d
02122 #define MPP_TYPE_ integer(INT_KIND)
02123 #define SHMEM_SUM_ SHMEM_INT4_SUM_TO_ALL
02124 #ifdef use_LAM_MPI
02125 #define MPI_TYPE_ MPI_INTEGER
02126 #else
02127 #define MPI_TYPE_ MPI_INTEGER4
02128 #endif
02129 #define MPP_TYPE_BYTELEN_ 4
02130 #include <mpp_sum.h>
02131 
02132 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
02133 !                                                                             !
02134 !           SYNCHRONIZATION ROUTINES: mpp_sync, mpp_sync_self                 !
02135 !                                                                             !
02136 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
02137 
02138     subroutine mpp_sync( pelist )
02139 !synchronize PEs in list
02140       integer, intent(in), optional :: pelist(:)
02141       integer :: n
02142 
02143       call mpp_sync_self(pelist)
02144 
02145       n = get_peset(pelist); if( peset(n)%count.EQ.1 )return
02146 
02147       if( current_clock.NE.0 )call SYSTEM_CLOCK(start_tick)
02148 #ifdef use_libSMA
02149       if( n.EQ.world_peset_num )then
02150           call SHMEM_BARRIER_ALL() !special call is faster
02151       else
02152           call SHMEM_BARRIER( peset(n)%start, peset(n)%log2stride, peset(n)%count, sync )
02153       end if
02154 #endif
02155 #ifdef use_libMPI
02156       call MPI_BARRIER( peset(n)%id, error )
02157 #endif
02158       if( current_clock.NE.0 )call increment_current_clock(EVENT_WAIT)
02159 
02160       return
02161     end subroutine mpp_sync
02162 
02163     subroutine mpp_sync_self( pelist )
02164 !this is to check if current PE's outstanding puts are complete
02165 !but we can't use shmem_fence because we are actually waiting for
02166 !a remote PE to complete its get
02167       integer, intent(in), optional :: pelist(:)
02168       integer :: i, m, n, stride
02169 
02170       n = get_peset(pelist); if( peset(n)%count.EQ.1 )return
02171 
02172       if( current_clock.NE.0 )call SYSTEM_CLOCK(start_tick)
02173 #ifdef use_libSMA
02174 #ifdef _CRAYT90
02175       call SHMEM_UDCFLUSH !invalidate data cache
02176 #endif
02177 #endif
02178       do m = 1,peset(n)%count
02179          i = peset(n)%list(m)
02180 #ifdef use_libSMA
02181          call SHMEM_INT8_WAIT( status(i), MPP_WAIT ) !wait for status.NE.MPP_WAIT
02182 #endif
02183 #ifdef use_libMPI
02184          if( mpp_request(i).NE.MPI_REQUEST_NULL )call MPI_WAIT( mpp_request(i), stat, error )
02185 #endif
02186       end do
02187       if( current_clock.NE.0 )call increment_current_clock(EVENT_WAIT)
02188       return
02189     end subroutine mpp_sync_self
02190 
02191 #ifdef use_libSMA
02192 !these local versions are written for grouping into shmem_integer_wait
02193     subroutine shmem_int4_wait_local( ivar, cmp_value )
02194 !dir$ INLINEALWAYS shmem_int4_wait_local
02195       integer(INT_KIND), intent(in) :: cmp_value
02196       integer(INT_KIND), intent(inout) :: ivar
02197       call SHMEM_INT4_WAIT( ivar, cmp_value )
02198       return
02199     end subroutine shmem_int4_wait_local
02200     subroutine shmem_int8_wait_local( ivar, cmp_value )
02201 !dir$ INLINEALWAYS shmem_int8_wait_local
02202       integer(LONG_KIND), intent(in) :: cmp_value
02203       integer(LONG_KIND), intent(inout) :: ivar
02204       call SHMEM_INT8_WAIT( ivar, cmp_value )
02205       return
02206     end subroutine shmem_int8_wait_local
02207 #endif
02208       
02209 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
02210 !                                                                             !
02211 !         MISCELLANEOUS UTILITIES: mpp_error, mpp_chksum, mpp_malloc          !
02212 !                                                                             !
02213 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
02214 
02215     subroutine mpp_error_basic( errortype, errormsg )
02216 !a very basic error handler
02217 !uses ABORT and FLUSH calls, may need to use cpp to rename
02218       integer, intent(in) :: errortype
02219       character(len=*), intent(in), optional :: errormsg
02220       character(len=128) :: text
02221       logical :: opened
02222       
02223       if( .NOT.module_is_initialized )call ABORT()
02224 
02225       select case( errortype )
02226       case(NOTE)
02227           text = 'NOTE'         !just FYI
02228       case(WARNING)
02229           text = 'WARNING'      !probable error
02230       case(FATAL)
02231           text = 'FATAL'        !fatal error
02232       case default
02233           text = 'WARNING: non-existent errortype (must be NOTE|WARNING|FATAL)'
02234       end select
02235 
02236       if( npes.GT.1 )write( text,'(a,i5)' )trim(text)//' from PE', pe !this is the mpp part
02237       if( PRESENT(errormsg) )text = trim(text)//': '//trim(errormsg)
02238 
02239       select case( errortype )
02240       case(NOTE)
02241           write( nulprt,'(a)' )trim(text)
02242       case default
02243           write( stderr(),'(/a/)' )trim(text)
02244           if( errortype.EQ.FATAL .OR. warnings_are_fatal )then
02245               call FLUSH(nulprt)
02246 #ifdef sgi_mipspro
02247 !              call TRACE_BACK_STACK_AND_PRINT()
02248 #endif
02249 #ifdef use_libMPI
02250 #ifndef sgi_mipspro
02251 !the call to MPI_ABORT is not trapped by TotalView on sgi
02252               call MPI_ABORT(  peset(0)%id, 1, error )
02253 #endif
02254 #endif
02255               call ABORT()  !automatically calls traceback on Cray systems
02256           end if
02257       end select
02258 
02259       error_state = errortype
02260       return
02261     end subroutine mpp_error_basic
02262 !overloads to mpp_error_basic
02263     subroutine mpp_error_mesg( routine, errormsg, errortype )
02264 !support for error_mesg routine in FMS
02265       character(len=*), intent(in) :: routine, errormsg
02266       integer, intent(in) :: errortype
02267       call mpp_error( errortype, trim(routine)//': '//trim(errormsg) )
02268       return
02269     end subroutine mpp_error_mesg
02270     subroutine mpp_error_noargs()
02271       call mpp_error(FATAL)
02272     end subroutine mpp_error_noargs
02273       
02274     subroutine mpp_set_warn_level(flag)
02275       integer, intent(in) :: flag
02276 
02277       if( flag.EQ.WARNING )then
02278           warnings_are_fatal = .FALSE.
02279       else if( flag.EQ.FATAL )then
02280           warnings_are_fatal = .TRUE.
02281       else
02282           call mpp_error( FATAL, 'MPP_SET_WARN_LEVEL: warning flag must be set to WARNING or FATAL.' )
02283       end if
02284       return
02285     end subroutine mpp_set_warn_level
02286 
02287     function mpp_error_state()
02288       integer :: mpp_error_state
02289       mpp_error_state = error_state
02290       return
02291     end function mpp_error_state
02292 
02293 #ifdef use_shmalloc
02294     subroutine mpp_malloc( ptr, newlen, len )
02295 !routine to perform symmetric allocation:
02296 !this is required on the t3e/O2k for variables that will be non-local arguments
02297 !to a shmem call (see man intro_shmem(3F)).
02298 !newlen is the required allocation length for the pointer ptr
02299 !   len is the current allocation (0 if unallocated)
02300       integer, intent(in) :: newlen
02301       integer, intent(inout) :: len
02302       real :: dummy
02303       integer :: words_per_long
02304       integer(LONG_KIND) :: long
02305 !argument ptr is a cray pointer, points to a dummy argument in this routine
02306       pointer( ptr, dummy )
02307 !      integer(LONG_KIND) :: error_8
02308 
02309       if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_MALLOC: You must first call mpp_init.' )
02310 !use existing allocation if it is enough
02311       if( newlen.LE.len )return
02312 
02313       call SHMEM_BARRIER_ALL()
02314 !if the pointer is already allocated, deallocate
02315 !      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)
02316       if( len.NE.0 )call SHPDEALLC( ptr, error, -1 )
02317 !allocate new length: assume that the array is KIND=8
02318       words_per_long = size(transfer(long,word))
02319       call SHPALLOC( ptr, newlen*words_per_long, error, -1 )
02320       len = newlen
02321       call SHMEM_BARRIER_ALL()
02322 
02323       if( debug )then
02324           call SYSTEM_CLOCK(tick)
02325           write( nulprt,'(a,i18,a,i5,a,2i8,i16)' )'T=', tick, ' PE=', pe, ' MPP_MALLOC: len, newlen, ptr=', len, newlen, ptr
02326       end if
02327       return
02328     end subroutine mpp_malloc
02329 #endif use_shmalloc
02330 
02331     subroutine mpp_set_stack_size(n)
02332 !set the mpp_stack variable to be at least n LONG words long
02333       integer, intent(in) :: n
02334       character(len=8) :: text
02335 #ifdef use_shmalloc
02336       call mpp_malloc( ptr_stack, n, mpp_stack_size )
02337 #else
02338       if( n.GT.mpp_stack_size .AND. allocated(mpp_stack) )deallocate(mpp_stack)
02339       if( .NOT.allocated(mpp_stack) )then
02340           allocate( mpp_stack(n) )
02341           mpp_stack_size = n
02342       end if
02343 #endif
02344       write( text,'(i8)' )n
02345       if( pe.EQ.root_pe )call mpp_error( NOTE, 'MPP_SET_STACK_SIZE: stack size set to '//text//'.' )
02346 
02347       return
02348     end subroutine mpp_set_stack_size
02349 
02350 #ifndef no_8byte_integers
02351 #define MPP_CHKSUM_INT_ mpp_chksum_i8_1d
02352 #define MPP_TYPE_ integer(LONG_KIND)
02353 #define MPP_RANK_  (:)
02354 #include <mpp_chksum_int.h>
02355 
02356 #define MPP_CHKSUM_INT_ mpp_chksum_i8_2d
02357 #define MPP_TYPE_ integer(LONG_KIND)
02358 #define MPP_RANK_  (:,:)
02359 #include <mpp_chksum_int.h>
02360 
02361 #define MPP_CHKSUM_INT_ mpp_chksum_i8_3d
02362 #define MPP_TYPE_ integer(LONG_KIND)
02363 #define MPP_RANK_  (:,:,:)
02364 #include <mpp_chksum_int.h>
02365 
02366 #define MPP_CHKSUM_INT_ mpp_chksum_i8_4d
02367 #define MPP_TYPE_ integer(LONG_KIND)
02368 #define MPP_RANK_  (:,:,:,:)
02369 #include <mpp_chksum_int.h>
02370 
02371 #define MPP_CHKSUM_INT_ mpp_chksum_i8_5d
02372 #define MPP_TYPE_ integer(LONG_KIND)
02373 #define MPP_RANK_  (:,:,:,:,:)
02374 #include <mpp_chksum_int.h>
02375 #endif
02376 
02377 #define MPP_CHKSUM_INT_ mpp_chksum_i4_1d
02378 #define MPP_TYPE_ integer(INT_KIND)
02379 #define MPP_RANK_  (:)
02380 #include <mpp_chksum_int.h>
02381 
02382 #define MPP_CHKSUM_INT_ mpp_chksum_i4_2d
02383 #define MPP_TYPE_ integer(INT_KIND)
02384 #define MPP_RANK_  (:,:)
02385 #include <mpp_chksum_int.h>
02386 
02387 #define MPP_CHKSUM_INT_ mpp_chksum_i4_3d
02388 #define MPP_TYPE_ integer(INT_KIND)
02389 #define MPP_RANK_  (:,:,:)
02390 #include <mpp_chksum_int.h>
02391 
02392 #define MPP_CHKSUM_INT_ mpp_chksum_i4_4d
02393 #define MPP_TYPE_ integer(INT_KIND)
02394 #define MPP_RANK_  (:,:,:,:)
02395 #include <mpp_chksum_int.h>
02396 
02397 #define MPP_CHKSUM_INT_ mpp_chksum_i4_5d
02398 #define MPP_TYPE_ integer(INT_KIND)
02399 #define MPP_RANK_  (:,:,:,:,:)
02400 #include <mpp_chksum_int.h>
02401 
02402 #define MPP_CHKSUM_ mpp_chksum_r8_0d
02403 #define MPP_TYPE_ real(DOUBLE_KIND)
02404 #define MPP_RANK_ !
02405 #include <mpp_chksum.h>
02406 
02407 #define MPP_CHKSUM_ mpp_chksum_r8_1d
02408 #define MPP_TYPE_ real(DOUBLE_KIND)
02409 #define MPP_RANK_ (:)
02410 #include <mpp_chksum.h>
02411 
02412 #define MPP_CHKSUM_ mpp_chksum_r8_2d
02413 #define MPP_TYPE_ real(DOUBLE_KIND)
02414 #define MPP_RANK_ (:,:)
02415 #include <mpp_chksum.h>
02416 
02417 #define MPP_CHKSUM_ mpp_chksum_r8_3d
02418 #define MPP_TYPE_ real(DOUBLE_KIND)
02419 #define MPP_RANK_ (:,:,:)
02420 #include <mpp_chksum.h>
02421 
02422 #define MPP_CHKSUM_ mpp_chksum_r8_4d
02423 #define MPP_TYPE_ real(DOUBLE_KIND)
02424 #define MPP_RANK_ (:,:,:,:)
02425 #include <mpp_chksum.h>
02426 
02427 #define MPP_CHKSUM_ mpp_chksum_r8_5d
02428 #define MPP_TYPE_ real(DOUBLE_KIND)
02429 #define MPP_RANK_ (:,:,:,:,:)
02430 #include <mpp_chksum.h>
02431 
02432 #define MPP_CHKSUM_ mpp_chksum_c8_0d
02433 #define MPP_TYPE_ complex(DOUBLE_KIND)
02434 #define MPP_RANK_ !
02435 #include <mpp_chksum.h>
02436 
02437 #define MPP_CHKSUM_ mpp_chksum_c8_1d
02438 #define MPP_TYPE_ complex(DOUBLE_KIND)
02439 #define MPP_RANK_ (:)
02440 #include <mpp_chksum.h>
02441 
02442 #define MPP_CHKSUM_ mpp_chksum_c8_2d
02443 #define MPP_TYPE_ complex(DOUBLE_KIND)
02444 #define MPP_RANK_ (:,:)
02445 #include <mpp_chksum.h>
02446 
02447 #define MPP_CHKSUM_ mpp_chksum_c8_3d
02448 #define MPP_TYPE_ complex(DOUBLE_KIND)
02449 #define MPP_RANK_ (:,:,:)
02450 #include <mpp_chksum.h>
02451 
02452 #define MPP_CHKSUM_ mpp_chksum_c8_4d
02453 #define MPP_TYPE_ complex(DOUBLE_KIND)
02454 #define MPP_RANK_ (:,:,:,:)
02455 #include <mpp_chksum.h>
02456 
02457 #define MPP_CHKSUM_ mpp_chksum_c8_5d
02458 #define MPP_TYPE_ complex(DOUBLE_KIND)
02459 #define MPP_RANK_ (:,:,:,:,:)
02460 #include <mpp_chksum.h>
02461 
02462 #ifndef no_4byte_reals
02463 !CAUTION: the r4 versions of these may produce
02464 !unpredictable results: I'm not sure what the result
02465 !of the TRANSFER() to integer(8) is from an odd number of real(4)s?
02466 !However the complex(4) will work, since it is guaranteed even.
02467 #define MPP_CHKSUM_ mpp_chksum_r4_0d
02468 #define MPP_TYPE_ real(FLOAT_KIND)
02469 #define MPP_RANK_ !
02470 #include <mpp_chksum.h>
02471 
02472 #define MPP_CHKSUM_ mpp_chksum_r4_1d
02473 #define MPP_TYPE_ real(FLOAT_KIND)
02474 #define MPP_RANK_ (:)
02475 #include <mpp_chksum.h>
02476 
02477 #define MPP_CHKSUM_ mpp_chksum_r4_2d
02478 #define MPP_TYPE_ real(FLOAT_KIND)
02479 #define MPP_RANK_ (:,:)
02480 #include <mpp_chksum.h>
02481 
02482 #define MPP_CHKSUM_ mpp_chksum_r4_3d
02483 #define MPP_TYPE_ real(FLOAT_KIND)
02484 #define MPP_RANK_ (:,:,:)
02485 #include <mpp_chksum.h>
02486 
02487 #define MPP_CHKSUM_ mpp_chksum_r4_4d
02488 #define MPP_TYPE_ real(FLOAT_KIND)
02489 #define MPP_RANK_ (:,:,:,:)
02490 #include <mpp_chksum.h>
02491 
02492 #define MPP_CHKSUM_ mpp_chksum_r4_5d
02493 #define MPP_TYPE_ real(FLOAT_KIND)
02494 #define MPP_RANK_ (:,:,:,:,:)
02495 #include <mpp_chksum.h>
02496 #endif
02497 
02498 #ifndef no_4byte_cmplx
02499 #define MPP_CHKSUM_ mpp_chksum_c4_0d
02500 #define MPP_TYPE_ complex(FLOAT_KIND)
02501 #define MPP_RANK_ !
02502 #include <mpp_chksum.h>
02503 
02504 #define MPP_CHKSUM_ mpp_chksum_c4_1d
02505 #define MPP_TYPE_ complex(FLOAT_KIND)
02506 #define MPP_RANK_ (:)
02507 #include <mpp_chksum.h>
02508 
02509 #define MPP_CHKSUM_ mpp_chksum_c4_2d
02510 #define MPP_TYPE_ complex(FLOAT_KIND)
02511 #define MPP_RANK_ (:,:)
02512 #include <mpp_chksum.h>
02513 
02514 #define MPP_CHKSUM_ mpp_chksum_c4_3d
02515 #define MPP_TYPE_ complex(FLOAT_KIND)
02516 #define MPP_RANK_ (:,:,:)
02517 #include <mpp_chksum.h>
02518 
02519 #define MPP_CHKSUM_ mpp_chksum_c4_4d
02520 #define MPP_TYPE_ complex(FLOAT_KIND)
02521 #define MPP_RANK_ (:,:,:,:)
02522 #include <mpp_chksum.h>
02523 
02524 #define MPP_CHKSUM_ mpp_chksum_c4_5d
02525 #define MPP_TYPE_ complex(FLOAT_KIND)
02526 #define MPP_RANK_ (:,:,:,:,:)
02527 #include <mpp_chksum.h>
02528 #endif
02529 
02530   end module mpp_mod_oa
02531 
02532 #endif
 All Data Structures Namespaces Files Functions Variables Defines