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.
lib_mpp.F90 in NEMO/trunk/src/OCE/LBC – NEMO

source: NEMO/trunk/src/OCE/LBC/lib_mpp.F90 @ 10440

Last change on this file since 10440 was 10437, checked in by smasson, 5 years ago

trunk: improve communication_report.txt for delayed global comm

  • Property svn:keywords set to Id
File size: 89.6 KB
RevLine 
[3]1MODULE lib_mpp
[13]2   !!======================================================================
3   !!                       ***  MODULE  lib_mpp  ***
[1344]4   !! Ocean numerics:  massively parallel processing library
[13]5   !!=====================================================================
[1344]6   !! History :  OPA  !  1994  (M. Guyon, J. Escobar, M. Imbard)  Original code
7   !!            7.0  !  1997  (A.M. Treguier)  SHMEM additions
8   !!            8.0  !  1998  (M. Imbard, J. Escobar, L. Colombet ) SHMEM and MPI
9   !!                 !  1998  (J.M. Molines) Open boundary conditions
[9019]10   !!   NEMO     1.0  !  2003  (J.M. Molines, G. Madec)  F90, free form
[1344]11   !!                 !  2003  (J.M. Molines) add mpp_ini_north(_3d,_2d)
12   !!             -   !  2004  (R. Bourdalle Badie)  isend option in mpi
13   !!                 !  2004  (J.M. Molines) minloc, maxloc
14   !!             -   !  2005  (G. Madec, S. Masson)  npolj=5,6 F-point & ice cases
15   !!             -   !  2005  (R. Redler) Replacement of MPI_COMM_WORLD except for MPI_Abort
16   !!             -   !  2005  (R. Benshila, G. Madec)  add extra halo case
17   !!             -   !  2008  (R. Benshila) add mpp_ini_ice
18   !!            3.2  !  2009  (R. Benshila) SHMEM suppression, north fold in lbc_nfd
[3764]19   !!            3.2  !  2009  (O. Marti)    add mpp_ini_znl
[2715]20   !!            4.0  !  2011  (G. Madec)  move ctl_ routines from in_out_manager
[9019]21   !!            3.5  !  2012  (S.Mocavero, I. Epicoco) Add mpp_lnk_bdy_3d/2d routines to optimize the BDY comm.
22   !!            3.5  !  2013  (C. Ethe, G. Madec)  message passing arrays as local variables
[6140]23   !!            3.5  !  2013  (S.Mocavero, I.Epicoco - CMCC) north fold optimizations
[9019]24   !!            3.6  !  2015  (O. Tintó and M. Castrillo - BSC) Added '_multiple' case for 2D lbc and max
25   !!            4.0  !  2017  (G. Madec) automatique allocation of array argument (use any 3rd dimension)
26   !!             -   !  2017  (G. Madec) create generic.h90 files to generate all lbc and north fold routines
[13]27   !!----------------------------------------------------------------------
[2715]28
29   !!----------------------------------------------------------------------
[6140]30   !!   ctl_stop      : update momentum and tracer Kz from a tke scheme
31   !!   ctl_warn      : initialization, namelist read, and parameters control
32   !!   ctl_opn       : Open file and check if required file is available.
33   !!   ctl_nam       : Prints informations when an error occurs while reading a namelist
34   !!   get_unit      : give the index of an unused logical unit
[2715]35   !!----------------------------------------------------------------------
[3764]36#if   defined key_mpp_mpi
[13]37   !!----------------------------------------------------------------------
[1344]38   !!   'key_mpp_mpi'             MPI massively parallel processing library
39   !!----------------------------------------------------------------------
[2715]40   !!   lib_mpp_alloc : allocate mpp arrays
41   !!   mynode        : indentify the processor unit
42   !!   mpp_lnk       : interface (defined in lbclnk) for message passing of 2d or 3d arrays (mpp_lnk_2d, mpp_lnk_3d)
[4990]43   !!   mpp_lnk_icb   : interface for message passing of 2d arrays with extra halo for icebergs (mpp_lnk_2d_icb)
[6140]44   !!   mpprecv       :
[9019]45   !!   mppsend       :
[2715]46   !!   mppscatter    :
47   !!   mppgather     :
48   !!   mpp_min       : generic interface for mppmin_int , mppmin_a_int , mppmin_real, mppmin_a_real
49   !!   mpp_max       : generic interface for mppmax_int , mppmax_a_int , mppmax_real, mppmax_a_real
50   !!   mpp_sum       : generic interface for mppsum_int , mppsum_a_int , mppsum_real, mppsum_a_real
51   !!   mpp_minloc    :
52   !!   mpp_maxloc    :
53   !!   mppsync       :
54   !!   mppstop       :
[1344]55   !!   mpp_ini_north : initialisation of north fold
[9019]56   !!   mpp_lbc_north_icb : alternative to mpp_nfd for extra outer halo with icebergs
[13]57   !!----------------------------------------------------------------------
[3764]58   USE dom_oce        ! ocean space and time domain
[2715]59   USE lbcnfd         ! north fold treatment
60   USE in_out_manager ! I/O manager
[3]61
[13]62   IMPLICIT NONE
[415]63   PRIVATE
[9019]64
65   INTERFACE mpp_nfd
[10425]66      MODULE PROCEDURE   mpp_nfd_2d    , mpp_nfd_3d    , mpp_nfd_4d
[9019]67      MODULE PROCEDURE   mpp_nfd_2d_ptr, mpp_nfd_3d_ptr, mpp_nfd_4d_ptr
68   END INTERFACE
69
70   ! Interface associated to the mpp_lnk_... routines is defined in lbclnk
[10425]71   PUBLIC   mpp_lnk_2d    , mpp_lnk_3d    , mpp_lnk_4d
[9019]72   PUBLIC   mpp_lnk_2d_ptr, mpp_lnk_3d_ptr, mpp_lnk_4d_ptr
73   !
74!!gm  this should be useless
75   PUBLIC   mpp_nfd_2d    , mpp_nfd_3d    , mpp_nfd_4d
76   PUBLIC   mpp_nfd_2d_ptr, mpp_nfd_3d_ptr, mpp_nfd_4d_ptr
77!!gm end
78   !
[4147]79   PUBLIC   ctl_stop, ctl_warn, get_unit, ctl_opn, ctl_nam
[1344]80   PUBLIC   mynode, mppstop, mppsync, mpp_comm_free
[9019]81   PUBLIC   mpp_ini_north
82   PUBLIC   mpp_lnk_2d_icb
83   PUBLIC   mpp_lbc_north_icb
[1344]84   PUBLIC   mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc
[10425]85   PUBLIC   mpp_delay_max, mpp_delay_sum, mpp_delay_rcv
[3294]86   PUBLIC   mppscatter, mppgather
[10425]87   PUBLIC   mpp_ini_znl
[3764]88   PUBLIC   mppsend, mpprecv                          ! needed by TAM and ICB routines
[9890]89   PUBLIC   mpp_lnk_bdy_2d, mpp_lnk_bdy_3d, mpp_lnk_bdy_4d
[5429]90   
[13]91   !! * Interfaces
92   !! define generic interface for these routine as they are called sometimes
[1344]93   !! with scalar arguments instead of array arguments, which causes problems
94   !! for the compilation on AIX system as well as NEC and SGI. Ok on COMPACQ
[13]95   INTERFACE mpp_min
96      MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real
97   END INTERFACE
98   INTERFACE mpp_max
[681]99      MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real
[13]100   END INTERFACE
101   INTERFACE mpp_sum
[6140]102      MODULE PROCEDURE mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real,   &
[9019]103         &             mppsum_realdd, mppsum_a_realdd
[13]104   END INTERFACE
[1344]105   INTERFACE mpp_minloc
106      MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d
107   END INTERFACE
108   INTERFACE mpp_maxloc
109      MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d
110   END INTERFACE
[6490]111
[51]112   !! ========================= !!
113   !!  MPI  variable definition !!
114   !! ========================= !!
[1629]115!$AGRIF_DO_NOT_TREAT
[2004]116   INCLUDE 'mpif.h'
[1629]117!$AGRIF_END_DO_NOT_TREAT
[3764]118
[1344]119   LOGICAL, PUBLIC, PARAMETER ::   lk_mpp = .TRUE.    !: mpp flag
[3]120
[1344]121   INTEGER, PARAMETER         ::   nprocmax = 2**10   ! maximun dimension (required to be a power of 2)
[3764]122
[10425]123   INTEGER, PUBLIC ::   mppsize        ! number of process
124   INTEGER, PUBLIC ::   mpprank        ! process number  [ 0 - size-1 ]
[2363]125!$AGRIF_DO_NOT_TREAT
[9570]126   INTEGER, PUBLIC ::   mpi_comm_oce   ! opa local communicator
[2363]127!$AGRIF_END_DO_NOT_TREAT
[3]128
[2480]129   INTEGER :: MPI_SUMDD
[1976]130
[1345]131   ! variables used for zonal integration
132   INTEGER, PUBLIC ::   ncomm_znl       !: communicator made by the processors on the same zonal average
[9019]133   LOGICAL, PUBLIC ::   l_znl_root      !: True on the 'left'most processor on the same row
134   INTEGER         ::   ngrp_znl        !  group ID for the znl processors
135   INTEGER         ::   ndim_rank_znl   !  number of processors on the same zonal average
[2715]136   INTEGER, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_znl  ! dimension ndim_rank_znl, number of the procs into the same znl domain
[3]137
[3764]138   ! North fold condition in mpp_mpi with jpni > 1 (PUBLIC for TAM)
[9019]139   INTEGER, PUBLIC ::   ngrp_world        !: group ID for the world processors
140   INTEGER, PUBLIC ::   ngrp_opa          !: group ID for the opa processors
141   INTEGER, PUBLIC ::   ngrp_north        !: group ID for the northern processors (to be fold)
142   INTEGER, PUBLIC ::   ncomm_north       !: communicator made by the processors belonging to ngrp_north
143   INTEGER, PUBLIC ::   ndim_rank_north   !: number of 'sea' processor in the northern line (can be /= jpni !)
144   INTEGER, PUBLIC ::   njmppmax          !: value of njmpp for the processors of the northern line
145   INTEGER, PUBLIC ::   north_root        !: number (in the comm_opa) of proc 0 in the northern comm
146   INTEGER, PUBLIC, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_north   !: dimension ndim_rank_north
[3764]147
[1344]148   ! Type of send : standard, buffered, immediate
[9019]149   CHARACTER(len=1), PUBLIC ::   cn_mpi_send        !: type od mpi send/recieve (S=standard, B=bsend, I=isend)
150   LOGICAL         , PUBLIC ::   l_isend = .FALSE.  !: isend use indicator (T if cn_mpi_send='I')
151   INTEGER         , PUBLIC ::   nn_buffer          !: size of the buffer in case of mpi_bsend
[3764]152
[10425]153   ! Communications summary report
154   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_lbc                   !: names of lbc_lnk calling routines
[10437]155   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_glb                   !: names of global comm calling routines
156   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_dlg                   !: names of delayed global comm calling routines
[10425]157   INTEGER, PUBLIC                               ::   ncom_stp = 0                 !: copy of time step # istp
158   INTEGER, PUBLIC                               ::   ncom_fsbc = 1                !: copy of sbc time step # nn_fsbc
159   INTEGER, PUBLIC                               ::   ncom_dttrc = 1               !: copy of top time step # nn_dttrc
160   INTEGER, PUBLIC                               ::   ncom_freq                    !: frequency of comm diagnostic
161   INTEGER, PUBLIC , DIMENSION(:,:), ALLOCATABLE ::   ncomm_sequence               !: size of communicated arrays (halos)
162   INTEGER, PARAMETER, PUBLIC                    ::   ncom_rec_max = 2000          !: max number of communication record
163   INTEGER, PUBLIC                               ::   n_sequence_lbc = 0           !: # of communicated arraysvia lbc
164   INTEGER, PUBLIC                               ::   n_sequence_glb = 0           !: # of global communications
[10437]165   INTEGER, PUBLIC                               ::   n_sequence_dlg = 0           !: # of delayed global communications
[10425]166   INTEGER, PUBLIC                               ::   numcom = -1                  !: logical unit for communicaton report
167   LOGICAL, PUBLIC                               ::   l_full_nf_update = .TRUE.    !: logical for a full (2lines) update of bc at North fold report
168   INTEGER,                    PARAMETER, PUBLIC ::   nbdelay = 2       !: number of delayed operations
169   !: name (used as id) of allreduce-delayed operations
170   CHARACTER(len=32), DIMENSION(nbdelay), PUBLIC ::   c_delaylist = (/ 'cflice', 'fwb' /)
171   !: component name where the allreduce-delayed operation is performed
172   CHARACTER(len=3),  DIMENSION(nbdelay), PUBLIC ::   c_delaycpnt = (/ 'ICE'   , 'OCE' /)
173   TYPE, PUBLIC ::   DELAYARR
174      REAL(   wp), POINTER, DIMENSION(:) ::  z1d => NULL()
175      COMPLEX(wp), POINTER, DIMENSION(:) ::  y1d => NULL()
176   END TYPE DELAYARR
177   TYPE( DELAYARR ), DIMENSION(nbdelay), PUBLIC  ::   todelay             
178   INTEGER,          DIMENSION(nbdelay), PUBLIC  ::   ndelayid = -1     !: mpi request id of the delayed operations
179
180   ! timing summary report
181   REAL(wp), DIMENSION(2), PUBLIC ::  waiting_time = 0._wp
182   REAL(wp)              , PUBLIC ::  compute_time = 0._wp, elapsed_time = 0._wp
183   
[9019]184   REAL(wp), DIMENSION(:), ALLOCATABLE, SAVE ::   tampon   ! buffer in case of bsend
[3]185
[9019]186   LOGICAL, PUBLIC ::   ln_nnogather                !: namelist control of northfold comms
187   LOGICAL, PUBLIC ::   l_north_nogather = .FALSE.  !: internal control of northfold comms
188
[51]189   !!----------------------------------------------------------------------
[9598]190   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[888]191   !! $Id$
[10068]192   !! Software governed by the CeCILL license (see ./LICENSE)
[1344]193   !!----------------------------------------------------------------------
[3]194CONTAINS
195
[9019]196   FUNCTION mynode( ldtxt, ldname, kumnam_ref, kumnam_cfg, kumond, kstop, localComm )
[2715]197      !!----------------------------------------------------------------------
[51]198      !!                  ***  routine mynode  ***
[3764]199      !!
[51]200      !! ** Purpose :   Find processor unit
201      !!----------------------------------------------------------------------
[6140]202      CHARACTER(len=*),DIMENSION(:), INTENT(  out) ::   ldtxt        !
203      CHARACTER(len=*)             , INTENT(in   ) ::   ldname       !
204      INTEGER                      , INTENT(in   ) ::   kumnam_ref   ! logical unit for reference namelist
205      INTEGER                      , INTENT(in   ) ::   kumnam_cfg   ! logical unit for configuration namelist
206      INTEGER                      , INTENT(inout) ::   kumond       ! logical unit for namelist output
207      INTEGER                      , INTENT(inout) ::   kstop        ! stop indicator
208      INTEGER         , OPTIONAL   , INTENT(in   ) ::   localComm    !
[2715]209      !
[4147]210      INTEGER ::   mynode, ierr, code, ji, ii, ios
[532]211      LOGICAL ::   mpi_was_called
[2715]212      !
[10428]213      NAMELIST/nammpp/ cn_mpi_send, nn_buffer, jpni, jpnj, ln_nnogather
[51]214      !!----------------------------------------------------------------------
[1344]215      !
[2481]216      ii = 1
[6140]217      WRITE(ldtxt(ii),*)                                                                  ;   ii = ii + 1
218      WRITE(ldtxt(ii),*) 'mynode : mpi initialisation'                                    ;   ii = ii + 1
219      WRITE(ldtxt(ii),*) '~~~~~~ '                                                        ;   ii = ii + 1
[1344]220      !
[4147]221      REWIND( kumnam_ref )              ! Namelist nammpp in reference namelist: mpi variables
222      READ  ( kumnam_ref, nammpp, IOSTAT = ios, ERR = 901)
[9168]223901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nammpp in reference namelist', lwp )
[9019]224      !
[4147]225      REWIND( kumnam_cfg )              ! Namelist nammpp in configuration namelist: mpi variables
226      READ  ( kumnam_cfg, nammpp, IOSTAT = ios, ERR = 902 )
[9168]227902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nammpp in configuration namelist', lwp )
[9019]228      !
[1344]229      !                              ! control print
[6140]230      WRITE(ldtxt(ii),*) '   Namelist nammpp'                                             ;   ii = ii + 1
231      WRITE(ldtxt(ii),*) '      mpi send type          cn_mpi_send = ', cn_mpi_send       ;   ii = ii + 1
232      WRITE(ldtxt(ii),*) '      size exported buffer   nn_buffer   = ', nn_buffer,' bytes';   ii = ii + 1
[9019]233      !
234      IF( jpni < 1 .OR. jpnj < 1  ) THEN
[10428]235         WRITE(ldtxt(ii),*) '      jpni and jpnj will be calculated automatically' ;   ii = ii + 1
[2715]236      ELSE
[6140]237         WRITE(ldtxt(ii),*) '      processor grid extent in i         jpni = ',jpni       ;   ii = ii + 1
238         WRITE(ldtxt(ii),*) '      processor grid extent in j         jpnj = ',jpnj       ;   ii = ii + 1
[9019]239      ENDIF
[2715]240
[3294]241      WRITE(ldtxt(ii),*) '      avoid use of mpi_allgather at the north fold  ln_nnogather = ', ln_nnogather  ; ii = ii + 1
242
[2480]243      CALL mpi_initialized ( mpi_was_called, code )
244      IF( code /= MPI_SUCCESS ) THEN
[3764]245         DO ji = 1, SIZE(ldtxt)
[2481]246            IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
[3764]247         END DO
[2480]248         WRITE(*, cform_err)
249         WRITE(*, *) 'lib_mpp: Error in routine mpi_initialized'
250         CALL mpi_abort( mpi_comm_world, code, ierr )
251      ENDIF
[415]252
[2480]253      IF( mpi_was_called ) THEN
254         !
255         SELECT CASE ( cn_mpi_send )
256         CASE ( 'S' )                ! Standard mpi send (blocking)
[6140]257            WRITE(ldtxt(ii),*) '           Standard blocking mpi send (send)'             ;   ii = ii + 1
[2480]258         CASE ( 'B' )                ! Buffer mpi send (blocking)
[6140]259            WRITE(ldtxt(ii),*) '           Buffer blocking mpi send (bsend)'              ;   ii = ii + 1
[9570]260            IF( Agrif_Root() )   CALL mpi_init_oce( ldtxt, ii, ierr )
[2480]261         CASE ( 'I' )                ! Immediate mpi send (non-blocking send)
[6140]262            WRITE(ldtxt(ii),*) '           Immediate non-blocking send (isend)'           ;   ii = ii + 1
[2480]263            l_isend = .TRUE.
264         CASE DEFAULT
[6140]265            WRITE(ldtxt(ii),cform_err)                                                    ;   ii = ii + 1
266            WRITE(ldtxt(ii),*) '           bad value for cn_mpi_send = ', cn_mpi_send     ;   ii = ii + 1
[2715]267            kstop = kstop + 1
[2480]268         END SELECT
[9019]269         !
270      ELSEIF ( PRESENT(localComm) .AND. .NOT. mpi_was_called ) THEN
[10425]271         WRITE(ldtxt(ii),cform_err)                                                    ;   ii = ii + 1
[6140]272         WRITE(ldtxt(ii),*) ' lib_mpp: You cannot provide a local communicator '          ;   ii = ii + 1
273         WRITE(ldtxt(ii),*) '          without calling MPI_Init before ! '                ;   ii = ii + 1
[2715]274         kstop = kstop + 1
[532]275      ELSE
[1601]276         SELECT CASE ( cn_mpi_send )
[524]277         CASE ( 'S' )                ! Standard mpi send (blocking)
[6140]278            WRITE(ldtxt(ii),*) '           Standard blocking mpi send (send)'             ;   ii = ii + 1
[2480]279            CALL mpi_init( ierr )
[524]280         CASE ( 'B' )                ! Buffer mpi send (blocking)
[6140]281            WRITE(ldtxt(ii),*) '           Buffer blocking mpi send (bsend)'              ;   ii = ii + 1
[9570]282            IF( Agrif_Root() )   CALL mpi_init_oce( ldtxt, ii, ierr )
[524]283         CASE ( 'I' )                ! Immediate mpi send (non-blocking send)
[6140]284            WRITE(ldtxt(ii),*) '           Immediate non-blocking send (isend)'           ;   ii = ii + 1
[524]285            l_isend = .TRUE.
[2480]286            CALL mpi_init( ierr )
[524]287         CASE DEFAULT
[6140]288            WRITE(ldtxt(ii),cform_err)                                                    ;   ii = ii + 1
289            WRITE(ldtxt(ii),*) '           bad value for cn_mpi_send = ', cn_mpi_send     ;   ii = ii + 1
[2715]290            kstop = kstop + 1
[524]291         END SELECT
[2480]292         !
[415]293      ENDIF
[570]294
[3764]295      IF( PRESENT(localComm) ) THEN
[2480]296         IF( Agrif_Root() ) THEN
[9570]297            mpi_comm_oce = localComm
[2480]298         ENDIF
299      ELSE
[9570]300         CALL mpi_comm_dup( mpi_comm_world, mpi_comm_oce, code)
[2480]301         IF( code /= MPI_SUCCESS ) THEN
[3764]302            DO ji = 1, SIZE(ldtxt)
[2481]303               IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
304            END DO
[2480]305            WRITE(*, cform_err)
306            WRITE(*, *) ' lib_mpp: Error in routine mpi_comm_dup'
307            CALL mpi_abort( mpi_comm_world, code, ierr )
308         ENDIF
[3764]309      ENDIF
[2480]310
[5656]311#if defined key_agrif
[9019]312      IF( Agrif_Root() ) THEN
[9570]313         CALL Agrif_MPI_Init(mpi_comm_oce)
[5656]314      ELSE
[9570]315         CALL Agrif_MPI_set_grid_comm(mpi_comm_oce)
[5656]316      ENDIF
317#endif
318
[9570]319      CALL mpi_comm_rank( mpi_comm_oce, mpprank, ierr )
320      CALL mpi_comm_size( mpi_comm_oce, mppsize, ierr )
[629]321      mynode = mpprank
[4624]322
323      IF( mynode == 0 ) THEN
[5407]324         CALL ctl_opn( kumond, TRIM(ldname), 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. , 1 )
325         WRITE(kumond, nammpp)     
[4624]326      ENDIF
[3764]327      !
[1976]328      CALL MPI_OP_CREATE(DDPDD_MPI, .TRUE., MPI_SUMDD, ierr)
329      !
[51]330   END FUNCTION mynode
[3]331
[9019]332   !!----------------------------------------------------------------------
333   !!                   ***  routine mpp_lnk_(2,3,4)d  ***
334   !!
335   !!   * Argument : dummy argument use in mpp_lnk_... routines
336   !!                ptab   :   array or pointer of arrays on which the boundary condition is applied
337   !!                cd_nat :   nature of array grid-points
338   !!                psgn   :   sign used across the north fold boundary
339   !!                kfld   :   optional, number of pt3d arrays
340   !!                cd_mpp :   optional, fill the overlap area only
341   !!                pval   :   optional, background value (used at closed boundaries)
342   !!----------------------------------------------------------------------
343   !
344   !                       !==  2D array and array of 2D pointer  ==!
345   !
346#  define DIM_2d
347#     define ROUTINE_LNK           mpp_lnk_2d
348#     include "mpp_lnk_generic.h90"
349#     undef ROUTINE_LNK
350#     define MULTI
351#     define ROUTINE_LNK           mpp_lnk_2d_ptr
352#     include "mpp_lnk_generic.h90"
353#     undef ROUTINE_LNK
354#     undef MULTI
355#  undef DIM_2d
356   !
357   !                       !==  3D array and array of 3D pointer  ==!
358   !
359#  define DIM_3d
360#     define ROUTINE_LNK           mpp_lnk_3d
361#     include "mpp_lnk_generic.h90"
362#     undef ROUTINE_LNK
363#     define MULTI
364#     define ROUTINE_LNK           mpp_lnk_3d_ptr
365#     include "mpp_lnk_generic.h90"
366#     undef ROUTINE_LNK
367#     undef MULTI
368#  undef DIM_3d
369   !
370   !                       !==  4D array and array of 4D pointer  ==!
371   !
372#  define DIM_4d
373#     define ROUTINE_LNK           mpp_lnk_4d
374#     include "mpp_lnk_generic.h90"
375#     undef ROUTINE_LNK
376#     define MULTI
377#     define ROUTINE_LNK           mpp_lnk_4d_ptr
378#     include "mpp_lnk_generic.h90"
379#     undef ROUTINE_LNK
380#     undef MULTI
381#  undef DIM_4d
[6140]382
[9019]383   !!----------------------------------------------------------------------
384   !!                   ***  routine mpp_nfd_(2,3,4)d  ***
385   !!
386   !!   * Argument : dummy argument use in mpp_nfd_... routines
387   !!                ptab   :   array or pointer of arrays on which the boundary condition is applied
388   !!                cd_nat :   nature of array grid-points
389   !!                psgn   :   sign used across the north fold boundary
390   !!                kfld   :   optional, number of pt3d arrays
391   !!                cd_mpp :   optional, fill the overlap area only
392   !!                pval   :   optional, background value (used at closed boundaries)
393   !!----------------------------------------------------------------------
394   !
395   !                       !==  2D array and array of 2D pointer  ==!
396   !
397#  define DIM_2d
398#     define ROUTINE_NFD           mpp_nfd_2d
399#     include "mpp_nfd_generic.h90"
400#     undef ROUTINE_NFD
401#     define MULTI
402#     define ROUTINE_NFD           mpp_nfd_2d_ptr
403#     include "mpp_nfd_generic.h90"
404#     undef ROUTINE_NFD
405#     undef MULTI
406#  undef DIM_2d
407   !
408   !                       !==  3D array and array of 3D pointer  ==!
409   !
410#  define DIM_3d
411#     define ROUTINE_NFD           mpp_nfd_3d
412#     include "mpp_nfd_generic.h90"
413#     undef ROUTINE_NFD
414#     define MULTI
415#     define ROUTINE_NFD           mpp_nfd_3d_ptr
416#     include "mpp_nfd_generic.h90"
417#     undef ROUTINE_NFD
418#     undef MULTI
419#  undef DIM_3d
420   !
421   !                       !==  4D array and array of 4D pointer  ==!
422   !
423#  define DIM_4d
424#     define ROUTINE_NFD           mpp_nfd_4d
425#     include "mpp_nfd_generic.h90"
426#     undef ROUTINE_NFD
427#     define MULTI
428#     define ROUTINE_NFD           mpp_nfd_4d_ptr
429#     include "mpp_nfd_generic.h90"
430#     undef ROUTINE_NFD
431#     undef MULTI
432#  undef DIM_4d
[3]433
[1344]434
[9019]435   !!----------------------------------------------------------------------
436   !!                   ***  routine mpp_lnk_bdy_(2,3,4)d  ***
437   !!
438   !!   * Argument : dummy argument use in mpp_lnk_... routines
439   !!                ptab   :   array or pointer of arrays on which the boundary condition is applied
440   !!                cd_nat :   nature of array grid-points
441   !!                psgn   :   sign used across the north fold boundary
442   !!                kb_bdy :   BDY boundary set
443   !!                kfld   :   optional, number of pt3d arrays
444   !!----------------------------------------------------------------------
445   !
446   !                       !==  2D array and array of 2D pointer  ==!
447   !
448#  define DIM_2d
449#     define ROUTINE_BDY           mpp_lnk_bdy_2d
450#     include "mpp_bdy_generic.h90"
451#     undef ROUTINE_BDY
452#  undef DIM_2d
453   !
454   !                       !==  3D array and array of 3D pointer  ==!
455   !
456#  define DIM_3d
457#     define ROUTINE_BDY           mpp_lnk_bdy_3d
458#     include "mpp_bdy_generic.h90"
459#     undef ROUTINE_BDY
460#  undef DIM_3d
461   !
462   !                       !==  4D array and array of 4D pointer  ==!
463   !
[9890]464#  define DIM_4d
465#     define ROUTINE_BDY           mpp_lnk_bdy_4d
466#     include "mpp_bdy_generic.h90"
467#     undef ROUTINE_BDY
468#  undef DIM_4d
[3]469
[9019]470   !!----------------------------------------------------------------------
471   !!
472   !!   load_array  &   mpp_lnk_2d_9    à generaliser a 3D et 4D
[5429]473   
474   
[9019]475   !!    mpp_lnk_sum_2d et 3D   ====>>>>>>   à virer du code !!!!
[5429]476   
[9019]477   
478   !!----------------------------------------------------------------------
[5429]479
480
[888]481
[1344]482   SUBROUTINE mppsend( ktyp, pmess, kbytes, kdest, md_req )
[51]483      !!----------------------------------------------------------------------
484      !!                  ***  routine mppsend  ***
[3764]485      !!
[51]486      !! ** Purpose :   Send messag passing array
487      !!
488      !!----------------------------------------------------------------------
[1344]489      REAL(wp), INTENT(inout) ::   pmess(*)   ! array of real
490      INTEGER , INTENT(in   ) ::   kbytes     ! size of the array pmess
491      INTEGER , INTENT(in   ) ::   kdest      ! receive process number
492      INTEGER , INTENT(in   ) ::   ktyp       ! tag of the message
493      INTEGER , INTENT(in   ) ::   md_req     ! argument for isend
494      !!
495      INTEGER ::   iflag
[51]496      !!----------------------------------------------------------------------
[1344]497      !
[1601]498      SELECT CASE ( cn_mpi_send )
[300]499      CASE ( 'S' )                ! Standard mpi send (blocking)
[9570]500         CALL mpi_send ( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce        , iflag )
[300]501      CASE ( 'B' )                ! Buffer mpi send (blocking)
[9570]502         CALL mpi_bsend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce        , iflag )
[300]503      CASE ( 'I' )                ! Immediate mpi send (non-blocking send)
[1344]504         ! be carefull, one more argument here : the mpi request identifier..
[9570]505         CALL mpi_isend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce, md_req, iflag )
[300]506      END SELECT
[1344]507      !
[51]508   END SUBROUTINE mppsend
[3]509
510
[3294]511   SUBROUTINE mpprecv( ktyp, pmess, kbytes, ksource )
[51]512      !!----------------------------------------------------------------------
513      !!                  ***  routine mpprecv  ***
514      !!
515      !! ** Purpose :   Receive messag passing array
516      !!
517      !!----------------------------------------------------------------------
[1344]518      REAL(wp), INTENT(inout) ::   pmess(*)   ! array of real
519      INTEGER , INTENT(in   ) ::   kbytes     ! suze of the array pmess
520      INTEGER , INTENT(in   ) ::   ktyp       ! Tag of the recevied message
[3764]521      INTEGER, OPTIONAL, INTENT(in) :: ksource    ! source process number
[1344]522      !!
[51]523      INTEGER :: istatus(mpi_status_size)
524      INTEGER :: iflag
[3294]525      INTEGER :: use_source
[1344]526      !!----------------------------------------------------------------------
527      !
[3764]528      ! If a specific process number has been passed to the receive call,
[3294]529      ! use that one. Default is to use mpi_any_source
[6140]530      use_source = mpi_any_source
531      IF( PRESENT(ksource) )   use_source = ksource
532      !
[9570]533      CALL mpi_recv( pmess, kbytes, mpi_double_precision, use_source, ktyp, mpi_comm_oce, istatus, iflag )
[1344]534      !
[51]535   END SUBROUTINE mpprecv
[3]536
537
[51]538   SUBROUTINE mppgather( ptab, kp, pio )
539      !!----------------------------------------------------------------------
540      !!                   ***  routine mppgather  ***
[3764]541      !!
542      !! ** Purpose :   Transfert between a local subdomain array and a work
[51]543      !!     array which is distributed following the vertical level.
544      !!
[1344]545      !!----------------------------------------------------------------------
[6140]546      REAL(wp), DIMENSION(jpi,jpj)      , INTENT(in   ) ::   ptab   ! subdomain input array
547      INTEGER                           , INTENT(in   ) ::   kp     ! record length
[1344]548      REAL(wp), DIMENSION(jpi,jpj,jpnij), INTENT(  out) ::   pio    ! subdomain input array
[51]549      !!
[1344]550      INTEGER :: itaille, ierror   ! temporary integer
[51]551      !!---------------------------------------------------------------------
[1344]552      !
553      itaille = jpi * jpj
554      CALL mpi_gather( ptab, itaille, mpi_double_precision, pio, itaille     ,   &
[9570]555         &                            mpi_double_precision, kp , mpi_comm_oce, ierror )
[1344]556      !
[51]557   END SUBROUTINE mppgather
[3]558
559
[51]560   SUBROUTINE mppscatter( pio, kp, ptab )
561      !!----------------------------------------------------------------------
562      !!                  ***  routine mppscatter  ***
563      !!
[3764]564      !! ** Purpose :   Transfert between awork array which is distributed
[51]565      !!      following the vertical level and the local subdomain array.
566      !!
567      !!----------------------------------------------------------------------
[6140]568      REAL(wp), DIMENSION(jpi,jpj,jpnij)  ::   pio    ! output array
569      INTEGER                             ::   kp     ! Tag (not used with MPI
570      REAL(wp), DIMENSION(jpi,jpj)        ::   ptab   ! subdomain array input
[1344]571      !!
572      INTEGER :: itaille, ierror   ! temporary integer
[51]573      !!---------------------------------------------------------------------
[1344]574      !
[6140]575      itaille = jpi * jpj
[1344]576      !
577      CALL mpi_scatter( pio, itaille, mpi_double_precision, ptab, itaille     ,   &
[9570]578         &                            mpi_double_precision, kp  , mpi_comm_oce, ierror )
[1344]579      !
[51]580   END SUBROUTINE mppscatter
[3]581
[10425]582   
583   SUBROUTINE mpp_delay_sum( cdname, cdelay, y_in, pout, ldlast, kcom )
584     !!----------------------------------------------------------------------
585      !!                   ***  routine mpp_delay_sum  ***
[1344]586      !!
[10425]587      !! ** Purpose :   performed delayed mpp_sum, the result is received on next call
588      !!
[1344]589      !!----------------------------------------------------------------------
[10425]590      CHARACTER(len=*), INTENT(in   )               ::   cdname  ! name of the calling subroutine
591      CHARACTER(len=*), INTENT(in   )               ::   cdelay  ! name (used as id) of the delayed operation
592      COMPLEX(wp),      INTENT(in   ), DIMENSION(:) ::   y_in
593      REAL(wp),         INTENT(  out), DIMENSION(:) ::   pout
594      LOGICAL,          INTENT(in   )               ::   ldlast  ! true if this is the last time we call this routine
595      INTEGER,          INTENT(in   ), OPTIONAL     ::   kcom
[1344]596      !!
[10425]597      INTEGER ::   ji, isz
598      INTEGER ::   idvar
599      INTEGER ::   ierr, ilocalcomm
600      COMPLEX(wp), ALLOCATABLE, DIMENSION(:) ::   ytmp
[1344]601      !!----------------------------------------------------------------------
[9570]602      ilocalcomm = mpi_comm_oce
[9019]603      IF( PRESENT(kcom) )   ilocalcomm = kcom
[3]604
[10425]605      isz = SIZE(y_in)
606     
[10437]607      IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ld_dlg = .TRUE. )
[13]608
[10425]609      idvar = -1
610      DO ji = 1, nbdelay
611         IF( TRIM(cdelay) == TRIM(c_delaylist(ji)) ) idvar = ji
612      END DO
613      IF ( idvar == -1 )   CALL ctl_stop( 'STOP',' mpp_delay_sum : please add a new delayed exchange for '//TRIM(cdname) )
614
615      IF ( ndelayid(idvar) == 0 ) THEN         ! first call    with restart: %z1d defined in iom_delay_rst
616         !                                       --------------------------
617         IF ( SIZE(todelay(idvar)%z1d) /= isz ) THEN                  ! Check dimension coherence
618            IF(lwp) WRITE(numout,*) ' WARNING: the nb of delayed variables in restart file is not the model one'
619            DEALLOCATE(todelay(idvar)%z1d)
620            ndelayid(idvar) = -1                                      ! do as if we had no restart
621         ELSE
622            ALLOCATE(todelay(idvar)%y1d(isz))
623            todelay(idvar)%y1d(:) = CMPLX(todelay(idvar)%z1d(:), 0., wp)   ! create %y1d, complex variable needed by mpi_sumdd
624         END IF
625      ENDIF
626     
627      IF( ndelayid(idvar) == -1 ) THEN         ! first call without restart: define %y1d and %z1d from y_in with blocking allreduce
628         !                                       --------------------------
629         ALLOCATE(todelay(idvar)%z1d(isz), todelay(idvar)%y1d(isz))   ! allocate also %z1d as used for the restart
630         CALL mpi_allreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ierr )   ! get %y1d
631         todelay(idvar)%z1d(:) = REAL(todelay(idvar)%y1d(:), wp)      ! define %z1d from %y1d
632      ENDIF
633
634      IF( ndelayid(idvar) > 0 )   CALL mpp_delay_rcv( idvar )         ! make sure %z1d is received
635
636      ! send back pout from todelay(idvar)%z1d defined at previous call
637      pout(:) = todelay(idvar)%z1d(:)
638
639      ! send y_in into todelay(idvar)%y1d with a non-blocking communication
640      CALL mpi_iallreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ndelayid(idvar), ierr )
641
642   END SUBROUTINE mpp_delay_sum
643
[9019]644   
[10425]645   SUBROUTINE mpp_delay_max( cdname, cdelay, p_in, pout, ldlast, kcom )
[9019]646      !!----------------------------------------------------------------------
[10425]647      !!                   ***  routine mpp_delay_max  ***
[9019]648      !!
[10425]649      !! ** Purpose :   performed delayed mpp_max, the result is received on next call
[9019]650      !!
651      !!----------------------------------------------------------------------
[10425]652      CHARACTER(len=*), INTENT(in   )                 ::   cdname  ! name of the calling subroutine
653      CHARACTER(len=*), INTENT(in   )                 ::   cdelay  ! name (used as id) of the delayed operation
654      REAL(wp),         INTENT(in   ), DIMENSION(:)   ::   p_in    !
655      REAL(wp),         INTENT(  out), DIMENSION(:)   ::   pout    !
656      LOGICAL,          INTENT(in   )                 ::   ldlast  ! true if this is the last time we call this routine
657      INTEGER,          INTENT(in   ), OPTIONAL       ::   kcom
[9019]658      !!
[10425]659      INTEGER ::   ji, isz
660      INTEGER ::   idvar
661      INTEGER ::   ierr, ilocalcomm
[9019]662      !!----------------------------------------------------------------------
[9570]663      ilocalcomm = mpi_comm_oce
[9019]664      IF( PRESENT(kcom) )   ilocalcomm = kcom
[6140]665
[10425]666      isz = SIZE(p_in)
[9019]667
[10437]668      IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ld_dlg = .TRUE. )
[13]669
[10425]670      idvar = -1
671      DO ji = 1, nbdelay
672         IF( TRIM(cdelay) == TRIM(c_delaylist(ji)) ) idvar = ji
673      END DO
674      IF ( idvar == -1 )   CALL ctl_stop( 'STOP',' mpp_delay_max : please add a new delayed exchange for '//TRIM(cdname) )
[3]675
[10425]676      IF ( ndelayid(idvar) == 0 ) THEN         ! first call    with restart: %z1d defined in iom_delay_rst
677         !                                       --------------------------
678         IF ( SIZE(todelay(idvar)%z1d) /= isz ) THEN                  ! Check dimension coherence
679            IF(lwp) WRITE(numout,*) ' WARNING: the nb of delayed variables in restart file is not the model one'
680            DEALLOCATE(todelay(idvar)%z1d)
681            ndelayid(idvar) = -1                                      ! do as if we had no restart
682         END IF
683      ENDIF
[13]684
[10425]685      IF( ndelayid(idvar) == -1 ) THEN         ! first call without restart: define %z1d from p_in with a blocking allreduce
686         !                                       --------------------------
687         ALLOCATE(todelay(idvar)%z1d(isz))
688         CALL mpi_allreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ierr )   ! get %z1d
689      ENDIF
[3]690
[10425]691      IF( ndelayid(idvar) > 0 )   CALL mpp_delay_rcv( idvar )         ! make sure %z1d is received
[3]692
[10425]693      ! send back pout from todelay(idvar)%z1d defined at previous call
694      pout(:) = todelay(idvar)%z1d(:)
[13]695
[10425]696      ! send p_in into todelay(idvar)%z1d with a non-blocking communication
697      CALL mpi_iallreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ndelayid(idvar), ierr )
698
699   END SUBROUTINE mpp_delay_max
700
701   
702   SUBROUTINE mpp_delay_rcv( kid )
703      !!----------------------------------------------------------------------
704      !!                   ***  routine mpp_delay_rcv  ***
[1344]705      !!
[10425]706      !! ** Purpose :  force barrier for delayed mpp (needed for restart)
[1344]707      !!
[10425]708      !!----------------------------------------------------------------------
709      INTEGER,INTENT(in   )      ::  kid 
710      INTEGER ::   ierr
711      !!----------------------------------------------------------------------
712      IF( ndelayid(kid) /= -2 ) THEN 
713         IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.)
714         CALL mpi_wait( ndelayid(kid), MPI_STATUS_IGNORE, ierr )                        ! make sure todelay(kid) is received
715         IF( ln_timing ) CALL tic_tac(.FALSE., ld_global = .TRUE.)
716         IF( ASSOCIATED(todelay(kid)%y1d) )   todelay(kid)%z1d(:) = REAL(todelay(kid)%y1d(:), wp)  ! define %z1d from %y1d
717         ndelayid(kid) = -2   ! add flag to know that mpi_wait was already called on kid
718      ENDIF
719   END SUBROUTINE mpp_delay_rcv
[3]720
[10425]721   
722   !!----------------------------------------------------------------------
723   !!    ***  mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real  ***
724   !!   
725   !!----------------------------------------------------------------------
726   !!
727#  define OPERATION_MAX
728#  define INTEGER_TYPE
729#  define DIM_0d
730#     define ROUTINE_ALLREDUCE           mppmax_int
731#     include "mpp_allreduce_generic.h90"
732#     undef ROUTINE_ALLREDUCE
733#  undef DIM_0d
734#  define DIM_1d
735#     define ROUTINE_ALLREDUCE           mppmax_a_int
736#     include "mpp_allreduce_generic.h90"
737#     undef ROUTINE_ALLREDUCE
738#  undef DIM_1d
739#  undef INTEGER_TYPE
740!
741#  define REAL_TYPE
742#  define DIM_0d
743#     define ROUTINE_ALLREDUCE           mppmax_real
744#     include "mpp_allreduce_generic.h90"
745#     undef ROUTINE_ALLREDUCE
746#  undef DIM_0d
747#  define DIM_1d
748#     define ROUTINE_ALLREDUCE           mppmax_a_real
749#     include "mpp_allreduce_generic.h90"
750#     undef ROUTINE_ALLREDUCE
751#  undef DIM_1d
752#  undef REAL_TYPE
753#  undef OPERATION_MAX
754   !!----------------------------------------------------------------------
755   !!    ***  mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real  ***
756   !!   
757   !!----------------------------------------------------------------------
758   !!
759#  define OPERATION_MIN
760#  define INTEGER_TYPE
761#  define DIM_0d
762#     define ROUTINE_ALLREDUCE           mppmin_int
763#     include "mpp_allreduce_generic.h90"
764#     undef ROUTINE_ALLREDUCE
765#  undef DIM_0d
766#  define DIM_1d
767#     define ROUTINE_ALLREDUCE           mppmin_a_int
768#     include "mpp_allreduce_generic.h90"
769#     undef ROUTINE_ALLREDUCE
770#  undef DIM_1d
771#  undef INTEGER_TYPE
772!
773#  define REAL_TYPE
774#  define DIM_0d
775#     define ROUTINE_ALLREDUCE           mppmin_real
776#     include "mpp_allreduce_generic.h90"
777#     undef ROUTINE_ALLREDUCE
778#  undef DIM_0d
779#  define DIM_1d
780#     define ROUTINE_ALLREDUCE           mppmin_a_real
781#     include "mpp_allreduce_generic.h90"
782#     undef ROUTINE_ALLREDUCE
783#  undef DIM_1d
784#  undef REAL_TYPE
785#  undef OPERATION_MIN
[869]786
[10425]787   !!----------------------------------------------------------------------
788   !!    ***  mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real  ***
789   !!   
790   !!   Global sum of 1D array or a variable (integer, real or complex)
791   !!----------------------------------------------------------------------
792   !!
793#  define OPERATION_SUM
794#  define INTEGER_TYPE
795#  define DIM_0d
796#     define ROUTINE_ALLREDUCE           mppsum_int
797#     include "mpp_allreduce_generic.h90"
798#     undef ROUTINE_ALLREDUCE
799#  undef DIM_0d
800#  define DIM_1d
801#     define ROUTINE_ALLREDUCE           mppsum_a_int
802#     include "mpp_allreduce_generic.h90"
803#     undef ROUTINE_ALLREDUCE
804#  undef DIM_1d
805#  undef INTEGER_TYPE
806!
807#  define REAL_TYPE
808#  define DIM_0d
809#     define ROUTINE_ALLREDUCE           mppsum_real
810#     include "mpp_allreduce_generic.h90"
811#     undef ROUTINE_ALLREDUCE
812#  undef DIM_0d
813#  define DIM_1d
814#     define ROUTINE_ALLREDUCE           mppsum_a_real
815#     include "mpp_allreduce_generic.h90"
816#     undef ROUTINE_ALLREDUCE
817#  undef DIM_1d
818#  undef REAL_TYPE
819#  undef OPERATION_SUM
820
821#  define OPERATION_SUM_DD
822#  define COMPLEX_TYPE
823#  define DIM_0d
824#     define ROUTINE_ALLREDUCE           mppsum_realdd
825#     include "mpp_allreduce_generic.h90"
826#     undef ROUTINE_ALLREDUCE
827#  undef DIM_0d
828#  define DIM_1d
829#     define ROUTINE_ALLREDUCE           mppsum_a_realdd
830#     include "mpp_allreduce_generic.h90"
831#     undef ROUTINE_ALLREDUCE
832#  undef DIM_1d
833#  undef COMPLEX_TYPE
834#  undef OPERATION_SUM_DD
835
836   !!----------------------------------------------------------------------
837   !!    ***  mpp_minloc2d, mpp_minloc3d, mpp_maxloc2d, mpp_maxloc3d
838   !!   
839   !!----------------------------------------------------------------------
840   !!
841#  define OPERATION_MINLOC
842#  define DIM_2d
843#     define ROUTINE_LOC           mpp_minloc2d
844#     include "mpp_loc_generic.h90"
845#     undef ROUTINE_LOC
846#  undef DIM_2d
847#  define DIM_3d
848#     define ROUTINE_LOC           mpp_minloc3d
849#     include "mpp_loc_generic.h90"
850#     undef ROUTINE_LOC
851#  undef DIM_3d
852#  undef OPERATION_MINLOC
853
854#  define OPERATION_MAXLOC
855#  define DIM_2d
856#     define ROUTINE_LOC           mpp_maxloc2d
857#     include "mpp_loc_generic.h90"
858#     undef ROUTINE_LOC
859#  undef DIM_2d
860#  define DIM_3d
861#     define ROUTINE_LOC           mpp_maxloc3d
862#     include "mpp_loc_generic.h90"
863#     undef ROUTINE_LOC
864#  undef DIM_3d
865#  undef OPERATION_MAXLOC
866
[1344]867   SUBROUTINE mppsync()
868      !!----------------------------------------------------------------------
869      !!                  ***  routine mppsync  ***
[3764]870      !!
[1344]871      !! ** Purpose :   Massively parallel processors, synchroneous
872      !!
873      !!-----------------------------------------------------------------------
874      INTEGER :: ierror
875      !!-----------------------------------------------------------------------
876      !
[9570]877      CALL mpi_barrier( mpi_comm_oce, ierror )
[1344]878      !
879   END SUBROUTINE mppsync
[3]880
881
[10425]882   SUBROUTINE mppstop( ldfinal, ld_force_abort ) 
[1344]883      !!----------------------------------------------------------------------
884      !!                  ***  routine mppstop  ***
[3764]885      !!
[3294]886      !! ** purpose :   Stop massively parallel processors method
[1344]887      !!
888      !!----------------------------------------------------------------------
[10425]889      LOGICAL, OPTIONAL, INTENT(in) :: ldfinal    ! source process number
890      LOGICAL, OPTIONAL, INTENT(in) :: ld_force_abort    ! source process number
891      LOGICAL ::   llfinal, ll_force_abort
[1344]892      INTEGER ::   info
893      !!----------------------------------------------------------------------
[10425]894      llfinal = .FALSE.
895      IF( PRESENT(ldfinal) ) llfinal = ldfinal
896      ll_force_abort = .FALSE.
897      IF( PRESENT(ld_force_abort) ) ll_force_abort = ld_force_abort
[1344]898      !
[10425]899      IF(ll_force_abort) THEN
900         CALL mpi_abort( MPI_COMM_WORLD )
901      ELSE
902         CALL mppsync
903         CALL mpi_finalize( info )
904      ENDIF
905      IF( .NOT. llfinal ) STOP 123456
[1344]906      !
907   END SUBROUTINE mppstop
[3]908
909
[1344]910   SUBROUTINE mpp_comm_free( kcom )
911      !!----------------------------------------------------------------------
912      INTEGER, INTENT(in) ::   kcom
913      !!
914      INTEGER :: ierr
915      !!----------------------------------------------------------------------
916      !
917      CALL MPI_COMM_FREE(kcom, ierr)
918      !
919   END SUBROUTINE mpp_comm_free
[3]920
[869]921
[2715]922   SUBROUTINE mpp_ini_znl( kumout )
[1345]923      !!----------------------------------------------------------------------
924      !!               ***  routine mpp_ini_znl  ***
925      !!
926      !! ** Purpose :   Initialize special communicator for computing zonal sum
927      !!
928      !! ** Method  : - Look for processors in the same row
929      !!              - Put their number in nrank_znl
930      !!              - Create group for the znl processors
931      !!              - Create a communicator for znl processors
932      !!              - Determine if processor should write znl files
933      !!
934      !! ** output
935      !!      ndim_rank_znl = number of processors on the same row
936      !!      ngrp_znl = group ID for the znl processors
937      !!      ncomm_znl = communicator for the ice procs.
938      !!      n_znl_root = number (in the world) of proc 0 in the ice comm.
939      !!
940      !!----------------------------------------------------------------------
[2715]941      INTEGER, INTENT(in) ::   kumout   ! ocean.output logical units
[1345]942      !
[2715]943      INTEGER :: jproc      ! dummy loop integer
944      INTEGER :: ierr, ii   ! local integer
945      INTEGER, ALLOCATABLE, DIMENSION(:) ::   kwork
946      !!----------------------------------------------------------------------
[1345]947      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_world     : ', ngrp_world
948      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_world : ', mpi_comm_world
[9570]949      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_oce   : ', mpi_comm_oce
[1345]950      !
[2715]951      ALLOCATE( kwork(jpnij), STAT=ierr )
952      IF( ierr /= 0 ) THEN
953         WRITE(kumout, cform_err)
954         WRITE(kumout,*) 'mpp_ini_znl : failed to allocate 1D array of length jpnij'
955         CALL mppstop
956      ENDIF
957
958      IF( jpnj == 1 ) THEN
[1345]959         ngrp_znl  = ngrp_world
[9570]960         ncomm_znl = mpi_comm_oce
[1345]961      ELSE
962         !
[9570]963         CALL MPI_ALLGATHER ( njmpp, 1, mpi_integer, kwork, 1, mpi_integer, mpi_comm_oce, ierr )
[1345]964         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - kwork pour njmpp : ', kwork
965         !-$$        CALL flush(numout)
966         !
967         ! Count number of processors on the same row
968         ndim_rank_znl = 0
969         DO jproc=1,jpnij
970            IF ( kwork(jproc) == njmpp ) THEN
971               ndim_rank_znl = ndim_rank_znl + 1
972            ENDIF
973         END DO
974         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ndim_rank_znl : ', ndim_rank_znl
975         !-$$        CALL flush(numout)
976         ! Allocate the right size to nrank_znl
[1441]977         IF (ALLOCATED (nrank_znl)) DEALLOCATE(nrank_znl)
[1345]978         ALLOCATE(nrank_znl(ndim_rank_znl))
[3764]979         ii = 0
[1345]980         nrank_znl (:) = 0
981         DO jproc=1,jpnij
982            IF ( kwork(jproc) == njmpp) THEN
983               ii = ii + 1
[3764]984               nrank_znl(ii) = jproc -1
[1345]985            ENDIF
986         END DO
987         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - nrank_znl : ', nrank_znl
988         !-$$        CALL flush(numout)
989
990         ! Create the opa group
[9570]991         CALL MPI_COMM_GROUP(mpi_comm_oce,ngrp_opa,ierr)
[1345]992         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_opa : ', ngrp_opa
993         !-$$        CALL flush(numout)
994
995         ! Create the znl group from the opa group
996         CALL MPI_GROUP_INCL  ( ngrp_opa, ndim_rank_znl, nrank_znl, ngrp_znl, ierr )
997         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_znl ', ngrp_znl
998         !-$$        CALL flush(numout)
999
1000         ! Create the znl communicator from the opa communicator, ie the pool of procs in the same row
[9570]1001         CALL MPI_COMM_CREATE ( mpi_comm_oce, ngrp_znl, ncomm_znl, ierr )
[1345]1002         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ncomm_znl ', ncomm_znl
1003         !-$$        CALL flush(numout)
1004         !
1005      END IF
1006
1007      ! Determines if processor if the first (starting from i=1) on the row
[3764]1008      IF ( jpni == 1 ) THEN
[1345]1009         l_znl_root = .TRUE.
1010      ELSE
1011         l_znl_root = .FALSE.
1012         kwork (1) = nimpp
[10425]1013         CALL mpp_min ( 'lib_mpp', kwork(1), kcom = ncomm_znl)
[1345]1014         IF ( nimpp == kwork(1)) l_znl_root = .TRUE.
1015      END IF
1016
[2715]1017      DEALLOCATE(kwork)
1018
[1345]1019   END SUBROUTINE mpp_ini_znl
1020
1021
[1344]1022   SUBROUTINE mpp_ini_north
1023      !!----------------------------------------------------------------------
1024      !!               ***  routine mpp_ini_north  ***
1025      !!
[3764]1026      !! ** Purpose :   Initialize special communicator for north folding
[1344]1027      !!      condition together with global variables needed in the mpp folding
1028      !!
1029      !! ** Method  : - Look for northern processors
1030      !!              - Put their number in nrank_north
1031      !!              - Create groups for the world processors and the north processors
1032      !!              - Create a communicator for northern processors
1033      !!
1034      !! ** output
1035      !!      njmppmax = njmpp for northern procs
1036      !!      ndim_rank_north = number of processors in the northern line
1037      !!      nrank_north (ndim_rank_north) = number  of the northern procs.
1038      !!      ngrp_world = group ID for the world processors
1039      !!      ngrp_north = group ID for the northern processors
1040      !!      ncomm_north = communicator for the northern procs.
1041      !!      north_root = number (in the world) of proc 0 in the northern comm.
1042      !!
1043      !!----------------------------------------------------------------------
1044      INTEGER ::   ierr
1045      INTEGER ::   jjproc
1046      INTEGER ::   ii, ji
1047      !!----------------------------------------------------------------------
1048      !
1049      njmppmax = MAXVAL( njmppt )
1050      !
1051      ! Look for how many procs on the northern boundary
1052      ndim_rank_north = 0
1053      DO jjproc = 1, jpnij
1054         IF( njmppt(jjproc) == njmppmax )   ndim_rank_north = ndim_rank_north + 1
1055      END DO
1056      !
1057      ! Allocate the right size to nrank_north
[1441]1058      IF (ALLOCATED (nrank_north)) DEALLOCATE(nrank_north)
[1344]1059      ALLOCATE( nrank_north(ndim_rank_north) )
[869]1060
[1344]1061      ! Fill the nrank_north array with proc. number of northern procs.
1062      ! Note : the rank start at 0 in MPI
1063      ii = 0
1064      DO ji = 1, jpnij
1065         IF ( njmppt(ji) == njmppmax   ) THEN
1066            ii=ii+1
1067            nrank_north(ii)=ji-1
1068         END IF
1069      END DO
1070      !
1071      ! create the world group
[9570]1072      CALL MPI_COMM_GROUP( mpi_comm_oce, ngrp_world, ierr )
[1344]1073      !
1074      ! Create the North group from the world group
1075      CALL MPI_GROUP_INCL( ngrp_world, ndim_rank_north, nrank_north, ngrp_north, ierr )
1076      !
1077      ! Create the North communicator , ie the pool of procs in the north group
[9570]1078      CALL MPI_COMM_CREATE( mpi_comm_oce, ngrp_north, ncomm_north, ierr )
[1344]1079      !
1080   END SUBROUTINE mpp_ini_north
[869]1081
1082
[9570]1083   SUBROUTINE mpi_init_oce( ldtxt, ksft, code )
[1344]1084      !!---------------------------------------------------------------------
1085      !!                   ***  routine mpp_init.opa  ***
1086      !!
1087      !! ** Purpose :: export and attach a MPI buffer for bsend
1088      !!
1089      !! ** Method  :: define buffer size in namelist, if 0 no buffer attachment
1090      !!            but classical mpi_init
[3764]1091      !!
1092      !! History :: 01/11 :: IDRIS initial version for IBM only
[1344]1093      !!            08/04 :: R. Benshila, generalisation
1094      !!---------------------------------------------------------------------
[3764]1095      CHARACTER(len=*),DIMENSION(:), INTENT(  out) ::   ldtxt
[2481]1096      INTEGER                      , INTENT(inout) ::   ksft
1097      INTEGER                      , INTENT(  out) ::   code
1098      INTEGER                                      ::   ierr, ji
1099      LOGICAL                                      ::   mpi_was_called
[1344]1100      !!---------------------------------------------------------------------
1101      !
1102      CALL mpi_initialized( mpi_was_called, code )      ! MPI initialization
[532]1103      IF ( code /= MPI_SUCCESS ) THEN
[3764]1104         DO ji = 1, SIZE(ldtxt)
[2481]1105            IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
[3764]1106         END DO
[2481]1107         WRITE(*, cform_err)
1108         WRITE(*, *) ' lib_mpp: Error in routine mpi_initialized'
[1344]1109         CALL mpi_abort( mpi_comm_world, code, ierr )
[532]1110      ENDIF
[1344]1111      !
1112      IF( .NOT. mpi_was_called ) THEN
1113         CALL mpi_init( code )
[9570]1114         CALL mpi_comm_dup( mpi_comm_world, mpi_comm_oce, code )
[532]1115         IF ( code /= MPI_SUCCESS ) THEN
[3764]1116            DO ji = 1, SIZE(ldtxt)
[2481]1117               IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
1118            END DO
1119            WRITE(*, cform_err)
1120            WRITE(*, *) ' lib_mpp: Error in routine mpi_comm_dup'
[532]1121            CALL mpi_abort( mpi_comm_world, code, ierr )
1122         ENDIF
1123      ENDIF
[1344]1124      !
[897]1125      IF( nn_buffer > 0 ) THEN
[2481]1126         WRITE(ldtxt(ksft),*) 'mpi_bsend, buffer allocation of  : ', nn_buffer   ;   ksft = ksft + 1
[897]1127         ! Buffer allocation and attachment
[2481]1128         ALLOCATE( tampon(nn_buffer), stat = ierr )
[3764]1129         IF( ierr /= 0 ) THEN
1130            DO ji = 1, SIZE(ldtxt)
[2481]1131               IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
1132            END DO
1133            WRITE(*, cform_err)
1134            WRITE(*, *) ' lib_mpp: Error in ALLOCATE', ierr
1135            CALL mpi_abort( mpi_comm_world, code, ierr )
1136         END IF
1137         CALL mpi_buffer_attach( tampon, nn_buffer, code )
[897]1138      ENDIF
[1344]1139      !
[9570]1140   END SUBROUTINE mpi_init_oce
[3]1141
[9019]1142
1143   SUBROUTINE DDPDD_MPI( ydda, yddb, ilen, itype )
[1976]1144      !!---------------------------------------------------------------------
1145      !!   Routine DDPDD_MPI: used by reduction operator MPI_SUMDD
1146      !!
1147      !!   Modification of original codes written by David H. Bailey
1148      !!   This subroutine computes yddb(i) = ydda(i)+yddb(i)
1149      !!---------------------------------------------------------------------
[9019]1150      INTEGER                     , INTENT(in)    ::   ilen, itype
1151      COMPLEX(wp), DIMENSION(ilen), INTENT(in)    ::   ydda
1152      COMPLEX(wp), DIMENSION(ilen), INTENT(inout) ::   yddb
[1976]1153      !
1154      REAL(wp) :: zerr, zt1, zt2    ! local work variables
[9019]1155      INTEGER  :: ji, ztmp           ! local scalar
1156      !!---------------------------------------------------------------------
1157      !
[1976]1158      ztmp = itype   ! avoid compilation warning
[9019]1159      !
[1976]1160      DO ji=1,ilen
1161      ! Compute ydda + yddb using Knuth's trick.
1162         zt1  = real(ydda(ji)) + real(yddb(ji))
1163         zerr = zt1 - real(ydda(ji))
1164         zt2  = ((real(yddb(ji)) - zerr) + (real(ydda(ji)) - (zt1 - zerr))) &
1165                + aimag(ydda(ji)) + aimag(yddb(ji))
1166
1167         ! The result is zt1 + zt2, after normalization.
1168         yddb(ji) = cmplx ( zt1 + zt2, zt2 - ((zt1 + zt2) - zt1),wp )
1169      END DO
[9019]1170      !
[1976]1171   END SUBROUTINE DDPDD_MPI
1172
[6140]1173
[9019]1174   SUBROUTINE mpp_lbc_north_icb( pt2d, cd_type, psgn, kextj)
[4990]1175      !!---------------------------------------------------------------------
1176      !!                   ***  routine mpp_lbc_north_icb  ***
1177      !!
1178      !! ** Purpose :   Ensure proper north fold horizontal bondary condition
1179      !!              in mpp configuration in case of jpn1 > 1 and for 2d
1180      !!              array with outer extra halo
1181      !!
1182      !! ** Method  :   North fold condition and mpp with more than one proc
1183      !!              in i-direction require a specific treatment. We gather
[9019]1184      !!              the 4+kextj northern lines of the global domain on 1
[4990]1185      !!              processor and apply lbc north-fold on this sub array.
1186      !!              Then we scatter the north fold array back to the processors.
[9019]1187      !!              This routine accounts for an extra halo with icebergs
1188      !!              and assumes ghost rows and columns have been suppressed.
[4990]1189      !!
1190      !!----------------------------------------------------------------------
1191      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pt2d     ! 2D array with extra halo
1192      CHARACTER(len=1)        , INTENT(in   ) ::   cd_type  ! nature of pt3d grid-points
1193      !                                                     !   = T ,  U , V , F or W -points
1194      REAL(wp)                , INTENT(in   ) ::   psgn     ! = -1. the sign change across the
1195      !!                                                    ! north fold, =  1. otherwise
[9019]1196      INTEGER                 , INTENT(in   ) ::   kextj    ! Extra halo width at north fold
[6140]1197      !
[4990]1198      INTEGER ::   ji, jj, jr
1199      INTEGER ::   ierr, itaille, ildi, ilei, iilb
[9019]1200      INTEGER ::   ipj, ij, iproc
[4990]1201      !
1202      REAL(wp), DIMENSION(:,:)  , ALLOCATABLE  ::  ztab_e, znorthloc_e
1203      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE  ::  znorthgloio_e
1204      !!----------------------------------------------------------------------
1205      !
[9019]1206      ipj=4
[9467]1207      ALLOCATE(        ztab_e(jpiglo, 1-kextj:ipj+kextj)       ,       &
1208     &            znorthloc_e(jpimax, 1-kextj:ipj+kextj)       ,       &
1209     &          znorthgloio_e(jpimax, 1-kextj:ipj+kextj,jpni)    )
[4990]1210      !
[9019]1211      ztab_e(:,:)      = 0._wp
1212      znorthloc_e(:,:) = 0._wp
[6140]1213      !
[9467]1214      ij = 1 - kextj
1215      ! put the last ipj+2*kextj lines of pt2d into znorthloc_e
1216      DO jj = jpj - ipj + 1 - kextj , jpj + kextj
1217         znorthloc_e(1:jpi,ij)=pt2d(1:jpi,jj)
[4990]1218         ij = ij + 1
1219      END DO
1220      !
[9467]1221      itaille = jpimax * ( ipj + 2*kextj )
[10425]1222      !
1223      IF( ln_timing ) CALL tic_tac(.TRUE.)
[9467]1224      CALL MPI_ALLGATHER( znorthloc_e(1,1-kextj)    , itaille, MPI_DOUBLE_PRECISION,    &
1225         &                znorthgloio_e(1,1-kextj,1), itaille, MPI_DOUBLE_PRECISION,    &
1226         &                ncomm_north, ierr )
[4990]1227      !
[10425]1228      IF( ln_timing ) CALL tic_tac(.FALSE.)
1229      !
[4990]1230      DO jr = 1, ndim_rank_north            ! recover the global north array
1231         iproc = nrank_north(jr) + 1
1232         ildi = nldit (iproc)
1233         ilei = nleit (iproc)
1234         iilb = nimppt(iproc)
[9467]1235         DO jj = 1-kextj, ipj+kextj
[4990]1236            DO ji = ildi, ilei
1237               ztab_e(ji+iilb-1,jj) = znorthgloio_e(ji,jj,jr)
1238            END DO
1239         END DO
1240      END DO
1241
1242      ! 2. North-Fold boundary conditions
1243      ! ----------------------------------
[9467]1244      CALL lbc_nfd( ztab_e(:,1-kextj:ipj+kextj), cd_type, psgn, kextj )
[4990]1245
[9467]1246      ij = 1 - kextj
[4990]1247      !! Scatter back to pt2d
[9467]1248      DO jj = jpj - ipj + 1 - kextj , jpj + kextj
[9019]1249         DO ji= 1, jpi
[4990]1250            pt2d(ji,jj) = ztab_e(ji+nimpp-1,ij)
1251         END DO
[9467]1252         ij  = ij +1
[4990]1253      END DO
1254      !
1255      DEALLOCATE( ztab_e, znorthloc_e, znorthgloio_e )
1256      !
1257   END SUBROUTINE mpp_lbc_north_icb
1258
[6140]1259
[10425]1260   SUBROUTINE mpp_lnk_2d_icb( cdname, pt2d, cd_type, psgn, kexti, kextj )
[4990]1261      !!----------------------------------------------------------------------
1262      !!                  ***  routine mpp_lnk_2d_icb  ***
1263      !!
[9019]1264      !! ** Purpose :   Message passing management for 2d array (with extra halo for icebergs)
1265      !!                This routine receives a (1-kexti:jpi+kexti,1-kexti:jpj+kextj)
1266      !!                array (usually (0:jpi+1, 0:jpj+1)) from lbc_lnk_icb calls.
[4990]1267      !!
1268      !! ** Method  :   Use mppsend and mpprecv function for passing mask
1269      !!      between processors following neighboring subdomains.
1270      !!            domain parameters
[9019]1271      !!                    jpi    : first dimension of the local subdomain
1272      !!                    jpj    : second dimension of the local subdomain
1273      !!                    kexti  : number of columns for extra outer halo
1274      !!                    kextj  : number of rows for extra outer halo
[4990]1275      !!                    nbondi : mark for "east-west local boundary"
1276      !!                    nbondj : mark for "north-south local boundary"
1277      !!                    noea   : number for local neighboring processors
1278      !!                    nowe   : number for local neighboring processors
1279      !!                    noso   : number for local neighboring processors
1280      !!                    nono   : number for local neighboring processors
1281      !!----------------------------------------------------------------------
[10425]1282      CHARACTER(len=*)                                        , INTENT(in   ) ::   cdname      ! name of the calling subroutine
[9019]1283      REAL(wp), DIMENSION(1-kexti:jpi+kexti,1-kextj:jpj+kextj), INTENT(inout) ::   pt2d     ! 2D array with extra halo
1284      CHARACTER(len=1)                                        , INTENT(in   ) ::   cd_type  ! nature of ptab array grid-points
1285      REAL(wp)                                                , INTENT(in   ) ::   psgn     ! sign used across the north fold
1286      INTEGER                                                 , INTENT(in   ) ::   kexti    ! extra i-halo width
1287      INTEGER                                                 , INTENT(in   ) ::   kextj    ! extra j-halo width
1288      !
[4990]1289      INTEGER  ::   jl   ! dummy loop indices
[9019]1290      INTEGER  ::   imigr, iihom, ijhom        ! local integers
1291      INTEGER  ::   ipreci, iprecj             !   -       -
[4990]1292      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend
1293      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend
1294      !!
[9019]1295      REAL(wp), DIMENSION(1-kexti:jpi+kexti,nn_hls+kextj,2) ::   r2dns, r2dsn
1296      REAL(wp), DIMENSION(1-kextj:jpj+kextj,nn_hls+kexti,2) ::   r2dwe, r2dew
[4990]1297      !!----------------------------------------------------------------------
1298
[9019]1299      ipreci = nn_hls + kexti      ! take into account outer extra 2D overlap area
1300      iprecj = nn_hls + kextj
[4990]1301
[10425]1302      IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, 1, 1, 1, ld_lbc = .TRUE. )
[4990]1303
1304      ! 1. standard boundary treatment
1305      ! ------------------------------
1306      ! Order matters Here !!!!
1307      !
1308      !                                      ! East-West boundaries
1309      !                                           !* Cyclic east-west
[9667]1310      IF( l_Iperio ) THEN
1311         pt2d(1-kexti:     1   ,:) = pt2d(jpim1-kexti: jpim1 ,:)       ! east
1312         pt2d(  jpi  :jpi+kexti,:) = pt2d(     2     :2+kexti,:)       ! west
[4990]1313         !
1314      ELSE                                        !* closed
[9667]1315         IF( .NOT. cd_type == 'F' )   pt2d(  1-kexti   :nn_hls   ,:) = 0._wp    ! east except at F-point
1316                                      pt2d(jpi-nn_hls+1:jpi+kexti,:) = 0._wp    ! west
[4990]1317      ENDIF
[9667]1318      !                                      ! North-South boundaries
1319      IF( l_Jperio ) THEN                         !* cyclic (only with no mpp j-split)
1320         pt2d(:,1-kextj:     1   ) = pt2d(:,jpjm1-kextj:  jpjm1)       ! north
1321         pt2d(:,  jpj  :jpj+kextj) = pt2d(:,     2     :2+kextj)       ! south
1322      ELSE                                        !* closed
1323         IF( .NOT. cd_type == 'F' )   pt2d(:,  1-kextj   :nn_hls   ) = 0._wp    ! north except at F-point
1324                                      pt2d(:,jpj-nn_hls+1:jpj+kextj) = 0._wp    ! south
1325      ENDIF
[4990]1326      !
1327
1328      ! north fold treatment
1329      ! -----------------------
1330      IF( npolj /= 0 ) THEN
1331         !
1332         SELECT CASE ( jpni )
[9019]1333                   CASE ( 1 )     ;   CALL lbc_nfd          ( pt2d(1:jpi,1:jpj+kextj), cd_type, psgn, kextj )
1334                   CASE DEFAULT   ;   CALL mpp_lbc_north_icb( pt2d(1:jpi,1:jpj+kextj), cd_type, psgn, kextj )
[4990]1335         END SELECT
1336         !
1337      ENDIF
1338
1339      ! 2. East and west directions exchange
1340      ! ------------------------------------
1341      ! we play with the neigbours AND the row number because of the periodicity
1342      !
1343      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions
1344      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case)
[9019]1345         iihom = jpi-nreci-kexti
[4990]1346         DO jl = 1, ipreci
[9019]1347            r2dew(:,jl,1) = pt2d(nn_hls+jl,:)
[4990]1348            r2dwe(:,jl,1) = pt2d(iihom +jl,:)
1349         END DO
1350      END SELECT
1351      !
1352      !                           ! Migrations
[9019]1353      imigr = ipreci * ( jpj + 2*kextj )
[4990]1354      !
[10425]1355      IF( ln_timing ) CALL tic_tac(.TRUE.)
1356      !
[4990]1357      SELECT CASE ( nbondi )
1358      CASE ( -1 )
[9019]1359         CALL mppsend( 2, r2dwe(1-kextj,1,1), imigr, noea, ml_req1 )
1360         CALL mpprecv( 1, r2dew(1-kextj,1,2), imigr, noea )
[4990]1361         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1362      CASE ( 0 )
[9019]1363         CALL mppsend( 1, r2dew(1-kextj,1,1), imigr, nowe, ml_req1 )
1364         CALL mppsend( 2, r2dwe(1-kextj,1,1), imigr, noea, ml_req2 )
1365         CALL mpprecv( 1, r2dew(1-kextj,1,2), imigr, noea )
1366         CALL mpprecv( 2, r2dwe(1-kextj,1,2), imigr, nowe )
[4990]1367         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1368         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err)
1369      CASE ( 1 )
[9019]1370         CALL mppsend( 1, r2dew(1-kextj,1,1), imigr, nowe, ml_req1 )
1371         CALL mpprecv( 2, r2dwe(1-kextj,1,2), imigr, nowe )
[4990]1372         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1373      END SELECT
1374      !
[10425]1375      IF( ln_timing ) CALL tic_tac(.FALSE.)
1376      !
[4990]1377      !                           ! Write Dirichlet lateral conditions
[9019]1378      iihom = jpi - nn_hls
[4990]1379      !
1380      SELECT CASE ( nbondi )
1381      CASE ( -1 )
1382         DO jl = 1, ipreci
1383            pt2d(iihom+jl,:) = r2dew(:,jl,2)
1384         END DO
1385      CASE ( 0 )
1386         DO jl = 1, ipreci
[9019]1387            pt2d(jl-kexti,:) = r2dwe(:,jl,2)
1388            pt2d(iihom+jl,:) = r2dew(:,jl,2)
[4990]1389         END DO
1390      CASE ( 1 )
1391         DO jl = 1, ipreci
[9019]1392            pt2d(jl-kexti,:) = r2dwe(:,jl,2)
[4990]1393         END DO
1394      END SELECT
1395
1396
1397      ! 3. North and south directions
1398      ! -----------------------------
1399      ! always closed : we play only with the neigbours
1400      !
1401      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions
[9019]1402         ijhom = jpj-nrecj-kextj
[4990]1403         DO jl = 1, iprecj
1404            r2dsn(:,jl,1) = pt2d(:,ijhom +jl)
[9019]1405            r2dns(:,jl,1) = pt2d(:,nn_hls+jl)
[4990]1406         END DO
1407      ENDIF
1408      !
1409      !                           ! Migrations
[9019]1410      imigr = iprecj * ( jpi + 2*kexti )
[4990]1411      !
[10425]1412      IF( ln_timing ) CALL tic_tac(.TRUE.)
1413      !
[4990]1414      SELECT CASE ( nbondj )
1415      CASE ( -1 )
[9019]1416         CALL mppsend( 4, r2dsn(1-kexti,1,1), imigr, nono, ml_req1 )
1417         CALL mpprecv( 3, r2dns(1-kexti,1,2), imigr, nono )
[4990]1418         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1419      CASE ( 0 )
[9019]1420         CALL mppsend( 3, r2dns(1-kexti,1,1), imigr, noso, ml_req1 )
1421         CALL mppsend( 4, r2dsn(1-kexti,1,1), imigr, nono, ml_req2 )
1422         CALL mpprecv( 3, r2dns(1-kexti,1,2), imigr, nono )
1423         CALL mpprecv( 4, r2dsn(1-kexti,1,2), imigr, noso )
[4990]1424         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1425         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err)
1426      CASE ( 1 )
[9019]1427         CALL mppsend( 3, r2dns(1-kexti,1,1), imigr, noso, ml_req1 )
1428         CALL mpprecv( 4, r2dsn(1-kexti,1,2), imigr, noso )
[4990]1429         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1430      END SELECT
1431      !
[10425]1432      IF( ln_timing ) CALL tic_tac(.FALSE.)
1433      !
[4990]1434      !                           ! Write Dirichlet lateral conditions
[9019]1435      ijhom = jpj - nn_hls
[4990]1436      !
1437      SELECT CASE ( nbondj )
1438      CASE ( -1 )
1439         DO jl = 1, iprecj
1440            pt2d(:,ijhom+jl) = r2dns(:,jl,2)
1441         END DO
1442      CASE ( 0 )
1443         DO jl = 1, iprecj
[9019]1444            pt2d(:,jl-kextj) = r2dsn(:,jl,2)
1445            pt2d(:,ijhom+jl) = r2dns(:,jl,2)
[4990]1446         END DO
1447      CASE ( 1 )
1448         DO jl = 1, iprecj
[9019]1449            pt2d(:,jl-kextj) = r2dsn(:,jl,2)
[4990]1450         END DO
1451      END SELECT
[9019]1452      !
[4990]1453   END SUBROUTINE mpp_lnk_2d_icb
[10425]1454
1455
[10437]1456   SUBROUTINE mpp_report( cdname, kpk, kpl, kpf, ld_lbc, ld_glb, ld_dlg )
[10425]1457      !!----------------------------------------------------------------------
1458      !!                  ***  routine mpp_report  ***
1459      !!
1460      !! ** Purpose :   report use of mpp routines per time-setp
1461      !!
1462      !!----------------------------------------------------------------------
1463      CHARACTER(len=*),           INTENT(in   ) ::   cdname      ! name of the calling subroutine
1464      INTEGER         , OPTIONAL, INTENT(in   ) ::   kpk, kpl, kpf
[10437]1465      LOGICAL         , OPTIONAL, INTENT(in   ) ::   ld_lbc, ld_glb, ld_dlg
[10425]1466      !!
[10437]1467      LOGICAL ::   ll_lbc, ll_glb, ll_dlg
[10425]1468      INTEGER ::    ji,  jj,  jk,  jh, jf   ! dummy loop indices
1469      !!----------------------------------------------------------------------
1470      !
1471      ll_lbc = .FALSE.
1472      IF( PRESENT(ld_lbc) ) ll_lbc = ld_lbc
1473      ll_glb = .FALSE.
1474      IF( PRESENT(ld_glb) ) ll_glb = ld_glb
[10437]1475      ll_dlg = .FALSE.
1476      IF( PRESENT(ld_dlg) ) ll_dlg = ld_dlg
[10425]1477      !
1478      ! find the smallest common frequency: default = frequency product, if multiple, choose the larger of the 2 frequency
1479      IF( ncom_dttrc /= 1 )   CALL ctl_stop( 'STOP', 'mpp_report, ncom_dttrc /= 1 not coded...' ) 
1480      ncom_freq = ncom_fsbc
1481      !
1482      IF ( ncom_stp == nit000+ncom_freq ) THEN   ! avoid to count extra communications in potential initializations at nit000
1483         IF( ll_lbc ) THEN
1484            IF( .NOT. ALLOCATED(ncomm_sequence) ) ALLOCATE( ncomm_sequence(ncom_rec_max,2) )
1485            IF( .NOT. ALLOCATED(    crname_lbc) ) ALLOCATE(     crname_lbc(ncom_rec_max  ) )
1486            n_sequence_lbc = n_sequence_lbc + 1
1487            IF( n_sequence_lbc > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' )   ! deadlock
1488            crname_lbc(n_sequence_lbc) = cdname     ! keep the name of the calling routine
1489            ncomm_sequence(n_sequence_lbc,1) = kpk*kpl   ! size of 3rd and 4th dimensions
1490            ncomm_sequence(n_sequence_lbc,2) = kpf       ! number of arrays to be treated (multi)
1491         ENDIF
1492         IF( ll_glb ) THEN
1493            IF( .NOT. ALLOCATED(crname_glb) ) ALLOCATE( crname_glb(ncom_rec_max) )
1494            n_sequence_glb = n_sequence_glb + 1
1495            IF( n_sequence_glb > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' )   ! deadlock
1496            crname_glb(n_sequence_glb) = cdname     ! keep the name of the calling routine
1497         ENDIF
[10437]1498         IF( ll_dlg ) THEN
1499            IF( .NOT. ALLOCATED(crname_dlg) ) ALLOCATE( crname_dlg(ncom_rec_max) )
1500            n_sequence_dlg = n_sequence_dlg + 1
1501            IF( n_sequence_dlg > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' )   ! deadlock
1502            crname_dlg(n_sequence_dlg) = cdname     ! keep the name of the calling routine
1503         ENDIF
[10425]1504      ELSE IF ( ncom_stp == nit000+2*ncom_freq ) THEN
1505         CALL ctl_opn( numcom, 'communication_report.txt', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE., narea )
1506         WRITE(numcom,*) ' '
1507         WRITE(numcom,*) ' ------------------------------------------------------------'
1508         WRITE(numcom,*) ' Communication pattern report (second oce+sbc+top time step):'
1509         WRITE(numcom,*) ' ------------------------------------------------------------'
1510         WRITE(numcom,*) ' '
1511         WRITE(numcom,'(A,I4)') ' Exchanged halos : ', n_sequence_lbc
1512         jj = 0; jk = 0; jf = 0; jh = 0
1513         DO ji = 1, n_sequence_lbc
1514            IF ( ncomm_sequence(ji,1) .GT. 1 ) jk = jk + 1
1515            IF ( ncomm_sequence(ji,2) .GT. 1 ) jf = jf + 1
1516            IF ( ncomm_sequence(ji,1) .GT. 1 .AND. ncomm_sequence(ji,2) .GT. 1 ) jj = jj + 1
1517            jh = MAX (jh, ncomm_sequence(ji,1)*ncomm_sequence(ji,2))
1518         END DO
1519         WRITE(numcom,'(A,I3)') ' 3D Exchanged halos : ', jk
1520         WRITE(numcom,'(A,I3)') ' Multi arrays exchanged halos : ', jf
1521         WRITE(numcom,'(A,I3)') '   from which 3D : ', jj
1522         WRITE(numcom,'(A,I10)') ' Array max size : ', jh*jpi*jpj
1523         WRITE(numcom,*) ' '
1524         WRITE(numcom,*) ' lbc_lnk called'
1525         jj = 1
1526         DO ji = 2, n_sequence_lbc
1527            IF( crname_lbc(ji-1) /= crname_lbc(ji) ) THEN
1528               WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_lbc(ji-1))
1529               jj = 0
1530            END IF
1531            jj = jj + 1 
1532         END DO
1533         WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_lbc(n_sequence_lbc))
1534         WRITE(numcom,*) ' '
1535         IF ( n_sequence_glb > 0 ) THEN
1536            WRITE(numcom,'(A,I4)') ' Global communications : ', n_sequence_glb
1537            jj = 1
1538            DO ji = 2, n_sequence_glb
1539               IF( crname_glb(ji-1) /= crname_glb(ji) ) THEN
1540                  WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_glb(ji-1))
1541                  jj = 0
1542               END IF
1543               jj = jj + 1 
1544            END DO
1545            WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_glb(n_sequence_glb))
1546            DEALLOCATE(crname_glb)
1547         ELSE
1548            WRITE(numcom,*) ' No MPI global communication '
1549         ENDIF
1550         WRITE(numcom,*) ' '
[10437]1551         IF ( n_sequence_dlg > 0 ) THEN
1552            WRITE(numcom,'(A,I4)') ' Delayed global communications : ', n_sequence_dlg
1553            jj = 1
1554            DO ji = 2, n_sequence_dlg
1555               IF( crname_dlg(ji-1) /= crname_dlg(ji) ) THEN
1556                  WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_dlg(ji-1))
1557                  jj = 0
1558               END IF
1559               jj = jj + 1 
1560            END DO
1561            WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_dlg(n_sequence_dlg))
1562            DEALLOCATE(crname_dlg)
1563         ELSE
1564            WRITE(numcom,*) ' No MPI delayed global communication '
1565         ENDIF
1566         WRITE(numcom,*) ' '
[10425]1567         WRITE(numcom,*) ' -----------------------------------------------'
1568         WRITE(numcom,*) ' '
1569         DEALLOCATE(ncomm_sequence)
1570         DEALLOCATE(crname_lbc)
1571      ENDIF
1572   END SUBROUTINE mpp_report
1573
[6140]1574   
[10425]1575   SUBROUTINE tic_tac (ld_tic, ld_global)
1576
1577    LOGICAL,           INTENT(IN) :: ld_tic
1578    LOGICAL, OPTIONAL, INTENT(IN) :: ld_global
1579    REAL(wp), DIMENSION(2), SAVE :: tic_wt
1580    REAL(wp),               SAVE :: tic_ct = 0._wp
1581    INTEGER :: ii
1582
1583    IF( ncom_stp <= nit000 ) RETURN
1584    IF( ncom_stp == nitend ) RETURN
1585    ii = 1
1586    IF( PRESENT( ld_global ) ) THEN
1587       IF( ld_global ) ii = 2
1588    END IF
1589   
1590    IF ( ld_tic ) THEN
1591       tic_wt(ii) = MPI_Wtime()                                                    ! start count tic->tac (waiting time)
1592       IF ( tic_ct > 0.0_wp ) compute_time = compute_time + MPI_Wtime() - tic_ct   ! cumulate count tac->tic
1593    ELSE
1594       waiting_time(ii) = waiting_time(ii) + MPI_Wtime() - tic_wt(ii)              ! cumulate count tic->tac
1595       tic_ct = MPI_Wtime()                                                        ! start count tac->tic (waiting time)
1596    ENDIF
1597   
1598   END SUBROUTINE tic_tac
1599
1600   
[13]1601#else
1602   !!----------------------------------------------------------------------
1603   !!   Default case:            Dummy module        share memory computing
1604   !!----------------------------------------------------------------------
[2715]1605   USE in_out_manager
[1976]1606
[13]1607   INTERFACE mpp_sum
[10425]1608      MODULE PROCEDURE mppsum_int, mppsum_a_int, mppsum_real, mppsum_a_real, mppsum_realdd, mppsum_a_realdd
[13]1609   END INTERFACE
1610   INTERFACE mpp_max
[681]1611      MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real
[13]1612   END INTERFACE
1613   INTERFACE mpp_min
1614      MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real
1615   END INTERFACE
[1344]1616   INTERFACE mpp_minloc
1617      MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d
1618   END INTERFACE
1619   INTERFACE mpp_maxloc
1620      MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d
1621   END INTERFACE
[3]1622
[13]1623   LOGICAL, PUBLIC, PARAMETER ::   lk_mpp = .FALSE.      !: mpp flag
[4147]1624   LOGICAL, PUBLIC            ::   ln_nnogather          !: namelist control of northfold comms (needed here in case "key_mpp_mpi" is not used)
[9570]1625   INTEGER, PUBLIC            ::   mpi_comm_oce          ! opa local communicator
[10425]1626
1627   INTEGER, PARAMETER, PUBLIC               ::   nbdelay = 0   ! make sure we don't enter loops: DO ji = 1, nbdelay
1628   CHARACTER(len=32), DIMENSION(1), PUBLIC  ::   c_delaylist = 'empty'
1629   CHARACTER(len=32), DIMENSION(1), PUBLIC  ::   c_delaycpnt = 'empty'
1630   TYPE ::   DELAYARR
1631      REAL(   wp), POINTER, DIMENSION(:) ::  z1d => NULL()
1632      COMPLEX(wp), POINTER, DIMENSION(:) ::  y1d => NULL()
1633   END TYPE DELAYARR
1634   TYPE( DELAYARR ), DIMENSION(1), PUBLIC  ::   todelay             
1635   INTEGER,  PUBLIC, DIMENSION(1)           ::   ndelayid = -1
[2715]1636   !!----------------------------------------------------------------------
[13]1637CONTAINS
[3]1638
[2715]1639   INTEGER FUNCTION lib_mpp_alloc(kumout)          ! Dummy function
1640      INTEGER, INTENT(in) ::   kumout
1641      lib_mpp_alloc = 0
1642   END FUNCTION lib_mpp_alloc
1643
[5407]1644   FUNCTION mynode( ldtxt, ldname, kumnam_ref, knumnam_cfg,  kumond , kstop, localComm ) RESULT (function_value)
[1579]1645      INTEGER, OPTIONAL            , INTENT(in   ) ::   localComm
[3764]1646      CHARACTER(len=*),DIMENSION(:) ::   ldtxt
[5407]1647      CHARACTER(len=*) ::   ldname
[4314]1648      INTEGER ::   kumnam_ref, knumnam_cfg , kumond , kstop
[9570]1649      IF( PRESENT( localComm ) ) mpi_comm_oce = localComm
[5412]1650      function_value = 0
[1579]1651      IF( .FALSE. )   ldtxt(:) = 'never done'
[5407]1652      CALL ctl_opn( kumond, TRIM(ldname), 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. , 1 )
[13]1653   END FUNCTION mynode
[3]1654
[13]1655   SUBROUTINE mppsync                       ! Dummy routine
1656   END SUBROUTINE mppsync
[3]1657
[10425]1658   !!----------------------------------------------------------------------
1659   !!    ***  mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real  ***
1660   !!   
1661   !!----------------------------------------------------------------------
1662   !!
1663#  define OPERATION_MAX
1664#  define INTEGER_TYPE
1665#  define DIM_0d
1666#     define ROUTINE_ALLREDUCE           mppmax_int
1667#     include "mpp_allreduce_generic.h90"
1668#     undef ROUTINE_ALLREDUCE
1669#  undef DIM_0d
1670#  define DIM_1d
1671#     define ROUTINE_ALLREDUCE           mppmax_a_int
1672#     include "mpp_allreduce_generic.h90"
1673#     undef ROUTINE_ALLREDUCE
1674#  undef DIM_1d
1675#  undef INTEGER_TYPE
1676!
1677#  define REAL_TYPE
1678#  define DIM_0d
1679#     define ROUTINE_ALLREDUCE           mppmax_real
1680#     include "mpp_allreduce_generic.h90"
1681#     undef ROUTINE_ALLREDUCE
1682#  undef DIM_0d
1683#  define DIM_1d
1684#     define ROUTINE_ALLREDUCE           mppmax_a_real
1685#     include "mpp_allreduce_generic.h90"
1686#     undef ROUTINE_ALLREDUCE
1687#  undef DIM_1d
1688#  undef REAL_TYPE
1689#  undef OPERATION_MAX
1690   !!----------------------------------------------------------------------
1691   !!    ***  mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real  ***
1692   !!   
1693   !!----------------------------------------------------------------------
1694   !!
1695#  define OPERATION_MIN
1696#  define INTEGER_TYPE
1697#  define DIM_0d
1698#     define ROUTINE_ALLREDUCE           mppmin_int
1699#     include "mpp_allreduce_generic.h90"
1700#     undef ROUTINE_ALLREDUCE
1701#  undef DIM_0d
1702#  define DIM_1d
1703#     define ROUTINE_ALLREDUCE           mppmin_a_int
1704#     include "mpp_allreduce_generic.h90"
1705#     undef ROUTINE_ALLREDUCE
1706#  undef DIM_1d
1707#  undef INTEGER_TYPE
1708!
1709#  define REAL_TYPE
1710#  define DIM_0d
1711#     define ROUTINE_ALLREDUCE           mppmin_real
1712#     include "mpp_allreduce_generic.h90"
1713#     undef ROUTINE_ALLREDUCE
1714#  undef DIM_0d
1715#  define DIM_1d
1716#     define ROUTINE_ALLREDUCE           mppmin_a_real
1717#     include "mpp_allreduce_generic.h90"
1718#     undef ROUTINE_ALLREDUCE
1719#  undef DIM_1d
1720#  undef REAL_TYPE
1721#  undef OPERATION_MIN
[3]1722
[10425]1723   !!----------------------------------------------------------------------
1724   !!    ***  mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real  ***
1725   !!   
1726   !!   Global sum of 1D array or a variable (integer, real or complex)
1727   !!----------------------------------------------------------------------
1728   !!
1729#  define OPERATION_SUM
1730#  define INTEGER_TYPE
1731#  define DIM_0d
1732#     define ROUTINE_ALLREDUCE           mppsum_int
1733#     include "mpp_allreduce_generic.h90"
1734#     undef ROUTINE_ALLREDUCE
1735#  undef DIM_0d
1736#  define DIM_1d
1737#     define ROUTINE_ALLREDUCE           mppsum_a_int
1738#     include "mpp_allreduce_generic.h90"
1739#     undef ROUTINE_ALLREDUCE
1740#  undef DIM_1d
1741#  undef INTEGER_TYPE
1742!
1743#  define REAL_TYPE
1744#  define DIM_0d
1745#     define ROUTINE_ALLREDUCE           mppsum_real
1746#     include "mpp_allreduce_generic.h90"
1747#     undef ROUTINE_ALLREDUCE
1748#  undef DIM_0d
1749#  define DIM_1d
1750#     define ROUTINE_ALLREDUCE           mppsum_a_real
1751#     include "mpp_allreduce_generic.h90"
1752#     undef ROUTINE_ALLREDUCE
1753#  undef DIM_1d
1754#  undef REAL_TYPE
1755#  undef OPERATION_SUM
[3]1756
[10425]1757#  define OPERATION_SUM_DD
1758#  define COMPLEX_TYPE
1759#  define DIM_0d
1760#     define ROUTINE_ALLREDUCE           mppsum_realdd
1761#     include "mpp_allreduce_generic.h90"
1762#     undef ROUTINE_ALLREDUCE
1763#  undef DIM_0d
1764#  define DIM_1d
1765#     define ROUTINE_ALLREDUCE           mppsum_a_realdd
1766#     include "mpp_allreduce_generic.h90"
1767#     undef ROUTINE_ALLREDUCE
1768#  undef DIM_1d
1769#  undef COMPLEX_TYPE
1770#  undef OPERATION_SUM_DD
[3]1771
[10425]1772   !!----------------------------------------------------------------------
1773   !!    ***  mpp_minloc2d, mpp_minloc3d, mpp_maxloc2d, mpp_maxloc3d
1774   !!   
1775   !!----------------------------------------------------------------------
1776   !!
1777#  define OPERATION_MINLOC
1778#  define DIM_2d
1779#     define ROUTINE_LOC           mpp_minloc2d
1780#     include "mpp_loc_generic.h90"
1781#     undef ROUTINE_LOC
1782#  undef DIM_2d
1783#  define DIM_3d
1784#     define ROUTINE_LOC           mpp_minloc3d
1785#     include "mpp_loc_generic.h90"
1786#     undef ROUTINE_LOC
1787#  undef DIM_3d
1788#  undef OPERATION_MINLOC
[2480]1789
[10425]1790#  define OPERATION_MAXLOC
1791#  define DIM_2d
1792#     define ROUTINE_LOC           mpp_maxloc2d
1793#     include "mpp_loc_generic.h90"
1794#     undef ROUTINE_LOC
1795#  undef DIM_2d
1796#  define DIM_3d
1797#     define ROUTINE_LOC           mpp_maxloc3d
1798#     include "mpp_loc_generic.h90"
1799#     undef ROUTINE_LOC
1800#  undef DIM_3d
1801#  undef OPERATION_MAXLOC
[13]1802
[10425]1803   SUBROUTINE mpp_delay_sum( cdname, cdelay, y_in, pout, ldlast, kcom )
1804      CHARACTER(len=*), INTENT(in   )               ::   cdname  ! name of the calling subroutine
1805      CHARACTER(len=*), INTENT(in   )               ::   cdelay  ! name (used as id) of the delayed operation
1806      COMPLEX(wp),      INTENT(in   ), DIMENSION(:) ::   y_in
1807      REAL(wp),         INTENT(  out), DIMENSION(:) ::   pout
1808      LOGICAL,          INTENT(in   )               ::   ldlast  ! true if this is the last time we call this routine
1809      INTEGER,          INTENT(in   ), OPTIONAL     ::   kcom
1810      !
1811      pout(:) = REAL(y_in(:), wp)
1812   END SUBROUTINE mpp_delay_sum
[3764]1813
[10425]1814   SUBROUTINE mpp_delay_max( cdname, cdelay, p_in, pout, ldlast, kcom )
1815      CHARACTER(len=*), INTENT(in   )               ::   cdname  ! name of the calling subroutine
1816      CHARACTER(len=*), INTENT(in   )               ::   cdelay  ! name (used as id) of the delayed operation
1817      REAL(wp),         INTENT(in   ), DIMENSION(:) ::   p_in
1818      REAL(wp),         INTENT(  out), DIMENSION(:) ::   pout
1819      LOGICAL,          INTENT(in   )               ::   ldlast  ! true if this is the last time we call this routine
1820      INTEGER,          INTENT(in   ), OPTIONAL     ::   kcom
1821      !
1822      pout(:) = p_in(:)
1823   END SUBROUTINE mpp_delay_max
[3294]1824
[10425]1825   SUBROUTINE mpp_delay_rcv( kid )
1826      INTEGER,INTENT(in   )      ::  kid 
1827      WRITE(*,*) 'mpp_delay_rcv: You should not have seen this print! error?', kid
1828   END SUBROUTINE mpp_delay_rcv
1829   
1830   SUBROUTINE mppstop( ldfinal, ld_force_abort )
1831      LOGICAL, OPTIONAL, INTENT(in) :: ldfinal    ! source process number
1832      LOGICAL, OPTIONAL, INTENT(in) :: ld_force_abort    ! source process number
[3799]1833      STOP      ! non MPP case, just stop the run
[51]1834   END SUBROUTINE mppstop
1835
[2715]1836   SUBROUTINE mpp_ini_znl( knum )
1837      INTEGER :: knum
1838      WRITE(*,*) 'mpp_ini_znl: You should not have seen this print! error?', knum
[1345]1839   END SUBROUTINE mpp_ini_znl
1840
[1344]1841   SUBROUTINE mpp_comm_free( kcom )
[869]1842      INTEGER :: kcom
[1344]1843      WRITE(*,*) 'mpp_comm_free: You should not have seen this print! error?', kcom
[869]1844   END SUBROUTINE mpp_comm_free
[9019]1845   
[3]1846#endif
[2715]1847
[13]1848   !!----------------------------------------------------------------------
[4147]1849   !!   All cases:         ctl_stop, ctl_warn, get_unit, ctl_opn, ctl_nam   routines
[2715]1850   !!----------------------------------------------------------------------
1851
1852   SUBROUTINE ctl_stop( cd1, cd2, cd3, cd4, cd5 ,   &
1853      &                 cd6, cd7, cd8, cd9, cd10 )
1854      !!----------------------------------------------------------------------
1855      !!                  ***  ROUTINE  stop_opa  ***
1856      !!
[3764]1857      !! ** Purpose :   print in ocean.outpput file a error message and
[2715]1858      !!                increment the error number (nstop) by one.
1859      !!----------------------------------------------------------------------
1860      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd1, cd2, cd3, cd4, cd5
1861      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd6, cd7, cd8, cd9, cd10
1862      !!----------------------------------------------------------------------
1863      !
[3764]1864      nstop = nstop + 1
[10425]1865
1866      ! force to open ocean.output file
1867      IF( numout == 6 ) CALL ctl_opn( numout, 'ocean.output', 'APPEND', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. )
1868       
1869      WRITE(numout,cform_err)
1870      IF( PRESENT(cd1 ) )   WRITE(numout,*) TRIM(cd1)
1871      IF( PRESENT(cd2 ) )   WRITE(numout,*) TRIM(cd2)
1872      IF( PRESENT(cd3 ) )   WRITE(numout,*) TRIM(cd3)
1873      IF( PRESENT(cd4 ) )   WRITE(numout,*) TRIM(cd4)
1874      IF( PRESENT(cd5 ) )   WRITE(numout,*) TRIM(cd5)
1875      IF( PRESENT(cd6 ) )   WRITE(numout,*) TRIM(cd6)
1876      IF( PRESENT(cd7 ) )   WRITE(numout,*) TRIM(cd7)
1877      IF( PRESENT(cd8 ) )   WRITE(numout,*) TRIM(cd8)
1878      IF( PRESENT(cd9 ) )   WRITE(numout,*) TRIM(cd9)
1879      IF( PRESENT(cd10) )   WRITE(numout,*) TRIM(cd10)
1880
[2715]1881                               CALL FLUSH(numout    )
1882      IF( numstp     /= -1 )   CALL FLUSH(numstp    )
[9019]1883      IF( numrun     /= -1 )   CALL FLUSH(numrun    )
[2715]1884      IF( numevo_ice /= -1 )   CALL FLUSH(numevo_ice)
1885      !
1886      IF( cd1 == 'STOP' ) THEN
[10425]1887         WRITE(numout,*)  'huge E-R-R-O-R : immediate stop'
1888         CALL mppstop(ld_force_abort = .true.)
[2715]1889      ENDIF
1890      !
1891   END SUBROUTINE ctl_stop
1892
1893
1894   SUBROUTINE ctl_warn( cd1, cd2, cd3, cd4, cd5,   &
1895      &                 cd6, cd7, cd8, cd9, cd10 )
1896      !!----------------------------------------------------------------------
1897      !!                  ***  ROUTINE  stop_warn  ***
1898      !!
[3764]1899      !! ** Purpose :   print in ocean.outpput file a error message and
[2715]1900      !!                increment the warning number (nwarn) by one.
1901      !!----------------------------------------------------------------------
1902      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd1, cd2, cd3, cd4, cd5
1903      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd6, cd7, cd8, cd9, cd10
1904      !!----------------------------------------------------------------------
[3764]1905      !
1906      nwarn = nwarn + 1
[2715]1907      IF(lwp) THEN
1908         WRITE(numout,cform_war)
[10425]1909         IF( PRESENT(cd1 ) ) WRITE(numout,*) TRIM(cd1)
1910         IF( PRESENT(cd2 ) ) WRITE(numout,*) TRIM(cd2)
1911         IF( PRESENT(cd3 ) ) WRITE(numout,*) TRIM(cd3)
1912         IF( PRESENT(cd4 ) ) WRITE(numout,*) TRIM(cd4)
1913         IF( PRESENT(cd5 ) ) WRITE(numout,*) TRIM(cd5)
1914         IF( PRESENT(cd6 ) ) WRITE(numout,*) TRIM(cd6)
1915         IF( PRESENT(cd7 ) ) WRITE(numout,*) TRIM(cd7)
1916         IF( PRESENT(cd8 ) ) WRITE(numout,*) TRIM(cd8)
1917         IF( PRESENT(cd9 ) ) WRITE(numout,*) TRIM(cd9)
1918         IF( PRESENT(cd10) ) WRITE(numout,*) TRIM(cd10)
[2715]1919      ENDIF
1920      CALL FLUSH(numout)
1921      !
1922   END SUBROUTINE ctl_warn
1923
1924
1925   SUBROUTINE ctl_opn( knum, cdfile, cdstat, cdform, cdacce, klengh, kout, ldwp, karea )
1926      !!----------------------------------------------------------------------
1927      !!                  ***  ROUTINE ctl_opn  ***
1928      !!
1929      !! ** Purpose :   Open file and check if required file is available.
1930      !!
1931      !! ** Method  :   Fortan open
1932      !!----------------------------------------------------------------------
1933      INTEGER          , INTENT(  out) ::   knum      ! logical unit to open
1934      CHARACTER(len=*) , INTENT(in   ) ::   cdfile    ! file name to open
1935      CHARACTER(len=*) , INTENT(in   ) ::   cdstat    ! disposition specifier
1936      CHARACTER(len=*) , INTENT(in   ) ::   cdform    ! formatting specifier
1937      CHARACTER(len=*) , INTENT(in   ) ::   cdacce    ! access specifier
1938      INTEGER          , INTENT(in   ) ::   klengh    ! record length
1939      INTEGER          , INTENT(in   ) ::   kout      ! number of logical units for write
1940      LOGICAL          , INTENT(in   ) ::   ldwp      ! boolean term for print
1941      INTEGER, OPTIONAL, INTENT(in   ) ::   karea     ! proc number
[5836]1942      !
[2715]1943      CHARACTER(len=80) ::   clfile
1944      INTEGER           ::   iost
1945      !!----------------------------------------------------------------------
[5836]1946      !
[2715]1947      ! adapt filename
1948      ! ----------------
1949      clfile = TRIM(cdfile)
1950      IF( PRESENT( karea ) ) THEN
1951         IF( karea > 1 )   WRITE(clfile, "(a,'_',i4.4)") TRIM(clfile), karea-1
1952      ENDIF
1953#if defined key_agrif
1954      IF( .NOT. Agrif_Root() )   clfile = TRIM(Agrif_CFixed())//'_'//TRIM(clfile)
1955      knum=Agrif_Get_Unit()
1956#else
1957      knum=get_unit()
1958#endif
[10425]1959      IF( TRIM(cdfile) == '/dev/null' )   clfile = TRIM(cdfile)   ! force the use of /dev/null
[5836]1960      !
[2715]1961      iost=0
[10425]1962      IF( cdacce(1:6) == 'DIRECT' )  THEN         ! cdacce has always more than 6 characters
1963         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat, RECL=klengh         , ERR=100, IOSTAT=iost )
1964      ELSE IF( TRIM(cdstat) == 'APPEND' )  THEN   ! cdstat can have less than 6 characters
1965         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS='UNKNOWN', POSITION='APPEND', ERR=100, IOSTAT=iost )
[2715]1966      ELSE
[10425]1967         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat                      , ERR=100, IOSTAT=iost )
[2715]1968      ENDIF
[10425]1969      IF( iost /= 0 .AND. TRIM(clfile) == '/dev/null' ) &   ! for windows
1970         &  OPEN(UNIT=knum,FILE='NUL', FORM=cdform, ACCESS=cdacce, STATUS=cdstat                      , ERR=100, IOSTAT=iost )   
[2715]1971      IF( iost == 0 ) THEN
1972         IF(ldwp) THEN
[10425]1973            WRITE(kout,*) '     file   : ', TRIM(clfile),' open ok'
[2715]1974            WRITE(kout,*) '     unit   = ', knum
1975            WRITE(kout,*) '     status = ', cdstat
1976            WRITE(kout,*) '     form   = ', cdform
1977            WRITE(kout,*) '     access = ', cdacce
1978            WRITE(kout,*)
1979         ENDIF
1980      ENDIF
1981100   CONTINUE
1982      IF( iost /= 0 ) THEN
1983         IF(ldwp) THEN
1984            WRITE(kout,*)
[10425]1985            WRITE(kout,*) ' ===>>>> : bad opening file: ', TRIM(clfile)
[2715]1986            WRITE(kout,*) ' =======   ===  '
1987            WRITE(kout,*) '           unit   = ', knum
1988            WRITE(kout,*) '           status = ', cdstat
1989            WRITE(kout,*) '           form   = ', cdform
1990            WRITE(kout,*) '           access = ', cdacce
1991            WRITE(kout,*) '           iostat = ', iost
1992            WRITE(kout,*) '           we stop. verify the file '
1993            WRITE(kout,*)
[9438]1994         ELSE  !!! Force writing to make sure we get the information - at least once - in this violent STOP!!
1995            WRITE(*,*)
[10425]1996            WRITE(*,*) ' ===>>>> : bad opening file: ', TRIM(clfile)
[9438]1997            WRITE(*,*) ' =======   ===  '
1998            WRITE(*,*) '           unit   = ', knum
1999            WRITE(*,*) '           status = ', cdstat
2000            WRITE(*,*) '           form   = ', cdform
2001            WRITE(*,*) '           access = ', cdacce
2002            WRITE(*,*) '           iostat = ', iost
2003            WRITE(*,*) '           we stop. verify the file '
2004            WRITE(*,*)
[2715]2005         ENDIF
[9019]2006         CALL FLUSH( kout ) 
[2715]2007         STOP 'ctl_opn bad opening'
2008      ENDIF
[5836]2009      !
[2715]2010   END SUBROUTINE ctl_opn
2011
[5836]2012
[4147]2013   SUBROUTINE ctl_nam ( kios, cdnam, ldwp )
2014      !!----------------------------------------------------------------------
2015      !!                  ***  ROUTINE ctl_nam  ***
2016      !!
2017      !! ** Purpose :   Informations when error while reading a namelist
2018      !!
2019      !! ** Method  :   Fortan open
2020      !!----------------------------------------------------------------------
[5836]2021      INTEGER         , INTENT(inout) ::   kios    ! IO status after reading the namelist
2022      CHARACTER(len=*), INTENT(in   ) ::   cdnam   ! group name of namelist for which error occurs
[7646]2023      CHARACTER(len=5)                ::   clios   ! string to convert iostat in character for print
[5836]2024      LOGICAL         , INTENT(in   ) ::   ldwp    ! boolean term for print
[4147]2025      !!----------------------------------------------------------------------
[5836]2026      !
[7646]2027      WRITE (clios, '(I5.0)')   kios
[4147]2028      IF( kios < 0 ) THEN         
[5836]2029         CALL ctl_warn( 'end of record or file while reading namelist '   &
2030            &           // TRIM(cdnam) // ' iostat = ' // TRIM(clios) )
[4147]2031      ENDIF
[5836]2032      !
[4147]2033      IF( kios > 0 ) THEN
[5836]2034         CALL ctl_stop( 'misspelled variable in namelist '   &
2035            &           // TRIM(cdnam) // ' iostat = ' // TRIM(clios) )
[4147]2036      ENDIF
2037      kios = 0
2038      RETURN
[5836]2039      !
[4147]2040   END SUBROUTINE ctl_nam
2041
[5836]2042
[2715]2043   INTEGER FUNCTION get_unit()
2044      !!----------------------------------------------------------------------
2045      !!                  ***  FUNCTION  get_unit  ***
2046      !!
2047      !! ** Purpose :   return the index of an unused logical unit
2048      !!----------------------------------------------------------------------
[3764]2049      LOGICAL :: llopn
[2715]2050      !!----------------------------------------------------------------------
2051      !
2052      get_unit = 15   ! choose a unit that is big enough then it is not already used in NEMO
2053      llopn = .TRUE.
2054      DO WHILE( (get_unit < 998) .AND. llopn )
2055         get_unit = get_unit + 1
2056         INQUIRE( unit = get_unit, opened = llopn )
2057      END DO
2058      IF( (get_unit == 999) .AND. llopn ) THEN
2059         CALL ctl_stop( 'get_unit: All logical units until 999 are used...' )
2060         get_unit = -1
2061      ENDIF
2062      !
2063   END FUNCTION get_unit
2064
2065   !!----------------------------------------------------------------------
[3]2066END MODULE lib_mpp
Note: See TracBrowser for help on using the repository browser.