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

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

source: branches/DEV_1879_mpp_sca/NEMO/OPA_SRC/lib_mpp.F90 @ 1926

Last change on this file since 1926 was 1926, checked in by acc, 14 years ago

First implementation of mpp scalability modifications (branch:DEV_1879_mpp_sca

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