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

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

Last change on this file since 2481 was 2481, checked in by smasson, 13 years ago

v33b: additional cleaning in libmpp, see ticket #779

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