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 @ 10300

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

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 2b: add waiting time for mpp_min/max/sum, see #2133

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