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

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

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 4b: reduce communications in si3, see #2133

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