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 @ 2359

Last change on this file since 2359 was 2304, checked in by rblod, 14 years ago

Choose one option for mpp reproducibility, see ticket #743

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