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

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

source: trunk/NEMO/OPA_SRC/lib_mpp.F90 @ 1502

Last change on this file since 1502 was 1490, checked in by rblod, 15 years ago

Supress remaining coments in French, see ticket #379

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