New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
lib_mpp.F90 in NEMO/branches/2018/dev_r9759_HPC09_ESIWACE/src/OCE/LBC – NEMO

source: NEMO/branches/2018/dev_r9759_HPC09_ESIWACE/src/OCE/LBC/lib_mpp.F90 @ 9772

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

dev_r9759_HPC09_ESIWACE: add benchmark features

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