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

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/lib_mpp.F90 @ 2412

Last change on this file since 2412 was 2363, checked in by rblod, 14 years ago

Fixes to compile with key_agrif and key_mpp_mpi

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