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

Last change on this file since 10172 was 10172, checked in by smasson, 6 years ago

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 2b: improve of timing, add computing and waiting time, see #2133

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