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/UKMO/r14075_coupling_sequence/src/OCE/LBC – NEMO

source: NEMO/branches/UKMO/r14075_coupling_sequence/src/OCE/LBC/lib_mpp.F90 @ 15814

Last change on this file since 15814 was 15424, checked in by jcastill, 3 years ago

Compiling code

File size: 57.9 KB
Line 
1MODULE lib_mpp
2   !!======================================================================
3   !!                       ***  MODULE  lib_mpp  ***
4   !! Ocean numerics:  massively parallel processing library
5   !!=====================================================================
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
10   !!   NEMO     1.0  !  2003  (J.M. Molines, G. Madec)  F90, free form
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
19   !!            3.2  !  2009  (O. Marti)    add mpp_ini_znl
20   !!            4.0  !  2011  (G. Madec)  move ctl_ routines from in_out_manager
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
23   !!            3.5  !  2013  (S.Mocavero, I.Epicoco - CMCC) north fold optimizations
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
27   !!----------------------------------------------------------------------
28
29   !!----------------------------------------------------------------------
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   !!----------------------------------------------------------------------
35   !!----------------------------------------------------------------------
36   !!   mpp_start     : get local communicator its size and rank
37   !!   mpp_lnk       : interface (defined in lbclnk) for message passing of 2d or 3d arrays (mpp_lnk_2d, mpp_lnk_3d)
38   !!   mpp_lnk_icb   : interface for message passing of 2d arrays with extra halo for icebergs (mpp_lnk_2d_icb)
39   !!   mpprecv       :
40   !!   mppsend       :
41   !!   mppscatter    :
42   !!   mppgather     :
43   !!   mpp_min       : generic interface for mppmin_int , mppmin_a_int , mppmin_real, mppmin_a_real
44   !!   mpp_max       : generic interface for mppmax_int , mppmax_a_int , mppmax_real, mppmax_a_real
45   !!   mpp_sum       : generic interface for mppsum_int , mppsum_a_int , mppsum_real, mppsum_a_real
46   !!   mpp_minloc    :
47   !!   mpp_maxloc    :
48   !!   mppsync       :
49   !!   mppstop       :
50   !!   mpp_ini_north : initialisation of north fold
51   !!   mpp_lbc_north_icb : alternative to mpp_nfd for extra outer halo with icebergs
52   !!----------------------------------------------------------------------
53   USE dom_oce        ! ocean space and time domain
54   USE in_out_manager ! I/O manager
55
56   IMPLICIT NONE
57   PRIVATE
58   !
59   PUBLIC   ctl_stop, ctl_warn, ctl_opn, ctl_nam
60   PUBLIC   mpp_start, mppstop, mppsync, mpp_comm_free
61   PUBLIC   mpp_ini_north
62   PUBLIC   mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc
63   PUBLIC   mpp_delay_max, mpp_delay_sum, mpp_delay_rcv
64   PUBLIC   mppscatter, mppgather
65   PUBLIC   mpp_ini_znl
66   PUBLIC   mppsend, mpprecv                          ! needed by TAM and ICB routines
67   PUBLIC   mpp_report
68   PUBLIC   tic_tac
69#if ! defined key_mpp_mpi
70   PUBLIC MPI_Wtime
71#endif
72   
73   !! * Interfaces
74   !! define generic interface for these routine as they are called sometimes
75   !! with scalar arguments instead of array arguments, which causes problems
76   !! for the compilation on AIX system as well as NEC and SGI. Ok on COMPACQ
77   INTERFACE mpp_min
78      MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real
79   END INTERFACE
80   INTERFACE mpp_max
81      MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real
82   END INTERFACE
83   INTERFACE mpp_sum
84      MODULE PROCEDURE mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real,   &
85         &             mppsum_realdd, mppsum_a_realdd
86   END INTERFACE
87   INTERFACE mpp_minloc
88      MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d
89   END INTERFACE
90   INTERFACE mpp_maxloc
91      MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d
92   END INTERFACE
93
94   !! ========================= !!
95   !!  MPI  variable definition !!
96   !! ========================= !!
97#if   defined key_mpp_mpi
98!$AGRIF_DO_NOT_TREAT
99   INCLUDE 'mpif.h'
100!$AGRIF_END_DO_NOT_TREAT
101   LOGICAL, PUBLIC, PARAMETER ::   lk_mpp = .TRUE.    !: mpp flag
102#else   
103   INTEGER, PUBLIC, PARAMETER ::   MPI_STATUS_SIZE = 1
104   INTEGER, PUBLIC, PARAMETER ::   MPI_DOUBLE_PRECISION = 8
105   LOGICAL, PUBLIC, PARAMETER ::   lk_mpp = .FALSE.    !: mpp flag
106#endif
107
108   INTEGER, PARAMETER         ::   nprocmax = 2**10   ! maximun dimension (required to be a power of 2)
109
110   INTEGER, PUBLIC ::   mppsize        ! number of process
111   INTEGER, PUBLIC ::   mpprank        ! process number  [ 0 - size-1 ]
112!$AGRIF_DO_NOT_TREAT
113   INTEGER, PUBLIC ::   mpi_comm_oce   ! opa local communicator
114!$AGRIF_END_DO_NOT_TREAT
115
116   INTEGER :: MPI_SUMDD
117
118   ! variables used for zonal integration
119   INTEGER, PUBLIC ::   ncomm_znl       !: communicator made by the processors on the same zonal average
120   LOGICAL, PUBLIC ::   l_znl_root      !: True on the 'left'most processor on the same row
121   INTEGER         ::   ngrp_znl        !  group ID for the znl processors
122   INTEGER         ::   ndim_rank_znl   !  number of processors on the same zonal average
123   INTEGER, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_znl  ! dimension ndim_rank_znl, number of the procs into the same znl domain
124
125   ! North fold condition in mpp_mpi with jpni > 1 (PUBLIC for TAM)
126   INTEGER, PUBLIC ::   ngrp_world        !: group ID for the world processors
127   INTEGER, PUBLIC ::   ngrp_opa          !: group ID for the opa processors
128   INTEGER, PUBLIC ::   ngrp_north        !: group ID for the northern processors (to be fold)
129   INTEGER, PUBLIC ::   ncomm_north       !: communicator made by the processors belonging to ngrp_north
130   INTEGER, PUBLIC ::   ndim_rank_north   !: number of 'sea' processor in the northern line (can be /= jpni !)
131   INTEGER, PUBLIC ::   njmppmax          !: value of njmpp for the processors of the northern line
132   INTEGER, PUBLIC ::   north_root        !: number (in the comm_opa) of proc 0 in the northern comm
133   INTEGER, PUBLIC, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_north   !: dimension ndim_rank_north
134
135   ! Communications summary report
136   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_lbc                   !: names of lbc_lnk calling routines
137   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_glb                   !: names of global comm calling routines
138   CHARACTER(len=128), DIMENSION(:), ALLOCATABLE ::   crname_dlg                   !: names of delayed global comm calling routines
139   INTEGER, PUBLIC                               ::   ncom_stp = 0                 !: copy of time step # istp
140   INTEGER, PUBLIC                               ::   ncom_fsbc = 1                !: copy of sbc time step # nn_fsbc
141   INTEGER, PUBLIC                               ::   ncom_dttrc = 1               !: copy of top time step # nn_dttrc
142   INTEGER, PUBLIC                               ::   ncom_freq                    !: frequency of comm diagnostic
143   INTEGER, PUBLIC , DIMENSION(:,:), ALLOCATABLE ::   ncomm_sequence               !: size of communicated arrays (halos)
144   INTEGER, PARAMETER, PUBLIC                    ::   ncom_rec_max = 5000          !: max number of communication record
145   INTEGER, PUBLIC                               ::   n_sequence_lbc = 0           !: # of communicated arraysvia lbc
146   INTEGER, PUBLIC                               ::   n_sequence_glb = 0           !: # of global communications
147   INTEGER, PUBLIC                               ::   n_sequence_dlg = 0           !: # of delayed global communications
148   INTEGER, PUBLIC                               ::   numcom = -1                  !: logical unit for communicaton report
149   LOGICAL, PUBLIC                               ::   l_full_nf_update = .TRUE.    !: logical for a full (2lines) update of bc at North fold report
150   INTEGER,                    PARAMETER, PUBLIC ::   nbdelay = 2       !: number of delayed operations
151   !: name (used as id) of allreduce-delayed operations
152   ! Warning: we must use the same character length in an array constructor (at least for gcc compiler)
153   CHARACTER(len=32), DIMENSION(nbdelay), PUBLIC ::   c_delaylist = (/ 'cflice', 'fwb   ' /)
154   !: component name where the allreduce-delayed operation is performed
155   CHARACTER(len=3),  DIMENSION(nbdelay), PUBLIC ::   c_delaycpnt = (/ 'ICE'   , 'OCE' /)
156   TYPE, PUBLIC ::   DELAYARR
157      REAL(   wp), POINTER, DIMENSION(:) ::  z1d => NULL()
158      COMPLEX(wp), POINTER, DIMENSION(:) ::  y1d => NULL()
159   END TYPE DELAYARR
160   TYPE( DELAYARR ), DIMENSION(nbdelay), PUBLIC, SAVE  ::   todelay         !: must have SAVE for default initialization of DELAYARR
161   INTEGER,          DIMENSION(nbdelay), PUBLIC        ::   ndelayid = -1   !: mpi request id of the delayed operations
162
163   ! timing summary report
164   REAL(wp), DIMENSION(2), PUBLIC ::  waiting_time = 0._wp
165   REAL(wp)              , PUBLIC ::  compute_time = 0._wp, elapsed_time = 0._wp
166   
167   REAL(wp), DIMENSION(:), ALLOCATABLE, SAVE ::   tampon   ! buffer in case of bsend
168
169   LOGICAL, PUBLIC ::   ln_nnogather                !: namelist control of northfold comms
170   LOGICAL, PUBLIC ::   l_north_nogather = .FALSE.  !: internal control of northfold comms
171   
172   !!----------------------------------------------------------------------
173   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
174   !! $Id$
175   !! Software governed by the CeCILL license (see ./LICENSE)
176   !!----------------------------------------------------------------------
177CONTAINS
178
179   SUBROUTINE mpp_start( localComm )
180      !!----------------------------------------------------------------------
181      !!                  ***  routine mpp_start  ***
182      !!
183      !! ** Purpose :   get mpi_comm_oce, mpprank and mppsize
184      !!----------------------------------------------------------------------
185      INTEGER         , OPTIONAL   , INTENT(in   ) ::   localComm    !
186      !
187      INTEGER ::   ierr
188      LOGICAL ::   llmpi_init
189      !!----------------------------------------------------------------------
190#if defined key_mpp_mpi
191      !
192      CALL mpi_initialized ( llmpi_init, ierr )
193      IF( ierr /= MPI_SUCCESS ) CALL ctl_stop( 'STOP', ' lib_mpp: Error in routine mpi_initialized' )
194
195      IF( .NOT. llmpi_init ) THEN
196         IF( PRESENT(localComm) ) THEN
197            WRITE(ctmp1,*) ' lib_mpp: You cannot provide a local communicator '
198            WRITE(ctmp2,*) '          without calling MPI_Init before ! '
199            CALL ctl_stop( 'STOP', ctmp1, ctmp2 )
200         ENDIF
201         CALL mpi_init( ierr )
202         IF( ierr /= MPI_SUCCESS ) CALL ctl_stop( 'STOP', ' lib_mpp: Error in routine mpi_init' )
203      ENDIF
204       
205      IF( PRESENT(localComm) ) THEN
206         IF( Agrif_Root() ) THEN
207            mpi_comm_oce = localComm
208         ENDIF
209      ELSE
210         CALL mpi_comm_dup( mpi_comm_world, mpi_comm_oce, ierr)
211         IF( ierr /= MPI_SUCCESS ) CALL ctl_stop( 'STOP', ' lib_mpp: Error in routine mpi_comm_dup' )
212      ENDIF
213
214# if defined key_agrif
215      IF( Agrif_Root() ) THEN
216         CALL Agrif_MPI_Init(mpi_comm_oce)
217      ELSE
218         CALL Agrif_MPI_set_grid_comm(mpi_comm_oce)
219      ENDIF
220# endif
221
222      CALL mpi_comm_rank( mpi_comm_oce, mpprank, ierr )
223      CALL mpi_comm_size( mpi_comm_oce, mppsize, ierr )
224      !
225      CALL MPI_OP_CREATE(DDPDD_MPI, .TRUE., MPI_SUMDD, ierr)
226      !
227#else
228      IF( PRESENT( localComm ) ) mpi_comm_oce = localComm
229      mppsize = 1
230      mpprank = 0
231#endif
232   END SUBROUTINE mpp_start
233
234
235   SUBROUTINE mppsend( ktyp, pmess, kbytes, kdest, md_req )
236      !!----------------------------------------------------------------------
237      !!                  ***  routine mppsend  ***
238      !!
239      !! ** Purpose :   Send messag passing array
240      !!
241      !!----------------------------------------------------------------------
242      REAL(wp), INTENT(inout) ::   pmess(*)   ! array of real
243      INTEGER , INTENT(in   ) ::   kbytes     ! size of the array pmess
244      INTEGER , INTENT(in   ) ::   kdest      ! receive process number
245      INTEGER , INTENT(in   ) ::   ktyp       ! tag of the message
246      INTEGER , INTENT(in   ) ::   md_req     ! argument for isend
247      !!
248      INTEGER ::   iflag
249      !!----------------------------------------------------------------------
250      !
251#if defined key_mpp_mpi
252      CALL mpi_isend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_oce, md_req, iflag )
253#endif
254      !
255   END SUBROUTINE mppsend
256
257
258   SUBROUTINE mpprecv( ktyp, pmess, kbytes, ksource )
259      !!----------------------------------------------------------------------
260      !!                  ***  routine mpprecv  ***
261      !!
262      !! ** Purpose :   Receive messag passing array
263      !!
264      !!----------------------------------------------------------------------
265      REAL(wp), INTENT(inout) ::   pmess(*)   ! array of real
266      INTEGER , INTENT(in   ) ::   kbytes     ! suze of the array pmess
267      INTEGER , INTENT(in   ) ::   ktyp       ! Tag of the recevied message
268      INTEGER, OPTIONAL, INTENT(in) :: ksource    ! source process number
269      !!
270      INTEGER :: istatus(mpi_status_size)
271      INTEGER :: iflag
272      INTEGER :: use_source
273      !!----------------------------------------------------------------------
274      !
275#if defined key_mpp_mpi
276      ! If a specific process number has been passed to the receive call,
277      ! use that one. Default is to use mpi_any_source
278      use_source = mpi_any_source
279      IF( PRESENT(ksource) )   use_source = ksource
280      !
281      CALL mpi_recv( pmess, kbytes, mpi_double_precision, use_source, ktyp, mpi_comm_oce, istatus, iflag )
282#endif
283      !
284   END SUBROUTINE mpprecv
285
286
287   SUBROUTINE mppgather( ptab, kp, pio )
288      !!----------------------------------------------------------------------
289      !!                   ***  routine mppgather  ***
290      !!
291      !! ** Purpose :   Transfert between a local subdomain array and a work
292      !!     array which is distributed following the vertical level.
293      !!
294      !!----------------------------------------------------------------------
295      REAL(wp), DIMENSION(jpi,jpj)      , INTENT(in   ) ::   ptab   ! subdomain input array
296      INTEGER                           , INTENT(in   ) ::   kp     ! record length
297      REAL(wp), DIMENSION(jpi,jpj,jpnij), INTENT(  out) ::   pio    ! subdomain input array
298      !!
299      INTEGER :: itaille, ierror   ! temporary integer
300      !!---------------------------------------------------------------------
301      !
302      itaille = jpi * jpj
303#if defined key_mpp_mpi
304      CALL mpi_gather( ptab, itaille, mpi_double_precision, pio, itaille     ,   &
305         &                            mpi_double_precision, kp , mpi_comm_oce, ierror )
306#else
307      pio(:,:,1) = ptab(:,:)
308#endif
309      !
310   END SUBROUTINE mppgather
311
312
313   SUBROUTINE mppscatter( pio, kp, ptab )
314      !!----------------------------------------------------------------------
315      !!                  ***  routine mppscatter  ***
316      !!
317      !! ** Purpose :   Transfert between awork array which is distributed
318      !!      following the vertical level and the local subdomain array.
319      !!
320      !!----------------------------------------------------------------------
321      REAL(wp), DIMENSION(jpi,jpj,jpnij)  ::   pio    ! output array
322      INTEGER                             ::   kp     ! Tag (not used with MPI
323      REAL(wp), DIMENSION(jpi,jpj)        ::   ptab   ! subdomain array input
324      !!
325      INTEGER :: itaille, ierror   ! temporary integer
326      !!---------------------------------------------------------------------
327      !
328      itaille = jpi * jpj
329      !
330#if defined key_mpp_mpi
331      CALL mpi_scatter( pio, itaille, mpi_double_precision, ptab, itaille     ,   &
332         &                            mpi_double_precision, kp  , mpi_comm_oce, ierror )
333#else
334      ptab(:,:) = pio(:,:,1)
335#endif
336      !
337   END SUBROUTINE mppscatter
338
339   
340   SUBROUTINE mpp_delay_sum( cdname, cdelay, y_in, pout, ldlast, kcom )
341     !!----------------------------------------------------------------------
342      !!                   ***  routine mpp_delay_sum  ***
343      !!
344      !! ** Purpose :   performed delayed mpp_sum, the result is received on next call
345      !!
346      !!----------------------------------------------------------------------
347      CHARACTER(len=*), INTENT(in   )               ::   cdname  ! name of the calling subroutine
348      CHARACTER(len=*), INTENT(in   )               ::   cdelay  ! name (used as id) of the delayed operation
349      COMPLEX(wp),      INTENT(in   ), DIMENSION(:) ::   y_in
350      REAL(wp),         INTENT(  out), DIMENSION(:) ::   pout
351      LOGICAL,          INTENT(in   )               ::   ldlast  ! true if this is the last time we call this routine
352      INTEGER,          INTENT(in   ), OPTIONAL     ::   kcom
353      !!
354      INTEGER ::   ji, isz
355      INTEGER ::   idvar
356      INTEGER ::   ierr, ilocalcomm
357      COMPLEX(wp), ALLOCATABLE, DIMENSION(:) ::   ytmp
358      !!----------------------------------------------------------------------
359#if defined key_mpp_mpi
360      ilocalcomm = mpi_comm_oce
361      IF( PRESENT(kcom) )   ilocalcomm = kcom
362
363      isz = SIZE(y_in)
364     
365      IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ld_dlg = .TRUE. )
366
367      idvar = -1
368      DO ji = 1, nbdelay
369         IF( TRIM(cdelay) == TRIM(c_delaylist(ji)) ) idvar = ji
370      END DO
371      IF ( idvar == -1 )   CALL ctl_stop( 'STOP',' mpp_delay_sum : please add a new delayed exchange for '//TRIM(cdname) )
372
373      IF ( ndelayid(idvar) == 0 ) THEN         ! first call    with restart: %z1d defined in iom_delay_rst
374         !                                       --------------------------
375         IF ( SIZE(todelay(idvar)%z1d) /= isz ) THEN                  ! Check dimension coherence
376            IF(lwp) WRITE(numout,*) ' WARNING: the nb of delayed variables in restart file is not the model one'
377            DEALLOCATE(todelay(idvar)%z1d)
378            ndelayid(idvar) = -1                                      ! do as if we had no restart
379         ELSE
380            ALLOCATE(todelay(idvar)%y1d(isz))
381            todelay(idvar)%y1d(:) = CMPLX(todelay(idvar)%z1d(:), 0., wp)   ! create %y1d, complex variable needed by mpi_sumdd
382            ndelayid(idvar) = MPI_REQUEST_NULL                             ! initialised request to a valid value
383         END IF
384      ENDIF
385     
386      IF( ndelayid(idvar) == -1 ) THEN         ! first call without restart: define %y1d and %z1d from y_in with blocking allreduce
387         !                                       --------------------------
388         ALLOCATE(todelay(idvar)%z1d(isz), todelay(idvar)%y1d(isz))   ! allocate also %z1d as used for the restart
389         CALL mpi_allreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ierr )   ! get %y1d
390         ndelayid(idvar) = MPI_REQUEST_NULL
391      ENDIF
392
393      CALL mpp_delay_rcv( idvar )         ! make sure %z1d is received
394
395      ! send back pout from todelay(idvar)%z1d defined at previous call
396      pout(:) = todelay(idvar)%z1d(:)
397
398      ! send y_in into todelay(idvar)%y1d with a non-blocking communication
399# if defined key_mpi2
400      IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.)
401      CALL  mpi_allreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ierr )
402      ndelayid(idvar) = MPI_REQUEST_NULL
403      IF( ln_timing ) CALL tic_tac(.FALSE., ld_global = .TRUE.)
404# else
405      CALL mpi_iallreduce( y_in(:), todelay(idvar)%y1d(:), isz, MPI_DOUBLE_COMPLEX, mpi_sumdd, ilocalcomm, ndelayid(idvar), ierr )
406# endif
407#else
408      pout(:) = REAL(y_in(:), wp)
409#endif
410
411   END SUBROUTINE mpp_delay_sum
412
413   
414   SUBROUTINE mpp_delay_max( cdname, cdelay, p_in, pout, ldlast, kcom )
415      !!----------------------------------------------------------------------
416      !!                   ***  routine mpp_delay_max  ***
417      !!
418      !! ** Purpose :   performed delayed mpp_max, the result is received on next call
419      !!
420      !!----------------------------------------------------------------------
421      CHARACTER(len=*), INTENT(in   )                 ::   cdname  ! name of the calling subroutine
422      CHARACTER(len=*), INTENT(in   )                 ::   cdelay  ! name (used as id) of the delayed operation
423      REAL(wp),         INTENT(in   ), DIMENSION(:)   ::   p_in    !
424      REAL(wp),         INTENT(  out), DIMENSION(:)   ::   pout    !
425      LOGICAL,          INTENT(in   )                 ::   ldlast  ! true if this is the last time we call this routine
426      INTEGER,          INTENT(in   ), OPTIONAL       ::   kcom
427      !!
428      INTEGER ::   ji, isz
429      INTEGER ::   idvar
430      INTEGER ::   ierr, ilocalcomm
431      !!----------------------------------------------------------------------
432     
433#if defined key_mpp_mpi
434      ilocalcomm = mpi_comm_oce
435      IF( PRESENT(kcom) )   ilocalcomm = kcom
436
437      isz = SIZE(p_in)
438
439      IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ld_dlg = .TRUE. )
440
441      idvar = -1
442      DO ji = 1, nbdelay
443         IF( TRIM(cdelay) == TRIM(c_delaylist(ji)) ) idvar = ji
444      END DO
445      IF ( idvar == -1 )   CALL ctl_stop( 'STOP',' mpp_delay_max : please add a new delayed exchange for '//TRIM(cdname) )
446
447      IF ( ndelayid(idvar) == 0 ) THEN         ! first call    with restart: %z1d defined in iom_delay_rst
448         !                                       --------------------------
449         IF ( SIZE(todelay(idvar)%z1d) /= isz ) THEN                  ! Check dimension coherence
450            IF(lwp) WRITE(numout,*) ' WARNING: the nb of delayed variables in restart file is not the model one'
451            DEALLOCATE(todelay(idvar)%z1d)
452            ndelayid(idvar) = -1                                      ! do as if we had no restart
453         ELSE
454            ndelayid(idvar) = MPI_REQUEST_NULL
455         END IF
456      ENDIF
457
458      IF( ndelayid(idvar) == -1 ) THEN         ! first call without restart: define %z1d from p_in with a blocking allreduce
459         !                                       --------------------------
460         ALLOCATE(todelay(idvar)%z1d(isz))
461         CALL mpi_allreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ierr )   ! get %z1d
462         ndelayid(idvar) = MPI_REQUEST_NULL
463      ENDIF
464
465      CALL mpp_delay_rcv( idvar )         ! make sure %z1d is received
466
467      ! send back pout from todelay(idvar)%z1d defined at previous call
468      pout(:) = todelay(idvar)%z1d(:)
469
470      ! send p_in into todelay(idvar)%z1d with a non-blocking communication
471      ! (PM) Should we get rid of MPI2 option ? MPI3 was release in 2013. Who is still using MPI2 ?
472# if defined key_mpi2
473      IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.)
474      CALL  mpi_allreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ierr )
475      IF( ln_timing ) CALL tic_tac(.FALSE., ld_global = .TRUE.)
476# else
477      CALL mpi_iallreduce( p_in(:), todelay(idvar)%z1d(:), isz, MPI_DOUBLE_PRECISION, mpi_max, ilocalcomm, ndelayid(idvar), ierr )
478# endif
479#else
480      pout(:) = p_in(:)
481#endif
482
483   END SUBROUTINE mpp_delay_max
484
485   
486   SUBROUTINE mpp_delay_rcv( kid )
487      !!----------------------------------------------------------------------
488      !!                   ***  routine mpp_delay_rcv  ***
489      !!
490      !! ** Purpose :  force barrier for delayed mpp (needed for restart)
491      !!
492      !!----------------------------------------------------------------------
493      INTEGER,INTENT(in   )      ::  kid 
494      INTEGER ::   ierr
495      !!----------------------------------------------------------------------
496#if defined key_mpp_mpi
497      IF( ln_timing ) CALL tic_tac( .TRUE., ld_global = .TRUE.)
498      ! test on ndelayid(kid) useless as mpi_wait return immediatly if the request handle is MPI_REQUEST_NULL
499      CALL mpi_wait( ndelayid(kid), MPI_STATUS_IGNORE, ierr ) ! after this ndelayid(kid) = MPI_REQUEST_NULL
500      IF( ln_timing ) CALL tic_tac( .FALSE., ld_global = .TRUE.)
501      IF( ASSOCIATED(todelay(kid)%y1d) )   todelay(kid)%z1d(:) = REAL(todelay(kid)%y1d(:), wp)  ! define %z1d from %y1d
502#endif
503   END SUBROUTINE mpp_delay_rcv
504
505   
506   !!----------------------------------------------------------------------
507   !!    ***  mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real  ***
508   !!   
509   !!----------------------------------------------------------------------
510   !!
511#  define OPERATION_MAX
512#  define INTEGER_TYPE
513#  define DIM_0d
514#     define ROUTINE_ALLREDUCE           mppmax_int
515#     include "mpp_allreduce_generic.h90"
516#     undef ROUTINE_ALLREDUCE
517#  undef DIM_0d
518#  define DIM_1d
519#     define ROUTINE_ALLREDUCE           mppmax_a_int
520#     include "mpp_allreduce_generic.h90"
521#     undef ROUTINE_ALLREDUCE
522#  undef DIM_1d
523#  undef INTEGER_TYPE
524!
525#  define REAL_TYPE
526#  define DIM_0d
527#     define ROUTINE_ALLREDUCE           mppmax_real
528#     include "mpp_allreduce_generic.h90"
529#     undef ROUTINE_ALLREDUCE
530#  undef DIM_0d
531#  define DIM_1d
532#     define ROUTINE_ALLREDUCE           mppmax_a_real
533#     include "mpp_allreduce_generic.h90"
534#     undef ROUTINE_ALLREDUCE
535#  undef DIM_1d
536#  undef REAL_TYPE
537#  undef OPERATION_MAX
538   !!----------------------------------------------------------------------
539   !!    ***  mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real  ***
540   !!   
541   !!----------------------------------------------------------------------
542   !!
543#  define OPERATION_MIN
544#  define INTEGER_TYPE
545#  define DIM_0d
546#     define ROUTINE_ALLREDUCE           mppmin_int
547#     include "mpp_allreduce_generic.h90"
548#     undef ROUTINE_ALLREDUCE
549#  undef DIM_0d
550#  define DIM_1d
551#     define ROUTINE_ALLREDUCE           mppmin_a_int
552#     include "mpp_allreduce_generic.h90"
553#     undef ROUTINE_ALLREDUCE
554#  undef DIM_1d
555#  undef INTEGER_TYPE
556!
557#  define REAL_TYPE
558#  define DIM_0d
559#     define ROUTINE_ALLREDUCE           mppmin_real
560#     include "mpp_allreduce_generic.h90"
561#     undef ROUTINE_ALLREDUCE
562#  undef DIM_0d
563#  define DIM_1d
564#     define ROUTINE_ALLREDUCE           mppmin_a_real
565#     include "mpp_allreduce_generic.h90"
566#     undef ROUTINE_ALLREDUCE
567#  undef DIM_1d
568#  undef REAL_TYPE
569#  undef OPERATION_MIN
570
571   !!----------------------------------------------------------------------
572   !!    ***  mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real  ***
573   !!   
574   !!   Global sum of 1D array or a variable (integer, real or complex)
575   !!----------------------------------------------------------------------
576   !!
577#  define OPERATION_SUM
578#  define INTEGER_TYPE
579#  define DIM_0d
580#     define ROUTINE_ALLREDUCE           mppsum_int
581#     include "mpp_allreduce_generic.h90"
582#     undef ROUTINE_ALLREDUCE
583#  undef DIM_0d
584#  define DIM_1d
585#     define ROUTINE_ALLREDUCE           mppsum_a_int
586#     include "mpp_allreduce_generic.h90"
587#     undef ROUTINE_ALLREDUCE
588#  undef DIM_1d
589#  undef INTEGER_TYPE
590!
591#  define REAL_TYPE
592#  define DIM_0d
593#     define ROUTINE_ALLREDUCE           mppsum_real
594#     include "mpp_allreduce_generic.h90"
595#     undef ROUTINE_ALLREDUCE
596#  undef DIM_0d
597#  define DIM_1d
598#     define ROUTINE_ALLREDUCE           mppsum_a_real
599#     include "mpp_allreduce_generic.h90"
600#     undef ROUTINE_ALLREDUCE
601#  undef DIM_1d
602#  undef REAL_TYPE
603#  undef OPERATION_SUM
604
605#  define OPERATION_SUM_DD
606#  define COMPLEX_TYPE
607#  define DIM_0d
608#     define ROUTINE_ALLREDUCE           mppsum_realdd
609#     include "mpp_allreduce_generic.h90"
610#     undef ROUTINE_ALLREDUCE
611#  undef DIM_0d
612#  define DIM_1d
613#     define ROUTINE_ALLREDUCE           mppsum_a_realdd
614#     include "mpp_allreduce_generic.h90"
615#     undef ROUTINE_ALLREDUCE
616#  undef DIM_1d
617#  undef COMPLEX_TYPE
618#  undef OPERATION_SUM_DD
619
620   !!----------------------------------------------------------------------
621   !!    ***  mpp_minloc2d, mpp_minloc3d, mpp_maxloc2d, mpp_maxloc3d
622   !!   
623   !!----------------------------------------------------------------------
624   !!
625#  define OPERATION_MINLOC
626#  define DIM_2d
627#     define ROUTINE_LOC           mpp_minloc2d
628#     include "mpp_loc_generic.h90"
629#     undef ROUTINE_LOC
630#  undef DIM_2d
631#  define DIM_3d
632#     define ROUTINE_LOC           mpp_minloc3d
633#     include "mpp_loc_generic.h90"
634#     undef ROUTINE_LOC
635#  undef DIM_3d
636#  undef OPERATION_MINLOC
637
638#  define OPERATION_MAXLOC
639#  define DIM_2d
640#     define ROUTINE_LOC           mpp_maxloc2d
641#     include "mpp_loc_generic.h90"
642#     undef ROUTINE_LOC
643#  undef DIM_2d
644#  define DIM_3d
645#     define ROUTINE_LOC           mpp_maxloc3d
646#     include "mpp_loc_generic.h90"
647#     undef ROUTINE_LOC
648#  undef DIM_3d
649#  undef OPERATION_MAXLOC
650
651   SUBROUTINE mppsync()
652      !!----------------------------------------------------------------------
653      !!                  ***  routine mppsync  ***
654      !!
655      !! ** Purpose :   Massively parallel processors, synchroneous
656      !!
657      !!-----------------------------------------------------------------------
658      INTEGER :: ierror
659      !!-----------------------------------------------------------------------
660      !
661#if defined key_mpp_mpi
662      CALL mpi_barrier( mpi_comm_oce, ierror )
663#endif
664      !
665   END SUBROUTINE mppsync
666
667
668   SUBROUTINE mppstop( ld_abort ) 
669      !!----------------------------------------------------------------------
670      !!                  ***  routine mppstop  ***
671      !!
672      !! ** purpose :   Stop massively parallel processors method
673      !!
674      !!----------------------------------------------------------------------
675#if defined key_oasis3 
676      USE mod_oasis      ! coupling routines
677#endif
678
679      LOGICAL, OPTIONAL, INTENT(in) :: ld_abort    ! source process number
680      LOGICAL ::   ll_abort
681      INTEGER ::   info
682      !!----------------------------------------------------------------------
683      ll_abort = .FALSE.
684      IF( PRESENT(ld_abort) ) ll_abort = ld_abort
685      !
686#if defined key_mpp_mpi
687      IF(ll_abort) THEN
688         CALL mpi_abort( MPI_COMM_WORLD )
689      ELSE
690#if defined key_oasis3   
691         ! If we're trying to shut down cleanly then we need to consider the fact   
692         ! that this could be part of an MPMD configuration - we don't want to   
693         ! leave other components deadlocked.   
694         CALL oasis_abort(nproc,"mppstop","NEMO initiated abort")   
695#else
696         CALL mppsync
697         CALL mpi_finalize( info )
698#endif
699      ENDIF 
700#endif
701      IF( ll_abort ) CALL ctl_stop ('STOP', 'NEMO abort mppstop')
702      !
703   END SUBROUTINE mppstop
704
705
706   SUBROUTINE mpp_comm_free( kcom )
707      !!----------------------------------------------------------------------
708      INTEGER, INTENT(in) ::   kcom
709      !!
710      INTEGER :: ierr
711      !!----------------------------------------------------------------------
712      !
713#if defined key_mpp_mpi
714      CALL MPI_COMM_FREE(kcom, ierr)
715#endif
716      !
717   END SUBROUTINE mpp_comm_free
718
719
720   SUBROUTINE mpp_ini_znl( kumout )
721      !!----------------------------------------------------------------------
722      !!               ***  routine mpp_ini_znl  ***
723      !!
724      !! ** Purpose :   Initialize special communicator for computing zonal sum
725      !!
726      !! ** Method  : - Look for processors in the same row
727      !!              - Put their number in nrank_znl
728      !!              - Create group for the znl processors
729      !!              - Create a communicator for znl processors
730      !!              - Determine if processor should write znl files
731      !!
732      !! ** output
733      !!      ndim_rank_znl = number of processors on the same row
734      !!      ngrp_znl = group ID for the znl processors
735      !!      ncomm_znl = communicator for the ice procs.
736      !!      n_znl_root = number (in the world) of proc 0 in the ice comm.
737      !!
738      !!----------------------------------------------------------------------
739      INTEGER, INTENT(in) ::   kumout   ! ocean.output logical units
740      !
741      INTEGER :: jproc      ! dummy loop integer
742      INTEGER :: ierr, ii   ! local integer
743      INTEGER, ALLOCATABLE, DIMENSION(:) ::   kwork
744      !!----------------------------------------------------------------------
745#if defined key_mpp_mpi
746      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_world     : ', ngrp_world
747      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_world : ', mpi_comm_world
748      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_oce   : ', mpi_comm_oce
749      !
750      ALLOCATE( kwork(jpnij), STAT=ierr )
751      IF( ierr /= 0 ) CALL ctl_stop( 'STOP', 'mpp_ini_znl : failed to allocate 1D array of length jpnij')
752
753      IF( jpnj == 1 ) THEN
754         ngrp_znl  = ngrp_world
755         ncomm_znl = mpi_comm_oce
756      ELSE
757         !
758         CALL MPI_ALLGATHER ( njmpp, 1, mpi_integer, kwork, 1, mpi_integer, mpi_comm_oce, ierr )
759         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - kwork pour njmpp : ', kwork
760         !-$$        CALL flush(numout)
761         !
762         ! Count number of processors on the same row
763         ndim_rank_znl = 0
764         DO jproc=1,jpnij
765            IF ( kwork(jproc) == njmpp ) THEN
766               ndim_rank_znl = ndim_rank_znl + 1
767            ENDIF
768         END DO
769         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ndim_rank_znl : ', ndim_rank_znl
770         !-$$        CALL flush(numout)
771         ! Allocate the right size to nrank_znl
772         IF (ALLOCATED (nrank_znl)) DEALLOCATE(nrank_znl)
773         ALLOCATE(nrank_znl(ndim_rank_znl))
774         ii = 0
775         nrank_znl (:) = 0
776         DO jproc=1,jpnij
777            IF ( kwork(jproc) == njmpp) THEN
778               ii = ii + 1
779               nrank_znl(ii) = jproc -1
780            ENDIF
781         END DO
782         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - nrank_znl : ', nrank_znl
783         !-$$        CALL flush(numout)
784
785         ! Create the opa group
786         CALL MPI_COMM_GROUP(mpi_comm_oce,ngrp_opa,ierr)
787         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_opa : ', ngrp_opa
788         !-$$        CALL flush(numout)
789
790         ! Create the znl group from the opa group
791         CALL MPI_GROUP_INCL  ( ngrp_opa, ndim_rank_znl, nrank_znl, ngrp_znl, ierr )
792         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_znl ', ngrp_znl
793         !-$$        CALL flush(numout)
794
795         ! Create the znl communicator from the opa communicator, ie the pool of procs in the same row
796         CALL MPI_COMM_CREATE ( mpi_comm_oce, ngrp_znl, ncomm_znl, ierr )
797         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ncomm_znl ', ncomm_znl
798         !-$$        CALL flush(numout)
799         !
800      END IF
801
802      ! Determines if processor if the first (starting from i=1) on the row
803      IF ( jpni == 1 ) THEN
804         l_znl_root = .TRUE.
805      ELSE
806         l_znl_root = .FALSE.
807         kwork (1) = nimpp
808         CALL mpp_min ( 'lib_mpp', kwork(1), kcom = ncomm_znl)
809         IF ( nimpp == kwork(1)) l_znl_root = .TRUE.
810      END IF
811
812      DEALLOCATE(kwork)
813#endif
814
815   END SUBROUTINE mpp_ini_znl
816
817
818   SUBROUTINE mpp_ini_north
819      !!----------------------------------------------------------------------
820      !!               ***  routine mpp_ini_north  ***
821      !!
822      !! ** Purpose :   Initialize special communicator for north folding
823      !!      condition together with global variables needed in the mpp folding
824      !!
825      !! ** Method  : - Look for northern processors
826      !!              - Put their number in nrank_north
827      !!              - Create groups for the world processors and the north processors
828      !!              - Create a communicator for northern processors
829      !!
830      !! ** output
831      !!      njmppmax = njmpp for northern procs
832      !!      ndim_rank_north = number of processors in the northern line
833      !!      nrank_north (ndim_rank_north) = number  of the northern procs.
834      !!      ngrp_world = group ID for the world processors
835      !!      ngrp_north = group ID for the northern processors
836      !!      ncomm_north = communicator for the northern procs.
837      !!      north_root = number (in the world) of proc 0 in the northern comm.
838      !!
839      !!----------------------------------------------------------------------
840      INTEGER ::   ierr
841      INTEGER ::   jjproc
842      INTEGER ::   ii, ji
843      !!----------------------------------------------------------------------
844      !
845#if defined key_mpp_mpi
846      njmppmax = MAXVAL( njmppt )
847      !
848      ! Look for how many procs on the northern boundary
849      ndim_rank_north = 0
850      DO jjproc = 1, jpnij
851         IF( njmppt(jjproc) == njmppmax )   ndim_rank_north = ndim_rank_north + 1
852      END DO
853      !
854      ! Allocate the right size to nrank_north
855      IF (ALLOCATED (nrank_north)) DEALLOCATE(nrank_north)
856      ALLOCATE( nrank_north(ndim_rank_north) )
857
858      ! Fill the nrank_north array with proc. number of northern procs.
859      ! Note : the rank start at 0 in MPI
860      ii = 0
861      DO ji = 1, jpnij
862         IF ( njmppt(ji) == njmppmax   ) THEN
863            ii=ii+1
864            nrank_north(ii)=ji-1
865         END IF
866      END DO
867      !
868      ! create the world group
869      CALL MPI_COMM_GROUP( mpi_comm_oce, ngrp_world, ierr )
870      !
871      ! Create the North group from the world group
872      CALL MPI_GROUP_INCL( ngrp_world, ndim_rank_north, nrank_north, ngrp_north, ierr )
873      !
874      ! Create the North communicator , ie the pool of procs in the north group
875      CALL MPI_COMM_CREATE( mpi_comm_oce, ngrp_north, ncomm_north, ierr )
876      !
877#endif
878   END SUBROUTINE mpp_ini_north
879
880
881   SUBROUTINE DDPDD_MPI( ydda, yddb, ilen, itype )
882      !!---------------------------------------------------------------------
883      !!   Routine DDPDD_MPI: used by reduction operator MPI_SUMDD
884      !!
885      !!   Modification of original codes written by David H. Bailey
886      !!   This subroutine computes yddb(i) = ydda(i)+yddb(i)
887      !!---------------------------------------------------------------------
888      INTEGER                     , INTENT(in)    ::   ilen, itype
889      COMPLEX(wp), DIMENSION(ilen), INTENT(in)    ::   ydda
890      COMPLEX(wp), DIMENSION(ilen), INTENT(inout) ::   yddb
891      !
892      REAL(wp) :: zerr, zt1, zt2    ! local work variables
893      INTEGER  :: ji, ztmp           ! local scalar
894      !!---------------------------------------------------------------------
895      !
896      ztmp = itype   ! avoid compilation warning
897      !
898      DO ji=1,ilen
899      ! Compute ydda + yddb using Knuth's trick.
900         zt1  = real(ydda(ji)) + real(yddb(ji))
901         zerr = zt1 - real(ydda(ji))
902         zt2  = ((real(yddb(ji)) - zerr) + (real(ydda(ji)) - (zt1 - zerr))) &
903                + aimag(ydda(ji)) + aimag(yddb(ji))
904
905         ! The result is zt1 + zt2, after normalization.
906         yddb(ji) = cmplx ( zt1 + zt2, zt2 - ((zt1 + zt2) - zt1),wp )
907      END DO
908      !
909   END SUBROUTINE DDPDD_MPI
910
911
912   SUBROUTINE mpp_report( cdname, kpk, kpl, kpf, ld_lbc, ld_glb, ld_dlg )
913      !!----------------------------------------------------------------------
914      !!                  ***  routine mpp_report  ***
915      !!
916      !! ** Purpose :   report use of mpp routines per time-setp
917      !!
918      !!----------------------------------------------------------------------
919      CHARACTER(len=*),           INTENT(in   ) ::   cdname      ! name of the calling subroutine
920      INTEGER         , OPTIONAL, INTENT(in   ) ::   kpk, kpl, kpf
921      LOGICAL         , OPTIONAL, INTENT(in   ) ::   ld_lbc, ld_glb, ld_dlg
922      !!
923      CHARACTER(len=128)                        ::   ccountname  ! name of a subroutine to count communications
924      LOGICAL ::   ll_lbc, ll_glb, ll_dlg
925      INTEGER ::    ji,  jj,  jk,  jh, jf, jcount   ! dummy loop indices
926      !!----------------------------------------------------------------------
927#if defined key_mpp_mpi
928      !
929      ll_lbc = .FALSE.
930      IF( PRESENT(ld_lbc) ) ll_lbc = ld_lbc
931      ll_glb = .FALSE.
932      IF( PRESENT(ld_glb) ) ll_glb = ld_glb
933      ll_dlg = .FALSE.
934      IF( PRESENT(ld_dlg) ) ll_dlg = ld_dlg
935      !
936      ! find the smallest common frequency: default = frequency product, if multiple, choose the larger of the 2 frequency
937      IF( ncom_dttrc /= 1 )   CALL ctl_stop( 'STOP', 'mpp_report, ncom_dttrc /= 1 not coded...' ) 
938      ncom_freq = ncom_fsbc
939      !
940      IF ( ncom_stp == nit000+ncom_freq ) THEN   ! avoid to count extra communications in potential initializations at nit000
941         IF( ll_lbc ) THEN
942            IF( .NOT. ALLOCATED(ncomm_sequence) ) ALLOCATE( ncomm_sequence(ncom_rec_max,2) )
943            IF( .NOT. ALLOCATED(    crname_lbc) ) ALLOCATE(     crname_lbc(ncom_rec_max  ) )
944            n_sequence_lbc = n_sequence_lbc + 1
945            IF( n_sequence_lbc > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' )   ! deadlock
946            crname_lbc(n_sequence_lbc) = cdname     ! keep the name of the calling routine
947            ncomm_sequence(n_sequence_lbc,1) = kpk*kpl   ! size of 3rd and 4th dimensions
948            ncomm_sequence(n_sequence_lbc,2) = kpf       ! number of arrays to be treated (multi)
949         ENDIF
950         IF( ll_glb ) THEN
951            IF( .NOT. ALLOCATED(crname_glb) ) ALLOCATE( crname_glb(ncom_rec_max) )
952            n_sequence_glb = n_sequence_glb + 1
953            IF( n_sequence_glb > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' )   ! deadlock
954            crname_glb(n_sequence_glb) = cdname     ! keep the name of the calling routine
955         ENDIF
956         IF( ll_dlg ) THEN
957            IF( .NOT. ALLOCATED(crname_dlg) ) ALLOCATE( crname_dlg(ncom_rec_max) )
958            n_sequence_dlg = n_sequence_dlg + 1
959            IF( n_sequence_dlg > ncom_rec_max ) CALL ctl_stop( 'STOP', 'lib_mpp, increase ncom_rec_max' )   ! deadlock
960            crname_dlg(n_sequence_dlg) = cdname     ! keep the name of the calling routine
961         ENDIF
962      ELSE IF ( ncom_stp == nit000+2*ncom_freq ) THEN
963         CALL ctl_opn( numcom, 'communication_report.txt', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE., narea )
964         WRITE(numcom,*) ' '
965         WRITE(numcom,*) ' ------------------------------------------------------------'
966         WRITE(numcom,*) ' Communication pattern report (second oce+sbc+top time step):'
967         WRITE(numcom,*) ' ------------------------------------------------------------'
968         WRITE(numcom,*) ' '
969         WRITE(numcom,'(A,I4)') ' Exchanged halos : ', n_sequence_lbc
970         jj = 0; jk = 0; jf = 0; jh = 0
971         DO ji = 1, n_sequence_lbc
972            IF ( ncomm_sequence(ji,1) .GT. 1 ) jk = jk + 1
973            IF ( ncomm_sequence(ji,2) .GT. 1 ) jf = jf + 1
974            IF ( ncomm_sequence(ji,1) .GT. 1 .AND. ncomm_sequence(ji,2) .GT. 1 ) jj = jj + 1
975            jh = MAX (jh, ncomm_sequence(ji,1)*ncomm_sequence(ji,2))
976         END DO
977         WRITE(numcom,'(A,I3)') ' 3D Exchanged halos : ', jk
978         WRITE(numcom,'(A,I3)') ' Multi arrays exchanged halos : ', jf
979         WRITE(numcom,'(A,I3)') '   from which 3D : ', jj
980         WRITE(numcom,'(A,I10)') ' Array max size : ', jh*jpi*jpj
981         WRITE(numcom,*) ' '
982         WRITE(numcom,*) ' lbc_lnk called'
983         DO ji = 1, n_sequence_lbc - 1
984            IF ( crname_lbc(ji) /= 'already counted' ) THEN
985               ccountname = crname_lbc(ji)
986               crname_lbc(ji) = 'already counted'
987               jcount = 1
988               DO jj = ji + 1, n_sequence_lbc
989                  IF ( ccountname ==  crname_lbc(jj) ) THEN
990                     jcount = jcount + 1
991                     crname_lbc(jj) = 'already counted'
992                  END IF
993               END DO
994               WRITE(numcom,'(A, I4, A, A)') ' - ', jcount,' times by subroutine ', TRIM(ccountname)
995            END IF
996         END DO
997         IF ( crname_lbc(n_sequence_lbc) /= 'already counted' ) THEN
998            WRITE(numcom,'(A, I4, A, A)') ' - ', 1,' times by subroutine ', TRIM(crname_lbc(ncom_rec_max))
999         END IF
1000         WRITE(numcom,*) ' '
1001         IF ( n_sequence_glb > 0 ) THEN
1002            WRITE(numcom,'(A,I4)') ' Global communications : ', n_sequence_glb
1003            jj = 1
1004            DO ji = 2, n_sequence_glb
1005               IF( crname_glb(ji-1) /= crname_glb(ji) ) THEN
1006                  WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_glb(ji-1))
1007                  jj = 0
1008               END IF
1009               jj = jj + 1 
1010            END DO
1011            WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_glb(n_sequence_glb))
1012            DEALLOCATE(crname_glb)
1013         ELSE
1014            WRITE(numcom,*) ' No MPI global communication '
1015         ENDIF
1016         WRITE(numcom,*) ' '
1017         IF ( n_sequence_dlg > 0 ) THEN
1018            WRITE(numcom,'(A,I4)') ' Delayed global communications : ', n_sequence_dlg
1019            jj = 1
1020            DO ji = 2, n_sequence_dlg
1021               IF( crname_dlg(ji-1) /= crname_dlg(ji) ) THEN
1022                  WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_dlg(ji-1))
1023                  jj = 0
1024               END IF
1025               jj = jj + 1 
1026            END DO
1027            WRITE(numcom,'(A, I4, A, A)') ' - ', jj,' times by subroutine ', TRIM(crname_dlg(n_sequence_dlg))
1028            DEALLOCATE(crname_dlg)
1029         ELSE
1030            WRITE(numcom,*) ' No MPI delayed global communication '
1031         ENDIF
1032         WRITE(numcom,*) ' '
1033         WRITE(numcom,*) ' -----------------------------------------------'
1034         WRITE(numcom,*) ' '
1035         DEALLOCATE(ncomm_sequence)
1036         DEALLOCATE(crname_lbc)
1037      ENDIF
1038#endif
1039   END SUBROUTINE mpp_report
1040
1041   
1042   SUBROUTINE tic_tac (ld_tic, ld_global)
1043
1044    LOGICAL,           INTENT(IN) :: ld_tic
1045    LOGICAL, OPTIONAL, INTENT(IN) :: ld_global
1046    REAL(wp), DIMENSION(2), SAVE :: tic_wt
1047    REAL(wp),               SAVE :: tic_ct = 0._wp
1048    INTEGER :: ii
1049#if defined key_mpp_mpi
1050
1051    IF( ncom_stp <= nit000 ) RETURN
1052    IF( ncom_stp == nitend ) RETURN
1053    ii = 1
1054    IF( PRESENT( ld_global ) ) THEN
1055       IF( ld_global ) ii = 2
1056    END IF
1057   
1058    IF ( ld_tic ) THEN
1059       tic_wt(ii) = MPI_Wtime()                                                    ! start count tic->tac (waiting time)
1060       IF ( tic_ct > 0.0_wp ) compute_time = compute_time + MPI_Wtime() - tic_ct   ! cumulate count tac->tic
1061    ELSE
1062       waiting_time(ii) = waiting_time(ii) + MPI_Wtime() - tic_wt(ii)              ! cumulate count tic->tac
1063       tic_ct = MPI_Wtime()                                                        ! start count tac->tic (waiting time)
1064    ENDIF
1065#endif
1066   
1067   END SUBROUTINE tic_tac
1068
1069#if ! defined key_mpp_mpi
1070   SUBROUTINE mpi_wait(request, status, ierror)
1071      INTEGER                            , INTENT(in   ) ::   request
1072      INTEGER, DIMENSION(MPI_STATUS_SIZE), INTENT(  out) ::   status
1073      INTEGER                            , INTENT(  out) ::   ierror
1074   END SUBROUTINE mpi_wait
1075
1076   
1077   FUNCTION MPI_Wtime()
1078      REAL(wp) ::  MPI_Wtime
1079      MPI_Wtime = -1.
1080   END FUNCTION MPI_Wtime
1081#endif
1082
1083   !!----------------------------------------------------------------------
1084   !!   ctl_stop, ctl_warn, get_unit, ctl_opn, ctl_nam   routines
1085   !!----------------------------------------------------------------------
1086
1087   SUBROUTINE ctl_stop( cd1, cd2, cd3, cd4, cd5 ,   &
1088      &                 cd6, cd7, cd8, cd9, cd10 )
1089      !!----------------------------------------------------------------------
1090      !!                  ***  ROUTINE  stop_opa  ***
1091      !!
1092      !! ** Purpose :   print in ocean.outpput file a error message and
1093      !!                increment the error number (nstop) by one.
1094      !!----------------------------------------------------------------------
1095      CHARACTER(len=*), INTENT(in   )           ::   cd1
1096      CHARACTER(len=*), INTENT(in   ), OPTIONAL ::        cd2, cd3, cd4, cd5
1097      CHARACTER(len=*), INTENT(in   ), OPTIONAL ::   cd6, cd7, cd8, cd9, cd10
1098      !
1099      INTEGER ::   inum
1100      !!----------------------------------------------------------------------
1101      !
1102      nstop = nstop + 1
1103      !
1104      IF( cd1 == 'STOP' .AND. narea /= 1 ) THEN    ! Immediate stop: add an arror message in 'ocean.output' file
1105         CALL ctl_opn( inum, 'ocean.output', 'APPEND', 'FORMATTED', 'SEQUENTIAL', -1, 6, .FALSE. )
1106         WRITE(inum,*)
1107         WRITE(inum,*) ' ==>>>   Look for "E R R O R" messages in all existing *ocean.output* files'
1108         CLOSE(inum)
1109      ENDIF
1110      IF( numout == 6 ) THEN                       ! force to open ocean.output file if not already opened
1111         CALL ctl_opn( numout, 'ocean.output', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, -1, .FALSE., narea )
1112      ENDIF
1113      !
1114                            WRITE(numout,*)
1115                            WRITE(numout,*) ' ===>>> : E R R O R'
1116                            WRITE(numout,*)
1117                            WRITE(numout,*) '         ==========='
1118                            WRITE(numout,*)
1119                            WRITE(numout,*) TRIM(cd1)
1120      IF( PRESENT(cd2 ) )   WRITE(numout,*) TRIM(cd2)
1121      IF( PRESENT(cd3 ) )   WRITE(numout,*) TRIM(cd3)
1122      IF( PRESENT(cd4 ) )   WRITE(numout,*) TRIM(cd4)
1123      IF( PRESENT(cd5 ) )   WRITE(numout,*) TRIM(cd5)
1124      IF( PRESENT(cd6 ) )   WRITE(numout,*) TRIM(cd6)
1125      IF( PRESENT(cd7 ) )   WRITE(numout,*) TRIM(cd7)
1126      IF( PRESENT(cd8 ) )   WRITE(numout,*) TRIM(cd8)
1127      IF( PRESENT(cd9 ) )   WRITE(numout,*) TRIM(cd9)
1128      IF( PRESENT(cd10) )   WRITE(numout,*) TRIM(cd10)
1129                            WRITE(numout,*)
1130      !
1131                               CALL FLUSH(numout    )
1132      IF( numstp     /= -1 )   CALL FLUSH(numstp    )
1133      IF( numrun     /= -1 )   CALL FLUSH(numrun    )
1134      IF( numevo_ice /= -1 )   CALL FLUSH(numevo_ice)
1135      !
1136      IF( cd1 == 'STOP' ) THEN
1137         WRITE(numout,*) 
1138         WRITE(numout,*)  'huge E-R-R-O-R : immediate stop'
1139         WRITE(numout,*) 
1140         CALL FLUSH(numout)
1141         CALL SLEEP(60)   ! make sure that all output and abort files are written by all cores. 60s should be enough...
1142         CALL mppstop( ld_abort = .true. )
1143      ENDIF
1144      !
1145   END SUBROUTINE ctl_stop
1146
1147
1148   SUBROUTINE ctl_warn( cd1, cd2, cd3, cd4, cd5,   &
1149      &                 cd6, cd7, cd8, cd9, cd10 )
1150      !!----------------------------------------------------------------------
1151      !!                  ***  ROUTINE  stop_warn  ***
1152      !!
1153      !! ** Purpose :   print in ocean.outpput file a error message and
1154      !!                increment the warning number (nwarn) by one.
1155      !!----------------------------------------------------------------------
1156      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd1, cd2, cd3, cd4, cd5
1157      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd6, cd7, cd8, cd9, cd10
1158      !!----------------------------------------------------------------------
1159      !
1160      nwarn = nwarn + 1
1161      !
1162      IF(lwp) THEN
1163                               WRITE(numout,*)
1164                               WRITE(numout,*) ' ===>>> : W A R N I N G'
1165                               WRITE(numout,*)
1166                               WRITE(numout,*) '         ==============='
1167                               WRITE(numout,*)
1168         IF( PRESENT(cd1 ) )   WRITE(numout,*) TRIM(cd1)
1169         IF( PRESENT(cd2 ) )   WRITE(numout,*) TRIM(cd2)
1170         IF( PRESENT(cd3 ) )   WRITE(numout,*) TRIM(cd3)
1171         IF( PRESENT(cd4 ) )   WRITE(numout,*) TRIM(cd4)
1172         IF( PRESENT(cd5 ) )   WRITE(numout,*) TRIM(cd5)
1173         IF( PRESENT(cd6 ) )   WRITE(numout,*) TRIM(cd6)
1174         IF( PRESENT(cd7 ) )   WRITE(numout,*) TRIM(cd7)
1175         IF( PRESENT(cd8 ) )   WRITE(numout,*) TRIM(cd8)
1176         IF( PRESENT(cd9 ) )   WRITE(numout,*) TRIM(cd9)
1177         IF( PRESENT(cd10) )   WRITE(numout,*) TRIM(cd10)
1178                               WRITE(numout,*)
1179      ENDIF
1180      CALL FLUSH(numout)
1181      !
1182   END SUBROUTINE ctl_warn
1183
1184
1185   SUBROUTINE ctl_opn( knum, cdfile, cdstat, cdform, cdacce, klengh, kout, ldwp, karea )
1186      !!----------------------------------------------------------------------
1187      !!                  ***  ROUTINE ctl_opn  ***
1188      !!
1189      !! ** Purpose :   Open file and check if required file is available.
1190      !!
1191      !! ** Method  :   Fortan open
1192      !!----------------------------------------------------------------------
1193      INTEGER          , INTENT(  out) ::   knum      ! logical unit to open
1194      CHARACTER(len=*) , INTENT(in   ) ::   cdfile    ! file name to open
1195      CHARACTER(len=*) , INTENT(in   ) ::   cdstat    ! disposition specifier
1196      CHARACTER(len=*) , INTENT(in   ) ::   cdform    ! formatting specifier
1197      CHARACTER(len=*) , INTENT(in   ) ::   cdacce    ! access specifier
1198      INTEGER          , INTENT(in   ) ::   klengh    ! record length
1199      INTEGER          , INTENT(in   ) ::   kout      ! number of logical units for write
1200      LOGICAL          , INTENT(in   ) ::   ldwp      ! boolean term for print
1201      INTEGER, OPTIONAL, INTENT(in   ) ::   karea     ! proc number
1202      !
1203      CHARACTER(len=80) ::   clfile
1204      CHARACTER(LEN=10) ::   clfmt            ! writing format
1205      INTEGER           ::   iost
1206      INTEGER           ::   idg              ! number of digits
1207      !!----------------------------------------------------------------------
1208      !
1209      ! adapt filename
1210      ! ----------------
1211      clfile = TRIM(cdfile)
1212      IF( PRESENT( karea ) ) THEN
1213         IF( karea > 1 ) THEN
1214            ! Warning: jpnij is maybe not already defined when calling ctl_opn -> use mppsize instead of jpnij
1215            idg = MAX( INT(LOG10(REAL(MAX(1,mppsize-1),wp))) + 1, 4 )      ! how many digits to we need to write? min=4, max=9
1216            WRITE(clfmt, "('(a,a,i', i1, '.', i1, ')')") idg, idg          ! '(a,a,ix.x)'
1217            WRITE(clfile, clfmt) TRIM(clfile), '_', karea-1
1218         ENDIF
1219      ENDIF
1220#if defined key_agrif
1221      IF( .NOT. Agrif_Root() )   clfile = TRIM(Agrif_CFixed())//'_'//TRIM(clfile)
1222      knum=Agrif_Get_Unit()
1223#else
1224      knum=get_unit()
1225#endif
1226      IF( TRIM(cdfile) == '/dev/null' )   clfile = TRIM(cdfile)   ! force the use of /dev/null
1227      !
1228      IF(       cdacce(1:6) == 'DIRECT' )  THEN   ! cdacce has always more than 6 characters
1229         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat, RECL=klengh         , ERR=100, IOSTAT=iost )
1230      ELSE IF( TRIM(cdstat) == 'APPEND' )  THEN   ! cdstat can have less than 6 characters
1231         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS='UNKNOWN', POSITION='APPEND', ERR=100, IOSTAT=iost )
1232      ELSE
1233         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat                      , ERR=100, IOSTAT=iost )
1234      ENDIF
1235      IF( iost /= 0 .AND. TRIM(clfile) == '/dev/null' ) &   ! for windows
1236         &  OPEN(UNIT=knum,FILE='NUL', FORM=cdform, ACCESS=cdacce, STATUS=cdstat                      , ERR=100, IOSTAT=iost )   
1237      IF( iost == 0 ) THEN
1238         IF(ldwp .AND. kout > 0) THEN
1239            WRITE(kout,*) '     file   : ', TRIM(clfile),' open ok'
1240            WRITE(kout,*) '     unit   = ', knum
1241            WRITE(kout,*) '     status = ', cdstat
1242            WRITE(kout,*) '     form   = ', cdform
1243            WRITE(kout,*) '     access = ', cdacce
1244            WRITE(kout,*)
1245         ENDIF
1246      ENDIF
1247100   CONTINUE
1248      IF( iost /= 0 ) THEN
1249         WRITE(ctmp1,*) ' ===>>>> : bad opening file: ', TRIM(clfile)
1250         WRITE(ctmp2,*) ' =======   ===  '
1251         WRITE(ctmp3,*) '           unit   = ', knum
1252         WRITE(ctmp4,*) '           status = ', cdstat
1253         WRITE(ctmp5,*) '           form   = ', cdform
1254         WRITE(ctmp6,*) '           access = ', cdacce
1255         WRITE(ctmp7,*) '           iostat = ', iost
1256         WRITE(ctmp8,*) '           we stop. verify the file '
1257         CALL ctl_stop( 'STOP', ctmp1, ctmp2, ctmp3, ctmp4, ctmp5, ctmp6, ctmp7, ctmp8 )
1258      ENDIF
1259      !
1260   END SUBROUTINE ctl_opn
1261
1262
1263   SUBROUTINE ctl_nam ( kios, cdnam )
1264      !!----------------------------------------------------------------------
1265      !!                  ***  ROUTINE ctl_nam  ***
1266      !!
1267      !! ** Purpose :   Informations when error while reading a namelist
1268      !!
1269      !! ** Method  :   Fortan open
1270      !!----------------------------------------------------------------------
1271      INTEGER                                , INTENT(inout) ::   kios    ! IO status after reading the namelist
1272      CHARACTER(len=*)                       , INTENT(in   ) ::   cdnam   ! group name of namelist for which error occurs
1273      !
1274      CHARACTER(len=5) ::   clios   ! string to convert iostat in character for print
1275      !!----------------------------------------------------------------------
1276      !
1277      WRITE (clios, '(I5.0)')   kios
1278      IF( kios < 0 ) THEN         
1279         CALL ctl_warn( 'end of record or file while reading namelist '   &
1280            &           // TRIM(cdnam) // ' iostat = ' // TRIM(clios) )
1281      ENDIF
1282      !
1283      IF( kios > 0 ) THEN
1284         CALL ctl_stop( 'misspelled variable in namelist '   &
1285            &           // TRIM(cdnam) // ' iostat = ' // TRIM(clios) )
1286      ENDIF
1287      kios = 0
1288      !
1289   END SUBROUTINE ctl_nam
1290
1291
1292   INTEGER FUNCTION get_unit()
1293      !!----------------------------------------------------------------------
1294      !!                  ***  FUNCTION  get_unit  ***
1295      !!
1296      !! ** Purpose :   return the index of an unused logical unit
1297      !!----------------------------------------------------------------------
1298      LOGICAL :: llopn
1299      !!----------------------------------------------------------------------
1300      !
1301      get_unit = 15   ! choose a unit that is big enough then it is not already used in NEMO
1302      llopn = .TRUE.
1303      DO WHILE( (get_unit < 998) .AND. llopn )
1304         get_unit = get_unit + 1
1305         INQUIRE( unit = get_unit, opened = llopn )
1306      END DO
1307      IF( (get_unit == 999) .AND. llopn ) THEN
1308         CALL ctl_stop( 'STOP', 'get_unit: All logical units until 999 are used...' )
1309      ENDIF
1310      !
1311   END FUNCTION get_unit
1312
1313   !!----------------------------------------------------------------------
1314END MODULE lib_mpp
Note: See TracBrowser for help on using the repository browser.