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 branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/LBC – NEMO

source: branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90 @ 9012

Last change on this file since 9012 was 9012, checked in by acc, 6 years ago

Branch dev_CNRS_2017. Merge in no_ghost changes from dev_r8126_ROBUST08_no_ghost. These changes include lib_mpp refresh and rationalisation of mppini from dev_r8126_ROBUST10_MPPINI

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