New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 10297 for NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/LBC/lib_mpp.F90 – NEMO

Ignore:
Timestamp:
2018-11-12T16:20:57+01:00 (5 years ago)
Author:
smasson
Message:

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 2a: add report calls of mppmin/max/sum, see #2133

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/LBC/lib_mpp.F90

    r10292 r10297  
    158158 
    159159   ! Communications summary report 
    160    CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname                       !: names of calling routines 
     160   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_lbc                   !: names of lbc_lnk calling routines 
     161   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_glb                   !: names of global comm  calling routines 
    161162   INTEGER, PUBLIC                               ::   ncom_stp = 0                 !: copy of time step # istp 
     163   INTEGER, PUBLIC                               ::   ncom_fsbc = 1                !: copy of sbc time step # nn_fsbc 
     164   INTEGER, PUBLIC                               ::   ncom_dttrc = 1               !: copy of top time step # nn_dttrc 
    162165   INTEGER, PUBLIC , DIMENSION(:,:), ALLOCATABLE ::   ncomm_sequence               !: size of communicated arrays (halos) 
    163    INTEGER, PUBLIC                               ::   n_sequence = 0               !: # of communicated arrays 
    164    LOGICAL                                       ::   l_comm_report_done = .false. !: print report only once 
    165     
     166   INTEGER, PUBLIC                               ::   n_sequence_lbc = 0           !: # of communicated arraysvia lbc 
     167   INTEGER, PUBLIC                               ::   n_sequence_glb = 0           !: # of global communications 
     168   INTEGER, PUBLIC                               ::   numcom = -1                  !: logical unit for communicaton report 
     169 
    166170   ! timing summary report 
    167171   REAL(wp), PUBLIC ::  waiting_time = 0._wp, compute_time = 0._wp, elapsed_time = 0._wp 
     
    577581   END SUBROUTINE mppscatter 
    578582 
    579    !!---------------------------------------------------------------------- 
    580    !!    ***  mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real  *** 
    581    !!    
    582    !!---------------------------------------------------------------------- 
    583    !! 
    584    SUBROUTINE mppmax_a_int( ktab, kdim, kcom ) 
    585       !!---------------------------------------------------------------------- 
    586       INTEGER , INTENT(in   )                  ::   kdim   ! size of array 
    587       INTEGER , INTENT(inout), DIMENSION(kdim) ::   ktab   ! input array 
    588       INTEGER , INTENT(in   ), OPTIONAL        ::   kcom   ! 
    589       INTEGER :: ierror, ilocalcomm   ! temporary integer 
    590       INTEGER, DIMENSION(kdim) ::   iwork 
    591       !!---------------------------------------------------------------------- 
    592       ilocalcomm = mpi_comm_oce 
    593       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    594       CALL mpi_allreduce( ktab, iwork, kdim, mpi_integer, mpi_max, ilocalcomm, ierror ) 
    595       ktab(:) = iwork(:) 
    596    END SUBROUTINE mppmax_a_int 
    597    !! 
    598    SUBROUTINE mppmax_int( ktab, kcom ) 
    599       !!---------------------------------------------------------------------- 
    600       INTEGER, INTENT(inout)           ::   ktab   ! ??? 
    601       INTEGER, INTENT(in   ), OPTIONAL ::   kcom   ! ??? 
    602       INTEGER ::   ierror, iwork, ilocalcomm   ! temporary integer 
    603       !!---------------------------------------------------------------------- 
    604       ilocalcomm = mpi_comm_oce 
    605       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    606       CALL mpi_allreduce( ktab, iwork, 1, mpi_integer, mpi_max, ilocalcomm, ierror ) 
    607       ktab = iwork 
    608    END SUBROUTINE mppmax_int 
    609    !! 
    610    SUBROUTINE mppmax_a_real( ptab, kdim, kcom ) 
    611       !!---------------------------------------------------------------------- 
    612       REAL(wp), DIMENSION(kdim), INTENT(inout) ::   ptab 
    613       INTEGER                  , INTENT(in   ) ::   kdim 
    614       INTEGER , OPTIONAL       , INTENT(in   ) ::   kcom 
    615       INTEGER :: ierror, ilocalcomm 
    616       REAL(wp), DIMENSION(kdim) ::  zwork 
    617       !!---------------------------------------------------------------------- 
    618       ilocalcomm = mpi_comm_oce 
    619       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    620       CALL mpi_allreduce( ptab, zwork, kdim, mpi_double_precision, mpi_max, ilocalcomm, ierror ) 
    621       ptab(:) = zwork(:) 
    622    END SUBROUTINE mppmax_a_real 
    623    !! 
    624    SUBROUTINE mppmax_real( ptab, kcom ) 
    625       !!---------------------------------------------------------------------- 
    626       REAL(wp), INTENT(inout)           ::   ptab   ! ??? 
    627       INTEGER , INTENT(in   ), OPTIONAL ::   kcom   ! ??? 
    628       INTEGER  ::   ierror, ilocalcomm 
    629       REAL(wp) ::   zwork 
    630       !!---------------------------------------------------------------------- 
    631       ilocalcomm = mpi_comm_oce 
    632       IF( PRESENT(kcom) )   ilocalcomm = kcom! 
    633       CALL mpi_allreduce( ptab, zwork, 1, mpi_double_precision, mpi_max, ilocalcomm, ierror ) 
    634       ptab = zwork 
    635    END SUBROUTINE mppmax_real 
    636583   !! 
    637584   SUBROUTINE mpp_ilor( ld_switch, ldlast, kcom ) 
     
    656603 
    657604   END SUBROUTINE mpp_ilor 
    658  
    659  
     605    
     606   !!---------------------------------------------------------------------- 
     607   !!    ***  mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real  *** 
     608   !!    
     609   !!---------------------------------------------------------------------- 
     610   !! 
     611#  define OPERATION_MAX 
     612#  define INTEGER_TYPE 
     613#  define DIM_0d 
     614#     define ROUTINE_ALLREDUCE           mppmax_int 
     615#     include "mpp_allreduce_generic.h90" 
     616#     undef ROUTINE_ALLREDUCE 
     617#  undef DIM_0d 
     618#  define DIM_1d 
     619#     define ROUTINE_ALLREDUCE           mppmax_a_int 
     620#     include "mpp_allreduce_generic.h90" 
     621#     undef ROUTINE_ALLREDUCE 
     622#  undef DIM_1d 
     623#  undef INTEGER_TYPE 
     624! 
     625#  define REAL_TYPE 
     626#  define DIM_0d 
     627#     define ROUTINE_ALLREDUCE           mppmax_real 
     628#     include "mpp_allreduce_generic.h90" 
     629#     undef ROUTINE_ALLREDUCE 
     630#  undef DIM_0d 
     631#  define DIM_1d 
     632#     define ROUTINE_ALLREDUCE           mppmax_a_real 
     633#     include "mpp_allreduce_generic.h90" 
     634#     undef ROUTINE_ALLREDUCE 
     635#  undef DIM_1d 
     636#  undef REAL_TYPE 
     637#  undef OPERATION_MAX 
    660638   !!---------------------------------------------------------------------- 
    661639   !!    ***  mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real  *** 
     
    663641   !!---------------------------------------------------------------------- 
    664642   !! 
    665    SUBROUTINE mppmin_a_int( ktab, kdim, kcom ) 
    666       !!---------------------------------------------------------------------- 
    667       INTEGER , INTENT( in  )                  ::   kdim   ! size of array 
    668       INTEGER , INTENT(inout), DIMENSION(kdim) ::   ktab   ! input array 
    669       INTEGER , INTENT( in  ), OPTIONAL        ::   kcom   ! input array 
    670       !! 
    671       INTEGER ::   ierror, ilocalcomm   ! temporary integer 
    672       INTEGER, DIMENSION(kdim) ::   iwork 
    673       !!---------------------------------------------------------------------- 
    674       ilocalcomm = mpi_comm_oce 
    675       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    676       CALL mpi_allreduce( ktab, iwork, kdim, mpi_integer, mpi_min, ilocalcomm, ierror ) 
    677       ktab(:) = iwork(:) 
    678    END SUBROUTINE mppmin_a_int 
    679    !! 
    680    SUBROUTINE mppmin_int( ktab, kcom ) 
    681       !!---------------------------------------------------------------------- 
    682       INTEGER, INTENT(inout) ::   ktab      ! ??? 
    683       INTEGER , INTENT( in  ), OPTIONAL        ::   kcom        ! input array 
    684       !! 
    685       INTEGER ::  ierror, iwork, ilocalcomm 
    686       !!---------------------------------------------------------------------- 
    687       ilocalcomm = mpi_comm_oce 
    688       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    689       CALL mpi_allreduce( ktab, iwork, 1, mpi_integer, mpi_min, ilocalcomm, ierror ) 
    690       ktab = iwork 
    691    END SUBROUTINE mppmin_int 
    692    !! 
    693    SUBROUTINE mppmin_a_real( ptab, kdim, kcom ) 
    694       !!---------------------------------------------------------------------- 
    695       INTEGER , INTENT(in   )                  ::   kdim 
    696       REAL(wp), INTENT(inout), DIMENSION(kdim) ::   ptab 
    697       INTEGER , INTENT(in   ), OPTIONAL        ::   kcom 
    698       INTEGER :: ierror, ilocalcomm 
    699       REAL(wp), DIMENSION(kdim) ::   zwork 
    700       !!----------------------------------------------------------------------- 
    701       ilocalcomm = mpi_comm_oce 
    702       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    703       CALL mpi_allreduce( ptab, zwork, kdim, mpi_double_precision, mpi_min, ilocalcomm, ierror ) 
    704       ptab(:) = zwork(:) 
    705    END SUBROUTINE mppmin_a_real 
    706    !! 
    707    SUBROUTINE mppmin_real( ptab, kcom ) 
    708       !!----------------------------------------------------------------------- 
    709       REAL(wp), INTENT(inout)           ::   ptab        ! 
    710       INTEGER , INTENT(in   ), OPTIONAL :: kcom 
    711       INTEGER  ::   ierror, ilocalcomm 
    712       REAL(wp) ::   zwork 
    713       !!----------------------------------------------------------------------- 
    714       ilocalcomm = mpi_comm_oce 
    715       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    716       CALL mpi_allreduce( ptab, zwork, 1, mpi_double_precision, mpi_min, ilocalcomm, ierror ) 
    717       ptab = zwork 
    718    END SUBROUTINE mppmin_real 
    719  
     643#  define OPERATION_MIN 
     644#  define INTEGER_TYPE 
     645#  define DIM_0d 
     646#     define ROUTINE_ALLREDUCE           mppmin_int 
     647#     include "mpp_allreduce_generic.h90" 
     648#     undef ROUTINE_ALLREDUCE 
     649#  undef DIM_0d 
     650#  define DIM_1d 
     651#     define ROUTINE_ALLREDUCE           mppmin_a_int 
     652#     include "mpp_allreduce_generic.h90" 
     653#     undef ROUTINE_ALLREDUCE 
     654#  undef DIM_1d 
     655#  undef INTEGER_TYPE 
     656! 
     657#  define REAL_TYPE 
     658#  define DIM_0d 
     659#     define ROUTINE_ALLREDUCE           mppmin_real 
     660#     include "mpp_allreduce_generic.h90" 
     661#     undef ROUTINE_ALLREDUCE 
     662#  undef DIM_0d 
     663#  define DIM_1d 
     664#     define ROUTINE_ALLREDUCE           mppmin_a_real 
     665#     include "mpp_allreduce_generic.h90" 
     666#     undef ROUTINE_ALLREDUCE 
     667#  undef DIM_1d 
     668#  undef REAL_TYPE 
     669#  undef OPERATION_MIN 
    720670 
    721671   !!---------------------------------------------------------------------- 
     
    725675   !!---------------------------------------------------------------------- 
    726676   !! 
    727    SUBROUTINE mppsum_a_int( ktab, kdim ) 
    728       !!---------------------------------------------------------------------- 
    729       INTEGER, INTENT(in   )                   ::   kdim   ! ??? 
    730       INTEGER, INTENT(inout), DIMENSION (kdim) ::   ktab   ! ??? 
    731       INTEGER :: ierror 
    732       INTEGER, DIMENSION (kdim) ::  iwork 
    733       !!---------------------------------------------------------------------- 
    734       CALL mpi_allreduce( ktab, iwork, kdim, mpi_integer, mpi_sum, mpi_comm_oce, ierror ) 
    735       ktab(:) = iwork(:) 
    736    END SUBROUTINE mppsum_a_int 
    737    !! 
    738    SUBROUTINE mppsum_int( ktab ) 
    739       !!---------------------------------------------------------------------- 
    740       INTEGER, INTENT(inout) ::   ktab 
    741       INTEGER :: ierror, iwork 
    742       !!---------------------------------------------------------------------- 
    743       CALL mpi_allreduce( ktab, iwork, 1, mpi_integer, mpi_sum, mpi_comm_oce, ierror ) 
    744       ktab = iwork 
    745    END SUBROUTINE mppsum_int 
    746    !! 
    747    SUBROUTINE mppsum_a_real( ptab, kdim, kcom ) 
    748       !!----------------------------------------------------------------------- 
    749       INTEGER                  , INTENT(in   ) ::   kdim   ! size of ptab 
    750       REAL(wp), DIMENSION(kdim), INTENT(inout) ::   ptab   ! input array 
    751       INTEGER , OPTIONAL       , INTENT(in   ) ::   kcom   ! specific communicator 
    752       INTEGER  ::   ierror, ilocalcomm    ! local integer 
    753       REAL(wp) ::   zwork(kdim)           ! local workspace 
    754       !!----------------------------------------------------------------------- 
    755       ilocalcomm = mpi_comm_oce 
    756       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    757       CALL mpi_allreduce( ptab, zwork, kdim, mpi_double_precision, mpi_sum, ilocalcomm, ierror ) 
    758       ptab(:) = zwork(:) 
    759    END SUBROUTINE mppsum_a_real 
    760    !! 
    761    SUBROUTINE mppsum_real( ptab, kcom ) 
    762       !!----------------------------------------------------------------------- 
    763       REAL(wp)          , INTENT(inout)           ::   ptab   ! input scalar 
    764       INTEGER , OPTIONAL, INTENT(in   ) ::   kcom 
    765       INTEGER  ::   ierror, ilocalcomm 
    766       REAL(wp) ::   zwork 
    767       !!----------------------------------------------------------------------- 
    768       ilocalcomm = mpi_comm_oce 
    769       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    770       CALL mpi_allreduce( ptab, zwork, 1, mpi_double_precision, mpi_sum, ilocalcomm, ierror ) 
    771       ptab = zwork 
    772    END SUBROUTINE mppsum_real 
    773    !! 
    774    SUBROUTINE mppsum_realdd( ytab, kcom ) 
    775       !!----------------------------------------------------------------------- 
    776       COMPLEX(wp)          , INTENT(inout) ::   ytab    ! input scalar 
    777       INTEGER    , OPTIONAL, INTENT(in   ) ::   kcom 
    778       INTEGER     ::   ierror, ilocalcomm 
    779       COMPLEX(wp) ::   zwork 
    780       !!----------------------------------------------------------------------- 
    781       ilocalcomm = mpi_comm_oce 
    782       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    783       CALL MPI_ALLREDUCE( ytab, zwork, 1, MPI_DOUBLE_COMPLEX, MPI_SUMDD, ilocalcomm, ierror ) 
    784       ytab = zwork 
    785    END SUBROUTINE mppsum_realdd 
    786    !! 
    787    SUBROUTINE mppsum_a_realdd( ytab, kdim, kcom ) 
    788       !!---------------------------------------------------------------------- 
    789       INTEGER                     , INTENT(in   ) ::   kdim   ! size of ytab 
    790       COMPLEX(wp), DIMENSION(kdim), INTENT(inout) ::   ytab   ! input array 
    791       INTEGER    , OPTIONAL       , INTENT(in   ) ::   kcom 
    792       INTEGER:: ierror, ilocalcomm    ! local integer 
    793       COMPLEX(wp), DIMENSION(kdim) :: zwork     ! temporary workspace 
    794       !!----------------------------------------------------------------------- 
    795       ilocalcomm = mpi_comm_oce 
    796       IF( PRESENT(kcom) )   ilocalcomm = kcom 
    797       CALL MPI_ALLREDUCE( ytab, zwork, kdim, MPI_DOUBLE_COMPLEX, MPI_SUMDD, ilocalcomm, ierror ) 
    798       ytab(:) = zwork(:) 
    799    END SUBROUTINE mppsum_a_realdd 
    800     
     677#  define OPERATION_SUM 
     678#  define INTEGER_TYPE 
     679#  define DIM_0d 
     680#     define ROUTINE_ALLREDUCE           mppsum_int 
     681#     include "mpp_allreduce_generic.h90" 
     682#     undef ROUTINE_ALLREDUCE 
     683#  undef DIM_0d 
     684#  define DIM_1d 
     685#     define ROUTINE_ALLREDUCE           mppsum_a_int 
     686#     include "mpp_allreduce_generic.h90" 
     687#     undef ROUTINE_ALLREDUCE 
     688#  undef DIM_1d 
     689#  undef INTEGER_TYPE 
     690! 
     691#  define REAL_TYPE 
     692#  define DIM_0d 
     693#     define ROUTINE_ALLREDUCE           mppsum_real 
     694#     include "mpp_allreduce_generic.h90" 
     695#     undef ROUTINE_ALLREDUCE 
     696#  undef DIM_0d 
     697#  define DIM_1d 
     698#     define ROUTINE_ALLREDUCE           mppsum_a_real 
     699#     include "mpp_allreduce_generic.h90" 
     700#     undef ROUTINE_ALLREDUCE 
     701#  undef DIM_1d 
     702#  undef REAL_TYPE 
     703#  undef OPERATION_SUM 
     704 
     705#  define OPERATION_SUM_DD 
     706#  define COMPLEX_TYPE 
     707#  define DIM_0d 
     708#     define ROUTINE_ALLREDUCE           mppsum_realdd 
     709#     include "mpp_allreduce_generic.h90" 
     710#     undef ROUTINE_ALLREDUCE 
     711#  undef DIM_0d 
     712#  define DIM_1d 
     713#     define ROUTINE_ALLREDUCE           mppsum_a_realdd 
     714#     include "mpp_allreduce_generic.h90" 
     715#     undef ROUTINE_ALLREDUCE 
     716#  undef DIM_1d 
     717#  undef COMPLEX_TYPE 
     718#  undef OPERATION_SUM_DD 
    801719 
    802720   SUBROUTINE mppmax_real_multiple( pt1d, kdim, kcom  ) 
     
    11171035         l_znl_root = .FALSE. 
    11181036         kwork (1) = nimpp 
    1119          CALL mpp_min ( kwork(1), kcom = ncomm_znl) 
     1037         CALL mpp_min ( 'lib_mpp', kwork(1), kcom = ncomm_znl) 
    11201038         IF ( nimpp == kwork(1)) l_znl_root = .TRUE. 
    11211039      END IF 
Note: See TracChangeset for help on using the changeset viewer.