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

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

source: branches/2011/dev_UKM0_2011/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90 @ 3062

Last change on this file since 3062 was 3062, checked in by rfurner, 12 years ago

ticket #885. added in changes from branches/2011/UKMO_MERCATOR_obc_bdy_merge@2888

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