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

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

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

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