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

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

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90 @ 3243

Last change on this file since 3243 was 3243, checked in by acc, 12 years ago

Branch 2011/dev_NEMO_MERGE_2011. Minor correction in lib_mpp.F90 to ensure that ctl_stop terminates execution when running without key_mpp_mpi

  • Property svn:keywords set to Id
File size: 129.3 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   !!----------------------------------------------------------------------
22
23   !!----------------------------------------------------------------------
24   !!   ctl_stop   : update momentum and tracer Kz from a tke scheme
25   !!   ctl_warn   : initialization, namelist read, and parameters control
26   !!   ctl_opn    : Open file and check if required file is available.
27   !!   get_unit    : give the index of an unused logical unit
28   !!----------------------------------------------------------------------
29#if   defined key_mpp_mpi 
30   !!----------------------------------------------------------------------
31   !!   'key_mpp_mpi'             MPI massively parallel processing library
32   !!----------------------------------------------------------------------
33   !!   lib_mpp_alloc : allocate mpp arrays
34   !!   mynode        : indentify the processor unit
35   !!   mpp_lnk       : interface (defined in lbclnk) for message passing of 2d or 3d arrays (mpp_lnk_2d, mpp_lnk_3d)
36   !!   mpp_lnk_3d_gather :  Message passing manadgement for two 3D arrays
37   !!   mpp_lnk_e     : interface (defined in lbclnk) for message passing of 2d array with extra halo (mpp_lnk_2d_e)
38   !!   mpprecv         :
39   !!   mppsend       :   SUBROUTINE mpp_ini_znl
40   !!   mppscatter    :
41   !!   mppgather     :
42   !!   mpp_min       : generic interface for mppmin_int , mppmin_a_int , mppmin_real, mppmin_a_real
43   !!   mpp_max       : generic interface for mppmax_int , mppmax_a_int , mppmax_real, mppmax_a_real
44   !!   mpp_sum       : generic interface for mppsum_int , mppsum_a_int , mppsum_real, mppsum_a_real
45   !!   mpp_minloc    :
46   !!   mpp_maxloc    :
47   !!   mppsync       :
48   !!   mppstop       :
49   !!   mppobc        : variant of mpp_lnk for open boundary condition
50   !!   mpp_ini_north : initialisation of north fold
51   !!   mpp_lbc_north : north fold processors gathering
52   !!   mpp_lbc_north_e : variant of mpp_lbc_north for extra outer halo
53   !!----------------------------------------------------------------------
54   USE dom_oce        ! ocean space and time domain
55   USE lbcnfd         ! north fold treatment
56   USE in_out_manager ! I/O manager
57
58   IMPLICIT NONE
59   PRIVATE
60   
61   PUBLIC   ctl_stop, ctl_warn, get_unit, ctl_opn
62   PUBLIC   mynode, mppstop, mppsync, mpp_comm_free
63   PUBLIC   mpp_ini_north, mpp_lbc_north, mpp_lbc_north_e
64   PUBLIC   mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc
65   PUBLIC   mpp_lnk_3d, mpp_lnk_3d_gather, mpp_lnk_2d, mpp_lnk_2d_e
66   PUBLIC   mppscatter, mppgather
67   PUBLIC   mppobc, mpp_ini_ice, mpp_ini_znl
68   PUBLIC   mppsize
69   PUBLIC   lib_mpp_alloc   ! Called in nemogcm.F90
70
71   !! * Interfaces
72   !! define generic interface for these routine as they are called sometimes
73   !! with scalar arguments instead of array arguments, which causes problems
74   !! for the compilation on AIX system as well as NEC and SGI. Ok on COMPACQ
75   INTERFACE mpp_min
76      MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real
77   END INTERFACE
78   INTERFACE mpp_max
79      MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real
80   END INTERFACE
81   INTERFACE mpp_sum
82# if defined key_mpp_rep
83      MODULE PROCEDURE mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real, &
84                       mppsum_realdd, mppsum_a_realdd
85# else
86      MODULE PROCEDURE mppsum_a_int, mppsum_int, mppsum_a_real, mppsum_real
87# endif
88   END INTERFACE
89   INTERFACE mpp_lbc_north
90      MODULE PROCEDURE mpp_lbc_north_3d, mpp_lbc_north_2d 
91   END INTERFACE
92   INTERFACE mpp_minloc
93      MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d
94   END INTERFACE
95   INTERFACE mpp_maxloc
96      MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d
97   END INTERFACE
98   
99   !! ========================= !!
100   !!  MPI  variable definition !!
101   !! ========================= !!
102!$AGRIF_DO_NOT_TREAT
103   INCLUDE 'mpif.h'
104!$AGRIF_END_DO_NOT_TREAT
105   
106   LOGICAL, PUBLIC, PARAMETER ::   lk_mpp = .TRUE.    !: mpp flag
107
108   INTEGER, PARAMETER         ::   nprocmax = 2**10   ! maximun dimension (required to be a power of 2)
109   
110   INTEGER ::   mppsize        ! number of process
111   INTEGER ::   mpprank        ! process number  [ 0 - size-1 ]
112!$AGRIF_DO_NOT_TREAT
113   INTEGER, PUBLIC ::   mpi_comm_opa   ! opa local communicator
114!$AGRIF_END_DO_NOT_TREAT
115
116# if defined key_mpp_rep
117   INTEGER :: MPI_SUMDD
118# endif
119
120   ! variables used in case of sea-ice
121   INTEGER, PUBLIC ::   ncomm_ice       !: communicator made by the processors with sea-ice
122   INTEGER ::   ngrp_ice        !  group ID for the ice processors (for rheology)
123   INTEGER ::   ndim_rank_ice   !  number of 'ice' processors
124   INTEGER ::   n_ice_root      !  number (in the comm_ice) of proc 0 in the ice comm
125   INTEGER, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_ice     ! dimension ndim_rank_ice
126
127   ! variables used for zonal integration
128   INTEGER, PUBLIC ::   ncomm_znl       !: communicator made by the processors on the same zonal average
129   LOGICAL, PUBLIC ::   l_znl_root      ! True on the 'left'most processor on the same row
130   INTEGER ::   ngrp_znl        ! group ID for the znl processors
131   INTEGER ::   ndim_rank_znl   ! number of processors on the same zonal average
132   INTEGER, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_znl  ! dimension ndim_rank_znl, number of the procs into the same znl domain
133   
134   ! North fold condition in mpp_mpi with jpni > 1
135   INTEGER ::   ngrp_world        ! group ID for the world processors
136   INTEGER ::   ngrp_opa          ! group ID for the opa processors
137   INTEGER ::   ngrp_north        ! group ID for the northern processors (to be fold)
138   INTEGER ::   ncomm_north       ! communicator made by the processors belonging to ngrp_north
139   INTEGER ::   ndim_rank_north   ! number of 'sea' processor in the northern line (can be /= jpni !)
140   INTEGER ::   njmppmax          ! value of njmpp for the processors of the northern line
141   INTEGER ::   north_root        ! number (in the comm_opa) of proc 0 in the northern comm
142   INTEGER, DIMENSION(:), ALLOCATABLE, SAVE ::   nrank_north   ! dimension ndim_rank_north
143
144   ! Type of send : standard, buffered, immediate
145   CHARACTER(len=1) ::   cn_mpi_send = 'S'    ! type od mpi send/recieve (S=standard, B=bsend, I=isend)
146   LOGICAL          ::   l_isend = .FALSE.   ! isend use indicator (T if cn_mpi_send='I')
147   INTEGER          ::   nn_buffer = 0       ! size of the buffer in case of mpi_bsend
148     
149   REAL(wp), DIMENSION(:), ALLOCATABLE, SAVE :: tampon  ! buffer in case of bsend
150
151   ! message passing arrays
152   REAL(wp), DIMENSION(:,:,:,:,:), ALLOCATABLE, SAVE ::   t4ns, t4sn   ! 2 x 3d for north-south & south-north
153   REAL(wp), DIMENSION(:,:,:,:,:), ALLOCATABLE, SAVE ::   t4ew, t4we   ! 2 x 3d for east-west & west-east
154   REAL(wp), DIMENSION(:,:,:,:,:), ALLOCATABLE, SAVE ::   t4p1, t4p2   ! 2 x 3d for north fold
155   REAL(wp), DIMENSION(:,:,:,:)  , ALLOCATABLE, SAVE ::   t3ns, t3sn   ! 3d for north-south & south-north
156   REAL(wp), DIMENSION(:,:,:,:)  , ALLOCATABLE, SAVE ::   t3ew, t3we   ! 3d for east-west & west-east
157   REAL(wp), DIMENSION(:,:,:,:)  , ALLOCATABLE, SAVE ::   t3p1, t3p2   ! 3d for north fold
158   REAL(wp), DIMENSION(:,:,:)    , ALLOCATABLE, SAVE ::   t2ns, t2sn   ! 2d for north-south & south-north
159   REAL(wp), DIMENSION(:,:,:)    , ALLOCATABLE, SAVE ::   t2ew, t2we   ! 2d for east-west & west-east
160   REAL(wp), DIMENSION(:,:,:)    , ALLOCATABLE, SAVE ::   t2p1, t2p2   ! 2d for north fold
161   REAL(wp), DIMENSION(:,:,:)    , ALLOCATABLE, SAVE ::   tr2ns, tr2sn ! 2d for north-south & south-north + extra outer halo
162   REAL(wp), DIMENSION(:,:,:)    , ALLOCATABLE, SAVE ::   tr2ew, tr2we ! 2d for east-west   & west-east   + extra outer halo
163
164   ! Arrays used in mpp_lbc_north_3d()
165   REAL(wp), DIMENSION(:,:,:)  , ALLOCATABLE, SAVE   ::   ztab, znorthloc
166   REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE, SAVE   ::   znorthgloio
167   REAL(wp), DIMENSION(:,:,:)  , ALLOCATABLE, SAVE   ::   zfoldwk      ! Workspace for message transfers avoiding mpi_allgather
168
169   ! Arrays used in mpp_lbc_north_2d()
170   REAL(wp), DIMENSION(:,:)  , ALLOCATABLE, SAVE    ::   ztab_2d, znorthloc_2d
171   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, SAVE    ::   znorthgloio_2d
172   REAL(wp), DIMENSION(:,:)  , ALLOCATABLE, SAVE    ::   zfoldwk_2d    ! Workspace for message transfers avoiding mpi_allgather
173
174   ! Arrays used in mpp_lbc_north_e()
175   REAL(wp), DIMENSION(:,:)  , ALLOCATABLE, SAVE    ::   ztab_e, znorthloc_e
176   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, SAVE    ::   znorthgloio_e
177
178   ! North fold arrays used to minimise the use of allgather operations. Set in nemo_northcomms (nemogcm) so need to be public
179   INTEGER, PUBLIC,  PARAMETER :: jpmaxngh = 8                 ! Assumed maximum number of active neighbours
180   INTEGER, PUBLIC,  PARAMETER :: jptyps   = 5                 ! Number of different neighbour lists to be used for northfold exchanges
181   INTEGER, PUBLIC,  DIMENSION (jpmaxngh,jptyps)    ::   isendto
182   INTEGER, PUBLIC,  DIMENSION (jptyps)             ::   nsndto
183   LOGICAL, PUBLIC                                  ::   ln_nnogather     = .FALSE.  ! namelist control of northfold comms
184   LOGICAL, PUBLIC                                  ::   l_north_nogather = .FALSE.  ! internal control of northfold comms
185   INTEGER, PUBLIC                                  ::   ityp
186   !!----------------------------------------------------------------------
187   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
188   !! $Id$
189   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
190   !!----------------------------------------------------------------------
191CONTAINS
192
193   INTEGER FUNCTION lib_mpp_alloc( kumout )
194      !!----------------------------------------------------------------------
195      !!              ***  routine lib_mpp_alloc  ***
196      !!----------------------------------------------------------------------
197      INTEGER, INTENT(in) ::   kumout   ! ocean.output logical unit
198      !!----------------------------------------------------------------------
199      !
200      ALLOCATE( t4ns(jpi,jprecj,jpk,2,2) , t4sn(jpi,jprecj,jpk,2,2) ,                                            &
201         &      t4ew(jpj,jpreci,jpk,2,2) , t4we(jpj,jpreci,jpk,2,2) ,                                            &
202         &      t4p1(jpi,jprecj,jpk,2,2) , t4p2(jpi,jprecj,jpk,2,2) ,                                            &
203         &      t3ns(jpi,jprecj,jpk,2)   , t3sn(jpi,jprecj,jpk,2)   ,                                            &
204         &      t3ew(jpj,jpreci,jpk,2)   , t3we(jpj,jpreci,jpk,2)   ,                                            &
205         &      t3p1(jpi,jprecj,jpk,2)   , t3p2(jpi,jprecj,jpk,2)   ,                                            &
206         &      t2ns(jpi,jprecj    ,2)   , t2sn(jpi,jprecj    ,2)   ,                                            &
207         &      t2ew(jpj,jpreci    ,2)   , t2we(jpj,jpreci    ,2)   ,                                            &
208         &      t2p1(jpi,jprecj    ,2)   , t2p2(jpi,jprecj    ,2)   ,                                            &
209         !
210         &      tr2ns(1-jpr2di:jpi+jpr2di,jprecj+jpr2dj,2) ,                                                     &
211         &      tr2sn(1-jpr2di:jpi+jpr2di,jprecj+jpr2dj,2) ,                                                     &
212         &      tr2ew(1-jpr2dj:jpj+jpr2dj,jpreci+jpr2di,2) ,                                                     &
213         &      tr2we(1-jpr2dj:jpj+jpr2dj,jpreci+jpr2di,2) ,                                                     &
214         !
215         &      ztab(jpiglo,4,jpk) , znorthloc(jpi,4,jpk) , znorthgloio(jpi,4,jpk,jpni) ,                        &
216         &      zfoldwk(jpi,4,jpk) ,                                                                             &
217         !
218         &      ztab_2d(jpiglo,4)  , znorthloc_2d(jpi,4)  , znorthgloio_2d(jpi,4,jpni)  ,                        &
219         &      zfoldwk_2d(jpi,4)  ,                                                                             &
220         !
221         &      ztab_e(jpiglo,4+2*jpr2dj) , znorthloc_e(jpi,4+2*jpr2dj) , znorthgloio_e(jpi,4+2*jpr2dj,jpni) ,   &
222         !
223         &      STAT=lib_mpp_alloc )
224         !
225      IF( lib_mpp_alloc /= 0 ) THEN
226         WRITE(kumout,cform_war)
227         WRITE(kumout,*) 'lib_mpp_alloc : failed to allocate arrays'
228      ENDIF
229      !
230   END FUNCTION lib_mpp_alloc
231
232
233   FUNCTION mynode( ldtxt, kumnam, kstop, localComm )
234      !!----------------------------------------------------------------------
235      !!                  ***  routine mynode  ***
236      !!                   
237      !! ** Purpose :   Find processor unit
238      !!----------------------------------------------------------------------
239      CHARACTER(len=*),DIMENSION(:), INTENT(  out) ::   ldtxt 
240      INTEGER                      , INTENT(in   ) ::   kumnam       ! namelist logical unit
241      INTEGER                      , INTENT(inout) ::   kstop        ! stop indicator
242      INTEGER, OPTIONAL            , INTENT(in   ) ::   localComm
243      !
244      INTEGER ::   mynode, ierr, code, ji, ii
245      LOGICAL ::   mpi_was_called
246      !
247      NAMELIST/nammpp/ cn_mpi_send, nn_buffer, jpni, jpnj, jpnij, ln_nnogather
248      !!----------------------------------------------------------------------
249      !
250      ii = 1
251      WRITE(ldtxt(ii),*)                                                                          ;   ii = ii + 1
252      WRITE(ldtxt(ii),*) 'mynode : mpi initialisation'                                            ;   ii = ii + 1
253      WRITE(ldtxt(ii),*) '~~~~~~ '                                                                ;   ii = ii + 1
254      !
255      jpni = -1; jpnj = -1; jpnij = -1
256      REWIND( kumnam )               ! Namelist namrun : parameters of the run
257      READ  ( kumnam, nammpp )
258      !                              ! control print
259      WRITE(ldtxt(ii),*) '   Namelist nammpp'                                                     ;   ii = ii + 1
260      WRITE(ldtxt(ii),*) '      mpi send type                      cn_mpi_send = ', cn_mpi_send   ;   ii = ii + 1
261      WRITE(ldtxt(ii),*) '      size in bytes of exported buffer   nn_buffer   = ', nn_buffer     ;   ii = ii + 1
262
263#if defined key_agrif
264      IF( .NOT. Agrif_Root() ) THEN
265         jpni  = Agrif_Parent(jpni ) 
266         jpnj  = Agrif_Parent(jpnj )
267         jpnij = Agrif_Parent(jpnij)
268      ENDIF
269#endif
270
271      IF(jpnij < 1)THEN
272         ! If jpnij is not specified in namelist then we calculate it - this
273         ! means there will be no land cutting out.
274         jpnij = jpni * jpnj
275      END IF
276
277      IF( (jpni < 1) .OR. (jpnj < 1) )THEN
278         WRITE(ldtxt(ii),*) '      jpni, jpnj and jpnij will be calculated automatically'; ii = ii + 1
279      ELSE
280         WRITE(ldtxt(ii),*) '      processor grid extent in i         jpni = ',jpni; ii = ii + 1
281         WRITE(ldtxt(ii),*) '      processor grid extent in j         jpnj = ',jpnj; ii = ii + 1
282         WRITE(ldtxt(ii),*) '      number of local domains           jpnij = ',jpnij; ii = ii +1
283      END IF
284
285      WRITE(ldtxt(ii),*) '      avoid use of mpi_allgather at the north fold  ln_nnogather = ', ln_nnogather  ; ii = ii + 1
286
287      CALL mpi_initialized ( mpi_was_called, code )
288      IF( code /= MPI_SUCCESS ) THEN
289         DO ji = 1, SIZE(ldtxt) 
290            IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
291         END DO         
292         WRITE(*, cform_err)
293         WRITE(*, *) 'lib_mpp: Error in routine mpi_initialized'
294         CALL mpi_abort( mpi_comm_world, code, ierr )
295      ENDIF
296
297      IF( mpi_was_called ) THEN
298         !
299         SELECT CASE ( cn_mpi_send )
300         CASE ( 'S' )                ! Standard mpi send (blocking)
301            WRITE(ldtxt(ii),*) '           Standard blocking mpi send (send)'                     ;   ii = ii + 1
302         CASE ( 'B' )                ! Buffer mpi send (blocking)
303            WRITE(ldtxt(ii),*) '           Buffer blocking mpi send (bsend)'                      ;   ii = ii + 1
304            IF( Agrif_Root() )   CALL mpi_init_opa( ldtxt, ii, ierr ) 
305         CASE ( 'I' )                ! Immediate mpi send (non-blocking send)
306            WRITE(ldtxt(ii),*) '           Immediate non-blocking send (isend)'                   ;   ii = ii + 1
307            l_isend = .TRUE.
308         CASE DEFAULT
309            WRITE(ldtxt(ii),cform_err)                                                            ;   ii = ii + 1
310            WRITE(ldtxt(ii),*) '           bad value for cn_mpi_send = ', cn_mpi_send             ;   ii = ii + 1
311            kstop = kstop + 1
312         END SELECT
313      ELSE IF ( PRESENT(localComm) .and. .not. mpi_was_called ) THEN
314         WRITE(ldtxt(ii),*) ' lib_mpp: You cannot provide a local communicator '                  ;   ii = ii + 1
315         WRITE(ldtxt(ii),*) '          without calling MPI_Init before ! '                        ;   ii = ii + 1
316         kstop = kstop + 1
317      ELSE
318         SELECT CASE ( cn_mpi_send )
319         CASE ( 'S' )                ! Standard mpi send (blocking)
320            WRITE(ldtxt(ii),*) '           Standard blocking mpi send (send)'                     ;   ii = ii + 1
321            CALL mpi_init( ierr )
322         CASE ( 'B' )                ! Buffer mpi send (blocking)
323            WRITE(ldtxt(ii),*) '           Buffer blocking mpi send (bsend)'                      ;   ii = ii + 1
324            IF( Agrif_Root() )   CALL mpi_init_opa( ldtxt, ii, ierr )
325         CASE ( 'I' )                ! Immediate mpi send (non-blocking send)
326            WRITE(ldtxt(ii),*) '           Immediate non-blocking send (isend)'                   ;   ii = ii + 1
327            l_isend = .TRUE.
328            CALL mpi_init( ierr )
329         CASE DEFAULT
330            WRITE(ldtxt(ii),cform_err)                                                            ;   ii = ii + 1
331            WRITE(ldtxt(ii),*) '           bad value for cn_mpi_send = ', cn_mpi_send             ;   ii = ii + 1
332            kstop = kstop + 1
333         END SELECT
334         !
335      ENDIF
336
337      IF( PRESENT(localComm) ) THEN
338         IF( Agrif_Root() ) THEN
339            mpi_comm_opa = localComm
340         ENDIF
341      ELSE
342         CALL mpi_comm_dup( mpi_comm_world, mpi_comm_opa, code)
343         IF( code /= MPI_SUCCESS ) THEN
344            DO ji = 1, SIZE(ldtxt) 
345               IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
346            END DO
347            WRITE(*, cform_err)
348            WRITE(*, *) ' lib_mpp: Error in routine mpi_comm_dup'
349            CALL mpi_abort( mpi_comm_world, code, ierr )
350         ENDIF
351      ENDIF
352
353      CALL mpi_comm_rank( mpi_comm_opa, mpprank, ierr )
354      CALL mpi_comm_size( mpi_comm_opa, mppsize, ierr )
355      mynode = mpprank
356      !
357#if defined key_mpp_rep
358      CALL MPI_OP_CREATE(DDPDD_MPI, .TRUE., MPI_SUMDD, ierr)
359#endif
360      !
361   END FUNCTION mynode
362
363
364   SUBROUTINE mpp_lnk_3d( ptab, cd_type, psgn, cd_mpp, pval )
365      !!----------------------------------------------------------------------
366      !!                  ***  routine mpp_lnk_3d  ***
367      !!
368      !! ** Purpose :   Message passing manadgement
369      !!
370      !! ** Method  :   Use mppsend and mpprecv function for passing mask
371      !!      between processors following neighboring subdomains.
372      !!            domain parameters
373      !!                    nlci   : first dimension of the local subdomain
374      !!                    nlcj   : second dimension of the local subdomain
375      !!                    nbondi : mark for "east-west local boundary"
376      !!                    nbondj : mark for "north-south local boundary"
377      !!                    noea   : number for local neighboring processors
378      !!                    nowe   : number for local neighboring processors
379      !!                    noso   : number for local neighboring processors
380      !!                    nono   : number for local neighboring processors
381      !!
382      !! ** Action  :   ptab with update value at its periphery
383      !!
384      !!----------------------------------------------------------------------
385      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   ptab     ! 3D array on which the boundary condition is applied
386      CHARACTER(len=1)                , INTENT(in   ) ::   cd_type  ! define the nature of ptab array grid-points
387      !                                                             ! = T , U , V , F , W points
388      REAL(wp)                        , INTENT(in   ) ::   psgn     ! =-1 the sign change across the north fold boundary
389      !                                                             ! =  1. , the sign is kept
390      CHARACTER(len=3), OPTIONAL      , INTENT(in   ) ::   cd_mpp   ! fill the overlap area only
391      REAL(wp)        , OPTIONAL      , INTENT(in   ) ::   pval     ! background value (used at closed boundaries)
392      !!
393      INTEGER  ::   ji, jj, jk, jl             ! dummy loop indices
394      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers
395      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend
396      REAL(wp) ::   zland
397      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend
398      !!----------------------------------------------------------------------
399
400      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value
401      ELSE                         ;   zland = 0.e0      ! zero by default
402      ENDIF
403
404      ! 1. standard boundary treatment
405      ! ------------------------------
406      IF( PRESENT( cd_mpp ) ) THEN      ! only fill added line/raw with existing values
407         !
408         ! WARNING ptab is defined only between nld and nle
409         DO jk = 1, jpk
410            DO jj = nlcj+1, jpj                 ! added line(s)   (inner only)
411               ptab(nldi  :nlei  , jj          ,jk) = ptab(nldi:nlei,     nlej,jk)   
412               ptab(1     :nldi-1, jj          ,jk) = ptab(nldi     ,     nlej,jk)
413               ptab(nlei+1:nlci  , jj          ,jk) = ptab(     nlei,     nlej,jk)
414            END DO
415            DO ji = nlci+1, jpi                 ! added column(s) (full)
416               ptab(ji           ,nldj  :nlej  ,jk) = ptab(     nlei,nldj:nlej,jk)
417               ptab(ji           ,1     :nldj-1,jk) = ptab(     nlei,nldj     ,jk)
418               ptab(ji           ,nlej+1:jpj   ,jk) = ptab(     nlei,     nlej,jk)
419            END DO
420         END DO
421         !
422      ELSE                              ! standard close or cyclic treatment
423         !
424         !                                   ! East-West boundaries
425         !                                        !* Cyclic east-west
426         IF( nbondi == 2 .AND. (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN
427            ptab( 1 ,:,:) = ptab(jpim1,:,:)
428            ptab(jpi,:,:) = ptab(  2  ,:,:)
429         ELSE                                     !* closed
430            IF( .NOT. cd_type == 'F' )   ptab(     1       :jpreci,:,:) = zland    ! south except F-point
431                                         ptab(nlci-jpreci+1:jpi   ,:,:) = zland    ! north
432         ENDIF
433         !                                   ! North-South boundaries (always closed)
434         IF( .NOT. cd_type == 'F' )   ptab(:,     1       :jprecj,:) = zland       ! south except F-point
435                                      ptab(:,nlcj-jprecj+1:jpj   ,:) = zland       ! north
436         !
437      ENDIF
438
439      ! 2. East and west directions exchange
440      ! ------------------------------------
441      ! we play with the neigbours AND the row number because of the periodicity
442      !
443      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions
444      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case)
445         iihom = nlci-nreci
446         DO jl = 1, jpreci
447            t3ew(:,jl,:,1) = ptab(jpreci+jl,:,:)
448            t3we(:,jl,:,1) = ptab(iihom +jl,:,:)
449         END DO
450      END SELECT 
451      !
452      !                           ! Migrations
453      imigr = jpreci * jpj * jpk
454      !
455      SELECT CASE ( nbondi ) 
456      CASE ( -1 )
457         CALL mppsend( 2, t3we(1,1,1,1), imigr, noea, ml_req1 )
458         CALL mpprecv( 1, t3ew(1,1,1,2), imigr, noea )
459         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
460      CASE ( 0 )
461         CALL mppsend( 1, t3ew(1,1,1,1), imigr, nowe, ml_req1 )
462         CALL mppsend( 2, t3we(1,1,1,1), imigr, noea, ml_req2 )
463         CALL mpprecv( 1, t3ew(1,1,1,2), imigr, noea )
464         CALL mpprecv( 2, t3we(1,1,1,2), imigr, nowe )
465         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
466         IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err)
467      CASE ( 1 )
468         CALL mppsend( 1, t3ew(1,1,1,1), imigr, nowe, ml_req1 )
469         CALL mpprecv( 2, t3we(1,1,1,2), imigr, nowe )
470         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
471      END SELECT
472      !
473      !                           ! Write Dirichlet lateral conditions
474      iihom = nlci-jpreci
475      !
476      SELECT CASE ( nbondi )
477      CASE ( -1 )
478         DO jl = 1, jpreci
479            ptab(iihom+jl,:,:) = t3ew(:,jl,:,2)
480         END DO
481      CASE ( 0 ) 
482         DO jl = 1, jpreci
483            ptab(jl      ,:,:) = t3we(:,jl,:,2)
484            ptab(iihom+jl,:,:) = t3ew(:,jl,:,2)
485         END DO
486      CASE ( 1 )
487         DO jl = 1, jpreci
488            ptab(jl      ,:,:) = t3we(:,jl,:,2)
489         END DO
490      END SELECT
491
492
493      ! 3. North and south directions
494      ! -----------------------------
495      ! always closed : we play only with the neigbours
496      !
497      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions
498         ijhom = nlcj-nrecj
499         DO jl = 1, jprecj
500            t3sn(:,jl,:,1) = ptab(:,ijhom +jl,:)
501            t3ns(:,jl,:,1) = ptab(:,jprecj+jl,:)
502         END DO
503      ENDIF
504      !
505      !                           ! Migrations
506      imigr = jprecj * jpi * jpk
507      !
508      SELECT CASE ( nbondj )     
509      CASE ( -1 )
510         CALL mppsend( 4, t3sn(1,1,1,1), imigr, nono, ml_req1 )
511         CALL mpprecv( 3, t3ns(1,1,1,2), imigr, nono )
512         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
513      CASE ( 0 )
514         CALL mppsend( 3, t3ns(1,1,1,1), imigr, noso, ml_req1 )
515         CALL mppsend( 4, t3sn(1,1,1,1), imigr, nono, ml_req2 )
516         CALL mpprecv( 3, t3ns(1,1,1,2), imigr, nono )
517         CALL mpprecv( 4, t3sn(1,1,1,2), imigr, noso )
518         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
519         IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err)
520      CASE ( 1 ) 
521         CALL mppsend( 3, t3ns(1,1,1,1), imigr, noso, ml_req1 )
522         CALL mpprecv( 4, t3sn(1,1,1,2), imigr, noso )
523         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
524      END SELECT
525      !
526      !                           ! Write Dirichlet lateral conditions
527      ijhom = nlcj-jprecj
528      !
529      SELECT CASE ( nbondj )
530      CASE ( -1 )
531         DO jl = 1, jprecj
532            ptab(:,ijhom+jl,:) = t3ns(:,jl,:,2)
533         END DO
534      CASE ( 0 ) 
535         DO jl = 1, jprecj
536            ptab(:,jl      ,:) = t3sn(:,jl,:,2)
537            ptab(:,ijhom+jl,:) = t3ns(:,jl,:,2)
538         END DO
539      CASE ( 1 )
540         DO jl = 1, jprecj
541            ptab(:,jl,:) = t3sn(:,jl,:,2)
542         END DO
543      END SELECT
544
545
546      ! 4. north fold treatment
547      ! -----------------------
548      !
549      IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN
550         !
551         SELECT CASE ( jpni )
552         CASE ( 1 )     ;   CALL lbc_nfd      ( ptab, cd_type, psgn )   ! only 1 northern proc, no mpp
553         CASE DEFAULT   ;   CALL mpp_lbc_north( ptab, cd_type, psgn )   ! for all northern procs.
554         END SELECT
555         !
556      ENDIF
557      !
558   END SUBROUTINE mpp_lnk_3d
559
560
561   SUBROUTINE mpp_lnk_2d( pt2d, cd_type, psgn, cd_mpp, pval )
562      !!----------------------------------------------------------------------
563      !!                  ***  routine mpp_lnk_2d  ***
564      !!                 
565      !! ** Purpose :   Message passing manadgement for 2d array
566      !!
567      !! ** Method  :   Use mppsend and mpprecv function for passing mask
568      !!      between processors following neighboring subdomains.
569      !!            domain parameters
570      !!                    nlci   : first dimension of the local subdomain
571      !!                    nlcj   : second dimension of the local subdomain
572      !!                    nbondi : mark for "east-west local boundary"
573      !!                    nbondj : mark for "north-south local boundary"
574      !!                    noea   : number for local neighboring processors
575      !!                    nowe   : number for local neighboring processors
576      !!                    noso   : number for local neighboring processors
577      !!                    nono   : number for local neighboring processors
578      !!
579      !!----------------------------------------------------------------------
580      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   pt2d     ! 2D array on which the boundary condition is applied
581      CHARACTER(len=1)            , INTENT(in   ) ::   cd_type  ! define the nature of ptab array grid-points
582      !                                                         ! = T , U , V , F , W and I points
583      REAL(wp)                    , INTENT(in   ) ::   psgn     ! =-1 the sign change across the north fold boundary
584      !                                                         ! =  1. , the sign is kept
585      CHARACTER(len=3), OPTIONAL  , INTENT(in   ) ::   cd_mpp   ! fill the overlap area only
586      REAL(wp)        , OPTIONAL  , INTENT(in   ) ::   pval     ! background value (used at closed boundaries)
587      !!
588      INTEGER  ::   ji, jj, jl   ! dummy loop indices
589      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers
590      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend
591      REAL(wp) ::   zland
592      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend
593      !!----------------------------------------------------------------------
594
595      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value
596      ELSE                         ;   zland = 0.e0      ! zero by default
597      ENDIF
598
599      ! 1. standard boundary treatment
600      ! ------------------------------
601      !
602      IF( PRESENT( cd_mpp ) ) THEN      ! only fill added line/raw with existing values
603         !
604         ! WARNING pt2d is defined only between nld and nle
605         DO jj = nlcj+1, jpj                 ! added line(s)   (inner only)
606            pt2d(nldi  :nlei  , jj          ) = pt2d(nldi:nlei,     nlej)   
607            pt2d(1     :nldi-1, jj          ) = pt2d(nldi     ,     nlej)
608            pt2d(nlei+1:nlci  , jj          ) = pt2d(     nlei,     nlej)
609         END DO
610         DO ji = nlci+1, jpi                 ! added column(s) (full)
611            pt2d(ji           ,nldj  :nlej  ) = pt2d(     nlei,nldj:nlej)
612            pt2d(ji           ,1     :nldj-1) = pt2d(     nlei,nldj     )
613            pt2d(ji           ,nlej+1:jpj   ) = pt2d(     nlei,     nlej)
614         END DO
615         !
616      ELSE                              ! standard close or cyclic treatment
617         !
618         !                                   ! East-West boundaries
619         IF( nbondi == 2 .AND.   &                ! Cyclic east-west
620            &    (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN
621            pt2d( 1 ,:) = pt2d(jpim1,:)                                    ! west
622            pt2d(jpi,:) = pt2d(  2  ,:)                                    ! east
623         ELSE                                     ! closed
624            IF( .NOT. cd_type == 'F' )   pt2d(     1       :jpreci,:) = zland    ! south except F-point
625                                         pt2d(nlci-jpreci+1:jpi   ,:) = zland    ! north
626         ENDIF
627         !                                   ! North-South boundaries (always closed)
628            IF( .NOT. cd_type == 'F' )   pt2d(:,     1       :jprecj) = zland    !south except F-point
629                                         pt2d(:,nlcj-jprecj+1:jpj   ) = zland    ! north
630         !
631      ENDIF
632
633      ! 2. East and west directions exchange
634      ! ------------------------------------
635      ! we play with the neigbours AND the row number because of the periodicity
636      !
637      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions
638      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case)
639         iihom = nlci-nreci
640         DO jl = 1, jpreci
641            t2ew(:,jl,1) = pt2d(jpreci+jl,:)
642            t2we(:,jl,1) = pt2d(iihom +jl,:)
643         END DO
644      END SELECT
645      !
646      !                           ! Migrations
647      imigr = jpreci * jpj
648      !
649      SELECT CASE ( nbondi )
650      CASE ( -1 )
651         CALL mppsend( 2, t2we(1,1,1), imigr, noea, ml_req1 )
652         CALL mpprecv( 1, t2ew(1,1,2), imigr, noea )
653         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
654      CASE ( 0 )
655         CALL mppsend( 1, t2ew(1,1,1), imigr, nowe, ml_req1 )
656         CALL mppsend( 2, t2we(1,1,1), imigr, noea, ml_req2 )
657         CALL mpprecv( 1, t2ew(1,1,2), imigr, noea )
658         CALL mpprecv( 2, t2we(1,1,2), imigr, nowe )
659         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
660         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err)
661      CASE ( 1 )
662         CALL mppsend( 1, t2ew(1,1,1), imigr, nowe, ml_req1 )
663         CALL mpprecv( 2, t2we(1,1,2), imigr, nowe )
664         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
665      END SELECT
666      !
667      !                           ! Write Dirichlet lateral conditions
668      iihom = nlci - jpreci
669      !
670      SELECT CASE ( nbondi )
671      CASE ( -1 )
672         DO jl = 1, jpreci
673            pt2d(iihom+jl,:) = t2ew(:,jl,2)
674         END DO
675      CASE ( 0 )
676         DO jl = 1, jpreci
677            pt2d(jl      ,:) = t2we(:,jl,2)
678            pt2d(iihom+jl,:) = t2ew(:,jl,2)
679         END DO
680      CASE ( 1 )
681         DO jl = 1, jpreci
682            pt2d(jl      ,:) = t2we(:,jl,2)
683         END DO
684      END SELECT
685
686
687      ! 3. North and south directions
688      ! -----------------------------
689      ! always closed : we play only with the neigbours
690      !
691      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions
692         ijhom = nlcj-nrecj
693         DO jl = 1, jprecj
694            t2sn(:,jl,1) = pt2d(:,ijhom +jl)
695            t2ns(:,jl,1) = pt2d(:,jprecj+jl)
696         END DO
697      ENDIF
698      !
699      !                           ! Migrations
700      imigr = jprecj * jpi
701      !
702      SELECT CASE ( nbondj )
703      CASE ( -1 )
704         CALL mppsend( 4, t2sn(1,1,1), imigr, nono, ml_req1 )
705         CALL mpprecv( 3, t2ns(1,1,2), imigr, nono )
706         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
707      CASE ( 0 )
708         CALL mppsend( 3, t2ns(1,1,1), imigr, noso, ml_req1 )
709         CALL mppsend( 4, t2sn(1,1,1), imigr, nono, ml_req2 )
710         CALL mpprecv( 3, t2ns(1,1,2), imigr, nono )
711         CALL mpprecv( 4, t2sn(1,1,2), imigr, noso )
712         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
713         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err)
714      CASE ( 1 )
715         CALL mppsend( 3, t2ns(1,1,1), imigr, noso, ml_req1 )
716         CALL mpprecv( 4, t2sn(1,1,2), imigr, noso )
717         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
718      END SELECT
719      !
720      !                           ! Write Dirichlet lateral conditions
721      ijhom = nlcj - jprecj
722      !
723      SELECT CASE ( nbondj )
724      CASE ( -1 )
725         DO jl = 1, jprecj
726            pt2d(:,ijhom+jl) = t2ns(:,jl,2)
727         END DO
728      CASE ( 0 )
729         DO jl = 1, jprecj
730            pt2d(:,jl      ) = t2sn(:,jl,2)
731            pt2d(:,ijhom+jl) = t2ns(:,jl,2)
732         END DO
733      CASE ( 1 ) 
734         DO jl = 1, jprecj
735            pt2d(:,jl      ) = t2sn(:,jl,2)
736         END DO
737      END SELECT
738
739
740      ! 4. north fold treatment
741      ! -----------------------
742      !
743      IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN
744         !
745         SELECT CASE ( jpni )
746         CASE ( 1 )     ;   CALL lbc_nfd      ( pt2d, cd_type, psgn )   ! only 1 northern proc, no mpp
747         CASE DEFAULT   ;   CALL mpp_lbc_north( pt2d, cd_type, psgn )   ! for all northern procs.
748         END SELECT
749         !
750      ENDIF
751      !
752   END SUBROUTINE mpp_lnk_2d
753
754
755   SUBROUTINE mpp_lnk_3d_gather( ptab1, cd_type1, ptab2, cd_type2, psgn )
756      !!----------------------------------------------------------------------
757      !!                  ***  routine mpp_lnk_3d_gather  ***
758      !!
759      !! ** Purpose :   Message passing manadgement for two 3D arrays
760      !!
761      !! ** Method  :   Use mppsend and mpprecv function for passing mask
762      !!      between processors following neighboring subdomains.
763      !!            domain parameters
764      !!                    nlci   : first dimension of the local subdomain
765      !!                    nlcj   : second dimension of the local subdomain
766      !!                    nbondi : mark for "east-west local boundary"
767      !!                    nbondj : mark for "north-south local boundary"
768      !!                    noea   : number for local neighboring processors
769      !!                    nowe   : number for local neighboring processors
770      !!                    noso   : number for local neighboring processors
771      !!                    nono   : number for local neighboring processors
772      !!
773      !! ** Action  :   ptab1 and ptab2  with update value at its periphery
774      !!
775      !!----------------------------------------------------------------------
776      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   ptab1     ! first and second 3D array on which
777      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   ptab2     ! the boundary condition is applied
778      CHARACTER(len=1)                , INTENT(in   ) ::   cd_type1  ! nature of ptab1 and ptab2 arrays
779      CHARACTER(len=1)                , INTENT(in   ) ::   cd_type2  ! i.e. grid-points = T , U , V , F or W points
780      REAL(wp)                        , INTENT(in   ) ::   psgn      ! =-1 the sign change across the north fold boundary
781      !!                                                             ! =  1. , the sign is kept
782      INTEGER  ::   jl   ! dummy loop indices
783      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers
784      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend
785      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend
786      !!----------------------------------------------------------------------
787
788      ! 1. standard boundary treatment
789      ! ------------------------------
790      !                                      ! East-West boundaries
791      !                                           !* Cyclic east-west
792      IF( nbondi == 2 .AND. (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN
793         ptab1( 1 ,:,:) = ptab1(jpim1,:,:)
794         ptab1(jpi,:,:) = ptab1(  2  ,:,:)
795         ptab2( 1 ,:,:) = ptab2(jpim1,:,:)
796         ptab2(jpi,:,:) = ptab2(  2  ,:,:)
797      ELSE                                        !* closed
798         IF( .NOT. cd_type1 == 'F' )   ptab1(     1       :jpreci,:,:) = 0.e0    ! south except at F-point
799         IF( .NOT. cd_type2 == 'F' )   ptab2(     1       :jpreci,:,:) = 0.e0
800                                       ptab1(nlci-jpreci+1:jpi   ,:,:) = 0.e0    ! north
801                                       ptab2(nlci-jpreci+1:jpi   ,:,:) = 0.e0
802      ENDIF
803
804     
805      !                                      ! North-South boundaries
806      IF( .NOT. cd_type1 == 'F' )   ptab1(:,     1       :jprecj,:) = 0.e0    ! south except at F-point
807      IF( .NOT. cd_type2 == 'F' )   ptab2(:,     1       :jprecj,:) = 0.e0
808                                    ptab1(:,nlcj-jprecj+1:jpj   ,:) = 0.e0    ! north
809                                    ptab2(:,nlcj-jprecj+1:jpj   ,:) = 0.e0
810
811
812      ! 2. East and west directions exchange
813      ! ------------------------------------
814      ! we play with the neigbours AND the row number because of the periodicity
815      !
816      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions
817      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case)
818         iihom = nlci-nreci
819         DO jl = 1, jpreci
820            t4ew(:,jl,:,1,1) = ptab1(jpreci+jl,:,:)
821            t4we(:,jl,:,1,1) = ptab1(iihom +jl,:,:)
822            t4ew(:,jl,:,2,1) = ptab2(jpreci+jl,:,:)
823            t4we(:,jl,:,2,1) = ptab2(iihom +jl,:,:)
824         END DO
825      END SELECT
826      !
827      !                           ! Migrations
828      imigr = jpreci * jpj * jpk *2
829      !
830      SELECT CASE ( nbondi ) 
831      CASE ( -1 )
832         CALL mppsend( 2, t4we(1,1,1,1,1), imigr, noea, ml_req1 )
833         CALL mpprecv( 1, t4ew(1,1,1,1,2), imigr, noea )
834         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
835      CASE ( 0 )
836         CALL mppsend( 1, t4ew(1,1,1,1,1), imigr, nowe, ml_req1 )
837         CALL mppsend( 2, t4we(1,1,1,1,1), imigr, noea, ml_req2 )
838         CALL mpprecv( 1, t4ew(1,1,1,1,2), imigr, noea )
839         CALL mpprecv( 2, t4we(1,1,1,1,2), imigr, nowe )
840         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
841         IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err)
842      CASE ( 1 )
843         CALL mppsend( 1, t4ew(1,1,1,1,1), imigr, nowe, ml_req1 )
844         CALL mpprecv( 2, t4we(1,1,1,1,2), imigr, nowe )
845         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
846      END SELECT
847      !
848      !                           ! Write Dirichlet lateral conditions
849      iihom = nlci - jpreci
850      !
851      SELECT CASE ( nbondi )
852      CASE ( -1 )
853         DO jl = 1, jpreci
854            ptab1(iihom+jl,:,:) = t4ew(:,jl,:,1,2)
855            ptab2(iihom+jl,:,:) = t4ew(:,jl,:,2,2)
856         END DO
857      CASE ( 0 ) 
858         DO jl = 1, jpreci
859            ptab1(jl      ,:,:) = t4we(:,jl,:,1,2)
860            ptab1(iihom+jl,:,:) = t4ew(:,jl,:,1,2)
861            ptab2(jl      ,:,:) = t4we(:,jl,:,2,2)
862            ptab2(iihom+jl,:,:) = t4ew(:,jl,:,2,2)
863         END DO
864      CASE ( 1 )
865         DO jl = 1, jpreci
866            ptab1(jl      ,:,:) = t4we(:,jl,:,1,2)
867            ptab2(jl      ,:,:) = t4we(:,jl,:,2,2)
868         END DO
869      END SELECT
870
871
872      ! 3. North and south directions
873      ! -----------------------------
874      ! always closed : we play only with the neigbours
875      !
876      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions
877         ijhom = nlcj - nrecj
878         DO jl = 1, jprecj
879            t4sn(:,jl,:,1,1) = ptab1(:,ijhom +jl,:)
880            t4ns(:,jl,:,1,1) = ptab1(:,jprecj+jl,:)
881            t4sn(:,jl,:,2,1) = ptab2(:,ijhom +jl,:)
882            t4ns(:,jl,:,2,1) = ptab2(:,jprecj+jl,:)
883         END DO
884      ENDIF
885      !
886      !                           ! Migrations
887      imigr = jprecj * jpi * jpk * 2
888      !
889      SELECT CASE ( nbondj )     
890      CASE ( -1 )
891         CALL mppsend( 4, t4sn(1,1,1,1,1), imigr, nono, ml_req1 )
892         CALL mpprecv( 3, t4ns(1,1,1,1,2), imigr, nono )
893         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
894      CASE ( 0 )
895         CALL mppsend( 3, t4ns(1,1,1,1,1), imigr, noso, ml_req1 )
896         CALL mppsend( 4, t4sn(1,1,1,1,1), imigr, nono, ml_req2 )
897         CALL mpprecv( 3, t4ns(1,1,1,1,2), imigr, nono )
898         CALL mpprecv( 4, t4sn(1,1,1,1,2), imigr, noso )
899         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
900         IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err)
901      CASE ( 1 ) 
902         CALL mppsend( 3, t4ns(1,1,1,1,1), imigr, noso, ml_req1 )
903         CALL mpprecv( 4, t4sn(1,1,1,1,2), imigr, noso )
904         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err)
905      END SELECT
906      !
907      !                           ! Write Dirichlet lateral conditions
908      ijhom = nlcj - jprecj
909      !
910      SELECT CASE ( nbondj )
911      CASE ( -1 )
912         DO jl = 1, jprecj
913            ptab1(:,ijhom+jl,:) = t4ns(:,jl,:,1,2)
914            ptab2(:,ijhom+jl,:) = t4ns(:,jl,:,2,2)
915         END DO
916      CASE ( 0 ) 
917         DO jl = 1, jprecj
918            ptab1(:,jl      ,:) = t4sn(:,jl,:,1,2)
919            ptab1(:,ijhom+jl,:) = t4ns(:,jl,:,1,2)
920            ptab2(:,jl      ,:) = t4sn(:,jl,:,2,2)
921            ptab2(:,ijhom+jl,:) = t4ns(:,jl,:,2,2)
922         END DO
923      CASE ( 1 )
924         DO jl = 1, jprecj
925            ptab1(:,jl,:) = t4sn(:,jl,:,1,2)
926            ptab2(:,jl,:) = t4sn(:,jl,:,2,2)
927         END DO
928      END SELECT
929
930
931      ! 4. north fold treatment
932      ! -----------------------
933      IF( npolj /= 0 ) THEN
934         !
935         SELECT CASE ( jpni )
936         CASE ( 1 )                                           
937            CALL lbc_nfd      ( ptab1, cd_type1, psgn )   ! only for northern procs.
938            CALL lbc_nfd      ( ptab2, cd_type2, psgn )
939         CASE DEFAULT
940            CALL mpp_lbc_north( ptab1, cd_type1, psgn )   ! for all northern procs.
941            CALL mpp_lbc_north (ptab2, cd_type2, psgn)
942         END SELECT 
943         !
944      ENDIF
945      !
946   END SUBROUTINE mpp_lnk_3d_gather
947
948
949   SUBROUTINE mpp_lnk_2d_e( pt2d, cd_type, psgn )
950      !!----------------------------------------------------------------------
951      !!                  ***  routine mpp_lnk_2d_e  ***
952      !!                 
953      !! ** Purpose :   Message passing manadgement for 2d array (with halo)
954      !!
955      !! ** Method  :   Use mppsend and mpprecv function for passing mask
956      !!      between processors following neighboring subdomains.
957      !!            domain parameters
958      !!                    nlci   : first dimension of the local subdomain
959      !!                    nlcj   : second dimension of the local subdomain
960      !!                    jpr2di : number of rows for extra outer halo
961      !!                    jpr2dj : number of columns for extra outer halo
962      !!                    nbondi : mark for "east-west local boundary"
963      !!                    nbondj : mark for "north-south local boundary"
964      !!                    noea   : number for local neighboring processors
965      !!                    nowe   : number for local neighboring processors
966      !!                    noso   : number for local neighboring processors
967      !!                    nono   : number for local neighboring processors
968      !!
969      !!----------------------------------------------------------------------
970      REAL(wp), DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj), INTENT(inout) ::   pt2d     ! 2D array with extra halo
971      CHARACTER(len=1)                                            , INTENT(in   ) ::   cd_type  ! nature of ptab array grid-points
972      !                                                                                         ! = T , U , V , F , W and I points
973      REAL(wp)                                                    , INTENT(in   ) ::   psgn     ! =-1 the sign change across the
974      !!                                                                                        ! north boundary, =  1. otherwise
975      INTEGER  ::   jl   ! dummy loop indices
976      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers
977      INTEGER  ::   ipreci, iprecj             ! temporary integers
978      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend
979      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend
980      !!----------------------------------------------------------------------
981
982      ipreci = jpreci + jpr2di      ! take into account outer extra 2D overlap area
983      iprecj = jprecj + jpr2dj
984
985
986      ! 1. standard boundary treatment
987      ! ------------------------------
988      ! Order matters Here !!!!
989      !
990      !                                      !* North-South boundaries (always colsed)
991      IF( .NOT. cd_type == 'F' )   pt2d(:,  1-jpr2dj   :  jprecj  ) = 0.e0    ! south except at F-point
992                                   pt2d(:,nlcj-jprecj+1:jpj+jpr2dj) = 0.e0    ! north
993                               
994      !                                      ! East-West boundaries
995      !                                           !* Cyclic east-west
996      IF( nbondi == 2 .AND. (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN
997         pt2d(1-jpr2di:     1    ,:) = pt2d(jpim1-jpr2di:  jpim1 ,:)       ! east
998         pt2d(   jpi  :jpi+jpr2di,:) = pt2d(     2      :2+jpr2di,:)       ! west
999         !
1000      ELSE                                        !* closed
1001         IF( .NOT. cd_type == 'F' )   pt2d(  1-jpr2di   :jpreci    ,:) = 0.e0    ! south except at F-point
1002                                      pt2d(nlci-jpreci+1:jpi+jpr2di,:) = 0.e0    ! north
1003      ENDIF
1004      !
1005
1006      ! north fold treatment
1007      ! -----------------------
1008      IF( npolj /= 0 ) THEN
1009         !
1010         SELECT CASE ( jpni )
1011         CASE ( 1 )     ;   CALL lbc_nfd        ( pt2d(1:jpi,1:jpj+jpr2dj), cd_type, psgn, pr2dj=jpr2dj )
1012         CASE DEFAULT   ;   CALL mpp_lbc_north_e( pt2d                    , cd_type, psgn               )
1013         END SELECT 
1014         !
1015      ENDIF
1016
1017      ! 2. East and west directions exchange
1018      ! ------------------------------------
1019      ! we play with the neigbours AND the row number because of the periodicity
1020      !
1021      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions
1022      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case)
1023         iihom = nlci-nreci-jpr2di
1024         DO jl = 1, ipreci
1025            tr2ew(:,jl,1) = pt2d(jpreci+jl,:)
1026            tr2we(:,jl,1) = pt2d(iihom +jl,:)
1027         END DO
1028      END SELECT
1029      !
1030      !                           ! Migrations
1031      imigr = ipreci * ( jpj + 2*jpr2dj)
1032      !
1033      SELECT CASE ( nbondi )
1034      CASE ( -1 )
1035         CALL mppsend( 2, tr2we(1-jpr2dj,1,1), imigr, noea, ml_req1 )
1036         CALL mpprecv( 1, tr2ew(1-jpr2dj,1,2), imigr, noea )
1037         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1038      CASE ( 0 )
1039         CALL mppsend( 1, tr2ew(1-jpr2dj,1,1), imigr, nowe, ml_req1 )
1040         CALL mppsend( 2, tr2we(1-jpr2dj,1,1), imigr, noea, ml_req2 )
1041         CALL mpprecv( 1, tr2ew(1-jpr2dj,1,2), imigr, noea )
1042         CALL mpprecv( 2, tr2we(1-jpr2dj,1,2), imigr, nowe )
1043         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1044         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err)
1045      CASE ( 1 )
1046         CALL mppsend( 1, tr2ew(1-jpr2dj,1,1), imigr, nowe, ml_req1 )
1047         CALL mpprecv( 2, tr2we(1-jpr2dj,1,2), imigr, nowe )
1048         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1049      END SELECT
1050      !
1051      !                           ! Write Dirichlet lateral conditions
1052      iihom = nlci - jpreci
1053      !
1054      SELECT CASE ( nbondi )
1055      CASE ( -1 )
1056         DO jl = 1, ipreci
1057            pt2d(iihom+jl,:) = tr2ew(:,jl,2)
1058         END DO
1059      CASE ( 0 )
1060         DO jl = 1, ipreci
1061            pt2d(jl-jpr2di,:) = tr2we(:,jl,2)
1062            pt2d( iihom+jl,:) = tr2ew(:,jl,2)
1063         END DO
1064      CASE ( 1 )
1065         DO jl = 1, ipreci
1066            pt2d(jl-jpr2di,:) = tr2we(:,jl,2)
1067         END DO
1068      END SELECT
1069
1070
1071      ! 3. North and south directions
1072      ! -----------------------------
1073      ! always closed : we play only with the neigbours
1074      !
1075      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions
1076         ijhom = nlcj-nrecj-jpr2dj
1077         DO jl = 1, iprecj
1078            tr2sn(:,jl,1) = pt2d(:,ijhom +jl)
1079            tr2ns(:,jl,1) = pt2d(:,jprecj+jl)
1080         END DO
1081      ENDIF
1082      !
1083      !                           ! Migrations
1084      imigr = iprecj * ( jpi + 2*jpr2di )
1085      !
1086      SELECT CASE ( nbondj )
1087      CASE ( -1 )
1088         CALL mppsend( 4, tr2sn(1-jpr2di,1,1), imigr, nono, ml_req1 )
1089         CALL mpprecv( 3, tr2ns(1-jpr2di,1,2), imigr, nono )
1090         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1091      CASE ( 0 )
1092         CALL mppsend( 3, tr2ns(1-jpr2di,1,1), imigr, noso, ml_req1 )
1093         CALL mppsend( 4, tr2sn(1-jpr2di,1,1), imigr, nono, ml_req2 )
1094         CALL mpprecv( 3, tr2ns(1-jpr2di,1,2), imigr, nono )
1095         CALL mpprecv( 4, tr2sn(1-jpr2di,1,2), imigr, noso )
1096         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1097         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err)
1098      CASE ( 1 )
1099         CALL mppsend( 3, tr2ns(1-jpr2di,1,1), imigr, noso, ml_req1 )
1100         CALL mpprecv( 4, tr2sn(1-jpr2di,1,2), imigr, noso )
1101         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err)
1102      END SELECT
1103      !
1104      !                           ! Write Dirichlet lateral conditions
1105      ijhom = nlcj - jprecj 
1106      !
1107      SELECT CASE ( nbondj )
1108      CASE ( -1 )
1109         DO jl = 1, iprecj
1110            pt2d(:,ijhom+jl) = tr2ns(:,jl,2)
1111         END DO
1112      CASE ( 0 )
1113         DO jl = 1, iprecj
1114            pt2d(:,jl-jpr2dj) = tr2sn(:,jl,2)
1115            pt2d(:,ijhom+jl ) = tr2ns(:,jl,2)
1116         END DO
1117      CASE ( 1 ) 
1118         DO jl = 1, iprecj
1119            pt2d(:,jl-jpr2dj) = tr2sn(:,jl,2)
1120         END DO
1121      END SELECT
1122
1123   END SUBROUTINE mpp_lnk_2d_e
1124
1125
1126   SUBROUTINE mppsend( ktyp, pmess, kbytes, kdest, md_req )
1127      !!----------------------------------------------------------------------
1128      !!                  ***  routine mppsend  ***
1129      !!                   
1130      !! ** Purpose :   Send messag passing array
1131      !!
1132      !!----------------------------------------------------------------------
1133      REAL(wp), INTENT(inout) ::   pmess(*)   ! array of real
1134      INTEGER , INTENT(in   ) ::   kbytes     ! size of the array pmess
1135      INTEGER , INTENT(in   ) ::   kdest      ! receive process number
1136      INTEGER , INTENT(in   ) ::   ktyp       ! tag of the message
1137      INTEGER , INTENT(in   ) ::   md_req     ! argument for isend
1138      !!
1139      INTEGER ::   iflag
1140      !!----------------------------------------------------------------------
1141      !
1142      SELECT CASE ( cn_mpi_send )
1143      CASE ( 'S' )                ! Standard mpi send (blocking)
1144         CALL mpi_send ( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_opa        , iflag )
1145      CASE ( 'B' )                ! Buffer mpi send (blocking)
1146         CALL mpi_bsend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_opa        , iflag )
1147      CASE ( 'I' )                ! Immediate mpi send (non-blocking send)
1148         ! be carefull, one more argument here : the mpi request identifier..
1149         CALL mpi_isend( pmess, kbytes, mpi_double_precision, kdest , ktyp, mpi_comm_opa, md_req, iflag )
1150      END SELECT
1151      !
1152   END SUBROUTINE mppsend
1153
1154
1155   SUBROUTINE mpprecv( ktyp, pmess, kbytes, ksource )
1156      !!----------------------------------------------------------------------
1157      !!                  ***  routine mpprecv  ***
1158      !!
1159      !! ** Purpose :   Receive messag passing array
1160      !!
1161      !!----------------------------------------------------------------------
1162      REAL(wp), INTENT(inout) ::   pmess(*)   ! array of real
1163      INTEGER , INTENT(in   ) ::   kbytes     ! suze of the array pmess
1164      INTEGER , INTENT(in   ) ::   ktyp       ! Tag of the recevied message
1165      INTEGER, OPTIONAL, INTENT(in) :: ksource    ! source process number
1166      !!
1167      INTEGER :: istatus(mpi_status_size)
1168      INTEGER :: iflag
1169      INTEGER :: use_source
1170      !!----------------------------------------------------------------------
1171      !
1172
1173      ! If a specific process number has been passed to the receive call,
1174      ! use that one. Default is to use mpi_any_source
1175      use_source=mpi_any_source
1176      if(present(ksource)) then
1177         use_source=ksource
1178      end if
1179
1180      CALL mpi_recv( pmess, kbytes, mpi_double_precision, use_source, ktyp, mpi_comm_opa, istatus, iflag )
1181      !
1182   END SUBROUTINE mpprecv
1183
1184
1185   SUBROUTINE mppgather( ptab, kp, pio )
1186      !!----------------------------------------------------------------------
1187      !!                   ***  routine mppgather  ***
1188      !!                   
1189      !! ** Purpose :   Transfert between a local subdomain array and a work
1190      !!     array which is distributed following the vertical level.
1191      !!
1192      !!----------------------------------------------------------------------
1193      REAL(wp), DIMENSION(jpi,jpj),       INTENT(in   ) ::   ptab   ! subdomain input array
1194      INTEGER ,                           INTENT(in   ) ::   kp     ! record length
1195      REAL(wp), DIMENSION(jpi,jpj,jpnij), INTENT(  out) ::   pio    ! subdomain input array
1196      !!
1197      INTEGER :: itaille, ierror   ! temporary integer
1198      !!---------------------------------------------------------------------
1199      !
1200      itaille = jpi * jpj
1201      CALL mpi_gather( ptab, itaille, mpi_double_precision, pio, itaille     ,   &
1202         &                            mpi_double_precision, kp , mpi_comm_opa, ierror ) 
1203      !
1204   END SUBROUTINE mppgather
1205
1206
1207   SUBROUTINE mppscatter( pio, kp, ptab )
1208      !!----------------------------------------------------------------------
1209      !!                  ***  routine mppscatter  ***
1210      !!
1211      !! ** Purpose :   Transfert between awork array which is distributed
1212      !!      following the vertical level and the local subdomain array.
1213      !!
1214      !!----------------------------------------------------------------------
1215      REAL(wp), DIMENSION(jpi,jpj,jpnij)  ::  pio        ! output array
1216      INTEGER                             ::   kp        ! Tag (not used with MPI
1217      REAL(wp), DIMENSION(jpi,jpj)        ::  ptab       ! subdomain array input
1218      !!
1219      INTEGER :: itaille, ierror   ! temporary integer
1220      !!---------------------------------------------------------------------
1221      !
1222      itaille=jpi*jpj
1223      !
1224      CALL mpi_scatter( pio, itaille, mpi_double_precision, ptab, itaille     ,   &
1225         &                            mpi_double_precision, kp  , mpi_comm_opa, ierror )
1226      !
1227   END SUBROUTINE mppscatter
1228
1229
1230   SUBROUTINE mppmax_a_int( ktab, kdim, kcom )
1231      !!----------------------------------------------------------------------
1232      !!                  ***  routine mppmax_a_int  ***
1233      !!
1234      !! ** Purpose :   Find maximum value in an integer layout array
1235      !!
1236      !!----------------------------------------------------------------------
1237      INTEGER , INTENT(in   )                  ::   kdim   ! size of array
1238      INTEGER , INTENT(inout), DIMENSION(kdim) ::   ktab   ! input array
1239      INTEGER , INTENT(in   ), OPTIONAL        ::   kcom   !
1240      !!
1241      INTEGER :: ierror, localcomm   ! temporary integer
1242      INTEGER, DIMENSION(kdim) ::   iwork
1243      !!----------------------------------------------------------------------
1244      !
1245      localcomm = mpi_comm_opa
1246      IF( PRESENT(kcom) )   localcomm = kcom
1247      !
1248      CALL mpi_allreduce( ktab, iwork, kdim, mpi_integer, mpi_max, localcomm, ierror )
1249      !
1250      ktab(:) = iwork(:)
1251      !
1252   END SUBROUTINE mppmax_a_int
1253
1254
1255   SUBROUTINE mppmax_int( ktab, kcom )
1256      !!----------------------------------------------------------------------
1257      !!                  ***  routine mppmax_int  ***
1258      !!
1259      !! ** Purpose :   Find maximum value in an integer layout array
1260      !!
1261      !!----------------------------------------------------------------------
1262      INTEGER, INTENT(inout)           ::   ktab      ! ???
1263      INTEGER, INTENT(in   ), OPTIONAL ::   kcom      ! ???
1264      !!
1265      INTEGER ::   ierror, iwork, localcomm   ! temporary integer
1266      !!----------------------------------------------------------------------
1267      !
1268      localcomm = mpi_comm_opa 
1269      IF( PRESENT(kcom) )   localcomm = kcom
1270      !
1271      CALL mpi_allreduce( ktab, iwork, 1, mpi_integer, mpi_max, localcomm, ierror)
1272      !
1273      ktab = iwork
1274      !
1275   END SUBROUTINE mppmax_int
1276
1277
1278   SUBROUTINE mppmin_a_int( ktab, kdim, kcom )
1279      !!----------------------------------------------------------------------
1280      !!                  ***  routine mppmin_a_int  ***
1281      !!
1282      !! ** Purpose :   Find minimum value in an integer layout array
1283      !!
1284      !!----------------------------------------------------------------------
1285      INTEGER , INTENT( in  )                  ::   kdim        ! size of array
1286      INTEGER , INTENT(inout), DIMENSION(kdim) ::   ktab        ! input array
1287      INTEGER , INTENT( in  ), OPTIONAL        ::   kcom        ! input array
1288      !!
1289      INTEGER ::   ierror, localcomm   ! temporary integer
1290      INTEGER, DIMENSION(kdim) ::   iwork
1291      !!----------------------------------------------------------------------
1292      !
1293      localcomm = mpi_comm_opa
1294      IF( PRESENT(kcom) )   localcomm = kcom
1295      !
1296      CALL mpi_allreduce( ktab, iwork, kdim, mpi_integer, mpi_min, localcomm, ierror )
1297      !
1298      ktab(:) = iwork(:)
1299      !
1300   END SUBROUTINE mppmin_a_int
1301
1302
1303   SUBROUTINE mppmin_int( ktab, kcom )
1304      !!----------------------------------------------------------------------
1305      !!                  ***  routine mppmin_int  ***
1306      !!
1307      !! ** Purpose :   Find minimum value in an integer layout array
1308      !!
1309      !!----------------------------------------------------------------------
1310      INTEGER, INTENT(inout) ::   ktab      ! ???
1311      INTEGER , INTENT( in  ), OPTIONAL        ::   kcom        ! input array
1312      !!
1313      INTEGER ::  ierror, iwork, localcomm
1314      !!----------------------------------------------------------------------
1315      !
1316      localcomm = mpi_comm_opa
1317      IF( PRESENT(kcom) )   localcomm = kcom
1318      !
1319     CALL mpi_allreduce( ktab, iwork, 1, mpi_integer, mpi_min, localcomm, ierror )
1320      !
1321      ktab = iwork
1322      !
1323   END SUBROUTINE mppmin_int
1324
1325
1326   SUBROUTINE mppsum_a_int( ktab, kdim )
1327      !!----------------------------------------------------------------------
1328      !!                  ***  routine mppsum_a_int  ***
1329      !!                   
1330      !! ** Purpose :   Global integer sum, 1D array case
1331      !!
1332      !!----------------------------------------------------------------------
1333      INTEGER, INTENT(in   )                   ::   kdim      ! ???
1334      INTEGER, INTENT(inout), DIMENSION (kdim) ::   ktab      ! ???
1335      !!
1336      INTEGER :: ierror
1337      INTEGER, DIMENSION (kdim) ::  iwork
1338      !!----------------------------------------------------------------------
1339      !
1340      CALL mpi_allreduce( ktab, iwork, kdim, mpi_integer, mpi_sum, mpi_comm_opa, ierror )
1341      !
1342      ktab(:) = iwork(:)
1343      !
1344   END SUBROUTINE mppsum_a_int
1345
1346
1347   SUBROUTINE mppsum_int( ktab )
1348      !!----------------------------------------------------------------------
1349      !!                 ***  routine mppsum_int  ***
1350      !!                 
1351      !! ** Purpose :   Global integer sum
1352      !!
1353      !!----------------------------------------------------------------------
1354      INTEGER, INTENT(inout) ::   ktab
1355      !!
1356      INTEGER :: ierror, iwork
1357      !!----------------------------------------------------------------------
1358      !
1359      CALL mpi_allreduce( ktab, iwork, 1, mpi_integer, mpi_sum, mpi_comm_opa, ierror )
1360      !
1361      ktab = iwork
1362      !
1363   END SUBROUTINE mppsum_int
1364
1365
1366   SUBROUTINE mppmax_a_real( ptab, kdim, kcom )
1367      !!----------------------------------------------------------------------
1368      !!                 ***  routine mppmax_a_real  ***
1369      !!                 
1370      !! ** Purpose :   Maximum
1371      !!
1372      !!----------------------------------------------------------------------
1373      INTEGER , INTENT(in   )                  ::   kdim
1374      REAL(wp), INTENT(inout), DIMENSION(kdim) ::   ptab
1375      INTEGER , INTENT(in   ), OPTIONAL        ::   kcom
1376      !!
1377      INTEGER :: ierror, localcomm
1378      REAL(wp), DIMENSION(kdim) ::  zwork
1379      !!----------------------------------------------------------------------
1380      !
1381      localcomm = mpi_comm_opa
1382      IF( PRESENT(kcom) ) localcomm = kcom
1383      !
1384      CALL mpi_allreduce( ptab, zwork, kdim, mpi_double_precision, mpi_max, localcomm, ierror )
1385      ptab(:) = zwork(:)
1386      !
1387   END SUBROUTINE mppmax_a_real
1388
1389
1390   SUBROUTINE mppmax_real( ptab, kcom )
1391      !!----------------------------------------------------------------------
1392      !!                  ***  routine mppmax_real  ***
1393      !!                   
1394      !! ** Purpose :   Maximum
1395      !!
1396      !!----------------------------------------------------------------------
1397      REAL(wp), INTENT(inout)           ::   ptab   ! ???
1398      INTEGER , INTENT(in   ), OPTIONAL ::   kcom   ! ???
1399      !!
1400      INTEGER  ::   ierror, localcomm
1401      REAL(wp) ::   zwork
1402      !!----------------------------------------------------------------------
1403      !
1404      localcomm = mpi_comm_opa 
1405      IF( PRESENT(kcom) )   localcomm = kcom
1406      !
1407      CALL mpi_allreduce( ptab, zwork, 1, mpi_double_precision, mpi_max, localcomm, ierror )
1408      ptab = zwork
1409      !
1410   END SUBROUTINE mppmax_real
1411
1412
1413   SUBROUTINE mppmin_a_real( ptab, kdim, kcom )
1414      !!----------------------------------------------------------------------
1415      !!                 ***  routine mppmin_a_real  ***
1416      !!                 
1417      !! ** Purpose :   Minimum of REAL, array case
1418      !!
1419      !!-----------------------------------------------------------------------
1420      INTEGER , INTENT(in   )                  ::   kdim
1421      REAL(wp), INTENT(inout), DIMENSION(kdim) ::   ptab
1422      INTEGER , INTENT(in   ), OPTIONAL        ::   kcom
1423      !!
1424      INTEGER :: ierror, localcomm
1425      REAL(wp), DIMENSION(kdim) ::   zwork
1426      !!-----------------------------------------------------------------------
1427      !
1428      localcomm = mpi_comm_opa 
1429      IF( PRESENT(kcom) ) localcomm = kcom
1430      !
1431      CALL mpi_allreduce( ptab, zwork, kdim, mpi_double_precision, mpi_min, localcomm, ierror )
1432      ptab(:) = zwork(:)
1433      !
1434   END SUBROUTINE mppmin_a_real
1435
1436
1437   SUBROUTINE mppmin_real( ptab, kcom )
1438      !!----------------------------------------------------------------------
1439      !!                  ***  routine mppmin_real  ***
1440      !!
1441      !! ** Purpose :   minimum of REAL, scalar case
1442      !!
1443      !!-----------------------------------------------------------------------
1444      REAL(wp), INTENT(inout)           ::   ptab        !
1445      INTEGER , INTENT(in   ), OPTIONAL :: kcom
1446      !!
1447      INTEGER  ::   ierror
1448      REAL(wp) ::   zwork
1449      INTEGER :: localcomm
1450      !!-----------------------------------------------------------------------
1451      !
1452      localcomm = mpi_comm_opa 
1453      IF( PRESENT(kcom) )   localcomm = kcom
1454      !
1455      CALL mpi_allreduce( ptab, zwork, 1, mpi_double_precision, mpi_min, localcomm, ierror )
1456      ptab = zwork
1457      !
1458   END SUBROUTINE mppmin_real
1459
1460
1461   SUBROUTINE mppsum_a_real( ptab, kdim, kcom )
1462      !!----------------------------------------------------------------------
1463      !!                  ***  routine mppsum_a_real  ***
1464      !!
1465      !! ** Purpose :   global sum, REAL ARRAY argument case
1466      !!
1467      !!-----------------------------------------------------------------------
1468      INTEGER , INTENT( in )                     ::   kdim      ! size of ptab
1469      REAL(wp), DIMENSION(kdim), INTENT( inout ) ::   ptab      ! input array
1470      INTEGER , INTENT( in ), OPTIONAL           :: kcom
1471      !!
1472      INTEGER                   ::   ierror    ! temporary integer
1473      INTEGER                   ::   localcomm 
1474      REAL(wp), DIMENSION(kdim) ::   zwork     ! temporary workspace
1475      !!-----------------------------------------------------------------------
1476      !
1477      localcomm = mpi_comm_opa 
1478      IF( PRESENT(kcom) )   localcomm = kcom
1479      !
1480      CALL mpi_allreduce( ptab, zwork, kdim, mpi_double_precision, mpi_sum, localcomm, ierror )
1481      ptab(:) = zwork(:)
1482      !
1483   END SUBROUTINE mppsum_a_real
1484
1485
1486   SUBROUTINE mppsum_real( ptab, kcom )
1487      !!----------------------------------------------------------------------
1488      !!                  ***  routine mppsum_real  ***
1489      !!             
1490      !! ** Purpose :   global sum, SCALAR argument case
1491      !!
1492      !!-----------------------------------------------------------------------
1493      REAL(wp), INTENT(inout)           ::   ptab   ! input scalar
1494      INTEGER , INTENT(in   ), OPTIONAL ::   kcom
1495      !!
1496      INTEGER  ::   ierror, localcomm 
1497      REAL(wp) ::   zwork
1498      !!-----------------------------------------------------------------------
1499      !
1500      localcomm = mpi_comm_opa 
1501      IF( PRESENT(kcom) ) localcomm = kcom
1502      !
1503      CALL mpi_allreduce( ptab, zwork, 1, mpi_double_precision, mpi_sum, localcomm, ierror )
1504      ptab = zwork
1505      !
1506   END SUBROUTINE mppsum_real
1507
1508# if defined key_mpp_rep
1509   SUBROUTINE mppsum_realdd( ytab, kcom )
1510      !!----------------------------------------------------------------------
1511      !!                  ***  routine mppsum_realdd ***
1512      !!
1513      !! ** Purpose :   global sum in Massively Parallel Processing
1514      !!                SCALAR argument case for double-double precision
1515      !!
1516      !!-----------------------------------------------------------------------
1517      COMPLEX(wp), INTENT(inout)         :: ytab    ! input scalar
1518      INTEGER , INTENT( in  ), OPTIONAL :: kcom
1519
1520      !! * Local variables   (MPI version)
1521      INTEGER  ::    ierror
1522      INTEGER  ::   localcomm
1523      COMPLEX(wp) :: zwork
1524
1525      localcomm = mpi_comm_opa
1526      IF( PRESENT(kcom) ) localcomm = kcom
1527
1528      ! reduce local sums into global sum
1529      CALL MPI_ALLREDUCE (ytab, zwork, 1, MPI_DOUBLE_COMPLEX, &
1530                       MPI_SUMDD,localcomm,ierror)
1531      ytab = zwork
1532
1533   END SUBROUTINE mppsum_realdd
1534 
1535 
1536   SUBROUTINE mppsum_a_realdd( ytab, kdim, kcom )
1537      !!----------------------------------------------------------------------
1538      !!                  ***  routine mppsum_a_realdd  ***
1539      !!
1540      !! ** Purpose :   global sum in Massively Parallel Processing
1541      !!                COMPLEX ARRAY case for double-double precision
1542      !!
1543      !!-----------------------------------------------------------------------
1544      INTEGER , INTENT( in )                        ::   kdim      ! size of ytab
1545      COMPLEX(wp), DIMENSION(kdim), INTENT( inout ) ::   ytab      ! input array
1546      INTEGER , INTENT( in  ), OPTIONAL :: kcom
1547
1548      !! * Local variables   (MPI version)
1549      INTEGER                      :: ierror    ! temporary integer
1550      INTEGER                      ::   localcomm
1551      COMPLEX(wp), DIMENSION(kdim) :: zwork     ! temporary workspace
1552
1553      localcomm = mpi_comm_opa
1554      IF( PRESENT(kcom) ) localcomm = kcom
1555
1556      CALL MPI_ALLREDUCE (ytab, zwork, kdim, MPI_DOUBLE_COMPLEX, &
1557                       MPI_SUMDD,localcomm,ierror)
1558      ytab(:) = zwork(:)
1559
1560   END SUBROUTINE mppsum_a_realdd
1561# endif   
1562   
1563   SUBROUTINE mpp_minloc2d( ptab, pmask, pmin, ki,kj )
1564      !!------------------------------------------------------------------------
1565      !!             ***  routine mpp_minloc  ***
1566      !!
1567      !! ** Purpose :   Compute the global minimum of an array ptab
1568      !!              and also give its global position
1569      !!
1570      !! ** Method  :   Use MPI_ALLREDUCE with MPI_MINLOC
1571      !!
1572      !!--------------------------------------------------------------------------
1573      REAL(wp), DIMENSION (jpi,jpj), INTENT(in   ) ::   ptab    ! Local 2D array
1574      REAL(wp), DIMENSION (jpi,jpj), INTENT(in   ) ::   pmask   ! Local mask
1575      REAL(wp)                     , INTENT(  out) ::   pmin    ! Global minimum of ptab
1576      INTEGER                      , INTENT(  out) ::   ki, kj   ! index of minimum in global frame
1577      !!
1578      INTEGER , DIMENSION(2)   ::   ilocs
1579      INTEGER :: ierror
1580      REAL(wp) ::   zmin   ! local minimum
1581      REAL(wp), DIMENSION(2,1) ::   zain, zaout
1582      !!-----------------------------------------------------------------------
1583      !
1584      zmin  = MINVAL( ptab(:,:) , mask= pmask == 1.e0 )
1585      ilocs = MINLOC( ptab(:,:) , mask= pmask == 1.e0 )
1586      !
1587      ki = ilocs(1) + nimpp - 1
1588      kj = ilocs(2) + njmpp - 1
1589      !
1590      zain(1,:)=zmin
1591      zain(2,:)=ki+10000.*kj
1592      !
1593      CALL MPI_ALLREDUCE( zain,zaout, 1, MPI_2DOUBLE_PRECISION,MPI_MINLOC,MPI_COMM_OPA,ierror)
1594      !
1595      pmin = zaout(1,1)
1596      kj = INT(zaout(2,1)/10000.)
1597      ki = INT(zaout(2,1) - 10000.*kj )
1598      !
1599   END SUBROUTINE mpp_minloc2d
1600
1601
1602   SUBROUTINE mpp_minloc3d( ptab, pmask, pmin, ki, kj ,kk)
1603      !!------------------------------------------------------------------------
1604      !!             ***  routine mpp_minloc  ***
1605      !!
1606      !! ** Purpose :   Compute the global minimum of an array ptab
1607      !!              and also give its global position
1608      !!
1609      !! ** Method  :   Use MPI_ALLREDUCE with MPI_MINLOC
1610      !!
1611      !!--------------------------------------------------------------------------
1612      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   ptab         ! Local 2D array
1613      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   pmask        ! Local mask
1614      REAL(wp)                         , INTENT(  out) ::   pmin         ! Global minimum of ptab
1615      INTEGER                          , INTENT(  out) ::   ki, kj, kk   ! index of minimum in global frame
1616      !!
1617      INTEGER  ::   ierror
1618      REAL(wp) ::   zmin     ! local minimum
1619      INTEGER , DIMENSION(3)   ::   ilocs
1620      REAL(wp), DIMENSION(2,1) ::   zain, zaout
1621      !!-----------------------------------------------------------------------
1622      !
1623      zmin  = MINVAL( ptab(:,:,:) , mask= pmask == 1.e0 )
1624      ilocs = MINLOC( ptab(:,:,:) , mask= pmask == 1.e0 )
1625      !
1626      ki = ilocs(1) + nimpp - 1
1627      kj = ilocs(2) + njmpp - 1
1628      kk = ilocs(3)
1629      !
1630      zain(1,:)=zmin
1631      zain(2,:)=ki+10000.*kj+100000000.*kk
1632      !
1633      CALL MPI_ALLREDUCE( zain,zaout, 1, MPI_2DOUBLE_PRECISION,MPI_MINLOC,MPI_COMM_OPA,ierror)
1634      !
1635      pmin = zaout(1,1)
1636      kk   = INT( zaout(2,1) / 100000000. )
1637      kj   = INT( zaout(2,1) - kk * 100000000. ) / 10000
1638      ki   = INT( zaout(2,1) - kk * 100000000. -kj * 10000. )
1639      !
1640   END SUBROUTINE mpp_minloc3d
1641
1642
1643   SUBROUTINE mpp_maxloc2d( ptab, pmask, pmax, ki, kj )
1644      !!------------------------------------------------------------------------
1645      !!             ***  routine mpp_maxloc  ***
1646      !!
1647      !! ** Purpose :   Compute the global maximum of an array ptab
1648      !!              and also give its global position
1649      !!
1650      !! ** Method  :   Use MPI_ALLREDUCE with MPI_MINLOC
1651      !!
1652      !!--------------------------------------------------------------------------
1653      REAL(wp), DIMENSION (jpi,jpj), INTENT(in   ) ::   ptab     ! Local 2D array
1654      REAL(wp), DIMENSION (jpi,jpj), INTENT(in   ) ::   pmask    ! Local mask
1655      REAL(wp)                     , INTENT(  out) ::   pmax     ! Global maximum of ptab
1656      INTEGER                      , INTENT(  out) ::   ki, kj   ! index of maximum in global frame
1657      !! 
1658      INTEGER  :: ierror
1659      INTEGER, DIMENSION (2)   ::   ilocs
1660      REAL(wp) :: zmax   ! local maximum
1661      REAL(wp), DIMENSION(2,1) ::   zain, zaout
1662      !!-----------------------------------------------------------------------
1663      !
1664      zmax  = MAXVAL( ptab(:,:) , mask= pmask == 1.e0 )
1665      ilocs = MAXLOC( ptab(:,:) , mask= pmask == 1.e0 )
1666      !
1667      ki = ilocs(1) + nimpp - 1
1668      kj = ilocs(2) + njmpp - 1
1669      !
1670      zain(1,:) = zmax
1671      zain(2,:) = ki + 10000. * kj
1672      !
1673      CALL MPI_ALLREDUCE( zain,zaout, 1, MPI_2DOUBLE_PRECISION,MPI_MAXLOC,MPI_COMM_OPA,ierror)
1674      !
1675      pmax = zaout(1,1)
1676      kj   = INT( zaout(2,1) / 10000.     )
1677      ki   = INT( zaout(2,1) - 10000.* kj )
1678      !
1679   END SUBROUTINE mpp_maxloc2d
1680
1681
1682   SUBROUTINE mpp_maxloc3d( ptab, pmask, pmax, ki, kj, kk )
1683      !!------------------------------------------------------------------------
1684      !!             ***  routine mpp_maxloc  ***
1685      !!
1686      !! ** Purpose :  Compute the global maximum of an array ptab
1687      !!              and also give its global position
1688      !!
1689      !! ** Method : Use MPI_ALLREDUCE with MPI_MINLOC
1690      !!
1691      !!--------------------------------------------------------------------------
1692      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   ptab         ! Local 2D array
1693      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   pmask        ! Local mask
1694      REAL(wp)                         , INTENT(  out) ::   pmax         ! Global maximum of ptab
1695      INTEGER                          , INTENT(  out) ::   ki, kj, kk   ! index of maximum in global frame
1696      !!   
1697      REAL(wp) :: zmax   ! local maximum
1698      REAL(wp), DIMENSION(2,1) ::   zain, zaout
1699      INTEGER , DIMENSION(3)   ::   ilocs
1700      INTEGER :: ierror
1701      !!-----------------------------------------------------------------------
1702      !
1703      zmax  = MAXVAL( ptab(:,:,:) , mask= pmask == 1.e0 )
1704      ilocs = MAXLOC( ptab(:,:,:) , mask= pmask == 1.e0 )
1705      !
1706      ki = ilocs(1) + nimpp - 1
1707      kj = ilocs(2) + njmpp - 1
1708      kk = ilocs(3)
1709      !
1710      zain(1,:)=zmax
1711      zain(2,:)=ki+10000.*kj+100000000.*kk
1712      !
1713      CALL MPI_ALLREDUCE( zain,zaout, 1, MPI_2DOUBLE_PRECISION,MPI_MAXLOC,MPI_COMM_OPA,ierror)
1714      !
1715      pmax = zaout(1,1)
1716      kk   = INT( zaout(2,1) / 100000000. )
1717      kj   = INT( zaout(2,1) - kk * 100000000. ) / 10000
1718      ki   = INT( zaout(2,1) - kk * 100000000. -kj * 10000. )
1719      !
1720   END SUBROUTINE mpp_maxloc3d
1721
1722
1723   SUBROUTINE mppsync()
1724      !!----------------------------------------------------------------------
1725      !!                  ***  routine mppsync  ***
1726      !!                   
1727      !! ** Purpose :   Massively parallel processors, synchroneous
1728      !!
1729      !!-----------------------------------------------------------------------
1730      INTEGER :: ierror
1731      !!-----------------------------------------------------------------------
1732      !
1733      CALL mpi_barrier( mpi_comm_opa, ierror )
1734      !
1735   END SUBROUTINE mppsync
1736
1737
1738   SUBROUTINE mppstop
1739      !!----------------------------------------------------------------------
1740      !!                  ***  routine mppstop  ***
1741      !!                   
1742      !! ** purpose :   Stop massively parallel processors method
1743      !!
1744      !!----------------------------------------------------------------------
1745      INTEGER ::   info
1746      !!----------------------------------------------------------------------
1747      !
1748      CALL mppsync
1749      CALL mpi_finalize( info )
1750      !
1751   END SUBROUTINE mppstop
1752
1753
1754   SUBROUTINE mppobc( ptab, kd1, kd2, kl, kk, ktype, kij , kumout)
1755      !!----------------------------------------------------------------------
1756      !!                  ***  routine mppobc  ***
1757      !!
1758      !! ** Purpose :   Message passing manadgement for open boundary
1759      !!     conditions array
1760      !!
1761      !! ** Method  :   Use mppsend and mpprecv function for passing mask
1762      !!       between processors following neighboring subdomains.
1763      !!       domain parameters
1764      !!                    nlci   : first dimension of the local subdomain
1765      !!                    nlcj   : second dimension of the local subdomain
1766      !!                    nbondi : mark for "east-west local boundary"
1767      !!                    nbondj : mark for "north-south local boundary"
1768      !!                    noea   : number for local neighboring processors
1769      !!                    nowe   : number for local neighboring processors
1770      !!                    noso   : number for local neighboring processors
1771      !!                    nono   : number for local neighboring processors
1772      !!
1773      !!----------------------------------------------------------------------
1774      USE wrk_nemo        ! Memory allocation
1775      !
1776      INTEGER , INTENT(in   )                     ::   kd1, kd2   ! starting and ending indices
1777      INTEGER , INTENT(in   )                     ::   kl         ! index of open boundary
1778      INTEGER , INTENT(in   )                     ::   kk         ! vertical dimension
1779      INTEGER , INTENT(in   )                     ::   ktype      ! define north/south or east/west cdt
1780      !                                                           !  = 1  north/south  ;  = 2  east/west
1781      INTEGER , INTENT(in   )                     ::   kij        ! horizontal dimension
1782      INTEGER , INTENT(in   )                     ::   kumout     ! ocean.output logical unit
1783      REAL(wp), INTENT(inout), DIMENSION(kij,kk)  ::   ptab       ! variable array
1784      !
1785      INTEGER ::   ji, jj, jk, jl        ! dummy loop indices
1786      INTEGER ::   iipt0, iipt1, ilpt1   ! local integers
1787      INTEGER ::   ijpt0, ijpt1          !   -       -
1788      INTEGER ::   imigr, iihom, ijhom   !   -       -
1789      INTEGER ::   ml_req1, ml_req2, ml_err    ! for key_mpi_isend
1790      INTEGER ::   ml_stat(MPI_STATUS_SIZE)    ! for key_mpi_isend
1791      REAL(wp), POINTER, DIMENSION(:,:) ::   ztab   ! temporary workspace
1792      !!----------------------------------------------------------------------
1793
1794      CALL wrk_alloc( jpi,jpj, ztab )
1795
1796      ! boundary condition initialization
1797      ! ---------------------------------
1798      ztab(:,:) = 0.e0
1799      !
1800      IF( ktype==1 ) THEN                                  ! north/south boundaries
1801         iipt0 = MAX( 1, MIN(kd1 - nimpp+1, nlci     ) )
1802         iipt1 = MAX( 0, MIN(kd2 - nimpp+1, nlci - 1 ) )
1803         ilpt1 = MAX( 1, MIN(kd2 - nimpp+1, nlci     ) )
1804         ijpt0 = MAX( 1, MIN(kl  - njmpp+1, nlcj     ) )
1805         ijpt1 = MAX( 0, MIN(kl  - njmpp+1, nlcj - 1 ) )
1806      ELSEIF( ktype==2 ) THEN                              ! east/west boundaries
1807         iipt0 = MAX( 1, MIN(kl  - nimpp+1, nlci     ) )
1808         iipt1 = MAX( 0, MIN(kl  - nimpp+1, nlci - 1 ) )
1809         ijpt0 = MAX( 1, MIN(kd1 - njmpp+1, nlcj     ) )
1810         ijpt1 = MAX( 0, MIN(kd2 - njmpp+1, nlcj - 1 ) )
1811         ilpt1 = MAX( 1, MIN(kd2 - njmpp+1, nlcj     ) )
1812      ELSE
1813         WRITE(kumout, cform_err)
1814         WRITE(kumout,*) 'mppobc : bad ktype'
1815         CALL mppstop
1816      ENDIF
1817     
1818      ! Communication level by level
1819      ! ----------------------------
1820!!gm Remark : this is very time consumming!!!
1821      !                                         ! ------------------------ !
1822      DO jk = 1, kk                             !   Loop over the levels   !
1823         !                                      ! ------------------------ !
1824         !
1825         IF( ktype == 1 ) THEN                               ! north/south boundaries
1826            DO jj = ijpt0, ijpt1
1827               DO ji = iipt0, iipt1
1828                  ztab(ji,jj) = ptab(ji,jk)
1829               END DO
1830            END DO
1831         ELSEIF( ktype == 2 ) THEN                           ! east/west boundaries
1832            DO jj = ijpt0, ijpt1
1833               DO ji = iipt0, iipt1
1834                  ztab(ji,jj) = ptab(jj,jk)
1835               END DO
1836            END DO
1837         ENDIF
1838
1839
1840         ! 1. East and west directions
1841         ! ---------------------------
1842         !
1843         IF( nbondi /= 2 ) THEN         ! Read Dirichlet lateral conditions
1844            iihom = nlci-nreci
1845            DO jl = 1, jpreci
1846               t2ew(:,jl,1) = ztab(jpreci+jl,:)
1847               t2we(:,jl,1) = ztab(iihom +jl,:)
1848            END DO
1849         ENDIF
1850         !
1851         !                              ! Migrations
1852         imigr=jpreci*jpj
1853         !
1854         IF( nbondi == -1 ) THEN
1855            CALL mppsend( 2, t2we(1,1,1), imigr, noea, ml_req1 )
1856            CALL mpprecv( 1, t2ew(1,1,2), imigr, noea )
1857            IF(l_isend)   CALL mpi_wait( ml_req1, ml_stat, ml_err )
1858         ELSEIF( nbondi == 0 ) THEN
1859            CALL mppsend( 1, t2ew(1,1,1), imigr, nowe, ml_req1 )
1860            CALL mppsend( 2, t2we(1,1,1), imigr, noea, ml_req2 )
1861            CALL mpprecv( 1, t2ew(1,1,2), imigr, noea )
1862            CALL mpprecv( 2, t2we(1,1,2), imigr, nowe )
1863            IF(l_isend)   CALL mpi_wait( ml_req1, ml_stat, ml_err )
1864            IF(l_isend)   CALL mpi_wait( ml_req2, ml_stat, ml_err )
1865         ELSEIF( nbondi == 1 ) THEN
1866            CALL mppsend( 1, t2ew(1,1,1), imigr, nowe, ml_req1 )
1867            CALL mpprecv( 2, t2we(1,1,2), imigr, nowe )
1868            IF(l_isend) CALL mpi_wait( ml_req1, ml_stat, ml_err )
1869         ENDIF
1870         !
1871         !                              ! Write Dirichlet lateral conditions
1872         iihom = nlci-jpreci
1873         !
1874         IF( nbondi == 0 .OR. nbondi == 1 ) THEN
1875            DO jl = 1, jpreci
1876               ztab(jl,:) = t2we(:,jl,2)
1877            END DO
1878         ENDIF
1879         IF( nbondi == -1 .OR. nbondi == 0 ) THEN
1880            DO jl = 1, jpreci
1881               ztab(iihom+jl,:) = t2ew(:,jl,2)
1882            END DO
1883         ENDIF
1884
1885
1886         ! 2. North and south directions
1887         ! -----------------------------
1888         !
1889         IF( nbondj /= 2 ) THEN         ! Read Dirichlet lateral conditions
1890            ijhom = nlcj-nrecj
1891            DO jl = 1, jprecj
1892               t2sn(:,jl,1) = ztab(:,ijhom +jl)
1893               t2ns(:,jl,1) = ztab(:,jprecj+jl)
1894            END DO
1895         ENDIF
1896         !
1897         !                              ! Migrations
1898         imigr = jprecj * jpi
1899         !
1900         IF( nbondj == -1 ) THEN
1901            CALL mppsend( 4, t2sn(1,1,1), imigr, nono, ml_req1 )
1902            CALL mpprecv( 3, t2ns(1,1,2), imigr, nono )
1903            IF(l_isend) CALL mpi_wait( ml_req1, ml_stat, ml_err )
1904         ELSEIF( nbondj == 0 ) THEN
1905            CALL mppsend( 3, t2ns(1,1,1), imigr, noso, ml_req1 )
1906            CALL mppsend( 4, t2sn(1,1,1), imigr, nono, ml_req2 )
1907            CALL mpprecv( 3, t2ns(1,1,2), imigr, nono )
1908            CALL mpprecv( 4, t2sn(1,1,2), imigr, noso )
1909            IF( l_isend )   CALL mpi_wait( ml_req1, ml_stat, ml_err )
1910            IF( l_isend )   CALL mpi_wait( ml_req2, ml_stat, ml_err )
1911         ELSEIF( nbondj == 1 ) THEN
1912            CALL mppsend( 3, t2ns(1,1,1), imigr, noso, ml_req1 )
1913            CALL mpprecv( 4, t2sn(1,1,2), imigr, noso)
1914            IF( l_isend )   CALL mpi_wait( ml_req1, ml_stat, ml_err )
1915         ENDIF
1916         !
1917         !                              ! Write Dirichlet lateral conditions
1918         ijhom = nlcj - jprecj
1919         IF( nbondj == 0 .OR. nbondj == 1 ) THEN
1920            DO jl = 1, jprecj
1921               ztab(:,jl) = t2sn(:,jl,2)
1922            END DO
1923         ENDIF
1924         IF( nbondj == 0 .OR. nbondj == -1 ) THEN
1925            DO jl = 1, jprecj
1926               ztab(:,ijhom+jl) = t2ns(:,jl,2)
1927            END DO
1928         ENDIF
1929         IF( ktype==1 .AND. kd1 <= jpi+nimpp-1 .AND. nimpp <= kd2 ) THEN
1930            DO jj = ijpt0, ijpt1            ! north/south boundaries
1931               DO ji = iipt0,ilpt1
1932                  ptab(ji,jk) = ztab(ji,jj) 
1933               END DO
1934            END DO
1935         ELSEIF( ktype==2 .AND. kd1 <= jpj+njmpp-1 .AND. njmpp <= kd2 ) THEN
1936            DO jj = ijpt0, ilpt1            ! east/west boundaries
1937               DO ji = iipt0,iipt1
1938                  ptab(jj,jk) = ztab(ji,jj) 
1939               END DO
1940            END DO
1941         ENDIF
1942         !
1943      END DO
1944      !
1945      CALL wrk_dealloc( jpi,jpj, ztab )
1946      !
1947   END SUBROUTINE mppobc
1948   
1949
1950   SUBROUTINE mpp_comm_free( kcom )
1951      !!----------------------------------------------------------------------
1952      !!----------------------------------------------------------------------
1953      INTEGER, INTENT(in) ::   kcom
1954      !!
1955      INTEGER :: ierr
1956      !!----------------------------------------------------------------------
1957      !
1958      CALL MPI_COMM_FREE(kcom, ierr)
1959      !
1960   END SUBROUTINE mpp_comm_free
1961
1962
1963   SUBROUTINE mpp_ini_ice( pindic, kumout )
1964      !!----------------------------------------------------------------------
1965      !!               ***  routine mpp_ini_ice  ***
1966      !!
1967      !! ** Purpose :   Initialize special communicator for ice areas
1968      !!      condition together with global variables needed in the ddmpp folding
1969      !!
1970      !! ** Method  : - Look for ice processors in ice routines
1971      !!              - Put their number in nrank_ice
1972      !!              - Create groups for the world processors and the ice processors
1973      !!              - Create a communicator for ice processors
1974      !!
1975      !! ** output
1976      !!      njmppmax = njmpp for northern procs
1977      !!      ndim_rank_ice = number of processors with ice
1978      !!      nrank_ice (ndim_rank_ice) = ice processors
1979      !!      ngrp_world = group ID for the world processors
1980      !!      ngrp_ice = group ID for the ice processors
1981      !!      ncomm_ice = communicator for the ice procs.
1982      !!      n_ice_root = number (in the world) of proc 0 in the ice comm.
1983      !!
1984      !!----------------------------------------------------------------------
1985      INTEGER, INTENT(in) ::   pindic
1986      INTEGER, INTENT(in) ::   kumout   ! ocean.output logical unit
1987      !!
1988      INTEGER :: jjproc
1989      INTEGER :: ii, ierr
1990      INTEGER, ALLOCATABLE, DIMENSION(:) ::   kice
1991      INTEGER, ALLOCATABLE, DIMENSION(:) ::   zwork
1992      !!----------------------------------------------------------------------
1993      !
1994      ! Since this is just an init routine and these arrays are of length jpnij
1995      ! then don't use wrk_nemo module - just allocate and deallocate.
1996      ALLOCATE( kice(jpnij), zwork(jpnij), STAT=ierr )
1997      IF( ierr /= 0 ) THEN
1998         WRITE(kumout, cform_err)
1999         WRITE(kumout,*) 'mpp_ini_ice : failed to allocate 2, 1D arrays (jpnij in length)'
2000         CALL mppstop
2001      ENDIF
2002
2003      ! Look for how many procs with sea-ice
2004      !
2005      kice = 0
2006      DO jjproc = 1, jpnij
2007         IF( jjproc == narea .AND. pindic .GT. 0 )   kice(jjproc) = 1   
2008      END DO
2009      !
2010      zwork = 0
2011      CALL MPI_ALLREDUCE( kice, zwork, jpnij, mpi_integer, mpi_sum, mpi_comm_opa, ierr )
2012      ndim_rank_ice = SUM( zwork )         
2013
2014      ! Allocate the right size to nrank_north
2015      IF( ALLOCATED ( nrank_ice ) )   DEALLOCATE( nrank_ice )
2016      ALLOCATE( nrank_ice(ndim_rank_ice) )
2017      !
2018      ii = 0     
2019      nrank_ice = 0
2020      DO jjproc = 1, jpnij
2021         IF( zwork(jjproc) == 1) THEN
2022            ii = ii + 1
2023            nrank_ice(ii) = jjproc -1 
2024         ENDIF
2025      END DO
2026
2027      ! Create the world group
2028      CALL MPI_COMM_GROUP( mpi_comm_opa, ngrp_world, ierr )
2029
2030      ! Create the ice group from the world group
2031      CALL MPI_GROUP_INCL( ngrp_world, ndim_rank_ice, nrank_ice, ngrp_ice, ierr )
2032
2033      ! Create the ice communicator , ie the pool of procs with sea-ice
2034      CALL MPI_COMM_CREATE( mpi_comm_opa, ngrp_ice, ncomm_ice, ierr )
2035
2036      ! Find proc number in the world of proc 0 in the north
2037      ! The following line seems to be useless, we just comment & keep it as reminder
2038      ! CALL MPI_GROUP_TRANSLATE_RANKS(ngrp_ice,1,0,ngrp_world,n_ice_root,ierr)
2039      !
2040      DEALLOCATE(kice, zwork)
2041      !
2042   END SUBROUTINE mpp_ini_ice
2043
2044
2045   SUBROUTINE mpp_ini_znl( kumout )
2046      !!----------------------------------------------------------------------
2047      !!               ***  routine mpp_ini_znl  ***
2048      !!
2049      !! ** Purpose :   Initialize special communicator for computing zonal sum
2050      !!
2051      !! ** Method  : - Look for processors in the same row
2052      !!              - Put their number in nrank_znl
2053      !!              - Create group for the znl processors
2054      !!              - Create a communicator for znl processors
2055      !!              - Determine if processor should write znl files
2056      !!
2057      !! ** output
2058      !!      ndim_rank_znl = number of processors on the same row
2059      !!      ngrp_znl = group ID for the znl processors
2060      !!      ncomm_znl = communicator for the ice procs.
2061      !!      n_znl_root = number (in the world) of proc 0 in the ice comm.
2062      !!
2063      !!----------------------------------------------------------------------
2064      INTEGER, INTENT(in) ::   kumout   ! ocean.output logical units
2065      !
2066      INTEGER :: jproc      ! dummy loop integer
2067      INTEGER :: ierr, ii   ! local integer
2068      INTEGER, ALLOCATABLE, DIMENSION(:) ::   kwork
2069      !!----------------------------------------------------------------------
2070      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_world     : ', ngrp_world
2071      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_world : ', mpi_comm_world
2072      !-$$     WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - mpi_comm_opa   : ', mpi_comm_opa
2073      !
2074      ALLOCATE( kwork(jpnij), STAT=ierr )
2075      IF( ierr /= 0 ) THEN
2076         WRITE(kumout, cform_err)
2077         WRITE(kumout,*) 'mpp_ini_znl : failed to allocate 1D array of length jpnij'
2078         CALL mppstop
2079      ENDIF
2080
2081      IF( jpnj == 1 ) THEN
2082         ngrp_znl  = ngrp_world
2083         ncomm_znl = mpi_comm_opa
2084      ELSE
2085         !
2086         CALL MPI_ALLGATHER ( njmpp, 1, mpi_integer, kwork, 1, mpi_integer, mpi_comm_opa, ierr )
2087         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - kwork pour njmpp : ', kwork
2088         !-$$        CALL flush(numout)
2089         !
2090         ! Count number of processors on the same row
2091         ndim_rank_znl = 0
2092         DO jproc=1,jpnij
2093            IF ( kwork(jproc) == njmpp ) THEN
2094               ndim_rank_znl = ndim_rank_znl + 1
2095            ENDIF
2096         END DO
2097         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ndim_rank_znl : ', ndim_rank_znl
2098         !-$$        CALL flush(numout)
2099         ! Allocate the right size to nrank_znl
2100         IF (ALLOCATED (nrank_znl)) DEALLOCATE(nrank_znl)
2101         ALLOCATE(nrank_znl(ndim_rank_znl))
2102         ii = 0     
2103         nrank_znl (:) = 0
2104         DO jproc=1,jpnij
2105            IF ( kwork(jproc) == njmpp) THEN
2106               ii = ii + 1
2107               nrank_znl(ii) = jproc -1 
2108            ENDIF
2109         END DO
2110         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - nrank_znl : ', nrank_znl
2111         !-$$        CALL flush(numout)
2112
2113         ! Create the opa group
2114         CALL MPI_COMM_GROUP(mpi_comm_opa,ngrp_opa,ierr)
2115         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_opa : ', ngrp_opa
2116         !-$$        CALL flush(numout)
2117
2118         ! Create the znl group from the opa group
2119         CALL MPI_GROUP_INCL  ( ngrp_opa, ndim_rank_znl, nrank_znl, ngrp_znl, ierr )
2120         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ngrp_znl ', ngrp_znl
2121         !-$$        CALL flush(numout)
2122
2123         ! Create the znl communicator from the opa communicator, ie the pool of procs in the same row
2124         CALL MPI_COMM_CREATE ( mpi_comm_opa, ngrp_znl, ncomm_znl, ierr )
2125         !-$$        WRITE (numout,*) 'mpp_ini_znl ', nproc, ' - ncomm_znl ', ncomm_znl
2126         !-$$        CALL flush(numout)
2127         !
2128      END IF
2129
2130      ! Determines if processor if the first (starting from i=1) on the row
2131      IF ( jpni == 1 ) THEN
2132         l_znl_root = .TRUE.
2133      ELSE
2134         l_znl_root = .FALSE.
2135         kwork (1) = nimpp
2136         CALL mpp_min ( kwork(1), kcom = ncomm_znl)
2137         IF ( nimpp == kwork(1)) l_znl_root = .TRUE.
2138      END IF
2139
2140      DEALLOCATE(kwork)
2141
2142   END SUBROUTINE mpp_ini_znl
2143
2144
2145   SUBROUTINE mpp_ini_north
2146      !!----------------------------------------------------------------------
2147      !!               ***  routine mpp_ini_north  ***
2148      !!
2149      !! ** Purpose :   Initialize special communicator for north folding
2150      !!      condition together with global variables needed in the mpp folding
2151      !!
2152      !! ** Method  : - Look for northern processors
2153      !!              - Put their number in nrank_north
2154      !!              - Create groups for the world processors and the north processors
2155      !!              - Create a communicator for northern processors
2156      !!
2157      !! ** output
2158      !!      njmppmax = njmpp for northern procs
2159      !!      ndim_rank_north = number of processors in the northern line
2160      !!      nrank_north (ndim_rank_north) = number  of the northern procs.
2161      !!      ngrp_world = group ID for the world processors
2162      !!      ngrp_north = group ID for the northern processors
2163      !!      ncomm_north = communicator for the northern procs.
2164      !!      north_root = number (in the world) of proc 0 in the northern comm.
2165      !!
2166      !!----------------------------------------------------------------------
2167      INTEGER ::   ierr
2168      INTEGER ::   jjproc
2169      INTEGER ::   ii, ji
2170      !!----------------------------------------------------------------------
2171      !
2172      njmppmax = MAXVAL( njmppt )
2173      !
2174      ! Look for how many procs on the northern boundary
2175      ndim_rank_north = 0
2176      DO jjproc = 1, jpnij
2177         IF( njmppt(jjproc) == njmppmax )   ndim_rank_north = ndim_rank_north + 1
2178      END DO
2179      !
2180      ! Allocate the right size to nrank_north
2181      IF (ALLOCATED (nrank_north)) DEALLOCATE(nrank_north)
2182      ALLOCATE( nrank_north(ndim_rank_north) )
2183
2184      ! Fill the nrank_north array with proc. number of northern procs.
2185      ! Note : the rank start at 0 in MPI
2186      ii = 0
2187      DO ji = 1, jpnij
2188         IF ( njmppt(ji) == njmppmax   ) THEN
2189            ii=ii+1
2190            nrank_north(ii)=ji-1
2191         END IF
2192      END DO
2193      !
2194      ! create the world group
2195      CALL MPI_COMM_GROUP( mpi_comm_opa, ngrp_world, ierr )
2196      !
2197      ! Create the North group from the world group
2198      CALL MPI_GROUP_INCL( ngrp_world, ndim_rank_north, nrank_north, ngrp_north, ierr )
2199      !
2200      ! Create the North communicator , ie the pool of procs in the north group
2201      CALL MPI_COMM_CREATE( mpi_comm_opa, ngrp_north, ncomm_north, ierr )
2202      !
2203   END SUBROUTINE mpp_ini_north
2204
2205
2206   SUBROUTINE mpp_lbc_north_3d( pt3d, cd_type, psgn )
2207      !!---------------------------------------------------------------------
2208      !!                   ***  routine mpp_lbc_north_3d  ***
2209      !!
2210      !! ** Purpose :   Ensure proper north fold horizontal bondary condition
2211      !!              in mpp configuration in case of jpn1 > 1
2212      !!
2213      !! ** Method  :   North fold condition and mpp with more than one proc
2214      !!              in i-direction require a specific treatment. We gather
2215      !!              the 4 northern lines of the global domain on 1 processor
2216      !!              and apply lbc north-fold on this sub array. Then we
2217      !!              scatter the north fold array back to the processors.
2218      !!
2219      !!----------------------------------------------------------------------
2220      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pt3d      ! 3D array on which the b.c. is applied
2221      CHARACTER(len=1)                , INTENT(in   ) ::   cd_type   ! nature of pt3d grid-points
2222      !                                                              !   = T ,  U , V , F or W  gridpoints
2223      REAL(wp)                        , INTENT(in   ) ::   psgn      ! = -1. the sign change across the north fold
2224      !!                                                             ! =  1. , the sign is kept
2225      INTEGER ::   ji, jj, jr
2226      INTEGER ::   ierr, itaille, ildi, ilei, iilb
2227      INTEGER ::   ijpj, ijpjm1, ij, iproc
2228      INTEGER, DIMENSION (jpmaxngh)          ::   ml_req_nf          ! for mpi_isend when avoiding mpi_allgather
2229      INTEGER                                ::   ml_err             ! for mpi_isend when avoiding mpi_allgather
2230      INTEGER, DIMENSION(MPI_STATUS_SIZE)    ::   ml_stat            ! for mpi_isend when avoiding mpi_allgather
2231      !!----------------------------------------------------------------------
2232      !   
2233      ijpj   = 4
2234      ityp = -1
2235      ijpjm1 = 3
2236      ztab(:,:,:) = 0.e0
2237      !
2238      DO jj = nlcj - ijpj +1, nlcj          ! put in znorthloc the last 4 jlines of pt3d
2239         ij = jj - nlcj + ijpj
2240         znorthloc(:,ij,:) = pt3d(:,jj,:)
2241      END DO
2242      !
2243      !                                     ! Build in procs of ncomm_north the znorthgloio
2244      itaille = jpi * jpk * ijpj
2245      IF ( l_north_nogather ) THEN
2246         !
2247         ! Avoid the use of mpi_allgather by exchanging only with the processes already identified
2248         ! (in nemo_northcomms) as being  involved in this process' northern boundary exchange
2249         !
2250         DO jj = nlcj-ijpj+1, nlcj          ! First put local values into the global array
2251            ij = jj - nlcj + ijpj
2252            DO ji = 1, nlci
2253               ztab(ji+nimpp-1,ij,:) = pt3d(ji,jj,:)
2254            END DO
2255         END DO
2256
2257         !
2258         ! Set the exchange type in order to access the correct list of active neighbours
2259         !
2260         SELECT CASE ( cd_type )
2261            CASE ( 'T' , 'W' )
2262               ityp = 1
2263            CASE ( 'U' )
2264               ityp = 2
2265            CASE ( 'V' )
2266               ityp = 3
2267            CASE ( 'F' )
2268               ityp = 4
2269            CASE ( 'I' )
2270               ityp = 5
2271            CASE DEFAULT
2272               ityp = -1                    ! Set a default value for unsupported types which
2273                                            ! will cause a fallback to the mpi_allgather method
2274         END SELECT
2275         IF ( ityp .gt. 0 ) THEN
2276
2277            DO jr = 1,nsndto(ityp)
2278               CALL mppsend(5, znorthloc, itaille, isendto(jr,ityp), ml_req_nf(jr) )
2279            END DO
2280            DO jr = 1,nsndto(ityp)
2281               CALL mpprecv(5, zfoldwk, itaille, isendto(jr,ityp))
2282               iproc = isendto(jr,ityp) + 1
2283               ildi = nldit (iproc)
2284               ilei = nleit (iproc)
2285               iilb = nimppt(iproc)
2286               DO jj = 1, ijpj
2287                  DO ji = ildi, ilei
2288                     ztab(ji+iilb-1,jj,:) = zfoldwk(ji,jj,:)
2289                  END DO
2290               END DO
2291            END DO
2292            IF (l_isend) THEN
2293               DO jr = 1,nsndto(ityp)
2294                  CALL mpi_wait(ml_req_nf(jr), ml_stat, ml_err)
2295               END DO
2296            ENDIF
2297
2298         ENDIF
2299
2300      ENDIF
2301
2302      IF ( ityp .lt. 0 ) THEN
2303         CALL MPI_ALLGATHER( znorthloc  , itaille, MPI_DOUBLE_PRECISION,                &
2304            &                znorthgloio, itaille, MPI_DOUBLE_PRECISION, ncomm_north, ierr )
2305         !
2306         DO jr = 1, ndim_rank_north         ! recover the global north array
2307            iproc = nrank_north(jr) + 1
2308            ildi  = nldit (iproc)
2309            ilei  = nleit (iproc)
2310            iilb  = nimppt(iproc)
2311            DO jj = 1, ijpj
2312               DO ji = ildi, ilei
2313                  ztab(ji+iilb-1,jj,:) = znorthgloio(ji,jj,:,jr)
2314               END DO
2315            END DO
2316         END DO
2317      ENDIF
2318      !
2319      ! The ztab array has been either:
2320      !  a. Fully populated by the mpi_allgather operation or
2321      !  b. Had the active points for this domain and northern neighbours populated
2322      !     by peer to peer exchanges
2323      ! Either way the array may be folded by lbc_nfd and the result for the span of
2324      ! this domain will be identical.
2325      !
2326      CALL lbc_nfd( ztab, cd_type, psgn )   ! North fold boundary condition
2327      !
2328      DO jj = nlcj-ijpj+1, nlcj             ! Scatter back to pt3d
2329         ij = jj - nlcj + ijpj
2330         DO ji= 1, nlci
2331            pt3d(ji,jj,:) = ztab(ji+nimpp-1,ij,:)
2332         END DO
2333      END DO
2334      !
2335   END SUBROUTINE mpp_lbc_north_3d
2336
2337
2338   SUBROUTINE mpp_lbc_north_2d( pt2d, cd_type, psgn)
2339      !!---------------------------------------------------------------------
2340      !!                   ***  routine mpp_lbc_north_2d  ***
2341      !!
2342      !! ** Purpose :   Ensure proper north fold horizontal bondary condition
2343      !!              in mpp configuration in case of jpn1 > 1 (for 2d array )
2344      !!
2345      !! ** Method  :   North fold condition and mpp with more than one proc
2346      !!              in i-direction require a specific treatment. We gather
2347      !!              the 4 northern lines of the global domain on 1 processor
2348      !!              and apply lbc north-fold on this sub array. Then we
2349      !!              scatter the north fold array back to the processors.
2350      !!
2351      !!----------------------------------------------------------------------
2352      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   pt2d      ! 3D array on which the b.c. is applied
2353      CHARACTER(len=1)            , INTENT(in   ) ::   cd_type   ! nature of pt3d grid-points
2354      !                                                          !   = T ,  U , V , F or W  gridpoints
2355      REAL(wp)                    , INTENT(in   ) ::   psgn      ! = -1. the sign change across the north fold
2356      !!                                                             ! =  1. , the sign is kept
2357      INTEGER ::   ji, jj, jr
2358      INTEGER ::   ierr, itaille, ildi, ilei, iilb
2359      INTEGER ::   ijpj, ijpjm1, ij, iproc
2360      INTEGER, DIMENSION (jpmaxngh)      ::   ml_req_nf          ! for mpi_isend when avoiding mpi_allgather
2361      INTEGER                            ::   ml_err             ! for mpi_isend when avoiding mpi_allgather
2362      INTEGER, DIMENSION(MPI_STATUS_SIZE)::   ml_stat            ! for mpi_isend when avoiding mpi_allgather
2363      !!----------------------------------------------------------------------
2364      !
2365      ijpj   = 4
2366      ityp = -1
2367      ijpjm1 = 3
2368      ztab_2d(:,:) = 0.e0
2369      !
2370      DO jj = nlcj-ijpj+1, nlcj             ! put in znorthloc_2d the last 4 jlines of pt2d
2371         ij = jj - nlcj + ijpj
2372         znorthloc_2d(:,ij) = pt2d(:,jj)
2373      END DO
2374
2375      !                                     ! Build in procs of ncomm_north the znorthgloio_2d
2376      itaille = jpi * ijpj
2377      IF ( l_north_nogather ) THEN
2378         !
2379         ! Avoid the use of mpi_allgather by exchanging only with the processes already identified
2380         ! (in nemo_northcomms) as being  involved in this process' northern boundary exchange
2381         !
2382         DO jj = nlcj-ijpj+1, nlcj          ! First put local values into the global array
2383            ij = jj - nlcj + ijpj
2384            DO ji = 1, nlci
2385               ztab_2d(ji+nimpp-1,ij) = pt2d(ji,jj)
2386            END DO
2387         END DO
2388
2389         !
2390         ! Set the exchange type in order to access the correct list of active neighbours
2391         !
2392         SELECT CASE ( cd_type )
2393            CASE ( 'T' , 'W' )
2394               ityp = 1
2395            CASE ( 'U' )
2396               ityp = 2
2397            CASE ( 'V' )
2398               ityp = 3
2399            CASE ( 'F' )
2400               ityp = 4
2401            CASE ( 'I' )
2402               ityp = 5
2403            CASE DEFAULT
2404               ityp = -1                    ! Set a default value for unsupported types which
2405                                            ! will cause a fallback to the mpi_allgather method
2406         END SELECT
2407
2408         IF ( ityp .gt. 0 ) THEN
2409
2410            DO jr = 1,nsndto(ityp)
2411               CALL mppsend(5, znorthloc_2d, itaille, isendto(jr,ityp), ml_req_nf(jr) )
2412            END DO
2413            DO jr = 1,nsndto(ityp)
2414               CALL mpprecv(5, zfoldwk_2d, itaille, isendto(jr,ityp))
2415               iproc = isendto(jr,ityp) + 1
2416               ildi = nldit (iproc)
2417               ilei = nleit (iproc)
2418               iilb = nimppt(iproc)
2419               DO jj = 1, ijpj
2420                  DO ji = ildi, ilei
2421                     ztab_2d(ji+iilb-1,jj) = zfoldwk_2d(ji,jj)
2422                  END DO
2423               END DO
2424            END DO
2425            IF (l_isend) THEN
2426               DO jr = 1,nsndto(ityp)
2427                  CALL mpi_wait(ml_req_nf(jr), ml_stat, ml_err)
2428               END DO
2429            ENDIF
2430
2431         ENDIF
2432
2433      ENDIF
2434
2435      IF ( ityp .lt. 0 ) THEN
2436         CALL MPI_ALLGATHER( znorthloc_2d  , itaille, MPI_DOUBLE_PRECISION,        &
2437            &                znorthgloio_2d, itaille, MPI_DOUBLE_PRECISION, ncomm_north, ierr )
2438         !
2439         DO jr = 1, ndim_rank_north            ! recover the global north array
2440            iproc = nrank_north(jr) + 1
2441            ildi = nldit (iproc)
2442            ilei = nleit (iproc)
2443            iilb = nimppt(iproc)
2444            DO jj = 1, ijpj
2445               DO ji = ildi, ilei
2446                  ztab_2d(ji+iilb-1,jj) = znorthgloio_2d(ji,jj,jr)
2447               END DO
2448            END DO
2449         END DO
2450      ENDIF
2451      !
2452      ! The ztab array has been either:
2453      !  a. Fully populated by the mpi_allgather operation or
2454      !  b. Had the active points for this domain and northern neighbours populated
2455      !     by peer to peer exchanges
2456      ! Either way the array may be folded by lbc_nfd and the result for the span of
2457      ! this domain will be identical.
2458      !
2459      CALL lbc_nfd( ztab_2d, cd_type, psgn )   ! North fold boundary condition
2460      !
2461      !
2462      DO jj = nlcj-ijpj+1, nlcj             ! Scatter back to pt2d
2463         ij = jj - nlcj + ijpj
2464         DO ji = 1, nlci
2465            pt2d(ji,jj) = ztab_2d(ji+nimpp-1,ij)
2466         END DO
2467      END DO
2468      !
2469   END SUBROUTINE mpp_lbc_north_2d
2470
2471
2472   SUBROUTINE mpp_lbc_north_e( pt2d, cd_type, psgn)
2473      !!---------------------------------------------------------------------
2474      !!                   ***  routine mpp_lbc_north_2d  ***
2475      !!
2476      !! ** Purpose :   Ensure proper north fold horizontal bondary condition
2477      !!              in mpp configuration in case of jpn1 > 1 and for 2d
2478      !!              array with outer extra halo
2479      !!
2480      !! ** Method  :   North fold condition and mpp with more than one proc
2481      !!              in i-direction require a specific treatment. We gather
2482      !!              the 4+2*jpr2dj northern lines of the global domain on 1
2483      !!              processor and apply lbc north-fold on this sub array.
2484      !!              Then we scatter the north fold array back to the processors.
2485      !!
2486      !!----------------------------------------------------------------------
2487      REAL(wp), DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj), INTENT(inout) ::   pt2d     ! 2D array with extra halo
2488      CHARACTER(len=1)                                            , INTENT(in   ) ::   cd_type  ! nature of pt3d grid-points
2489      !                                                                                         !   = T ,  U , V , F or W -points
2490      REAL(wp)                                                    , INTENT(in   ) ::   psgn     ! = -1. the sign change across the 
2491      !!                                                                                        ! north fold, =  1. otherwise
2492      INTEGER ::   ji, jj, jr
2493      INTEGER ::   ierr, itaille, ildi, ilei, iilb
2494      INTEGER ::   ijpj, ij, iproc
2495      !!----------------------------------------------------------------------
2496      !
2497      ijpj=4
2498      ztab_e(:,:) = 0.e0
2499
2500      ij=0
2501      ! put in znorthloc_e the last 4 jlines of pt2d
2502      DO jj = nlcj - ijpj + 1 - jpr2dj, nlcj +jpr2dj
2503         ij = ij + 1
2504         DO ji = 1, jpi
2505            znorthloc_e(ji,ij)=pt2d(ji,jj)
2506         END DO
2507      END DO
2508      !
2509      itaille = jpi * ( ijpj + 2 * jpr2dj )
2510      CALL MPI_ALLGATHER( znorthloc_e(1,1)  , itaille, MPI_DOUBLE_PRECISION,    &
2511         &                znorthgloio_e(1,1,1), itaille, MPI_DOUBLE_PRECISION, ncomm_north, ierr )
2512      !
2513      DO jr = 1, ndim_rank_north            ! recover the global north array
2514         iproc = nrank_north(jr) + 1
2515         ildi = nldit (iproc)
2516         ilei = nleit (iproc)
2517         iilb = nimppt(iproc)
2518         DO jj = 1, ijpj+2*jpr2dj
2519            DO ji = ildi, ilei
2520               ztab_e(ji+iilb-1,jj) = znorthgloio_e(ji,jj,jr)
2521            END DO
2522         END DO
2523      END DO
2524
2525
2526      ! 2. North-Fold boundary conditions
2527      ! ----------------------------------
2528      CALL lbc_nfd( ztab_e(:,:), cd_type, psgn, pr2dj = jpr2dj )
2529
2530      ij = jpr2dj
2531      !! Scatter back to pt2d
2532      DO jj = nlcj - ijpj + 1 , nlcj +jpr2dj
2533      ij  = ij +1 
2534         DO ji= 1, nlci
2535            pt2d(ji,jj) = ztab_e(ji+nimpp-1,ij)
2536         END DO
2537      END DO
2538      !
2539   END SUBROUTINE mpp_lbc_north_e
2540
2541
2542   SUBROUTINE mpi_init_opa( ldtxt, ksft, code )
2543      !!---------------------------------------------------------------------
2544      !!                   ***  routine mpp_init.opa  ***
2545      !!
2546      !! ** Purpose :: export and attach a MPI buffer for bsend
2547      !!
2548      !! ** Method  :: define buffer size in namelist, if 0 no buffer attachment
2549      !!            but classical mpi_init
2550      !!
2551      !! History :: 01/11 :: IDRIS initial version for IBM only 
2552      !!            08/04 :: R. Benshila, generalisation
2553      !!---------------------------------------------------------------------
2554      CHARACTER(len=*),DIMENSION(:), INTENT(  out) ::   ldtxt 
2555      INTEGER                      , INTENT(inout) ::   ksft
2556      INTEGER                      , INTENT(  out) ::   code
2557      INTEGER                                      ::   ierr, ji
2558      LOGICAL                                      ::   mpi_was_called
2559      !!---------------------------------------------------------------------
2560      !
2561      CALL mpi_initialized( mpi_was_called, code )      ! MPI initialization
2562      IF ( code /= MPI_SUCCESS ) THEN
2563         DO ji = 1, SIZE(ldtxt) 
2564            IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
2565         END DO         
2566         WRITE(*, cform_err)
2567         WRITE(*, *) ' lib_mpp: Error in routine mpi_initialized'
2568         CALL mpi_abort( mpi_comm_world, code, ierr )
2569      ENDIF
2570      !
2571      IF( .NOT. mpi_was_called ) THEN
2572         CALL mpi_init( code )
2573         CALL mpi_comm_dup( mpi_comm_world, mpi_comm_opa, code )
2574         IF ( code /= MPI_SUCCESS ) THEN
2575            DO ji = 1, SIZE(ldtxt) 
2576               IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
2577            END DO
2578            WRITE(*, cform_err)
2579            WRITE(*, *) ' lib_mpp: Error in routine mpi_comm_dup'
2580            CALL mpi_abort( mpi_comm_world, code, ierr )
2581         ENDIF
2582      ENDIF
2583      !
2584      IF( nn_buffer > 0 ) THEN
2585         WRITE(ldtxt(ksft),*) 'mpi_bsend, buffer allocation of  : ', nn_buffer   ;   ksft = ksft + 1
2586         ! Buffer allocation and attachment
2587         ALLOCATE( tampon(nn_buffer), stat = ierr )
2588         IF( ierr /= 0 ) THEN
2589            DO ji = 1, SIZE(ldtxt) 
2590               IF( TRIM(ldtxt(ji)) /= '' )   WRITE(*,*) ldtxt(ji)      ! control print of mynode
2591            END DO
2592            WRITE(*, cform_err)
2593            WRITE(*, *) ' lib_mpp: Error in ALLOCATE', ierr
2594            CALL mpi_abort( mpi_comm_world, code, ierr )
2595         END IF
2596         CALL mpi_buffer_attach( tampon, nn_buffer, code )
2597      ENDIF
2598      !
2599   END SUBROUTINE mpi_init_opa
2600
2601#if defined key_mpp_rep
2602   SUBROUTINE DDPDD_MPI (ydda, yddb, ilen, itype)
2603      !!---------------------------------------------------------------------
2604      !!   Routine DDPDD_MPI: used by reduction operator MPI_SUMDD
2605      !!
2606      !!   Modification of original codes written by David H. Bailey
2607      !!   This subroutine computes yddb(i) = ydda(i)+yddb(i)
2608      !!---------------------------------------------------------------------
2609      INTEGER, INTENT(in)                         :: ilen, itype
2610      COMPLEX(wp), DIMENSION(ilen), INTENT(in)     :: ydda
2611      COMPLEX(wp), DIMENSION(ilen), INTENT(inout)  :: yddb
2612      !
2613      REAL(wp) :: zerr, zt1, zt2    ! local work variables
2614      INTEGER :: ji, ztmp           ! local scalar
2615
2616      ztmp = itype   ! avoid compilation warning
2617
2618      DO ji=1,ilen
2619      ! Compute ydda + yddb using Knuth's trick.
2620         zt1  = real(ydda(ji)) + real(yddb(ji))
2621         zerr = zt1 - real(ydda(ji))
2622         zt2  = ((real(yddb(ji)) - zerr) + (real(ydda(ji)) - (zt1 - zerr))) &
2623                + aimag(ydda(ji)) + aimag(yddb(ji))
2624
2625         ! The result is zt1 + zt2, after normalization.
2626         yddb(ji) = cmplx ( zt1 + zt2, zt2 - ((zt1 + zt2) - zt1),wp )
2627      END DO
2628
2629   END SUBROUTINE DDPDD_MPI
2630#endif
2631
2632#else
2633   !!----------------------------------------------------------------------
2634   !!   Default case:            Dummy module        share memory computing
2635   !!----------------------------------------------------------------------
2636   USE in_out_manager
2637
2638   INTERFACE mpp_sum
2639      MODULE PROCEDURE mpp_sum_a2s, mpp_sum_as, mpp_sum_ai, mpp_sum_s, mpp_sum_i
2640   END INTERFACE
2641   INTERFACE mpp_max
2642      MODULE PROCEDURE mppmax_a_int, mppmax_int, mppmax_a_real, mppmax_real
2643   END INTERFACE
2644   INTERFACE mpp_min
2645      MODULE PROCEDURE mppmin_a_int, mppmin_int, mppmin_a_real, mppmin_real
2646   END INTERFACE
2647   INTERFACE mppobc
2648      MODULE PROCEDURE mppobc_1d, mppobc_2d, mppobc_3d, mppobc_4d
2649   END INTERFACE
2650   INTERFACE mpp_minloc
2651      MODULE PROCEDURE mpp_minloc2d ,mpp_minloc3d
2652   END INTERFACE
2653   INTERFACE mpp_maxloc
2654      MODULE PROCEDURE mpp_maxloc2d ,mpp_maxloc3d
2655   END INTERFACE
2656
2657   LOGICAL, PUBLIC, PARAMETER ::   lk_mpp = .FALSE.      !: mpp flag
2658   LOGICAL, PUBLIC            ::   ln_nnogather  = .FALSE.  !: namelist control of northfold comms (needed here in case "key_mpp_mpi" is not used)
2659   INTEGER :: ncomm_ice
2660   !!----------------------------------------------------------------------
2661CONTAINS
2662
2663   INTEGER FUNCTION lib_mpp_alloc(kumout)          ! Dummy function
2664      INTEGER, INTENT(in) ::   kumout
2665      lib_mpp_alloc = 0
2666   END FUNCTION lib_mpp_alloc
2667
2668   FUNCTION mynode( ldtxt, kumnam, kstop, localComm ) RESULT (function_value)
2669      INTEGER, OPTIONAL            , INTENT(in   ) ::   localComm
2670      CHARACTER(len=*),DIMENSION(:) ::   ldtxt 
2671      INTEGER ::   kumnam, kstop
2672      IF( PRESENT( localComm ) .OR. .NOT.PRESENT( localComm ) )   function_value = 0
2673      IF( .FALSE. )   ldtxt(:) = 'never done'
2674   END FUNCTION mynode
2675
2676   SUBROUTINE mppsync                       ! Dummy routine
2677   END SUBROUTINE mppsync
2678
2679   SUBROUTINE mpp_sum_as( parr, kdim, kcom )      ! Dummy routine
2680      REAL   , DIMENSION(:) :: parr
2681      INTEGER               :: kdim
2682      INTEGER, OPTIONAL     :: kcom 
2683      WRITE(*,*) 'mpp_sum_as: You should not have seen this print! error?', kdim, parr(1), kcom
2684   END SUBROUTINE mpp_sum_as
2685
2686   SUBROUTINE mpp_sum_a2s( parr, kdim, kcom )      ! Dummy routine
2687      REAL   , DIMENSION(:,:) :: parr
2688      INTEGER               :: kdim
2689      INTEGER, OPTIONAL     :: kcom 
2690      WRITE(*,*) 'mpp_sum_a2s: You should not have seen this print! error?', kdim, parr(1,1), kcom
2691   END SUBROUTINE mpp_sum_a2s
2692
2693   SUBROUTINE mpp_sum_ai( karr, kdim, kcom )      ! Dummy routine
2694      INTEGER, DIMENSION(:) :: karr
2695      INTEGER               :: kdim
2696      INTEGER, OPTIONAL     :: kcom 
2697      WRITE(*,*) 'mpp_sum_ai: You should not have seen this print! error?', kdim, karr(1), kcom
2698   END SUBROUTINE mpp_sum_ai
2699
2700   SUBROUTINE mpp_sum_s( psca, kcom )            ! Dummy routine
2701      REAL                  :: psca
2702      INTEGER, OPTIONAL     :: kcom 
2703      WRITE(*,*) 'mpp_sum_s: You should not have seen this print! error?', psca, kcom
2704   END SUBROUTINE mpp_sum_s
2705
2706   SUBROUTINE mpp_sum_i( kint, kcom )            ! Dummy routine
2707      integer               :: kint
2708      INTEGER, OPTIONAL     :: kcom 
2709      WRITE(*,*) 'mpp_sum_i: You should not have seen this print! error?', kint, kcom
2710   END SUBROUTINE mpp_sum_i
2711
2712   SUBROUTINE mppmax_a_real( parr, kdim, kcom )
2713      REAL   , DIMENSION(:) :: parr
2714      INTEGER               :: kdim
2715      INTEGER, OPTIONAL     :: kcom 
2716      WRITE(*,*) 'mppmax_a_real: You should not have seen this print! error?', kdim, parr(1), kcom
2717   END SUBROUTINE mppmax_a_real
2718
2719   SUBROUTINE mppmax_real( psca, kcom )
2720      REAL                  :: psca
2721      INTEGER, OPTIONAL     :: kcom 
2722      WRITE(*,*) 'mppmax_real: You should not have seen this print! error?', psca, kcom
2723   END SUBROUTINE mppmax_real
2724
2725   SUBROUTINE mppmin_a_real( parr, kdim, kcom )
2726      REAL   , DIMENSION(:) :: parr
2727      INTEGER               :: kdim
2728      INTEGER, OPTIONAL     :: kcom 
2729      WRITE(*,*) 'mppmin_a_real: You should not have seen this print! error?', kdim, parr(1), kcom
2730   END SUBROUTINE mppmin_a_real
2731
2732   SUBROUTINE mppmin_real( psca, kcom )
2733      REAL                  :: psca
2734      INTEGER, OPTIONAL     :: kcom 
2735      WRITE(*,*) 'mppmin_real: You should not have seen this print! error?', psca, kcom
2736   END SUBROUTINE mppmin_real
2737
2738   SUBROUTINE mppmax_a_int( karr, kdim ,kcom)
2739      INTEGER, DIMENSION(:) :: karr
2740      INTEGER               :: kdim
2741      INTEGER, OPTIONAL     :: kcom 
2742      WRITE(*,*) 'mppmax_a_int: You should not have seen this print! error?', kdim, karr(1), kcom
2743   END SUBROUTINE mppmax_a_int
2744
2745   SUBROUTINE mppmax_int( kint, kcom)
2746      INTEGER               :: kint
2747      INTEGER, OPTIONAL     :: kcom 
2748      WRITE(*,*) 'mppmax_int: You should not have seen this print! error?', kint, kcom
2749   END SUBROUTINE mppmax_int
2750
2751   SUBROUTINE mppmin_a_int( karr, kdim, kcom )
2752      INTEGER, DIMENSION(:) :: karr
2753      INTEGER               :: kdim
2754      INTEGER, OPTIONAL     :: kcom 
2755      WRITE(*,*) 'mppmin_a_int: You should not have seen this print! error?', kdim, karr(1), kcom
2756   END SUBROUTINE mppmin_a_int
2757
2758   SUBROUTINE mppmin_int( kint, kcom )
2759      INTEGER               :: kint
2760      INTEGER, OPTIONAL     :: kcom 
2761      WRITE(*,*) 'mppmin_int: You should not have seen this print! error?', kint, kcom
2762   END SUBROUTINE mppmin_int
2763
2764   SUBROUTINE mppobc_1d( parr, kd1, kd2, kl, kk, ktype, kij, knum )
2765      INTEGER  ::   kd1, kd2, kl , kk, ktype, kij, knum
2766      REAL, DIMENSION(:) ::   parr           ! variable array
2767      WRITE(*,*) 'mppobc: You should not have seen this print! error?', parr(1), kd1, kd2, kl, kk, ktype, kij, knum
2768   END SUBROUTINE mppobc_1d
2769
2770   SUBROUTINE mppobc_2d( parr, kd1, kd2, kl, kk, ktype, kij, knum )
2771      INTEGER  ::   kd1, kd2, kl , kk, ktype, kij, knum
2772      REAL, DIMENSION(:,:) ::   parr           ! variable array
2773      WRITE(*,*) 'mppobc: You should not have seen this print! error?', parr(1,1), kd1, kd2, kl, kk, ktype, kij, knum
2774   END SUBROUTINE mppobc_2d
2775
2776   SUBROUTINE mppobc_3d( parr, kd1, kd2, kl, kk, ktype, kij, knum )
2777      INTEGER  ::   kd1, kd2, kl , kk, ktype, kij, knum
2778      REAL, DIMENSION(:,:,:) ::   parr           ! variable array
2779      WRITE(*,*) 'mppobc: You should not have seen this print! error?', parr(1,1,1), kd1, kd2, kl, kk, ktype, kij, knum
2780   END SUBROUTINE mppobc_3d
2781
2782   SUBROUTINE mppobc_4d( parr, kd1, kd2, kl, kk, ktype, kij, knum )
2783      INTEGER  ::   kd1, kd2, kl , kk, ktype, kij, knum
2784      REAL, DIMENSION(:,:,:,:) ::   parr           ! variable array
2785      WRITE(*,*) 'mppobc: You should not have seen this print! error?', parr(1,1,1,1), kd1, kd2, kl, kk, ktype, kij, knum
2786   END SUBROUTINE mppobc_4d
2787
2788   SUBROUTINE mpp_minloc2d( ptab, pmask, pmin, ki, kj )
2789      REAL                   :: pmin
2790      REAL , DIMENSION (:,:) :: ptab, pmask
2791      INTEGER :: ki, kj
2792      WRITE(*,*) 'mpp_minloc2d: You should not have seen this print! error?', pmin, ki, kj, ptab(1,1), pmask(1,1)
2793   END SUBROUTINE mpp_minloc2d
2794
2795   SUBROUTINE mpp_minloc3d( ptab, pmask, pmin, ki, kj, kk )
2796      REAL                     :: pmin
2797      REAL , DIMENSION (:,:,:) :: ptab, pmask
2798      INTEGER :: ki, kj, kk
2799      WRITE(*,*) 'mpp_minloc3d: You should not have seen this print! error?', pmin, ki, kj, kk, ptab(1,1,1), pmask(1,1,1)
2800   END SUBROUTINE mpp_minloc3d
2801
2802   SUBROUTINE mpp_maxloc2d( ptab, pmask, pmax, ki, kj )
2803      REAL                   :: pmax
2804      REAL , DIMENSION (:,:) :: ptab, pmask
2805      INTEGER :: ki, kj
2806      WRITE(*,*) 'mpp_maxloc2d: You should not have seen this print! error?', pmax, ki, kj, ptab(1,1), pmask(1,1)
2807   END SUBROUTINE mpp_maxloc2d
2808
2809   SUBROUTINE mpp_maxloc3d( ptab, pmask, pmax, ki, kj, kk )
2810      REAL                     :: pmax
2811      REAL , DIMENSION (:,:,:) :: ptab, pmask
2812      INTEGER :: ki, kj, kk
2813      WRITE(*,*) 'mpp_maxloc3d: You should not have seen this print! error?', pmax, ki, kj, kk, ptab(1,1,1), pmask(1,1,1)
2814   END SUBROUTINE mpp_maxloc3d
2815
2816   SUBROUTINE mppstop
2817      WRITE(*,*) 'mppstop: You should not have seen this print if running in mpp mode! error?...'
2818      WRITE(*,*) 'mppstop: ..otherwise this is a stop condition raised by ctl_stop in single processor mode'
2819      STOP
2820   END SUBROUTINE mppstop
2821
2822   SUBROUTINE mpp_ini_ice( kcom, knum )
2823      INTEGER :: kcom, knum
2824      WRITE(*,*) 'mpp_ini_ice: You should not have seen this print! error?', kcom, knum
2825   END SUBROUTINE mpp_ini_ice
2826
2827   SUBROUTINE mpp_ini_znl( knum )
2828      INTEGER :: knum
2829      WRITE(*,*) 'mpp_ini_znl: You should not have seen this print! error?', knum
2830   END SUBROUTINE mpp_ini_znl
2831
2832   SUBROUTINE mpp_comm_free( kcom )
2833      INTEGER :: kcom
2834      WRITE(*,*) 'mpp_comm_free: You should not have seen this print! error?', kcom
2835   END SUBROUTINE mpp_comm_free
2836#endif
2837
2838   !!----------------------------------------------------------------------
2839   !!   All cases:         ctl_stop, ctl_warn, get_unit, ctl_opn   routines
2840   !!----------------------------------------------------------------------
2841
2842   SUBROUTINE ctl_stop( cd1, cd2, cd3, cd4, cd5 ,   &
2843      &                 cd6, cd7, cd8, cd9, cd10 )
2844      !!----------------------------------------------------------------------
2845      !!                  ***  ROUTINE  stop_opa  ***
2846      !!
2847      !! ** Purpose :   print in ocean.outpput file a error message and
2848      !!                increment the error number (nstop) by one.
2849      !!----------------------------------------------------------------------
2850      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd1, cd2, cd3, cd4, cd5
2851      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd6, cd7, cd8, cd9, cd10
2852      !!----------------------------------------------------------------------
2853      !
2854      nstop = nstop + 1 
2855      IF(lwp) THEN
2856         WRITE(numout,cform_err)
2857         IF( PRESENT(cd1 ) )   WRITE(numout,*) cd1
2858         IF( PRESENT(cd2 ) )   WRITE(numout,*) cd2
2859         IF( PRESENT(cd3 ) )   WRITE(numout,*) cd3
2860         IF( PRESENT(cd4 ) )   WRITE(numout,*) cd4
2861         IF( PRESENT(cd5 ) )   WRITE(numout,*) cd5
2862         IF( PRESENT(cd6 ) )   WRITE(numout,*) cd6
2863         IF( PRESENT(cd7 ) )   WRITE(numout,*) cd7
2864         IF( PRESENT(cd8 ) )   WRITE(numout,*) cd8
2865         IF( PRESENT(cd9 ) )   WRITE(numout,*) cd9
2866         IF( PRESENT(cd10) )   WRITE(numout,*) cd10
2867      ENDIF
2868                               CALL FLUSH(numout    )
2869      IF( numstp     /= -1 )   CALL FLUSH(numstp    )
2870      IF( numsol     /= -1 )   CALL FLUSH(numsol    )
2871      IF( numevo_ice /= -1 )   CALL FLUSH(numevo_ice)
2872      !
2873      IF( cd1 == 'STOP' ) THEN
2874         IF(lwp) WRITE(numout,*)  'huge E-R-R-O-R : immediate stop'
2875         CALL mppstop()
2876      ENDIF
2877      !
2878   END SUBROUTINE ctl_stop
2879
2880
2881   SUBROUTINE ctl_warn( cd1, cd2, cd3, cd4, cd5,   &
2882      &                 cd6, cd7, cd8, cd9, cd10 )
2883      !!----------------------------------------------------------------------
2884      !!                  ***  ROUTINE  stop_warn  ***
2885      !!
2886      !! ** Purpose :   print in ocean.outpput file a error message and
2887      !!                increment the warning number (nwarn) by one.
2888      !!----------------------------------------------------------------------
2889      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd1, cd2, cd3, cd4, cd5
2890      CHARACTER(len=*), INTENT(in), OPTIONAL ::  cd6, cd7, cd8, cd9, cd10
2891      !!----------------------------------------------------------------------
2892      !
2893      nwarn = nwarn + 1 
2894      IF(lwp) THEN
2895         WRITE(numout,cform_war)
2896         IF( PRESENT(cd1 ) ) WRITE(numout,*) cd1
2897         IF( PRESENT(cd2 ) ) WRITE(numout,*) cd2
2898         IF( PRESENT(cd3 ) ) WRITE(numout,*) cd3
2899         IF( PRESENT(cd4 ) ) WRITE(numout,*) cd4
2900         IF( PRESENT(cd5 ) ) WRITE(numout,*) cd5
2901         IF( PRESENT(cd6 ) ) WRITE(numout,*) cd6
2902         IF( PRESENT(cd7 ) ) WRITE(numout,*) cd7
2903         IF( PRESENT(cd8 ) ) WRITE(numout,*) cd8
2904         IF( PRESENT(cd9 ) ) WRITE(numout,*) cd9
2905         IF( PRESENT(cd10) ) WRITE(numout,*) cd10
2906      ENDIF
2907      CALL FLUSH(numout)
2908      !
2909   END SUBROUTINE ctl_warn
2910
2911
2912   SUBROUTINE ctl_opn( knum, cdfile, cdstat, cdform, cdacce, klengh, kout, ldwp, karea )
2913      !!----------------------------------------------------------------------
2914      !!                  ***  ROUTINE ctl_opn  ***
2915      !!
2916      !! ** Purpose :   Open file and check if required file is available.
2917      !!
2918      !! ** Method  :   Fortan open
2919      !!----------------------------------------------------------------------
2920      INTEGER          , INTENT(  out) ::   knum      ! logical unit to open
2921      CHARACTER(len=*) , INTENT(in   ) ::   cdfile    ! file name to open
2922      CHARACTER(len=*) , INTENT(in   ) ::   cdstat    ! disposition specifier
2923      CHARACTER(len=*) , INTENT(in   ) ::   cdform    ! formatting specifier
2924      CHARACTER(len=*) , INTENT(in   ) ::   cdacce    ! access specifier
2925      INTEGER          , INTENT(in   ) ::   klengh    ! record length
2926      INTEGER          , INTENT(in   ) ::   kout      ! number of logical units for write
2927      LOGICAL          , INTENT(in   ) ::   ldwp      ! boolean term for print
2928      INTEGER, OPTIONAL, INTENT(in   ) ::   karea     ! proc number
2929      !!
2930      CHARACTER(len=80) ::   clfile
2931      INTEGER           ::   iost
2932      !!----------------------------------------------------------------------
2933
2934      ! adapt filename
2935      ! ----------------
2936      clfile = TRIM(cdfile)
2937      IF( PRESENT( karea ) ) THEN
2938         IF( karea > 1 )   WRITE(clfile, "(a,'_',i4.4)") TRIM(clfile), karea-1
2939      ENDIF
2940#if defined key_agrif
2941      IF( .NOT. Agrif_Root() )   clfile = TRIM(Agrif_CFixed())//'_'//TRIM(clfile)
2942      knum=Agrif_Get_Unit()
2943#else
2944      knum=get_unit()
2945#endif
2946
2947      iost=0
2948      IF( cdacce(1:6) == 'DIRECT' )  THEN
2949         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat, RECL=klengh, ERR=100, IOSTAT=iost )
2950      ELSE
2951         OPEN( UNIT=knum, FILE=clfile, FORM=cdform, ACCESS=cdacce, STATUS=cdstat             , ERR=100, IOSTAT=iost )
2952      ENDIF
2953      IF( iost == 0 ) THEN
2954         IF(ldwp) THEN
2955            WRITE(kout,*) '     file   : ', clfile,' open ok'
2956            WRITE(kout,*) '     unit   = ', knum
2957            WRITE(kout,*) '     status = ', cdstat
2958            WRITE(kout,*) '     form   = ', cdform
2959            WRITE(kout,*) '     access = ', cdacce
2960            WRITE(kout,*)
2961         ENDIF
2962      ENDIF
2963100   CONTINUE
2964      IF( iost /= 0 ) THEN
2965         IF(ldwp) THEN
2966            WRITE(kout,*)
2967            WRITE(kout,*) ' ===>>>> : bad opening file: ', clfile
2968            WRITE(kout,*) ' =======   ===  '
2969            WRITE(kout,*) '           unit   = ', knum
2970            WRITE(kout,*) '           status = ', cdstat
2971            WRITE(kout,*) '           form   = ', cdform
2972            WRITE(kout,*) '           access = ', cdacce
2973            WRITE(kout,*) '           iostat = ', iost
2974            WRITE(kout,*) '           we stop. verify the file '
2975            WRITE(kout,*)
2976         ENDIF
2977         STOP 'ctl_opn bad opening'
2978      ENDIF
2979     
2980   END SUBROUTINE ctl_opn
2981
2982
2983   INTEGER FUNCTION get_unit()
2984      !!----------------------------------------------------------------------
2985      !!                  ***  FUNCTION  get_unit  ***
2986      !!
2987      !! ** Purpose :   return the index of an unused logical unit
2988      !!----------------------------------------------------------------------
2989      LOGICAL :: llopn 
2990      !!----------------------------------------------------------------------
2991      !
2992      get_unit = 15   ! choose a unit that is big enough then it is not already used in NEMO
2993      llopn = .TRUE.
2994      DO WHILE( (get_unit < 998) .AND. llopn )
2995         get_unit = get_unit + 1
2996         INQUIRE( unit = get_unit, opened = llopn )
2997      END DO
2998      IF( (get_unit == 999) .AND. llopn ) THEN
2999         CALL ctl_stop( 'get_unit: All logical units until 999 are used...' )
3000         get_unit = -1
3001      ENDIF
3002      !
3003   END FUNCTION get_unit
3004
3005   !!----------------------------------------------------------------------
3006END MODULE lib_mpp
Note: See TracBrowser for help on using the repository browser.