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/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/LBC – NEMO

source: NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/LBC/lib_mpp.F90 @ 10330

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

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 6: find the best mpi domain decomposition, see #2133

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