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.
mppini.F90 in NEMO/branches/2021/dev_r14447_HPC-07_Irrmann_try_new_pt2pt/src/OCE/LBC – NEMO

source: NEMO/branches/2021/dev_r14447_HPC-07_Irrmann_try_new_pt2pt/src/OCE/LBC/mppini.F90 @ 15072

Last change on this file since 15072 was 15072, checked in by girrmann, 3 years ago

add structure for new RMA communications, small asynchronous communications debugging

  • Property svn:keywords set to Id
File size: 65.6 KB
Line 
1MODULE mppini
2   !!======================================================================
3   !!                       ***  MODULE mppini   ***
4   !! Ocean initialization : distributed memory computing initialization
5   !!======================================================================
6   !! History :  6.0  !  1994-11  (M. Guyon)  Original code
7   !!  OPA       7.0  !  1995-04  (J. Escobar, M. Imbard)
8   !!            8.0  !  1998-05  (M. Imbard, J. Escobar, L. Colombet )  SHMEM and MPI versions
9   !!  NEMO      1.0  !  2004-01  (G. Madec, J.M Molines)  F90 : free form , north fold jpni > 1
10   !!            3.4  !  2011-10  (A. C. Coward, NOCS & J. Donners, PRACE)  add init_nfdcom
11   !!            3.   !  2013-06  (I. Epicoco, S. Mocavero, CMCC)  init_nfdcom: setup avoiding MPI communication
12   !!            4.0  !  2016-06  (G. Madec)  use domain configuration file instead of bathymetry file
13   !!            4.0  !  2017-06  (J.M. Molines, T. Lovato) merge of mppini and mppini_2
14   !!----------------------------------------------------------------------
15
16   !!----------------------------------------------------------------------
17   !!  mpp_init       : Lay out the global domain over processors with/without land processor elimination
18   !!      init_ioipsl: IOIPSL initialization in mpp
19   !!      init_nfdcom: Setup for north fold exchanges with explicit point-to-point messaging
20   !!      init_doloop: set the starting/ending indices of DO-loop used in do_loop_substitute
21   !!----------------------------------------------------------------------
22   USE dom_oce        ! ocean space and time domain
23   USE bdy_oce        ! open BounDarY
24   !
25   USE lbcnfd  , ONLY : isendto, nsndto ! Setup of north fold exchanges
26   USE lib_mpp        ! distribued memory computing library
27   USE iom            ! nemo I/O library
28   USE ioipsl         ! I/O IPSL library
29   USE in_out_manager ! I/O Manager
30
31   IMPLICIT NONE
32   PRIVATE
33
34   PUBLIC   mpp_init       ! called by nemogcm.F90
35   PUBLIC   mpp_getnum     ! called by prtctl
36   PUBLIC   mpp_basesplit  ! called by prtctl
37   PUBLIC   mpp_is_ocean   ! called by prtctl
38
39   INTEGER ::   numbot = -1   ! 'bottom_level' local logical unit
40   INTEGER ::   numbdy = -1   ! 'bdy_msk'      local logical unit
41   !!----------------------------------------------------------------------
42   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
43   !! $Id$
44   !! Software governed by the CeCILL license (see ./LICENSE)
45   !!----------------------------------------------------------------------
46CONTAINS
47
48#if defined key_mpi_off
49   !!----------------------------------------------------------------------
50   !!   Default option :                            shared memory computing
51   !!----------------------------------------------------------------------
52
53   SUBROUTINE mpp_init
54      !!----------------------------------------------------------------------
55      !!                  ***  ROUTINE mpp_init  ***
56      !!
57      !! ** Purpose :   Lay out the global domain over processors.
58      !!
59      !! ** Method  :   Shared memory computing, set the local processor
60      !!              variables to the value of the global domain
61      !!----------------------------------------------------------------------
62      !
63      nn_hls = 1
64      jpiglo = Ni0glo + 2 * nn_hls
65      jpjglo = Nj0glo + 2 * nn_hls
66      jpimax = jpiglo
67      jpjmax = jpjglo
68      jpi    = jpiglo
69      jpj    = jpjglo
70      jpk    = MAX( 2, jpkglo )
71      jpij   = jpi*jpj
72      jpni   = 1
73      jpnj   = 1
74      jpnij  = jpni*jpnj
75      nimpp  = 1
76      njmpp  = 1
77      nidom  = FLIO_DOM_NONE
78      !
79      CALL init_doloop                       ! set start/end indices or do-loop depending on the halo width value (nn_hls)
80      !
81      IF(lwp) THEN
82         WRITE(numout,*)
83         WRITE(numout,*) 'mpp_init : NO massively parallel processing'
84         WRITE(numout,*) '~~~~~~~~ '
85         WRITE(numout,*) '   l_Iperio = ', l_Iperio, '    l_Jperio = ', l_Jperio
86         WRITE(numout,*) '     njmpp  = ', njmpp
87      ENDIF
88      !
89#if defined key_agrif
90    IF (.NOT.agrif_root()) THEN
91      call agrif_nemo_init()
92    ENDIF
93#endif
94   END SUBROUTINE mpp_init
95
96#else
97   !!----------------------------------------------------------------------
98   !!                   MPI massively parallel processing
99   !!----------------------------------------------------------------------
100
101
102   SUBROUTINE mpp_init
103      !!----------------------------------------------------------------------
104      !!                  ***  ROUTINE mpp_init  ***
105      !!
106      !! ** Purpose :   Lay out the global domain over processors.
107      !!      If land processors are to be eliminated, this program requires the
108      !!      presence of the domain configuration file. Land processors elimination
109      !!      is performed if jpni x jpnj /= jpnij. In this case, using the MPP_PREP
110      !!      preprocessing tool, help for defining the best cutting out.
111      !!
112      !! ** Method  :   Global domain is distributed in smaller local domains.
113      !!      Periodic condition is a function of the local domain position
114      !!      (global boundary or neighbouring domain) and of the global periodic
115      !!
116      !! ** Action : - set domain parameters
117      !!                    nimpp     : longitudinal index
118      !!                    njmpp     : latitudinal  index
119      !!                    narea     : number for local area
120      !!                    mpinei    : number of neighboring domains (starting at 0, -1 if no neighbourg)
121      !!----------------------------------------------------------------------
122      INTEGER ::   ji, jj, jn, jp, jh
123      INTEGER ::   ii, ij, ii2, ij2
124      INTEGER ::   inijmin   ! number of oce subdomains
125      INTEGER ::   inum, inum0
126      INTEGER ::   ifreq, il1, imil, il2, ijm1
127      INTEGER ::   ierr, ios
128      INTEGER ::   inbi, inbj, iimax, ijmax, icnt1, icnt2
129      INTEGER, DIMENSION(16*n_hlsmax) :: ichanged
130      INTEGER, ALLOCATABLE, DIMENSION(:    ) ::   iin, ijn
131      INTEGER, ALLOCATABLE, DIMENSION(:,:  ) ::   iimppt, ijpi, ipproc
132      INTEGER, ALLOCATABLE, DIMENSION(:,:  ) ::   ijmppt, ijpj
133      INTEGER, ALLOCATABLE, DIMENSION(:,:  ) ::   impi
134      INTEGER, ALLOCATABLE, DIMENSION(:,:,:) ::   inei
135      LOGICAL ::   llbest, llauto
136      LOGICAL ::   llwrtlay
137      LOGICAL ::   llmpi_Iperio, llmpi_Jperio, llmpiNFold
138      LOGICAL ::   ln_listonly
139      LOGICAL, ALLOCATABLE, DIMENSION(:,:  ) ::   llisOce  ! is not land-domain only?
140      LOGICAL, ALLOCATABLE, DIMENSION(:,:,:) ::   llnei    ! are neighbourgs existing?
141      NAMELIST/nambdy/ ln_bdy, nb_bdy, ln_coords_file, cn_coords_file,           &
142           &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
143           &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
144           &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
145           &             cn_ice, nn_ice_dta,                                     &
146           &             ln_vol, nn_volctl, nn_rimwidth
147      NAMELIST/nammpp/ jpni, jpnj, nn_hls, ln_nnogather, ln_listonly, nn_comm, ln_tspers, ln_async, ln_RMA
148      !!----------------------------------------------------------------------
149      !
150      llwrtlay = lwm .OR. sn_cfctl%l_layout
151      !
152      !  0. read namelists parameters
153      ! -----------------------------------
154      !
155      READ  ( numnam_ref, nammpp, IOSTAT = ios, ERR = 901 )
156901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nammpp in reference namelist' )
157      READ  ( numnam_cfg, nammpp, IOSTAT = ios, ERR = 902 )
158902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nammpp in configuration namelist' )
159      !
160      nn_hls = MAX(1, nn_hls)   ! nn_hls must be > 0
161      IF(lwp) THEN
162            WRITE(numout,*) '   Namelist nammpp'
163         IF( jpni < 1 .OR. jpnj < 1  ) THEN
164            WRITE(numout,*) '      jpni and jpnj will be calculated automatically'
165         ELSE
166            WRITE(numout,*) '      processor grid extent in i                            jpni = ', jpni
167            WRITE(numout,*) '      processor grid extent in j                            jpnj = ', jpnj
168         ENDIF
169            WRITE(numout,*) '      avoid use of mpi_allgather at the north fold  ln_nnogather = ', ln_nnogather
170            WRITE(numout,*) '      halo width (applies to both rows and columns)       nn_hls = ', nn_hls
171            WRITE(numout,*) '      communication type                                 nn_comm = ', nn_comm
172            WRITE(numout,*) '      switch to persistent calls in time-splitting    ln_tspers = ', ln_tspers
173            WRITE(numout,*) '      switch to asynchronous com in time-splitting    ln_async = ' , ln_async
174            WRITE(numout,*) '      switch to RMA communications in time-splitting    ln_RMA = ' , ln_RMA
175      ENDIF
176      !
177      IF(lwm)   WRITE( numond, nammpp )
178      !
179      jpiglo = Ni0glo + 2 * nn_hls
180      jpjglo = Nj0glo + 2 * nn_hls
181      !
182      ! do we need to take into account bdy_msk?
183      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 903)
184903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nambdy in reference namelist (mppini)' )
185      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 904 )
186904   IF( ios >  0 )   CALL ctl_nam ( ios , 'nambdy in configuration namelist (mppini)' )
187      !
188      IF(               ln_read_cfg ) CALL iom_open( cn_domcfg,    numbot )
189      IF( ln_bdy .AND. ln_mask_file ) CALL iom_open( cn_mask_file, numbdy )
190      !
191      IF( ln_listonly )   CALL bestpartition( MAX(mppsize,jpni*jpnj), ldlist = .TRUE. )   ! must be done by all core
192      !
193      !  1. Dimension arrays for subdomains
194      ! -----------------------------------
195      !
196      ! If dimensions of MPI processes grid weren't specified in the namelist file
197      ! then we calculate them here now that we have our communicator size
198      IF(lwp) THEN
199         WRITE(numout,*)
200         WRITE(numout,*) 'mpp_init:'
201         WRITE(numout,*) '~~~~~~~~ '
202      ENDIF
203      IF( jpni < 1 .OR. jpnj < 1 ) THEN
204         CALL bestpartition( mppsize, jpni, jpnj )           ! best mpi decomposition for mppsize mpi processes
205         llauto = .TRUE.
206         llbest = .TRUE.
207      ELSE
208         llauto = .FALSE.
209         CALL bestpartition( mppsize, inbi, inbj, icnt2 )    ! best mpi decomposition for mppsize mpi processes
210         ! largest subdomain size for mpi decoposition jpni*jpnj given in the namelist
211         CALL mpp_basesplit( jpiglo, jpjglo, nn_hls, jpni, jpnj, jpimax, jpjmax )
212         ! largest subdomain size for mpi decoposition inbi*inbj given by bestpartition
213         CALL mpp_basesplit( jpiglo, jpjglo, nn_hls, inbi, inbj,  iimax,  ijmax )
214         icnt1 = jpni*jpnj - mppsize   ! number of land subdomains that should be removed to use mppsize mpi processes
215         IF(lwp) THEN
216            WRITE(numout,9000) '   The chosen domain decomposition ', jpni, ' x ', jpnj, ' with ', icnt1, ' land subdomains'
217            WRITE(numout,9002) '      - uses a total of ',  mppsize,' mpi process'
218            WRITE(numout,9000) '      - has mpi subdomains with a maximum size of (jpi = ', jpimax, ', jpj = ', jpjmax,   &
219               &                                                                ', jpi*jpj = ', jpimax*jpjmax, ')'
220            WRITE(numout,9000) '   The best domain decompostion ', inbi, ' x ', inbj, ' with ', icnt2, ' land subdomains'
221            WRITE(numout,9002) '      - uses a total of ',  inbi*inbj-icnt2,' mpi process'
222            WRITE(numout,9000) '      - has mpi subdomains with a maximum size of (jpi = ',  iimax, ', jpj = ',  ijmax,   &
223               &                                                             ', jpi*jpj = ',  iimax* ijmax, ')'
224         ENDIF
225         IF( iimax*ijmax < jpimax*jpjmax ) THEN   ! chosen subdomain size is larger that the best subdomain size
226            llbest = .FALSE.
227            IF ( inbi*inbj-icnt2 < mppsize ) THEN
228               WRITE(ctmp1,*) '   ==> You could therefore have smaller mpi subdomains with less mpi processes'
229            ELSE
230               WRITE(ctmp1,*) '   ==> You could therefore have smaller mpi subdomains with the same number of mpi processes'
231            ENDIF
232            CALL ctl_warn( ' ', ctmp1, ' ', '    ---   YOU ARE WASTING CPU...   ---', ' ' )
233         ELSE IF ( iimax*ijmax == jpimax*jpjmax .AND. (inbi*inbj-icnt2) <  mppsize) THEN
234            llbest = .FALSE.
235            WRITE(ctmp1,*) '   ==> You could therefore have the same mpi subdomains size with less mpi processes'
236            CALL ctl_warn( ' ', ctmp1, ' ', '    ---   YOU ARE WASTING CPU...   ---', ' ' )
237         ELSE
238            llbest = .TRUE.
239         ENDIF
240      ENDIF
241
242      ! look for land mpi subdomains...
243      ALLOCATE( llisOce(jpni,jpnj) )
244      CALL mpp_is_ocean( llisOce )
245      inijmin = COUNT( llisOce )   ! number of oce subdomains
246
247      IF( mppsize < inijmin ) THEN   ! too many oce subdomains: can happen only if jpni and jpnj are prescribed...
248         WRITE(ctmp1,9001) '   With this specified domain decomposition: jpni = ', jpni, ' jpnj = ', jpnj
249         WRITE(ctmp2,9002) '   we can eliminate only ', jpni*jpnj - inijmin, ' land mpi subdomains therefore '
250         WRITE(ctmp3,9001) '   the number of ocean mpi subdomains (', inijmin,') exceed the number of MPI processes:', mppsize
251         WRITE(ctmp4,*) '   ==>>> There is the list of best domain decompositions you should use: '
252         CALL ctl_stop( ctmp1, ctmp2, ctmp3, ' ', ctmp4, ' ' )
253         CALL bestpartition( mppsize, ldlist = .TRUE. )   ! must be done by all core
254      ENDIF
255
256      IF( mppsize > jpni*jpnj ) THEN   ! not enough mpi subdomains for the total number of mpi processes
257         IF(lwp) THEN
258            WRITE(numout,9003) '   The number of mpi processes: ', mppsize
259            WRITE(numout,9003) '   exceeds the maximum number of subdomains (ocean+land) = ', jpni*jpnj
260            WRITE(numout,9001) '   defined by the following domain decomposition: jpni = ', jpni, ' jpnj = ', jpnj
261            WRITE(numout,   *) '   You should: '
262           IF( llauto ) THEN
263               WRITE(numout,*) '     - either prescribe your domain decomposition with the namelist variables'
264               WRITE(numout,*) '       jpni and jpnj to match the number of mpi process you want to use, '
265               WRITE(numout,*) '       even IF it not the best choice...'
266               WRITE(numout,*) '     - or keep the automatic and optimal domain decomposition by picking up one'
267               WRITE(numout,*) '       of the number of mpi process proposed in the list bellow'
268            ELSE
269               WRITE(numout,*) '     - either properly prescribe your domain decomposition with jpni and jpnj'
270               WRITE(numout,*) '       in order to be consistent with the number of mpi process you want to use'
271               WRITE(numout,*) '       even IF it not the best choice...'
272               WRITE(numout,*) '     - or use the automatic and optimal domain decomposition and pick up one of'
273               WRITE(numout,*) '       the domain decomposition proposed in the list bellow'
274            ENDIF
275            WRITE(numout,*)
276         ENDIF
277         CALL bestpartition( mppsize, ldlist = .TRUE. )   ! must be done by all core
278      ENDIF
279
280      jpnij = mppsize   ! force jpnij definition <-- remove as much land subdomains as needed to reach this condition
281      IF( mppsize > inijmin ) THEN
282         WRITE(ctmp1,9003) '   The number of mpi processes: ', mppsize
283         WRITE(ctmp2,9003) '   exceeds the maximum number of ocean subdomains = ', inijmin
284         WRITE(ctmp3,9002) '   we suppressed ', jpni*jpnj - mppsize, ' land subdomains '
285         WRITE(ctmp4,9002) '   BUT we had to keep ', mppsize - inijmin, ' land subdomains that are useless...'
286         CALL ctl_warn( ctmp1, ctmp2, ctmp3, ctmp4, ' ', '    --- YOU ARE WASTING CPU... ---', ' ' )
287      ELSE   ! mppsize = inijmin
288         IF(lwp) THEN
289            IF(llbest) WRITE(numout,*) '   ==> you use the best mpi decomposition'
290            WRITE(numout,*)
291            WRITE(numout,9003) '   Number of mpi processes: ', mppsize
292            WRITE(numout,9003) '   Number of ocean subdomains = ', inijmin
293            WRITE(numout,9003) '   Number of suppressed land subdomains = ', jpni*jpnj - inijmin
294            WRITE(numout,*)
295         ENDIF
296      ENDIF
2979000  FORMAT (a, i4, a, i4, a, i7, a)
2989001  FORMAT (a, i4, a, i4)
2999002  FORMAT (a, i4, a)
3009003  FORMAT (a, i5)
301
302      ALLOCATE( nfimpp(jpni), nfproc(jpni), nfjpi(jpni),   &
303         &      iin(jpnij), ijn(jpnij),   &
304         &      iimppt(jpni,jpnj), ijmppt(jpni,jpnj), ijpi(jpni,jpnj), ijpj(jpni,jpnj), ipproc(jpni,jpnj),   &
305         &      inei(8,jpni,jpnj), llnei(8,jpni,jpnj),   &
306         &      impi(8,jpnij),   &
307         &      STAT=ierr )
308      CALL mpp_sum( 'mppini', ierr )
309      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'mpp_init: unable to allocate standard ocean arrays' )
310
311#if defined key_agrif
312      IF( .NOT. Agrif_Root() ) THEN       ! AGRIF children: specific setting (cf. agrif_user.F90)
313         CALL agrif_nemo_init()
314      ENDIF
315#endif
316      !
317      !  2. Index arrays for subdomains
318      ! -----------------------------------
319      !
320      CALL mpp_basesplit( jpiglo, jpjglo, nn_hls, jpni, jpnj, jpimax, jpjmax, iimppt, ijmppt, ijpi, ijpj )
321      CALL mpp_getnum( llisOce, ipproc, iin, ijn )
322      !
323      ii = iin(narea)
324      ij = ijn(narea)
325      jpi   = ijpi(ii,ij)
326      jpj   = ijpj(ii,ij)
327      jpk   = MAX( 2, jpkglo )
328      jpij  = jpi*jpj
329      nimpp = iimppt(ii,ij)
330      njmpp = ijmppt(ii,ij)
331      !
332      CALL init_doloop                          ! set start/end indices of do-loop, depending on the halo width value (nn_hls)
333      !
334      IF(lwp) THEN
335         WRITE(numout,*)
336         WRITE(numout,*) 'MPI Message Passing MPI - domain lay out over processors'
337         WRITE(numout,*)
338         WRITE(numout,*) '   defines mpp subdomains'
339         WRITE(numout,*) '      jpni = ', jpni
340         WRITE(numout,*) '      jpnj = ', jpnj
341         WRITE(numout,*) '     jpnij = ', jpnij
342         WRITE(numout,*) '     nimpp = ', nimpp
343         WRITE(numout,*) '     njmpp = ', njmpp
344         WRITE(numout,*)
345         WRITE(numout,*) '      sum ijpi(i,1) = ', sum(ijpi(:,1)), ' jpiglo = ', jpiglo
346         WRITE(numout,*) '      sum ijpj(1,j) = ', SUM(ijpj(1,:)), ' jpjglo = ', jpjglo
347         
348         ! Subdomain grid print
349         ifreq = 4
350         il1 = 1
351         DO jn = 1, (jpni-1)/ifreq+1
352            il2 = MIN(jpni,il1+ifreq-1)
353            WRITE(numout,*)
354            WRITE(numout,9400) ('***',ji=il1,il2-1)
355            DO jj = jpnj, 1, -1
356               WRITE(numout,9403) ('   ',ji=il1,il2-1)
357               WRITE(numout,9402) jj, (ijpi(ji,jj),ijpj(ji,jj),ji=il1,il2)
358               WRITE(numout,9404) (ipproc(ji,jj),ji=il1,il2)
359               WRITE(numout,9403) ('   ',ji=il1,il2-1)
360               WRITE(numout,9400) ('***',ji=il1,il2-1)
361            END DO
362            WRITE(numout,9401) (ji,ji=il1,il2)
363            il1 = il1+ifreq
364         END DO
365 9400    FORMAT('           ***'   ,20('*************',a3)    )
366 9403    FORMAT('           *     ',20('         *   ',a3)    )
367 9401    FORMAT('              '   ,20('   ',i3,'          ') )
368 9402    FORMAT('       ',i3,' *  ',20(i3,'  x',i3,'   *   ') )
369 9404    FORMAT('           *  '   ,20('     ' ,i4,'   *   ') )
370      ENDIF
371      !
372      ! Store informations for the north pole folding communications
373      nfproc(:) = ipproc(:,jpnj)
374      nfimpp(:) = iimppt(:,jpnj)
375      nfjpi (:) =   ijpi(:,jpnj)
376      !
377      ! 3. Define Western, Eastern, Southern and Northern neighbors + corners in the subdomain grid reference
378      ! ------------------------------------------------------------------------------------------------------
379      !
380      ! note that North fold is has specific treatment for its MPI communications.
381      ! This must not be treated as a "usual" communication with a northern neighbor.
382      !    -> North fold processes have no Northern neighbor in the definition done bellow
383      !
384      llmpi_Iperio = jpni > 1 .AND. l_Iperio                         ! do i-periodicity with an MPI communication?
385      llmpi_Jperio = jpnj > 1 .AND. l_Jperio                         ! do j-periodicity with an MPI communication?
386      !
387      l_SelfPerio(1:2) = l_Iperio .AND. jpni == 1                    !  west,  east periodicity by itself
388      l_SelfPerio(3:4) = l_Jperio .AND. jpnj == 1                    ! south, north periodicity by itself
389      l_SelfPerio(5:8) = l_SelfPerio(jpwe) .AND. l_SelfPerio(jpso)   ! corners bi-periodicity by itself
390      !
391      ! define neighbors mapping (1/2): default definition: ignore if neighbours are land-only subdomains or not
392      DO jj = 1, jpnj
393         DO ji = 1, jpni
394            !
395            IF ( llisOce(ji,jj) ) THEN                     ! this subdomain has some ocean: it has neighbours
396               !
397               inum0 = ji - 1 + ( jj - 1 ) * jpni             ! index in the subdomains grid. start at 0
398               !
399               ! Is there a neighbor?
400               llnei(jpwe,ji,jj) = ji >   1  .OR. llmpi_Iperio           ! West  nei exists if not the first column or llmpi_Iperio
401               llnei(jpea,ji,jj) = ji < jpni .OR. llmpi_Iperio           ! East  nei exists if not the last  column or llmpi_Iperio
402               llnei(jpso,ji,jj) = jj >   1  .OR. llmpi_Jperio           ! South nei exists if not the first line   or llmpi_Jperio
403               llnei(jpno,ji,jj) = jj < jpnj .OR. llmpi_Jperio           ! North nei exists if not the last  line   or llmpi_Jperio
404               llnei(jpsw,ji,jj) = llnei(jpwe,ji,jj) .AND. llnei(jpso,ji,jj)   ! So-We nei exists if both South and West nei exist
405               llnei(jpse,ji,jj) = llnei(jpea,ji,jj) .AND. llnei(jpso,ji,jj)   ! So-Ea nei exists if both South and East nei exist
406               llnei(jpnw,ji,jj) = llnei(jpwe,ji,jj) .AND. llnei(jpno,ji,jj)   ! No-We nei exists if both North and West nei exist
407               llnei(jpne,ji,jj) = llnei(jpea,ji,jj) .AND. llnei(jpno,ji,jj)   ! No-Ea nei exists if both North and East nei exist
408               !
409               ! Which index (starting at 0) have neighbors in the subdomains grid?
410               IF( llnei(jpwe,ji,jj) )   inei(jpwe,ji,jj) =            inum0 -    1 + jpni        * COUNT( (/ ji ==    1 /) )
411               IF( llnei(jpea,ji,jj) )   inei(jpea,ji,jj) =            inum0 +    1 - jpni        * COUNT( (/ ji == jpni /) )
412               IF( llnei(jpso,ji,jj) )   inei(jpso,ji,jj) =            inum0 - jpni + jpni * jpnj * COUNT( (/ jj ==    1 /) )
413               IF( llnei(jpno,ji,jj) )   inei(jpno,ji,jj) =            inum0 + jpni - jpni * jpnj * COUNT( (/ jj == jpnj /) )
414               IF( llnei(jpsw,ji,jj) )   inei(jpsw,ji,jj) = inei(jpso,ji,jj) -    1 + jpni        * COUNT( (/ ji ==    1 /) )
415               IF( llnei(jpse,ji,jj) )   inei(jpse,ji,jj) = inei(jpso,ji,jj) +    1 - jpni        * COUNT( (/ ji == jpni /) )
416               IF( llnei(jpnw,ji,jj) )   inei(jpnw,ji,jj) = inei(jpno,ji,jj) -    1 + jpni        * COUNT( (/ ji ==    1 /) )
417               IF( llnei(jpne,ji,jj) )   inei(jpne,ji,jj) = inei(jpno,ji,jj) +    1 - jpni        * COUNT( (/ ji == jpni /) )
418               !
419            ELSE                                           ! land-only domain has no neighbour
420               llnei(:,ji,jj) = .FALSE.
421            ENDIF
422            !
423         END DO
424      END DO
425      !
426      ! define neighbors mapping (2/2): check if neighbours are not land-only subdomains
427      DO jj = 1, jpnj
428         DO ji = 1, jpni
429            DO jn = 1, 8
430               IF( llnei(jn,ji,jj) ) THEN   ! if a neighbour is existing -> this should not be a land-only domain
431                  ii = 1 + MOD( inei(jn,ji,jj) , jpni )
432                  ij = 1 +      inei(jn,ji,jj) / jpni
433                  llnei(jn,ji,jj) = llisOce( ii, ij )
434               ENDIF
435            END DO
436         END DO
437      END DO
438      !
439      ! update index of the neighbours in the subdomains grid
440      WHERE( .NOT. llnei )   inei = -1
441      !
442      ! Save processor layout in ascii file
443      IF (llwrtlay) THEN
444         CALL ctl_opn( inum, 'layout.dat', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE., narea )
445         WRITE(inum,'(a)') '  jpnij jpimax jpjmax    jpk jpiglo jpjglo ( local:   narea    jpi    jpj )'
446         WRITE(inum,'(6i7,a,3i7,a)') jpnij,jpimax,jpjmax,jpk,jpiglo,jpjglo,' ( local: ',narea,jpi,jpj,' )'
447         WRITE(inum,*) 
448         WRITE(inum,       *) '------------------------------------'
449         WRITE(inum,'(a,i2)') ' Mapping of the default neighnourgs '
450         WRITE(inum,       *) '------------------------------------'
451         WRITE(inum,*) 
452         WRITE(inum,'(a)') '  rank    ii    ij   jpi   jpj nimpp njmpp mpiwe mpiea mpiso mpino mpisw mpise mpinw mpine'
453         DO jp = 1, jpnij
454            ii = iin(jp)
455            ij = ijn(jp)
456            WRITE(inum,'(15i6)')  jp-1, ii, ij, ijpi(ii,ij), ijpj(ii,ij), iimppt(ii,ij), ijmppt(ii,ij), inei(:,ii,ij)
457         END DO
458      ENDIF
459
460      !
461      ! 4. Define Western, Eastern, Southern and Northern neighbors + corners for each mpi process
462      ! ------------------------------------------------------------------------------------------
463      !
464      ! rewrite information from "subdomain grid" to mpi process list
465      ! Warning, for example:
466      !    position of the northern neighbor in the "subdomain grid"
467      !    position of the northern neighbor in the "mpi process list"
468     
469      ! default definition: no neighbors
470      impi(:,:) = -1   ! (starting at 0, -1 if no neighbourg)
471     
472      DO jp = 1, jpnij
473         ii = iin(jp)
474         ij = ijn(jp)
475         DO jn = 1, 8
476            IF( llnei(jn,ii,ij) ) THEN   ! must be tested as some land-domain can be kept to fit mppsize
477               ii2 = 1 + MOD( inei(jn,ii,ij) , jpni )
478               ij2 = 1 +      inei(jn,ii,ij) / jpni
479               impi(jn,jp) = ipproc( ii2, ij2 )
480            ENDIF
481         END DO
482      END DO
483
484      !
485      ! 4. keep information for the local process
486      ! -----------------------------------------
487      !
488      ! set default neighbours
489      mpinei(:) = impi(:,narea)
490      DO jh = 1, n_hlsmax
491         mpiSnei(jh,:) = impi(:,narea)   ! default definition
492         mpiRnei(jh,:) = impi(:,narea)
493      END DO
494      !
495      IF(lwp) THEN
496         WRITE(numout,*)
497         WRITE(numout,*) '   resulting internal parameters : '
498         WRITE(numout,*) '      narea = ', narea
499         WRITE(numout,*) '      mpi nei  west = ', mpinei(jpwe)  , '   mpi nei  east = ', mpinei(jpea)
500         WRITE(numout,*) '      mpi nei south = ', mpinei(jpso)  , '   mpi nei north = ', mpinei(jpno)
501         WRITE(numout,*) '      mpi nei so-we = ', mpinei(jpsw)  , '   mpi nei so-ea = ', mpinei(jpse)
502         WRITE(numout,*) '      mpi nei no-we = ', mpinei(jpnw)  , '   mpi nei no-ea = ', mpinei(jpne)
503      ENDIF
504      !                          ! Prepare mpp north fold
505      !
506      llmpiNFold =          jpni  > 1 .AND. l_NFold   ! is the North fold done with an MPI communication?
507      l_IdoNFold = ijn(narea) == jpnj .AND. l_NFold   ! is this process doing North fold?
508      !
509      IF( llmpiNFold ) THEN
510         CALL mpp_ini_north
511         IF (lwp) THEN
512            WRITE(numout,*)
513            WRITE(numout,*) '   ==>>>   North fold boundary prepared for jpni >1'
514         ENDIF
515         IF (llwrtlay) THEN      ! additional prints in layout.dat
516            WRITE(inum,*)
517            WRITE(inum,*)
518            WRITE(inum,*) 'Number of subdomains located along the north fold : ', ndim_rank_north
519            WRITE(inum,*) 'Rank of the subdomains located along the north fold : ', ndim_rank_north
520            DO jp = 1, ndim_rank_north, 5
521               WRITE(inum,*) nrank_north( jp:MINVAL( (/jp+4,ndim_rank_north/) ) )
522            END DO
523         ENDIF
524         IF ( l_IdoNFold .AND. ln_nnogather ) THEN
525            CALL init_nfdcom     ! northfold neighbour lists
526            IF (llwrtlay) THEN
527               WRITE(inum,*)
528               WRITE(inum,*) 'north fold exchanges with explicit point-to-point messaging :'
529               WRITE(inum,*) '   nsndto  : ', nsndto
530               WRITE(inum,*) '   isendto : ', isendto(1:nsndto)
531            ENDIF
532         ENDIF
533      ENDIF
534      !
535      CALL mpp_ini_nc(nn_hls)    ! Initialize communicator for neighbourhood collective communications
536      DO jh = 1, n_hlsmax
537         mpi_nc_com4(jh) = mpi_nc_com4(nn_hls)   ! default definition
538         mpi_nc_com8(jh) = mpi_nc_com8(nn_hls)
539      END DO
540      !
541      CALL init_excl_landpt      ! exclude exchanges which contain only land points
542      !
543      ! Save processor layout changes in ascii file
544      DO jh = 1, n_hlsmax    ! different halo size
545         DO ji = 1, 8
546            ichanged(16*(jh-1)  +ji) = COUNT( mpinei(ji:ji) /= mpiSnei(jh,ji:ji) )
547            ichanged(16*(jh-1)+8+ji) = COUNT( mpinei(ji:ji) /= mpiRnei(jh,ji:ji) )
548         END DO
549      END DO
550      CALL mpp_sum( "mpp_init", ichanged )   ! must be called by all processes
551      IF (llwrtlay) THEN
552         WRITE(inum,*) 
553         WRITE(inum,       *) '----------------------------------------------------------------------'
554         WRITE(inum,'(a,i2)') ' Mapping of the neighnourgs once excluding comm with only land points '
555         WRITE(inum,       *) '----------------------------------------------------------------------'
556         DO jh = 1, n_hlsmax    ! different halo size
557            WRITE(inum,*) 
558            WRITE(inum,'(a,i2)') 'halo size: ', jh
559            WRITE(inum,       *) '---------'
560            WRITE(inum,'(a)') '  rank    ii    ij mpiwe mpiea mpiso mpino mpisw mpise mpinw mpine'
561            WRITE(inum,   '(11i6,a)')  narea-1, iin(narea), ijn(narea),   mpinei(:), ' <- Org'
562            WRITE(inum,'(18x,8i6,a,i1,a)')   mpiSnei(jh,:), ' <- Send ', COUNT( mpinei(:) /= mpiSnei(jh,:) ), ' modif'
563            WRITE(inum,'(18x,8i6,a,i1,a)')   mpiRnei(jh,:), ' <- Recv ', COUNT( mpinei(:) /= mpiRnei(jh,:) ), ' modif'
564            WRITE(inum,*) ' total changes among all mpi tasks:'
565            WRITE(inum,*) '       mpiwe mpiea mpiso mpino mpisw mpise mpinw mpine'
566            WRITE(inum,'(a,8i6)') ' Send: ', ichanged(jh*16-15:jh*16-8)
567            WRITE(inum,'(a,8i6)') ' Recv: ', ichanged(jh*16 -7:jh*16  )
568         END DO
569      ENDIF
570      !
571      IF( ln_tspers ) CALL mpp_ini_pers   ! initialize persistent call
572      !
573      CALL init_ioipsl           ! Prepare NetCDF output file (if necessary)
574      !
575      IF (llwrtlay) CLOSE(inum)
576      !
577      DEALLOCATE(iin, ijn, iimppt, ijmppt, ijpi, ijpj, ipproc, inei, llnei, impi, llisOce)
578      !
579    END SUBROUTINE mpp_init
580
581#endif
582
583    SUBROUTINE mpp_basesplit( kiglo, kjglo, khls, knbi, knbj, kimax, kjmax, kimppt, kjmppt, klci, klcj)
584      !!----------------------------------------------------------------------
585      !!                  ***  ROUTINE mpp_basesplit  ***
586      !!
587      !! ** Purpose :   Lay out the global domain over processors.
588      !!
589      !! ** Method  :   Global domain is distributed in smaller local domains.
590      !!
591      !! ** Action : - set for all knbi*knbj domains:
592      !!                    kimppt     : longitudinal index
593      !!                    kjmppt     : latitudinal  index
594      !!                    klci       : first dimension
595      !!                    klcj       : second dimension
596      !!----------------------------------------------------------------------
597      INTEGER,                                 INTENT(in   ) ::   kiglo, kjglo
598      INTEGER,                                 INTENT(in   ) ::   khls
599      INTEGER,                                 INTENT(in   ) ::   knbi, knbj
600      INTEGER,                                 INTENT(  out) ::   kimax, kjmax
601      INTEGER, DIMENSION(knbi,knbj), OPTIONAL, INTENT(  out) ::   kimppt, kjmppt
602      INTEGER, DIMENSION(knbi,knbj), OPTIONAL, INTENT(  out) ::   klci, klcj
603      !
604      INTEGER ::   ji, jj
605      INTEGER ::   i2hls
606      INTEGER ::   iresti, irestj, irm, ijpjmin
607      !!----------------------------------------------------------------------
608      i2hls = 2*khls
609      !
610#if defined key_nemocice_decomp
611      kimax = ( nx_global+2-i2hls + (knbi-1) ) / knbi + i2hls    ! first  dim.
612      kjmax = ( ny_global+2-i2hls + (knbj-1) ) / knbj + i2hls    ! second dim.
613#else
614      kimax = ( kiglo - i2hls + (knbi-1) ) / knbi + i2hls    ! first  dim.
615      kjmax = ( kjglo - i2hls + (knbj-1) ) / knbj + i2hls    ! second dim.
616#endif
617      IF( .NOT. PRESENT(kimppt) ) RETURN
618      !
619      !  1. Dimension arrays for subdomains
620      ! -----------------------------------
621      !  Computation of local domain sizes klci() klcj()
622      !  These dimensions depend on global sizes knbi,knbj and kiglo,kjglo
623      !  The subdomains are squares lesser than or equal to the global
624      !  dimensions divided by the number of processors minus the overlap array.
625      !
626      iresti = 1 + MOD( kiglo - i2hls - 1 , knbi )
627      irestj = 1 + MOD( kjglo - i2hls - 1 , knbj )
628      !
629      !  Need to use kimax and kjmax here since jpi and jpj not yet defined
630#if defined key_nemocice_decomp
631      ! Change padding to be consistent with CICE
632      klci(1:knbi-1,:       ) = kimax
633      klci(  knbi  ,:       ) = kiglo - (knbi - 1) * (kimax - i2hls)
634      klcj(:       ,1:knbj-1) = kjmax
635      klcj(:       ,  knbj  ) = kjglo - (knbj - 1) * (kjmax - i2hls)
636#else
637      klci(1:iresti      ,:) = kimax
638      klci(iresti+1:knbi ,:) = kimax-1
639      IF( MINVAL(klci) < 2*i2hls ) THEN
640         WRITE(ctmp1,*) '   mpp_basesplit: minimum value of jpi must be >= ', 2*i2hls
641         WRITE(ctmp2,*) '   We have ', MINVAL(klci)
642        CALL ctl_stop( 'STOP', ctmp1, ctmp2 )
643      ENDIF
644      IF( l_NFold ) THEN
645         ! minimize the size of the last row to compensate for the north pole folding coast
646         IF( c_NFtype == 'T' )   ijpjmin = 2+3*khls   ! V and F folding must be outside of southern halos
647         IF( c_NFtype == 'F' )   ijpjmin = 1+3*khls   ! V and F folding must be outside of southern halos
648         irm = knbj - irestj                          ! total number of lines to be removed
649         klcj(:,knbj) = MAX( ijpjmin, kjmax-irm )     ! we must have jpj >= ijpjmin in the last row
650         irm = irm - ( kjmax - klcj(1,knbj) )         ! remaining number of lines to remove
651         irestj = knbj - 1 - irm
652         klcj(:, irestj+1:knbj-1) = kjmax-1
653      ELSE
654         klcj(:, irestj+1:knbj  ) = kjmax-1
655      ENDIF
656      klcj(:,1:irestj) = kjmax
657      IF( MINVAL(klcj) < 2*i2hls ) THEN
658         WRITE(ctmp1,*) '   mpp_basesplit: minimum value of jpj must be >= ', 2*i2hls
659         WRITE(ctmp2,*) '   We have ', MINVAL(klcj)
660         CALL ctl_stop( 'STOP', ctmp1, ctmp2 )
661      ENDIF
662#endif
663
664      !  2. Index arrays for subdomains
665      ! -------------------------------
666      kimppt(:,:) = 1
667      kjmppt(:,:) = 1
668      !
669      IF( knbi > 1 ) THEN
670         DO jj = 1, knbj
671            DO ji = 2, knbi
672               kimppt(ji,jj) = kimppt(ji-1,jj) + klci(ji-1,jj) - i2hls
673            END DO
674         END DO
675      ENDIF
676      !
677      IF( knbj > 1 )THEN
678         DO jj = 2, knbj
679            DO ji = 1, knbi
680               kjmppt(ji,jj) = kjmppt(ji,jj-1) + klcj(ji,jj-1) - i2hls
681            END DO
682         END DO
683      ENDIF
684
685   END SUBROUTINE mpp_basesplit
686
687
688   SUBROUTINE bestpartition( knbij, knbi, knbj, knbcnt, ldlist )
689      !!----------------------------------------------------------------------
690      !!                 ***  ROUTINE bestpartition  ***
691      !!
692      !! ** Purpose :
693      !!
694      !! ** Method  :
695      !!----------------------------------------------------------------------
696      INTEGER,           INTENT(in   ) ::   knbij         ! total number of subdomains (knbi*knbj)
697      INTEGER, OPTIONAL, INTENT(  out) ::   knbi, knbj    ! number if subdomains along i and j (knbi and knbj)
698      INTEGER, OPTIONAL, INTENT(  out) ::   knbcnt        ! number of land subdomains
699      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldlist        ! .true.: print the list the best domain decompositions (with land)
700      !
701      INTEGER :: ji, jj, ii, iitarget
702      INTEGER :: iszitst, iszjtst
703      INTEGER :: isziref, iszjref
704      INTEGER :: iszimin, iszjmin
705      INTEGER :: inbij, iszij
706      INTEGER :: inbimax, inbjmax, inbijmax, inbijold
707      INTEGER :: isz0, isz1
708      INTEGER, DIMENSION(  :), ALLOCATABLE :: indexok
709      INTEGER, DIMENSION(  :), ALLOCATABLE :: inbi0, inbj0, inbij0   ! number of subdomains along i,j
710      INTEGER, DIMENSION(  :), ALLOCATABLE :: iszi0, iszj0, iszij0   ! max size of the subdomains along i,j
711      INTEGER, DIMENSION(  :), ALLOCATABLE :: inbi1, inbj1, inbij1   ! number of subdomains along i,j
712      INTEGER, DIMENSION(  :), ALLOCATABLE :: iszi1, iszj1, iszij1   ! max size of the subdomains along i,j
713      LOGICAL :: llist
714      LOGICAL, DIMENSION(:,:), ALLOCATABLE :: llmsk2d                 ! max size of the subdomains along i,j
715      LOGICAL, DIMENSION(:,:), ALLOCATABLE :: llisOce              !  -     -
716      REAL(wp)::   zpropland
717      !!----------------------------------------------------------------------
718      !
719      llist = .FALSE.
720      IF( PRESENT(ldlist) ) llist = ldlist
721
722      CALL mpp_init_landprop( zpropland )                      ! get the proportion of land point over the gloal domain
723      inbij = NINT( REAL(knbij, wp) / ( 1.0 - zpropland ) )    ! define the largest possible value for jpni*jpnj
724      !
725      IF( llist ) THEN   ;   inbijmax = inbij*2
726      ELSE               ;   inbijmax = inbij
727      ENDIF
728      !
729      ALLOCATE(inbi0(inbijmax),inbj0(inbijmax),iszi0(inbijmax),iszj0(inbijmax))
730      !
731      inbimax = 0
732      inbjmax = 0
733      isziref = jpiglo*jpjglo+1   ! define a value that is larger than the largest possible
734      iszjref = jpiglo*jpjglo+1
735      !
736      iszimin = 4*nn_hls          ! minimum size of the MPI subdomain so halos are always adressing neighbor inner domain
737      iszjmin = 4*nn_hls
738      IF( c_NFtype == 'T' )   iszjmin = MAX(iszjmin, 2+3*nn_hls)   ! V and F folding must be outside of southern halos
739      IF( c_NFtype == 'F' )   iszjmin = MAX(iszjmin, 1+3*nn_hls)   ! V and F folding must be outside of southern halos
740      !
741      ! get the list of knbi that gives a smaller jpimax than knbi-1
742      ! get the list of knbj that gives a smaller jpjmax than knbj-1
743      DO ji = 1, inbijmax
744#if defined key_nemocice_decomp
745         iszitst = ( nx_global+2-2*nn_hls + (ji-1) ) / ji + 2*nn_hls    ! first  dim.
746#else
747         iszitst = ( Ni0glo + (ji-1) ) / ji + 2*nn_hls   ! max subdomain i-size
748#endif
749         IF( iszitst < isziref .AND. iszitst >= iszimin ) THEN
750            isziref = iszitst
751            inbimax = inbimax + 1
752            inbi0(inbimax) = ji
753            iszi0(inbimax) = isziref
754         ENDIF
755#if defined key_nemocice_decomp
756         iszjtst = ( ny_global+2-2*nn_hls + (ji-1) ) / ji + 2*nn_hls    ! first  dim.
757#else
758         iszjtst = ( Nj0glo + (ji-1) ) / ji + 2*nn_hls   ! max subdomain j-size
759#endif
760         IF( iszjtst < iszjref .AND. iszjtst >= iszjmin ) THEN
761            iszjref = iszjtst
762            inbjmax = inbjmax + 1
763            inbj0(inbjmax) = ji
764            iszj0(inbjmax) = iszjref
765         ENDIF
766      END DO
767
768      ! combine these 2 lists to get all possible knbi*knbj <  inbijmax
769      ALLOCATE( llmsk2d(inbimax,inbjmax) )
770      DO jj = 1, inbjmax
771         DO ji = 1, inbimax
772            IF ( inbi0(ji) * inbj0(jj) <= inbijmax ) THEN   ;   llmsk2d(ji,jj) = .TRUE.
773            ELSE                                            ;   llmsk2d(ji,jj) = .FALSE.
774            ENDIF
775         END DO
776      END DO
777      isz1 = COUNT(llmsk2d)
778      ALLOCATE( inbi1(isz1), inbj1(isz1), iszi1(isz1), iszj1(isz1) )
779      ii = 0
780      DO jj = 1, inbjmax
781         DO ji = 1, inbimax
782            IF( llmsk2d(ji,jj) .EQV. .TRUE. ) THEN
783               ii = ii + 1
784               inbi1(ii) = inbi0(ji)
785               inbj1(ii) = inbj0(jj)
786               iszi1(ii) = iszi0(ji)
787               iszj1(ii) = iszj0(jj)
788            ENDIF
789         END DO
790      END DO
791      DEALLOCATE( inbi0, inbj0, iszi0, iszj0 )
792      DEALLOCATE( llmsk2d )
793
794      ALLOCATE( inbij1(isz1), iszij1(isz1) )
795      inbij1(:) = inbi1(:) * inbj1(:)
796      iszij1(:) = iszi1(:) * iszj1(:)
797
798      ! if there is no land and no print
799      IF( .NOT. llist .AND. numbot == -1 .AND. numbdy == -1 ) THEN
800         ! get the smaller partition which gives the smallest subdomain size
801         ii = MINLOC(inbij1, mask = iszij1 == MINVAL(iszij1), dim = 1)
802         knbi = inbi1(ii)
803         knbj = inbj1(ii)
804         IF(PRESENT(knbcnt))   knbcnt = 0
805         DEALLOCATE( inbi1, inbj1, inbij1, iszi1, iszj1, iszij1 )
806         RETURN
807      ENDIF
808
809      ! extract only the partitions which reduce the subdomain size in comparison with smaller partitions
810      ALLOCATE( indexok(isz1) )                                 ! to store indices of the best partitions
811      isz0 = 0                                                  ! number of best partitions
812      inbij = 1                                                 ! start with the min value of inbij1 => 1
813      iszij = jpiglo*jpjglo+1                                   ! default: larger than global domain
814      DO WHILE( inbij <= inbijmax )                             ! if we did not reach the max of inbij1
815         ii = MINLOC(iszij1, mask = inbij1 == inbij, dim = 1)   ! warning: send back the first occurence if multiple results
816         IF ( iszij1(ii) < iszij ) THEN
817            ii = MINLOC( iszi1+iszj1, mask = iszij1 == iszij1(ii) .AND. inbij1 == inbij, dim = 1)  ! select the smaller perimeter if multiple min
818            isz0 = isz0 + 1
819            indexok(isz0) = ii
820            iszij = iszij1(ii)
821         ENDIF
822         inbij = MINVAL(inbij1, mask = inbij1 > inbij)   ! warning: return largest integer value if mask = .false. everywhere
823      END DO
824      DEALLOCATE( inbij1, iszij1 )
825
826      ! keep only the best partitions (sorted by increasing order of subdomains number and decreassing subdomain size)
827      ALLOCATE( inbi0(isz0), inbj0(isz0), iszi0(isz0), iszj0(isz0) )
828      DO ji = 1, isz0
829         ii = indexok(ji)
830         inbi0(ji) = inbi1(ii)
831         inbj0(ji) = inbj1(ii)
832         iszi0(ji) = iszi1(ii)
833         iszj0(ji) = iszj1(ii)
834      END DO
835      DEALLOCATE( indexok, inbi1, inbj1, iszi1, iszj1 )
836
837      IF( llist ) THEN
838         IF(lwp) THEN
839            WRITE(numout,*)
840            WRITE(numout,*) '                  For your information:'
841            WRITE(numout,*) '  list of the best partitions including land supression'
842            WRITE(numout,*) '  -----------------------------------------------------'
843            WRITE(numout,*)
844         ENDIF
845         ji = isz0   ! initialization with the largest value
846         ALLOCATE( llisOce(inbi0(ji), inbj0(ji)) )
847         CALL mpp_is_ocean( llisOce )   ! Warning: must be call by all cores (call mpp_sum)
848         inbijold = COUNT(llisOce)
849         DEALLOCATE( llisOce )
850         DO ji =isz0-1,1,-1
851            ALLOCATE( llisOce(inbi0(ji), inbj0(ji)) )
852            CALL mpp_is_ocean( llisOce )   ! Warning: must be call by all cores (call mpp_sum)
853            inbij = COUNT(llisOce)
854            DEALLOCATE( llisOce )
855            IF(lwp .AND. inbij < inbijold) THEN
856               WRITE(numout,'(a, i6, a, i6, a, f4.1, a, i9, a, i6, a, i6, a)')                                 &
857                  &   'nb_cores oce: ', inbij, ', land domains excluded: ', inbi0(ji)*inbj0(ji) - inbij,       &
858                  &   ' (', REAL(inbi0(ji)*inbj0(ji) - inbij,wp) / REAL(inbi0(ji)*inbj0(ji),wp) *100.,         &
859                  &   '%), largest oce domain: ', iszi0(ji)*iszj0(ji), ' ( ', iszi0(ji),' x ', iszj0(ji), ' )'
860               inbijold = inbij
861            ENDIF
862         END DO
863         DEALLOCATE( inbi0, inbj0, iszi0, iszj0 )
864         IF(lwp) THEN
865            WRITE(numout,*)
866            WRITE(numout,*)  '  -----------------------------------------------------------'
867         ENDIF
868         CALL mppsync
869         CALL mppstop( ld_abort = .TRUE. )
870      ENDIF
871
872      DEALLOCATE( iszi0, iszj0 )
873      inbij = inbijmax + 1        ! default: larger than possible
874      ii = isz0+1                 ! start from the end of the list (smaller subdomains)
875      DO WHILE( inbij > knbij )   ! while the number of ocean subdomains exceed the number of procs
876         ii = ii -1
877         ALLOCATE( llisOce(inbi0(ii), inbj0(ii)) )
878         CALL mpp_is_ocean( llisOce )            ! must be done by all core
879         inbij = COUNT(llisOce)
880         DEALLOCATE( llisOce )
881      END DO
882      knbi = inbi0(ii)
883      knbj = inbj0(ii)
884      IF(PRESENT(knbcnt))   knbcnt = knbi * knbj - inbij
885      DEALLOCATE( inbi0, inbj0 )
886      !
887   END SUBROUTINE bestpartition
888
889
890   SUBROUTINE mpp_init_landprop( propland )
891      !!----------------------------------------------------------------------
892      !!                  ***  ROUTINE mpp_init_landprop  ***
893      !!
894      !! ** Purpose : the the proportion of land points in the surface land-sea mask
895      !!
896      !! ** Method  : read iproc strips (of length Ni0glo) of the land-sea mask
897      !!----------------------------------------------------------------------
898      REAL(wp), INTENT(  out) :: propland    ! proportion of land points in the global domain (between 0 and 1)
899      !
900      INTEGER, DIMENSION(jpni*jpnj) ::   kusedom_1d
901      INTEGER :: inboce, iarea
902      INTEGER :: iproc, idiv, ijsz
903      INTEGER :: ijstr
904      LOGICAL, ALLOCATABLE, DIMENSION(:,:) ::   lloce
905      !!----------------------------------------------------------------------
906      ! do nothing if there is no land-sea mask
907      IF( numbot == -1 .and. numbdy == -1 ) THEN
908         propland = 0.
909         RETURN
910      ENDIF
911
912      ! number of processes reading the bathymetry file
913      iproc = MINVAL( (/mppsize, Nj0glo/2, 100/) )  ! read a least 2 lines, no more that 100 processes reading at the same time
914
915      ! we want to read iproc strips of the land-sea mask. -> pick up iproc processes every idiv processes starting at 1
916      IF( iproc == 1 ) THEN   ;   idiv = mppsize
917      ELSE                    ;   idiv = ( mppsize - 1 ) / ( iproc - 1 )
918      ENDIF
919
920      iarea = (narea-1)/idiv   ! involed process number (starting counting at 0)
921      IF( MOD( narea-1, idiv ) == 0 .AND. iarea < iproc ) THEN   ! beware idiv can be = to 1
922         !
923         ijsz = Nj0glo / iproc                                               ! width of the stripe to read
924         IF( iarea < MOD(Nj0glo,iproc) ) ijsz = ijsz + 1
925         ijstr = iarea*(Nj0glo/iproc) + MIN(iarea, MOD(Nj0glo,iproc)) + 1    ! starting j position of the reading
926         !
927         ALLOCATE( lloce(Ni0glo, ijsz) )                                     ! allocate the strip
928         CALL read_mask( 1, ijstr, Ni0glo, ijsz, lloce )
929         inboce = COUNT(lloce)                                               ! number of ocean point in the stripe
930         DEALLOCATE(lloce)
931         !
932      ELSE
933         inboce = 0
934      ENDIF
935      CALL mpp_sum( 'mppini', inboce )   ! total number of ocean points over the global domain
936      !
937      propland = REAL( Ni0glo*Nj0glo - inboce, wp ) / REAL( Ni0glo*Nj0glo, wp )
938      !
939   END SUBROUTINE mpp_init_landprop
940
941
942   SUBROUTINE mpp_is_ocean( ldIsOce )
943      !!----------------------------------------------------------------------
944      !!                  ***  ROUTINE mpp_is_ocean  ***
945      !!
946      !! ** Purpose : Check for a mpi domain decomposition inbi x inbj which
947      !!              subdomains, including 1 halo (even if nn_hls>1), contain
948      !!              at least 1 ocean point.
949      !!              We must indeed ensure that each subdomain that is a neighbour
950      !!              of a land subdomain, has only land points on its boundary
951      !!              (inside the inner subdomain) with the land subdomain.
952      !!              This is needed to get the proper bondary conditions on
953      !!              a subdomain with a closed boundary.
954      !!
955      !! ** Method  : read inbj strips (of length Ni0glo) of the land-sea mask
956      !!----------------------------------------------------------------------
957      LOGICAL, DIMENSION(:,:), INTENT(  out) ::   ldIsOce        ! .true. if a sub domain constains 1 ocean point
958      !
959      INTEGER :: idiv, iimax, ijmax, iarea
960      INTEGER :: inbi, inbj, inx, iny, inry, isty
961      INTEGER :: ji, jn
962      INTEGER, ALLOCATABLE, DIMENSION(:,:) ::   inboce           ! number oce oce pint in each mpi subdomain
963      INTEGER, ALLOCATABLE, DIMENSION(:  ) ::   inboce_1d
964      INTEGER, ALLOCATABLE, DIMENSION(:,:) ::   iimppt, ijpi
965      INTEGER, ALLOCATABLE, DIMENSION(:,:) ::   ijmppt, ijpj
966      LOGICAL, ALLOCATABLE, DIMENSION(:,:) ::   lloce            ! lloce(i,j) = .true. if the point (i,j) is ocean
967      !!----------------------------------------------------------------------
968      ! do nothing if there is no land-sea mask
969      IF( numbot == -1 .AND. numbdy == -1 ) THEN
970         ldIsOce(:,:) = .TRUE.
971         RETURN
972      ENDIF
973      !
974      inbi = SIZE( ldIsOce, dim = 1 )
975      inbj = SIZE( ldIsOce, dim = 2 )
976      !
977      ! we want to read inbj strips of the land-sea mask. -> pick up inbj processes every idiv processes starting at 1
978      IF           ( inbj == 1 ) THEN   ;   idiv = mppsize
979      ELSE IF ( mppsize < inbj ) THEN   ;   idiv = 1
980      ELSE                              ;   idiv = ( mppsize - 1 ) / ( inbj - 1 )
981      ENDIF
982      !
983      ALLOCATE( inboce(inbi,inbj), inboce_1d(inbi*inbj) )
984      inboce(:,:) = 0          ! default no ocean point found
985      !
986      DO jn = 0, (inbj-1)/mppsize   ! if mppsize < inbj : more strips than mpi processes (because of potential land domains)
987         !
988         iarea = (narea-1)/idiv + jn * mppsize + 1                     ! involed process number (starting counting at 1)
989         IF( MOD( narea-1, idiv ) == 0 .AND. iarea <= inbj ) THEN      ! beware idiv can be = to 1
990            !
991            ALLOCATE( iimppt(inbi,inbj), ijmppt(inbi,inbj), ijpi(inbi,inbj), ijpj(inbi,inbj) )
992            CALL mpp_basesplit( Ni0glo, Nj0glo, 0, inbi, inbj, iimax, ijmax, iimppt, ijmppt, ijpi, ijpj )
993            !
994            inx = Ni0glo + 2   ;   iny = ijpj(1,iarea) + 2             ! strip size + 1 halo on each direction (even if nn_hls>1)
995            ALLOCATE( lloce(inx, iny) )                                ! allocate the strip
996            inry = iny - COUNT( (/ iarea == 1, iarea == inbj /) )      ! number of point to read in y-direction
997            isty = 1 + COUNT( (/ iarea == 1 /) )                       ! read from the first or the second line?
998            CALL read_mask( 1, ijmppt(1,iarea) - 2 + isty, Ni0glo, inry, lloce(2:inx-1, isty:inry+isty-1) )   ! read the strip
999            !
1000            IF( iarea == 1    ) THEN                                   ! the first line was not read
1001               IF( l_Jperio ) THEN                                     !   north-south periodocity
1002                  CALL read_mask( 1, Nj0glo, Ni0glo, 1, lloce(2:inx-1, 1) )   !   read the last line -> first line of lloce
1003               ELSE
1004                  lloce(2:inx-1,  1) = .FALSE.                         !   closed boundary
1005               ENDIF
1006            ENDIF
1007            IF( iarea == inbj ) THEN                                   ! the last line was not read
1008               IF( l_Jperio ) THEN                                     !   north-south periodocity
1009                  CALL read_mask( 1, 1, Ni0glo, 1, lloce(2:inx-1,iny) )   !      read the first line -> last line of lloce
1010               ELSEIF( c_NFtype == 'T' ) THEN                          !   north-pole folding T-pivot, T-point
1011                  lloce(2,iny) = lloce(2,iny-2)                        !      here we have 1 halo (even if nn_hls>1)
1012                  DO ji = 3,inx-1
1013                     lloce(ji,iny  ) = lloce(inx-ji+2,iny-2)           !      ok, we have at least 3 lines
1014                  END DO
1015                  DO ji = inx/2+2,inx-1
1016                     lloce(ji,iny-1) = lloce(inx-ji+2,iny-1)
1017                  END DO
1018               ELSEIF( c_NFtype == 'F' ) THEN                          !   north-pole folding F-pivot, T-point, 1 halo
1019                  lloce(inx/2+1,iny-1) = lloce(inx/2,iny-1)            !      here we have 1 halo (even if nn_hls>1)
1020                  lloce(inx  -1,iny-1) = lloce(2    ,iny-1)
1021                  DO ji = 2,inx-1
1022                     lloce(ji,iny) = lloce(inx-ji+1,iny-1)
1023                  END DO
1024               ELSE                                                    !   closed boundary
1025                  lloce(2:inx-1,iny) = .FALSE.
1026               ENDIF
1027            ENDIF
1028            !                                                          ! first and last column were not read
1029            IF( l_Iperio ) THEN
1030               lloce(1,:) = lloce(inx-1,:)   ;   lloce(inx,:) = lloce(2,:)   ! east-west periodocity
1031            ELSE
1032               lloce(1,:) = .FALSE.          ;   lloce(inx,:) = .FALSE.      ! closed boundary
1033            ENDIF
1034            !
1035            DO  ji = 1, inbi
1036               inboce(ji,iarea) = COUNT( lloce(iimppt(ji,1):iimppt(ji,1)+ijpi(ji,1)+1,:) )   ! lloce as 2 points more than Ni0glo
1037            END DO
1038            !
1039            DEALLOCATE(lloce)
1040            DEALLOCATE(iimppt, ijmppt, ijpi, ijpj)
1041            !
1042         ENDIF
1043      END DO
1044
1045      inboce_1d = RESHAPE(inboce, (/ inbi*inbj /))
1046      CALL mpp_sum( 'mppini', inboce_1d )
1047      inboce = RESHAPE(inboce_1d, (/inbi, inbj/))
1048      ldIsOce(:,:) = inboce(:,:) /= 0
1049      DEALLOCATE(inboce, inboce_1d)
1050      !
1051   END SUBROUTINE mpp_is_ocean
1052
1053
1054   SUBROUTINE read_mask( kistr, kjstr, kicnt, kjcnt, ldoce )
1055      !!----------------------------------------------------------------------
1056      !!                  ***  ROUTINE read_mask  ***
1057      !!
1058      !! ** Purpose : Read relevant bathymetric information in order to
1059      !!              provide a land/sea mask used for the elimination
1060      !!              of land domains, in an mpp computation.
1061      !!
1062      !! ** Method  : read stipe of size (Ni0glo,...)
1063      !!----------------------------------------------------------------------
1064      INTEGER                        , INTENT(in   ) ::   kistr, kjstr   ! starting i and j position of the reading
1065      INTEGER                        , INTENT(in   ) ::   kicnt, kjcnt   ! number of points to read in i and j directions
1066      LOGICAL, DIMENSION(kicnt,kjcnt), INTENT(  out) ::   ldoce          ! ldoce(i,j) = .true. if the point (i,j) is ocean
1067      !
1068      INTEGER                          ::   inumsave                     ! local logical unit
1069      REAL(wp), DIMENSION(kicnt,kjcnt) ::   zbot, zbdy
1070      !!----------------------------------------------------------------------
1071      !
1072      inumsave = numout   ;   numout = numnul   !   redirect all print to /dev/null
1073      !
1074      IF( numbot /= -1 ) THEN
1075         CALL iom_get( numbot, jpdom_unknown, 'bottom_level', zbot, kstart = (/kistr,kjstr/), kcount = (/kicnt, kjcnt/) )
1076      ELSE
1077         zbot(:,:) = 1._wp                      ! put a non-null value
1078      ENDIF
1079      !
1080      IF( numbdy /= -1 ) THEN                   ! Adjust with bdy_msk if it exists
1081         CALL iom_get ( numbdy, jpdom_unknown,     'bdy_msk', zbdy, kstart = (/kistr,kjstr/), kcount = (/kicnt, kjcnt/) )
1082         zbot(:,:) = zbot(:,:) * zbdy(:,:)
1083      ENDIF
1084      !
1085      ldoce(:,:) = NINT(zbot(:,:)) > 0
1086      numout = inumsave
1087      !
1088   END SUBROUTINE read_mask
1089
1090
1091   SUBROUTINE mpp_getnum( ldIsOce, kproc, kipos, kjpos )
1092      !!----------------------------------------------------------------------
1093      !!                  ***  ROUTINE mpp_getnum  ***
1094      !!
1095      !! ** Purpose : give a number to each MPI subdomains (starting at 0)
1096      !!
1097      !! ** Method  : start from bottom left. First skip land subdomain, and finally use them if needed
1098      !!----------------------------------------------------------------------
1099      LOGICAL, DIMENSION(:,:), INTENT(in   ) ::   ldIsOce     ! F if land process
1100      INTEGER, DIMENSION(:,:), INTENT(  out) ::   kproc       ! subdomain number (-1 if not existing, starting at 0)
1101      INTEGER, DIMENSION(  :), INTENT(  out) ::   kipos       ! i-position of the subdomain (from 1 to jpni)
1102      INTEGER, DIMENSION(  :), INTENT(  out) ::   kjpos       ! j-position of the subdomain (from 1 to jpnj)
1103      !
1104      INTEGER :: ii, ij, jarea, iarea0
1105      INTEGER :: icont, i2add , ini, inj, inij
1106      !!----------------------------------------------------------------------
1107      !
1108      ini = SIZE(ldIsOce, dim = 1)
1109      inj = SIZE(ldIsOce, dim = 2)
1110      inij = SIZE(kipos)
1111      !
1112      ! specify which subdomains are oce subdomains; other are land subdomains
1113      kproc(:,:) = -1
1114      icont = -1
1115      DO jarea = 1, ini*inj
1116         iarea0 = jarea - 1
1117         ii = 1 + MOD(iarea0,ini)
1118         ij = 1 +     iarea0/ini
1119         IF( ldIsOce(ii,ij) ) THEN
1120            icont = icont + 1
1121            kproc(ii,ij) = icont
1122            kipos(icont+1) = ii
1123            kjpos(icont+1) = ij
1124         ENDIF
1125      END DO
1126      ! if needed add some land subdomains to reach inij active subdomains
1127      i2add = inij - COUNT( ldIsOce )
1128      DO jarea = 1, ini*inj
1129         iarea0 = jarea - 1
1130         ii = 1 + MOD(iarea0,ini)
1131         ij = 1 +     iarea0/ini
1132         IF( .NOT. ldIsOce(ii,ij) .AND. i2add > 0 ) THEN
1133            icont = icont + 1
1134            kproc(ii,ij) = icont
1135            kipos(icont+1) = ii
1136            kjpos(icont+1) = ij
1137            i2add = i2add - 1
1138         ENDIF
1139      END DO
1140      !
1141   END SUBROUTINE mpp_getnum
1142
1143
1144   SUBROUTINE init_excl_landpt
1145      !!----------------------------------------------------------------------
1146      !!                  ***  ROUTINE   ***
1147      !!
1148      !! ** Purpose : exclude exchanges which contain only land points
1149      !!
1150      !! ** Method  : if a send or receive buffer constains only land point we
1151      !!              flag off the corresponding communication
1152      !!              Warning: this selection depend on the halo size -> loop on halo size
1153      !!
1154      !!----------------------------------------------------------------------
1155      INTEGER ::   inumsave
1156      INTEGER ::   jh
1157      INTEGER ::   ipi, ipj
1158      INTEGER ::   iiwe, iiea, iist, iisz 
1159      INTEGER ::   ijso, ijno, ijst, ijsz 
1160      LOGICAL ::   llsave
1161      REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   zmsk
1162      LOGICAL , DIMENSION(Ni_0,Nj_0,1)      ::   lloce
1163      !!----------------------------------------------------------------------
1164      !
1165      ! read the land-sea mask on the inner domain
1166      CALL read_mask( nimpp, njmpp, Ni_0, Nj_0, lloce(:,:,1) )
1167      !
1168      ! Here we look only at communications excluding the NP folding.
1169      ! As lbcnfd not validated if halo size /= nn_hls, we switch if off temporary...
1170      llsave = l_IdoNFold
1171      l_IdoNFold = .FALSE.
1172      !
1173      DO jh = 1, n_hlsmax    ! different halo size
1174         !
1175         ipi = Ni_0 + 2*jh   ! local domain size
1176         ipj = Nj_0 + 2*jh
1177         !
1178         ALLOCATE( zmsk(ipi,ipj) )
1179         zmsk(jh+1:jh+Ni_0,jh+1:jh+Nj_0) = REAL(COUNT(lloce, dim = 3), wp)   ! define inner domain -> need REAL to use lbclnk
1180         CALL lbc_lnk('mppini', zmsk, 'T', 1._wp, khls = jh)                 ! fill halos
1181         !       
1182         iiwe = jh   ;   iiea = Ni_0   ! bottom-left corfer - 1 of the sent data
1183         ijso = jh   ;   ijno = Nj_0
1184         IF( nn_comm == 1 ) THEN
1185            iist =  0   ;   iisz = ipi
1186            ijst =  0   ;   ijsz = ipj
1187         ELSE
1188            iist = jh   ;   iisz = Ni_0
1189            ijst = jh   ;   ijsz = Nj_0
1190         ENDIF
1191IF( nn_comm == 1 ) THEN       ! SM: NOT WORKING FOR NEIGHBOURHOOD COLLECTIVE COMMUNICATIONS, I DON'T KNOW WHY...
1192         ! do not send if we send only land points
1193         IF( NINT(SUM( zmsk(iiwe+1:iiwe+jh  ,ijst+1:ijst+ijsz) )) == 0 )   mpiSnei(jh,jpwe) = -1
1194         IF( NINT(SUM( zmsk(iiea+1:iiea+jh  ,ijst+1:ijst+ijsz) )) == 0 )   mpiSnei(jh,jpea) = -1
1195         IF( NINT(SUM( zmsk(iist+1:iist+iisz,ijso+1:ijso+jh  ) )) == 0 )   mpiSnei(jh,jpso) = -1
1196         IF( NINT(SUM( zmsk(iist+1:iist+iisz,ijno+1:ijno+jh  ) )) == 0 )   mpiSnei(jh,jpno) = -1
1197         IF( NINT(SUM( zmsk(iiwe+1:iiwe+jh  ,ijso+1:ijso+jh  ) )) == 0 )   mpiSnei(jh,jpsw) = -1
1198         IF( NINT(SUM( zmsk(iiea+1:iiea+jh  ,ijso+1:ijso+jh  ) )) == 0 )   mpiSnei(jh,jpse) = -1
1199         IF( NINT(SUM( zmsk(iiwe+1:iiwe+jh  ,ijno+1:ijno+jh  ) )) == 0 )   mpiSnei(jh,jpnw) = -1
1200         IF( NINT(SUM( zmsk(iiea+1:iiea+jh  ,ijno+1:ijno+jh  ) )) == 0 )   mpiSnei(jh,jpne) = -1
1201         !
1202         iiwe = iiwe-jh   ;   iiea = iiea+jh   ! bottom-left corfer - 1 of the received data
1203         ijso = ijso-jh   ;   ijno = ijno+jh
1204         ! do not send if we send only land points
1205         IF( NINT(SUM( zmsk(iiwe+1:iiwe+jh  ,ijst+1:ijst+ijsz) )) == 0 )   mpiRnei(jh,jpwe) = -1
1206         IF( NINT(SUM( zmsk(iiea+1:iiea+jh  ,ijst+1:ijst+ijsz) )) == 0 )   mpiRnei(jh,jpea) = -1
1207         IF( NINT(SUM( zmsk(iist+1:iist+iisz,ijso+1:ijso+jh  ) )) == 0 )   mpiRnei(jh,jpso) = -1
1208         IF( NINT(SUM( zmsk(iist+1:iist+iisz,ijno+1:ijno+jh  ) )) == 0 )   mpiRnei(jh,jpno) = -1
1209         IF( NINT(SUM( zmsk(iiwe+1:iiwe+jh  ,ijso+1:ijso+jh  ) )) == 0 )   mpiRnei(jh,jpsw) = -1
1210         IF( NINT(SUM( zmsk(iiea+1:iiea+jh  ,ijso+1:ijso+jh  ) )) == 0 )   mpiRnei(jh,jpse) = -1
1211         IF( NINT(SUM( zmsk(iiwe+1:iiwe+jh  ,ijno+1:ijno+jh  ) )) == 0 )   mpiRnei(jh,jpnw) = -1
1212         IF( NINT(SUM( zmsk(iiea+1:iiea+jh  ,ijno+1:ijno+jh  ) )) == 0 )   mpiRnei(jh,jpne) = -1
1213ENDIF
1214         !
1215         ! Specific (and rare) problem in corner treatment because we do 1st West-East comm, next South-North comm
1216         IF( nn_comm == 1 ) THEN
1217            IF( mpiSnei(jh,jpwe) > -1 )   mpiSnei(jh, (/jpsw,jpnw/) ) = -1   ! SW and NW corners already sent through West nei
1218            IF( mpiSnei(jh,jpea) > -1 )   mpiSnei(jh, (/jpse,jpne/) ) = -1   ! SE and NE corners already sent through East nei
1219            IF( mpiRnei(jh,jpso) > -1 )   mpiRnei(jh, (/jpsw,jpse/) ) = -1   ! SW and SE corners will be received through South nei
1220            IF( mpiRnei(jh,jpno) > -1 )   mpiRnei(jh, (/jpnw,jpne/) ) = -1   ! NW and NE corners will be received through North nei
1221        ENDIF
1222         !
1223         DEALLOCATE( zmsk )
1224         !
1225         CALL mpp_ini_nc(jh)    ! Initialize/Update communicator for neighbourhood collective communications
1226         !
1227      END DO
1228      l_IdoNFold = llsave
1229
1230   END SUBROUTINE init_excl_landpt
1231
1232
1233   SUBROUTINE init_ioipsl
1234      !!----------------------------------------------------------------------
1235      !!                  ***  ROUTINE init_ioipsl  ***
1236      !!
1237      !! ** Purpose :
1238      !!
1239      !! ** Method  :
1240      !!
1241      !! History :
1242      !!   9.0  !  04-03  (G. Madec )  MPP-IOIPSL
1243      !!   " "  !  08-12  (A. Coward)  addition in case of jpni*jpnj < jpnij
1244      !!----------------------------------------------------------------------
1245      INTEGER, DIMENSION(2) ::   iglo, iloc, iabsf, iabsl, ihals, ihale, idid
1246      !!----------------------------------------------------------------------
1247
1248      ! The domain is split only horizontally along i- or/and j- direction
1249      ! So we need at the most only 1D arrays with 2 elements.
1250      ! Set idompar values equivalent to the jpdom_local_noextra definition
1251      ! used in IOM. This works even if jpnij .ne. jpni*jpnj.
1252      iglo( :) = (/ Ni0glo, Nj0glo /)
1253      iloc( :) = (/ Ni_0  , Nj_0   /)
1254      iabsf(:) = (/ Nis0  , Njs0   /) + (/ nimpp, njmpp /) - 1 - nn_hls   ! corresponds to mig0(Nis0) but mig0 is not yet defined!
1255      iabsl(:) = iabsf(:) + iloc(:) - 1
1256      ihals(:) = (/ 0     , 0      /)
1257      ihale(:) = (/ 0     , 0      /)
1258      idid( :) = (/ 1     , 2      /)
1259
1260      IF(lwp) THEN
1261          WRITE(numout,*)
1262          WRITE(numout,*) 'mpp init_ioipsl :   iloc  = ', iloc
1263          WRITE(numout,*) '~~~~~~~~~~~~~~~     iabsf = ', iabsf
1264          WRITE(numout,*) '                    ihals = ', ihals
1265          WRITE(numout,*) '                    ihale = ', ihale
1266      ENDIF
1267      !
1268      CALL flio_dom_set ( jpnij, narea-1, idid, iglo, iloc, iabsf, iabsl, ihals, ihale, 'BOX', nidom)
1269      !
1270   END SUBROUTINE init_ioipsl
1271
1272
1273   SUBROUTINE init_nfdcom
1274      !!----------------------------------------------------------------------
1275      !!                     ***  ROUTINE  init_nfdcom  ***
1276      !! ** Purpose :   Setup for north fold exchanges with explicit
1277      !!                point-to-point messaging
1278      !!
1279      !! ** Method :   Initialization of the northern neighbours lists.
1280      !!----------------------------------------------------------------------
1281      !!    1.0  ! 2011-10  (A. C. Coward, NOCS & J. Donners, PRACE)
1282      !!    2.0  ! 2013-06 Setup avoiding MPI communication (I. Epicoco, S. Mocavero, CMCC)
1283      !!----------------------------------------------------------------------
1284      INTEGER  ::   sxM, dxM, sxT, dxT, jn
1285      !!----------------------------------------------------------------------
1286      !
1287      !sxM is the first point (in the global domain) needed to compute the north-fold for the current process
1288      sxM = jpiglo - nimpp - jpi + 1
1289      !dxM is the last point (in the global domain) needed to compute the north-fold for the current process
1290      dxM = jpiglo - nimpp + 2
1291      !
1292      ! loop over the other north-fold processes to find the processes
1293      ! managing the points belonging to the sxT-dxT range
1294      !
1295      nsndto = 0
1296      DO jn = 1, jpni
1297         !
1298         sxT = nfimpp(jn)                    ! sxT = 1st  point (in the global domain) of the jn process
1299         dxT = nfimpp(jn) + nfjpi(jn) - 1    ! dxT = last point (in the global domain) of the jn process
1300         !
1301         IF    ( sxT < sxM  .AND.  sxM < dxT ) THEN
1302            nsndto          = nsndto + 1
1303            isendto(nsndto) = jn
1304         ELSEIF( sxM <= sxT  .AND.  dxM >= dxT ) THEN
1305            nsndto          = nsndto + 1
1306            isendto(nsndto) = jn
1307         ELSEIF( dxM <  dxT  .AND.  sxT <  dxM ) THEN
1308            nsndto          = nsndto + 1
1309            isendto(nsndto) = jn
1310         ENDIF
1311         !
1312      END DO
1313      !
1314   END SUBROUTINE init_nfdcom
1315
1316
1317   SUBROUTINE init_doloop
1318      !!----------------------------------------------------------------------
1319      !!                  ***  ROUTINE init_doloop  ***
1320      !!
1321      !! ** Purpose :   set the starting/ending indices of DO-loop
1322      !!              These indices are used in do_loop_substitute.h90
1323      !!----------------------------------------------------------------------
1324      !
1325      Nis0 =   1+nn_hls
1326      Njs0 =   1+nn_hls
1327      Nie0 = jpi-nn_hls
1328      Nje0 = jpj-nn_hls
1329      !
1330      Ni_0 = Nie0 - Nis0 + 1
1331      Nj_0 = Nje0 - Njs0 + 1
1332      !
1333      ! old indices to be removed...
1334      jpim1 = jpi-1                             ! inner domain indices
1335      jpjm1 = jpj-1                             !   "           "
1336      jpkm1 = jpk-1                             !   "           "
1337      !
1338   END SUBROUTINE init_doloop
1339
1340   !!======================================================================
1341END MODULE mppini
Note: See TracBrowser for help on using the repository browser.