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/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/LBC – NEMO

source: NEMO/branches/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/LBC/mppini.F90 @ 13236

Last change on this file since 13236 was 13236, checked in by smasson, 4 years ago

dev_r12558_HPC-08_epico_Extra_Halo: fix merge with trunk@13218, see #2366

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