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/releases/r4.0/r4.0-HEAD/src/OCE/LBC – NEMO

source: NEMO/releases/r4.0/r4.0-HEAD/src/OCE/LBC/mppini.F90 @ 12737

Last change on this file since 12737 was 12737, checked in by jchanut, 4 years ago

Fixes AGRIF reproductibility with land processors removal, i.e. #2240. Trunk is not concerned by this problem since nbondi/nbondj variables are not used anymore.

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