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 trunk/NEMOGCM/NEMO/OPA_SRC/LBC – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90 @ 3918

Last change on this file since 3918 was 3918, checked in by smasson, 11 years ago

trunk: better fortran and error in the conv of agrif, see ticket #1111 and #1112

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