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

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

source: NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/LBC/lib_mpp.F90 @ 10297

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

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 2a: add report calls of mppmin/max/sum, see #2133

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