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.
Changeset 15267 – NEMO

Changeset 15267


Ignore:
Timestamp:
2021-09-17T11:04:34+02:00 (3 years ago)
Author:
smasson
Message:

trunk: new nogather nolding, #2724

Location:
NEMO/trunk
Files:
1 deleted
14 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/DOM/dom_oce.F90

    r15062 r15267  
    274274      ! 
    275275      ii = ii+1 
    276       ALLOCATE( mig(jpi), mjg(jpj), mig0(jpi), mjg0(jpj), STAT=ierr(ii) ) 
    277          ! 
    278       ii = ii+1 
    279       ALLOCATE( mi0(jpiglo), mi1(jpiglo), mj0(jpjglo), mj1(jpjglo), STAT=ierr(ii) ) 
    280          ! 
    281       ii = ii+1 
    282276      ALLOCATE( glamt(jpi,jpj) ,    glamu(jpi,jpj) ,  glamv(jpi,jpj) ,  glamf(jpi,jpj) ,     & 
    283277         &      gphit(jpi,jpj) ,    gphiu(jpi,jpj) ,  gphiv(jpi,jpj) ,  gphif(jpi,jpj) ,     & 
  • NEMO/trunk/src/OCE/DOM/domain.F90

    r15023 r15267  
    124124      !           !==  Reference coordinate system  ==! 
    125125      ! 
    126       CALL dom_glo                      ! global domain versus local domain 
     126!      CALL dom_glo                      ! global domain versus local domain 
    127127      CALL dom_nam                      ! read namelist ( namrun, namdom ) 
    128128      CALL dom_tile_init                ! Tile domain 
  • NEMO/trunk/src/OCE/DYN/dynldf_lap_blp_lf.F90

    r15033 r15267  
    102102         DO_3D( iij-1, iij-1, iij-1, iij-1, 1, jpkm1 )                           ! Horizontal slab 
    103103            !                                      ! ahm * e3 * curl  (warning: computed for ji-1,jj-1) 
    104             zcur     = ahmf(ji,jj,jk) * e3f(ji,jj,jk) * r1_e1e2f(ji,jj)               &   ! ahmf already * by fmask    
    105                &       * ( e2v(ji+1,jj) * pv(ji+1,jj,jk) - e2v(ji,jj) * pv(ji,jj,jk)  & 
    106                &       - e1u(ji,jj+1) * pu(ji,jj+1,jk) + e1u(ji,jj) * pu(ji,jj,jk) ) 
    107             zcur_jm1 = ahmf(ji  ,jj-1,jk) * e3f(ji,jj-1,jk) * r1_e1e2f(ji,jj-1)               &   ! ahmf already * by fmask 
     104            zcur     = ahmf(ji  ,jj  ,jk) * e3f(ji  ,jj  ,jk) * r1_e1e2f(ji  ,jj  )               &   ! ahmf already * by fmask    
     105               &       * ( e2v(ji+1,jj  ) * pv(ji+1,jj  ,jk) - e2v(ji,jj) * pv(ji,jj,jk)  & 
     106               &         - e1u(ji  ,jj+1) * pu(ji  ,jj+1,jk) + e1u(ji,jj) * pu(ji,jj,jk) ) 
     107            zcur_jm1 = ahmf(ji  ,jj-1,jk) * e3f(ji  ,jj-1,jk) * r1_e1e2f(ji  ,jj-1)               &   ! ahmf already * by fmask 
    108108               &       * ( e2v(ji+1,jj-1) * pv(ji+1,jj-1,jk) - e2v(ji,jj-1) * pv(ji,jj-1,jk)  & 
    109                &       - e1u(ji,jj) * pu(ji,jj,jk) + e1u(ji,jj-1) * pu(ji,jj-1,jk) ) 
    110             zcur_im1 = ahmf(ji-1,jj,jk) * e3f(ji-1,jj,jk) * r1_e1e2f(ji-1,jj)         &   ! ahmf already * by fmask 
    111                &       * ( e2v(ji,jj) * pv(ji,jj,jk) - e2v(ji-1,jj) * pv(ji-1,jj,jk)  & 
    112                &       - e1u(ji-1,jj+1) * pu(ji-1,jj+1,jk) + e1u(ji-1,jj) * pu(ji-1,jj,jk) ) 
     109               &         - e1u(ji  ,jj  ) * pu(ji  ,jj  ,jk) + e1u(ji,jj-1) * pu(ji,jj-1,jk) ) 
     110            zcur_im1 = ahmf(ji-1,jj  ,jk) * e3f(ji-1,jj  ,jk) * r1_e1e2f(ji-1,jj  )         &   ! ahmf already * by fmask 
     111               &       * ( e2v(ji  ,jj  ) * pv(ji  ,jj  ,jk) - e2v(ji-1,jj) * pv(ji-1,jj,jk)  & 
     112               &         - e1u(ji-1,jj+1) * pu(ji-1,jj+1,jk) + e1u(ji-1,jj) * pu(ji-1,jj,jk) ) 
    113113            !                                      ! ahm * div        (warning: computed for ji,jj) 
    114114            zdiv     = ahmt(ji,jj,jk) * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kbb)               &   ! ahmt already * by tmask 
  • NEMO/trunk/src/OCE/LBC/lbc_nfd_generic.h90

    r14433 r15267  
    2929            SELECT CASE ( cd_nat(jf) ) 
    3030            CASE ( 'T' , 'W' )                         ! T-, W-point 
    31                DO jl = 1, ipl; DO jk = 1, ipk 
     31               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    3232                  ! 
    3333                  ! last khls lines (from ipj to ipj-khls+1) : full 
     
    5353                     DO ji = 1, 1                 ! point ipi - khls + 1 
    5454                        ii1 = ipi - khls + ji 
    55                         ii2 =          khls + ji 
     55                        ii2 =       khls + ji 
    5656                        ptab(jf)%pt4d(ii1,ij1,jk,jl) = psgn(jf) * ptab(jf)%pt4d(ii2,ij2,jk,jl) 
    5757                     END DO 
     
    8484               END DO; END DO 
    8585            CASE ( 'U' )                               ! U-point 
    86                DO jl = 1, ipl; DO jk = 1, ipk 
     86               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    8787                  ! 
    8888                  ! last khls lines (from ipj to ipj-khls+1) : full 
     
    103103                     DO ji = 1, khls              ! last khls points 
    104104                        ii1 = ipi - khls + ji            ! ends at: ipi - khls + khls = ipi 
    105                         ii2 = ipi - khls + 1 - ji        ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
     105                        ii2 = ipi - khls - ji + 1        ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
    106106                        ptab(jf)%pt4d(ii1,ij1,jk,jl) = psgn(jf) * ptab(jf)%pt4d(ii2,ij2,jk,jl) 
    107107                     END DO 
     
    129129               END DO; END DO 
    130130            CASE ( 'V' )                               ! V-point 
    131                DO jl = 1, ipl; DO jk = 1, ipk 
     131               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    132132                  ! 
    133133                  ! last khls+1 lines (from ipj to ipj-khls) : full 
     
    153153                     DO ji = 1, 1                 ! point ipi - khls + 1 
    154154                        ii1 = ipi - khls + ji 
    155                         ii2 =          khls + ji 
     155                        ii2 =       khls + ji 
    156156                        ptab(jf)%pt4d(ii1,ij1,jk,jl) = psgn(jf) * ptab(jf)%pt4d(ii2,ij2,jk,jl) 
    157157                     END DO 
     
    165165               END DO; END DO 
    166166            CASE ( 'F' )                               ! F-point 
    167                DO jl = 1, ipl; DO jk = 1, ipk 
     167               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    168168                  ! 
    169169                  ! last khls+1 lines (from ipj to ipj-khls) : full 
     
    184184                     DO ji = 1, khls              ! last khls points 
    185185                        ii1 = ipi - khls + ji            ! ends at: ipi - khls + khls = ipi 
    186                         ii2 = ipi - khls + 1 - ji        ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
     186                        ii2 = ipi - khls - ji + 1        ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
    187187                        ptab(jf)%pt4d(ii1,ij1,jk,jl) = psgn(jf) * ptab(jf)%pt4d(ii2,ij2,jk,jl) 
    188188                     END DO 
     
    198198            SELECT CASE ( cd_nat(jf) ) 
    199199            CASE ( 'T' , 'W' )                         ! T-, W-point 
    200                DO jl = 1, ipl; DO jk = 1, ipk 
     200               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    201201                  ! 
    202202                  ! first: line number ipj-khls : 3 points 
     
    212212                     DO ji = 1, 1                 ! points ipi - khls 
    213213                        ii1 = ipi - khls + ji - 1 
    214                         ii2 =          khls + ji 
     214                        ii2 =       khls + ji 
    215215                        ptab(jf)%pt4d(ii1,ij1,jk,jl) =            ptab(jf)%pt4d(ii2,ij2,jk,jl)   ! Warning: pb with sign... 
    216216                     END DO 
     
    240240                     DO ji = 1, khls              ! last khls points 
    241241                        ii1 = ipi - khls + ji            ! ends at: ipi - khls + khls = ipi 
    242                         ii2 = ipi - khls + 1 - ji        ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
     242                        ii2 = ipi - khls - ji + 1        ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
    243243                        ptab(jf)%pt4d(ii1,ij1,jk,jl) = psgn(jf) * ptab(jf)%pt4d(ii2,ij2,jk,jl) 
    244244                     END DO 
     
    247247               END DO; END DO 
    248248            CASE ( 'U' )                               ! U-point 
    249                DO jl = 1, ipl; DO jk = 1, ipk 
     249               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    250250                  ! 
    251251                  ! last khls lines (from ipj to ipj-khls+1) : full 
     
    283283               END DO; END DO 
    284284            CASE ( 'V' )                               ! V-point 
    285                DO jl = 1, ipl; DO jk = 1, ipk 
     285               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    286286                  ! 
    287287                  ! last khls lines (from ipj to ipj-khls+1) : full 
     
    302302                     DO ji = 1, khls            ! last khls points 
    303303                        ii1 = ipi - khls + ji          ! ends at: ipi - khls + khls = ipi 
    304                         ii2 = ipi - khls + 1 - ji      ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
     304                        ii2 = ipi - khls - ji + 1      ! ends at: ipi - khls + 1 - khls = ipi - 2*khls + 1 
    305305                        ptab(jf)%pt4d(ii1,ij1,jk,jl) = psgn(jf) * ptab(jf)%pt4d(ii2,ij2,jk,jl) 
    306306                     END DO 
     
    328328               END DO; END DO 
    329329            CASE ( 'F' )                               ! F-point 
    330                DO jl = 1, ipl; DO jk = 1, ipk 
     330               DO jl = 1, ipl   ;  DO jk = 1, ipk 
    331331                  ! 
    332332                  ! last khls lines (from ipj to ipj-khls+1) : full 
  • NEMO/trunk/src/OCE/LBC/lbcnfd.F90

    r14433 r15267  
    1010 
    1111   !!---------------------------------------------------------------------- 
    12    !!   lbc_nfd       : generic interface for lbc_nfd_3d and lbc_nfd_2d routines 
    13    !!   lbc_nfd_3d    : lateral boundary condition: North fold treatment for a 3D arrays   (lbc_nfd) 
    14    !!   lbc_nfd_2d    : lateral boundary condition: North fold treatment for a 2D arrays   (lbc_nfd) 
    15    !!   lbc_nfd_nogather       : generic interface for lbc_nfd_nogather_3d and  
    16    !!                            lbc_nfd_nogather_2d routines (designed for use 
    17    !!                            with ln_nnogather to avoid global width arrays 
    18    !!                            mpi all gather operations) 
     12   !!   lbc_nfd       : generic interface for lbc_nfd_sp and lbc_nfd_dp routines that is doing the north fold in a non-mpi case  
     13   !!   mpp_nfd       : generic interface for mpp_nfd_sp and mpp_nfd_dp routines that will use lbc_nfd directly or indirectly 
    1914   !!---------------------------------------------------------------------- 
    2015   USE dom_oce        ! ocean space and time domain  
     
    3631      MODULE PROCEDURE   mpp_nfd_sp, mpp_nfd_dp 
    3732   END INTERFACE 
    38  
    39    INTERFACE lbc_nfd_nogather   ! called by mpp_nfd 
    40       MODULE PROCEDURE   lbc_nfd_nogather_sp, lbc_nfd_nogather_dp 
    41    END INTERFACE 
    4233    
    4334   PUBLIC   mpp_nfd            ! mpi north fold conditions 
    4435   PUBLIC   lbc_nfd            ! north fold conditions 
    45    PUBLIC   lbc_nfd_nogather   ! north fold conditions (no allgather case) 
    4636 
    47    INTEGER, PUBLIC, PARAMETER            ::   jpmaxngh = 3               !: 
    48    INTEGER, PUBLIC                       ::   nsndto                     !: 
    49    INTEGER, PUBLIC, DIMENSION (jpmaxngh) ::   isendto                    !: processes to which communicate 
    50    INTEGER, PUBLIC                       ::   ijpj 
    51  
     37   INTEGER, PUBLIC                               :: nfd_nbnei 
     38   INTEGER, PUBLIC, ALLOCATABLE, DIMENSION (:  ) :: nfd_rknei 
     39   INTEGER, PUBLIC, ALLOCATABLE, DIMENSION (:,:) :: nfd_rksnd 
     40   INTEGER, PUBLIC, ALLOCATABLE, DIMENSION (:,:) :: nfd_jisnd 
     41    
    5242   !!---------------------------------------------------------------------- 
    5343   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    5949   !!---------------------------------------------------------------------- 
    6050   !!                   ***  routine lbc_nfd_[sd]p  *** 
    61    !!               ***  routine lbc_nfd_nogather_[sd]p  *** 
    6251   !!                   ***  routine lbc_nfd_ext_[sd]p  *** 
    6352   !!---------------------------------------------------------------------- 
     
    7564#define PRECISION sp 
    7665#  include "lbc_nfd_generic.h90" 
    77 #  include "lbc_nfd_nogather_generic.h90" 
    7866#  include "lbc_nfd_ext_generic.h90" 
    7967#undef PRECISION 
     
    8371#define PRECISION dp 
    8472#  include "lbc_nfd_generic.h90" 
    85 #  include "lbc_nfd_nogather_generic.h90" 
    8673#  include "lbc_nfd_ext_generic.h90" 
    8774#undef PRECISION 
  • NEMO/trunk/src/OCE/LBC/lib_mpp.F90

    r15023 r15267  
    4949   !!   mppsync       : 
    5050   !!   mppstop       : 
    51    !!   mpp_ini_north : initialisation of north fold 
     51   !!   mpp_ini_northgather : initialisation of north fold with gathering of the communications 
    5252   !!   mpp_lbc_north_icb : alternative to mpp_nfd for extra outer halo with icebergs 
    5353   !!   mpp_bcast_nml : broadcast/receive namelist character buffer from reading process to all others 
     
    6464   PUBLIC   ctl_stop, ctl_warn, ctl_opn, ctl_nam, load_nml 
    6565   PUBLIC   mpp_start, mppstop, mppsync, mpp_comm_free 
    66    PUBLIC   mpp_ini_north 
     66   PUBLIC   mpp_ini_northgather 
    6767   PUBLIC   mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc 
    6868   PUBLIC   mpp_delay_max, mpp_delay_sum, mpp_delay_rcv 
     
    11521152 
    11531153 
    1154    SUBROUTINE mpp_ini_north 
    1155       !!---------------------------------------------------------------------- 
    1156       !!               ***  routine mpp_ini_north  *** 
     1154   SUBROUTINE mpp_ini_northgather 
     1155      !!---------------------------------------------------------------------- 
     1156      !!               ***  routine mpp_ini_northgather  *** 
    11571157      !! 
    11581158      !! ** Purpose :   Initialize special communicator for north folding 
     
    12101210      ! 
    12111211#endif 
    1212    END SUBROUTINE mpp_ini_north 
     1212   END SUBROUTINE mpp_ini_northgather 
    12131213 
    12141214 
  • NEMO/trunk/src/OCE/LBC/mpp_nfd_generic.h90

    r14433 r15267  
    1010      ! 
    1111      LOGICAL  ::   ll_add_line 
    12       INTEGER  ::   ji,  jj,  jk,  jl, jh, jf, jr   ! dummy loop indices 
     12      INTEGER  ::   ji,  jj,  jk,  jl, jf, jr, jg, jn   ! dummy loop indices 
    1313      INTEGER  ::   ipi, ipj, ipj2, ipk, ipl, ipf   ! dimension of the input array 
    14       INTEGER  ::   imigr, iihom, ijhom             ! local integers 
    1514      INTEGER  ::   ierr, ibuffsize, iis0, iie0, impp 
    16       INTEGER  ::   ii1, ii2, ij1, ij2 
    17       INTEGER  ::   ipimax, i0max 
     15      INTEGER  ::   ii1, ii2, ij1, ij2, iis, iie, iib, iig, iin 
     16      INTEGER  ::   i0max 
    1817      INTEGER  ::   ij, iproc, ipni, ijnr 
    19       INTEGER, DIMENSION (jpmaxngh)       ::   ml_req_nf   ! for mpi_isend when avoiding mpi_allgather 
    20       INTEGER                             ::   ml_err      ! for mpi_isend when avoiding mpi_allgather 
    21       !                                                    ! Workspace for message transfers avoiding mpi_allgather 
    22       INTEGER                             ::   ipj_b       ! sum of lines for all multi fields 
    23       INTEGER                             ::   i012        ! 0, 1 or 2 
    24       INTEGER , DIMENSION(:,:)        , ALLOCATABLE ::   jj_s  ! position of sent lines 
    25       INTEGER , DIMENSION(:,:)        , ALLOCATABLE ::   jj_b  ! position of buffer lines 
    26       INTEGER , DIMENSION(:)          , ALLOCATABLE ::   ipj_s ! number of sent lines 
    27       REAL(PRECISION), DIMENSION(:,:,:,:)    , ALLOCATABLE ::   ztabb, ztabr, ztabw  ! buffer, receive and work arrays 
     18      INTEGER, DIMENSION (:), ALLOCATABLE ::   ireq_s, ireq_r   ! for mpi_isend when avoiding mpi_allgather 
     19      INTEGER                             ::   ipjtot           ! sum of lines for all multi fields 
     20      INTEGER                             ::   i012             ! 0, 1 or 2 
     21      INTEGER , DIMENSION(:,:)        , ALLOCATABLE ::   ijsnd  ! j-position of sent lines for each field 
     22      INTEGER , DIMENSION(:,:)        , ALLOCATABLE ::   ijbuf  ! j-position of send buffer lines for each field 
     23      INTEGER , DIMENSION(:,:)        , ALLOCATABLE ::   ijrcv  ! j-position of recv buffer lines for each field 
     24      INTEGER , DIMENSION(:,:)        , ALLOCATABLE ::   ii1st, iiend 
     25      INTEGER , DIMENSION(:)          , ALLOCATABLE ::   ipjfld ! number of sent lines for each field 
     26      REAL(PRECISION), DIMENSION(:,:,:,:)    , ALLOCATABLE ::   zbufs  ! buffer, receive and work arrays 
     27      REAL(PRECISION), DIMENSION(:,:,:,:,:)  , ALLOCATABLE ::   zbufr  ! buffer, receive and work arrays 
    2828      REAL(PRECISION), DIMENSION(:,:,:,:,:)  , ALLOCATABLE ::   znorthloc 
    2929      REAL(PRECISION), DIMENSION(:,:,:,:,:,:), ALLOCATABLE ::   znorthglo 
     
    6262         ll_add_line = l_full_nf_update .OR. ( ncom_stp <= nit000+1 .AND. .NOT. ln_rstart ) 
    6363          
    64          ALLOCATE(ipj_s(ipf))                ! how many lines do we exchange? 
     64         ALLOCATE(ipjfld(ipf))                 ! how many lines do we exchange for each field? 
    6565         IF( ll_add_line ) THEN 
    66             DO jf = 1, ipf                      ! Loop over the number of arrays to be processed 
    67                ipj_s(jf) = khls + COUNT( (/ c_NFtype == 'T' .OR. cd_nat(jf) == 'V' .OR. cd_nat(jf) == 'F' /) )  
     66            DO jf = 1, ipf                     ! Loop over the number of arrays to be processed 
     67               ipjfld(jf) = khls + COUNT( (/ c_NFtype == 'T' .OR. cd_nat(jf) == 'V' .OR. cd_nat(jf) == 'F' /) ) 
    6868            END DO 
    6969         ELSE 
    70             ipj_s(:) = khls 
     70            ipjfld(:) = khls 
    7171         ENDIF 
    7272          
    73          ipj   = MAXVAL(ipj_s(:))            ! Max 2nd dimension of message transfers 
    74          ipj_b = SUM(   ipj_s(:))            ! Total number of lines to be exchanged 
    75          ALLOCATE( jj_s(ipj, ipf), jj_b(ipj, ipf) ) 
     73         ipj    = MAXVAL(ipjfld(:))            ! Max 2nd dimension of message transfers 
     74         ipjtot = SUM(   ipjfld(:))            ! Total number of lines to be exchanged 
    7675 
    7776         ! Index of modifying lines in input 
     77         ALLOCATE( ijsnd(ipj, ipf), ijbuf(ipj, ipf), ijrcv(ipj, ipf), ii1st(ipj, ipf), iiend(ipj, ipf) ) 
     78 
    7879         ij1 = 0 
    79          DO jf = 1, ipf                      ! Loop over the number of arrays to be processed 
    80             ! 
     80         DO jf = 1, ipf                        ! Loop over the number of arrays to be processed 
     81            ! 
     82            DO jj = 1, khls   ! first khls lines (starting from top) must be fully defined 
     83               ii1st(jj, jf) = 1 
     84               iiend(jj, jf) = jpi 
     85            END DO 
     86            ! 
     87            ! what do we do with line khls+1 (starting from top) 
    8188            IF( c_NFtype == 'T' ) THEN          ! *  North fold  T-point pivot 
    8289               SELECT CASE ( cd_nat(jf) ) 
    83                CASE ( 'T', 'W', 'U' )   ;   i012 = 1   ! T-, U-, W-point 
    84                CASE ( 'V', 'F'      )   ;   i012 = 2   ! V-, F-point 
     90               CASE ('T','W')   ;   i012 = 1   ;   ii1st(khls+1, jf) = mi0(jpiglo/2+2)   ;   iiend(khls+1, jf) = mi1(jpiglo-khls) 
     91               CASE ('U'    )   ;   i012 = 1   ;   ii1st(khls+1, jf) = mi0(jpiglo/2+1)   ;   iiend(khls+1, jf) = mi1(jpiglo-khls) 
     92               CASE ('V'    )   ;   i012 = 2   ;   ii1st(khls+1, jf) = 1                 ;   iiend(khls+1, jf) = jpi 
     93               CASE ('F'    )   ;   i012 = 2   ;   ii1st(khls+1, jf) = 1                 ;   iiend(khls+1, jf) = jpi 
    8594               END SELECT 
    8695            ENDIF 
    8796            IF( c_NFtype == 'F' ) THEN          ! *  North fold  F-point pivot 
    8897               SELECT CASE ( cd_nat(jf) ) 
    89                CASE ( 'T', 'W', 'U' )   ;   i012 = 0   ! T-, U-, W-point 
    90                CASE ( 'V', 'F'      )   ;   i012 = 1   ! V-, F-point 
     98               CASE ('T','W')   ;   i012 = 0   ! we don't touch line khls+1 
     99               CASE ('U'    )   ;   i012 = 0   ! we don't touch line khls+1 
     100               CASE ('V'    )   ;   i012 = 1   ;   ii1st(khls+1, jf) = mi0(jpiglo/2+1)   ;   iiend(khls+1, jf) = mi1(jpiglo-khls  ) 
     101               CASE ('F'    )   ;   i012 = 1   ;   ii1st(khls+1, jf) = mi0(jpiglo/2+1)   ;   iiend(khls+1, jf) = mi1(jpiglo-khls-1) 
    91102               END SELECT 
    92103            ENDIF 
    93                ! 
    94             DO jj = 1, ipj_s(jf) 
     104            ! 
     105            DO jj = 1, ipjfld(jf) 
    95106               ij1 = ij1 + 1 
    96                jj_b(jj,jf) = ij1 
    97                jj_s(jj,jf) = jpj - 2*khls + jj - i012 
     107               ijsnd(jj,jf) = jpj - 2*khls + jj - i012   ! sent lines (from bottom of sent lines) 
     108               ijbuf(jj,jf) = ij1                        ! gather all lines in the snd/rcv buffers 
     109               ijrcv(jj,jf) = jpj - jj + 1               ! recv lines (from the top -> reverse order for jj) 
    98110            END DO 
    99111            ! 
    100112         END DO 
    101113         ! 
    102          ALLOCATE( ztabb(jpimax,ipj_b,ipk,ipl) )   ! store all the data to be sent in a buffer array 
    103          ibuffsize = jpimax * ipj_b * ipk * ipl 
    104          ! 
     114         i0max = jpimax - 2 * khls                                    ! we are not sending the halos 
     115         ALLOCATE( zbufs(i0max,ipjtot,ipk,ipl), ireq_s(nfd_nbnei) )   ! store all the data to be sent in a buffer array 
     116         ibuffsize = i0max * ipjtot * ipk * ipl 
     117         ! 
     118         ! fill the send buffer with all the lines 
    105119         DO jf = 1, ipf   ;   DO jl = 1, ipl   ;   DO jk = 1, ipk 
    106             DO jj = 1, ipj_s(jf) 
    107                ij1 = jj_b(jj,jf) 
    108                ij2 = jj_s(jj,jf) 
    109                DO ji = 1, jpi 
    110                   ztabb(ji,ij1,jk,jl) = ptab(jf)%pt4d(ji,ij2,jk,jl) 
    111                END DO 
    112                DO ji = jpi+1, jpimax 
    113                   ztabb(ji,ij1,jk,jl) = HUGE(0._/**/PRECISION)   ! avoid sending uninitialized values (make sure we don't use it) 
     120            DO jj = 1, ipjfld(jf) 
     121               ij1 = ijbuf(jj,jf) 
     122               ij2 = ijsnd(jj,jf) 
     123               DO ji = Nis0, Nie0       ! should not use any other value 
     124                  iib = ji - Nis0 + 1 
     125                  zbufs(iib,ij1,jk,jl) = ptab(jf)%pt4d(ji,ij2,jk,jl) 
     126               END DO 
     127               DO ji = Ni_0+1, i0max    ! avoid sending uninitialized values (make sure we don't use it) 
     128                  zbufs(ji,ij1,jk,jl) = HUGE(0._/**/PRECISION)   ! make sure we don't use it... 
    114129               END DO 
    115130            END DO 
     
    119134         IF( ln_timing ) CALL tic_tac(.TRUE.) 
    120135         ! 
    121          ! send the data as soon as possible 
    122          DO jr = 1, nsndto 
    123             iproc = nfproc(isendto(jr)) 
     136         ! send the same buffer data to all neighbourgs as soon as possible 
     137         DO jn = 1, nfd_nbnei 
     138            iproc = nfd_rknei(jn) 
    124139            IF( iproc /= narea-1 .AND. iproc /= -1 ) THEN 
    125140#if ! defined key_mpi_off 
    126                CALL MPI_ISEND( ztabb, ibuffsize, MPI_TYPE, iproc, 5, mpi_comm_oce, ml_req_nf(jr), ierr ) 
     141               CALL MPI_Isend( zbufs, ibuffsize, MPI_TYPE, iproc, 5, mpi_comm_oce, ireq_s(jn), ierr ) 
    127142#endif 
     143            ELSE 
     144               ireq_s(jn) = MPI_REQUEST_NULL 
    128145            ENDIF 
    129146         END DO 
    130147         ! 
    131          ipimax = jpimax * jpmaxngh 
    132          ALLOCATE( ztabw(jpimax,ipj_b,ipk,ipl), ztabr(ipimax,ipj_b,ipk,ipl) )  
    133          ! 
    134          DO jr = 1, nsndto 
    135             ! 
    136             ipni  = isendto(jr) 
    137             iproc = nfproc(ipni) 
    138             ipi   = nfjpi (ipni) 
    139             ! 
    140             IF( ipni ==   1  ) THEN   ;   iis0 =   1          ! domain  left side: as e-w comm already done -> from 1st column 
    141             ELSE                      ;   iis0 =   1 + khls   ! default: -> from inner domain  
    142             ENDIF 
    143             IF( ipni == jpni ) THEN   ;   iie0 = ipi          ! domain right side: as e-w comm already done -> until last column 
    144             ELSE                      ;   iie0 = ipi - khls   ! default: -> until inner domain  
    145             ENDIF 
    146             impp = nfimpp(ipni) - nfimpp(isendto(1)) 
     148         ALLOCATE( zbufr(i0max,ipjtot,ipk,ipl,nfd_nbnei), ireq_r(nfd_nbnei) )  
     149         ! 
     150         DO jn = 1, nfd_nbnei 
     151            ! 
     152            iproc = nfd_rknei(jn) 
    147153            ! 
    148154            IF(           iproc == -1 ) THEN   ! No neighbour (land proc that was suppressed) 
    149155               ! 
     156               ireq_r(jn) = MPI_REQUEST_NULL                ! no message to be received 
     157               zbufr(:,:,:,:,jn) = HUGE(0._/**/PRECISION)   ! default: define it and make sure we don't use it... 
    150158               SELECT CASE ( kfillmode ) 
    151                CASE ( jpfillnothing )               ! no filling  
    152                CASE ( jpfillcopy    )               ! filling with inner domain values 
     159               CASE ( jpfillnothing )                       ! no filling  
     160               CASE ( jpfillcopy    )                       ! filling with inner domain values 
    153161                  DO jf = 1, ipf   ;   DO jl = 1, ipl   ;   DO jk = 1, ipk 
    154                      DO jj = 1, ipj_s(jf) 
    155                         ij1 = jj_b(jj,jf) 
    156                         ij2 = jj_s(jj,jf) 
    157                         DO ji = iis0, iie0 
    158                            ztabr(impp+ji,ij1,jk,jl) = ptab(jf)%pt4d(Nis0,ij2,jk,jl)   ! chose to take the 1st iner domain point 
    159                         END DO 
     162                     DO jj = 1, ipjfld(jf) 
     163                        ij1 = ijbuf(jj,jf) 
     164                        ij2 = ijsnd(jj,jf)                                      ! we will use only the first value, see init_nfdcom 
     165                        zbufr(1,ij1,jk,jl,jn) = ptab(jf)%pt4d(Nis0,ij2,jk,jl)   ! chose to take the 1st inner domain point 
    160166                     END DO 
    161167                  END DO   ;   END DO   ;   END DO 
    162                CASE ( jpfillcst     )               ! filling with constant value 
    163                   DO jl = 1, ipl   ;   DO jk = 1, ipk 
    164                      DO jj = 1, ipj_b 
    165                         DO ji = iis0, iie0 
    166                            ztabr(impp+ji,jj,jk,jl) = pfillval 
    167                         END DO 
    168                      END DO 
    169                   END DO   ;   END DO 
     168               CASE ( jpfillcst     )                       ! filling with constant value 
     169                  zbufr(1,:,:,:,jn) = pfillval              ! we will use only the first value, see init_nfdcom 
    170170               END SELECT 
    171171               ! 
    172172            ELSE IF( iproc == narea-1 ) THEN   ! get data from myself! 
    173173               ! 
     174               ireq_r(jn) = MPI_REQUEST_NULL                ! no message to be received 
    174175               DO jf = 1, ipf   ;   DO jl = 1, ipl  ;   DO jk = 1, ipk 
    175                   DO jj = 1, ipj_s(jf) 
    176                      ij1 = jj_b(jj,jf) 
    177                      ij2 = jj_s(jj,jf) 
    178                      DO ji = iis0, iie0 
    179                         ztabr(impp+ji,ij1,jk,jl) = ptab(jf)%pt4d(ji,ij2,jk,jl) 
     176                  DO jj = 1, ipjfld(jf) 
     177                     ij1 = ijbuf(jj,jf) 
     178                     ij2 = ijsnd(jj,jf) 
     179                     DO ji = Nis0, Nie0                     ! should not use any other value 
     180                        iib = ji - Nis0 + 1 
     181                        zbufr(iib,ij1,jk,jl,jn) = ptab(jf)%pt4d(ji,ij2,jk,jl) 
    180182                     END DO 
    181183                  END DO 
     
    183185               ! 
    184186            ELSE                               ! get data from a neighbour trough communication 
    185                !   
    186187#if ! defined key_mpi_off 
    187                CALL MPI_RECV( ztabw, ibuffsize, MPI_TYPE, iproc, 5, mpi_comm_oce, MPI_STATUS_IGNORE, ierr ) 
     188               CALL MPI_Irecv( zbufr(:,:,:,:,jn), ibuffsize, MPI_TYPE, iproc, 5, mpi_comm_oce, ireq_r(jn), ierr ) 
    188189#endif 
    189                DO jl = 1, ipl   ;   DO jk = 1, ipk 
    190                   DO jj = 1, ipj_b 
    191                      DO ji = iis0, iie0 
    192                         ztabr(impp+ji,jj,jk,jl) = ztabw(ji,jj,jk,jl) 
    193                      END DO 
    194                   END DO 
    195                END DO   ;   END DO 
    196                 
    197             ENDIF 
    198             ! 
    199          END DO   ! nsndto 
     190            ENDIF 
     191            ! 
     192         END DO   ! nfd_nbnei 
     193         ! 
     194         CALL mpi_waitall(nfd_nbnei, ireq_r, MPI_STATUSES_IGNORE, ierr)   ! wait for all Irecv 
    200195         ! 
    201196         IF( ln_timing ) CALL tic_tac(.FALSE.) 
     
    204199         ! 
    205200         DO jf = 1, ipf 
    206             ij1 = jj_b(       1 ,jf) 
    207             ij2 = jj_b(ipj_s(jf),jf) 
    208             CALL lbc_nfd_nogather( ptab(jf)%pt4d(:,:,:,:), ztabr(:,ij1:ij2,:,:), cd_nat(jf), psgn(jf), khls ) 
    209          END DO 
    210          ! 
    211          DEALLOCATE( ztabr, ztabw, jj_s, jj_b, ipj_s ) 
    212          ! 
    213          DO jr = 1,nsndto 
    214             iproc = nfproc(isendto(jr)) 
    215             IF( iproc /= narea-1 .AND. iproc /= -1 ) THEN 
    216                CALL mpi_wait( ml_req_nf(jr), MPI_STATUS_IGNORE, ml_err )   ! put the wait at the very end just before the deallocate 
    217             ENDIF 
    218          END DO 
    219          DEALLOCATE( ztabb ) 
     201            ! 
     202            SELECT CASE ( cd_nat(jf) )     ! which grid number? 
     203            CASE ('T','W')   ;   iig = 1   ! T-, W-point 
     204            CASE ('U')       ;   iig = 2   ! U-point 
     205            CASE ('V')       ;   iig = 3   ! V-point 
     206            CASE ('F')       ;   iig = 4   ! F-point 
     207            END SELECT 
     208            ! 
     209            DO jl = 1, ipl   ;   DO jk = 1, ipk 
     210               ! 
     211               ! if T point with F-point pivot : must be done first 
     212               !    --> specific correction of 3 points near the 2 pivots (to be clean, usually masked -> so useless)  
     213               IF( c_NFtype == 'F' .AND. iig == 1 ) THEN 
     214                  ij1 = jpj - khls     ! j-index in the receiving array 
     215                  ij2 = 1              ! only 1 line in the buffer 
     216                  DO ji = mi0(khls), mi1(khls) 
     217                     iib = nfd_jisnd(mi0(       khls),iig)   ! i-index in the buffer 
     218                     iin = nfd_rksnd(mi0(       khls),iig)   ! neigbhour-index in the buffer 
     219                     IF( nfd_rknei(iin) == -1 .AND. kfillmode == jpfillnothing )   CYCLE 
     220                     ptab(jf)%pt4d(ji,ij1,jk,jl) = zbufr(iib,ij2,jk,jl,iin)   ! no psgn(jf) 
     221                  END DO 
     222                  DO ji = mi0(jpiglo/2+1), mi1(jpiglo/2+1) 
     223                     iib = nfd_jisnd(mi0( jpiglo/2+1),iig)   ! i-index in the buffer 
     224                     iin = nfd_rksnd(mi0( jpiglo/2+1),iig)   ! neigbhour-index in the buffer 
     225                     IF( nfd_rknei(iin) == -1 .AND. kfillmode == jpfillnothing )   CYCLE 
     226                     ptab(jf)%pt4d(ji,ij1,jk,jl) = zbufr(iib,ij2,jk,jl,iin)   ! no psgn(jf) 
     227                  END DO 
     228                  DO ji = mi0(jpiglo-khls), mi1(jpiglo-khls) 
     229                     iib = nfd_jisnd(mi0(jpiglo-khls),iig)   ! i-index in the buffer 
     230                     iin = nfd_rksnd(mi0(jpiglo-khls),iig)   ! neigbhour-index in the buffer 
     231                     IF( nfd_rknei(iin) == -1 .AND. kfillmode == jpfillnothing )   CYCLE 
     232                     ptab(jf)%pt4d(ji,ij1,jk,jl) = zbufr(iib,ij2,jk,jl,iin)   ! no psgn(jf) 
     233                  END DO 
     234               ENDIF 
     235               ! 
     236               ! Apply the North pole folding. 
     237               DO jj = 1, ipjfld(jf)   ! for all lines to be exchanged for this field 
     238                  ij1 = ijrcv(jj,jf)   ! j-index in the receiving array 
     239                  ij2 = ijbuf(jj,jf)   ! j-index in the buffer 
     240                  iis = ii1st(jj,jf)   ! stating i-index in the receiving array 
     241                  iie = iiend(jj,jf)   !  ending i-index in the receiving array 
     242                  DO ji = iis, iie  
     243                     iib = nfd_jisnd(ji,iig)   ! i-index in the buffer 
     244                     iin = nfd_rksnd(ji,iig)   ! neigbhour-index in the buffer 
     245                     IF( nfd_rknei(iin) == -1 .AND. kfillmode == jpfillnothing )   CYCLE 
     246                     ptab(jf)%pt4d(ji,ij1,jk,jl) = psgn(jf) * zbufr(iib,ij2,jk,jl,iin) 
     247                  END DO 
     248               END DO 
     249               ! 
     250               ! re-apply periodocity when we modified the eastern side of the inner domain (and not the full line) 
     251               IF( c_NFtype == 'T' ) THEN          ! *  North fold  T-point pivot 
     252                  IF(     iig <= 2 ) THEN   ;   iis = mi0(1)   ;   iie = mi1(khls)   ! 'T','W','U': update west halo 
     253                  ELSE                      ;   iis = 1        ;   iie = 0           ! 'V','F'    : full line already exchanged 
     254                  ENDIF 
     255               ENDIF 
     256               IF( c_NFtype == 'F' ) THEN          ! *  North fold  F-point pivot 
     257                  IF(     iig <= 2 ) THEN   ;   iis = 1        ;   iie = 0           ! 'T','W','U': nothing to do 
     258                  ELSEIF( iig == 3 ) THEN   ;   iis = mi0(1)   ;   iie = mi1(khls)   ! 'V'        : update west halo 
     259                  ELSEIF( khls > 1 ) THEN   ;   iis = mi0(1)   ;   iie = mi1(khls-1) ! 'F' and khls > 1 
     260                  ELSE                      ;   iis = 1        ;   iie = 0           ! 'F' and khls == 1 : nothing to do 
     261                  ENDIF 
     262               ENDIF 
     263               jj  = ipjfld(jf)     ! only for the last line of this field 
     264               ij1 = ijrcv(jj,jf)   ! j-index in the receiving array 
     265               ij2 = ijbuf(jj,jf)   ! j-index in the buffer 
     266               DO ji = iis, iie 
     267                  iib = nfd_jisnd(ji,iig)   ! i-index in the buffer 
     268                  iin = nfd_rksnd(ji,iig)   ! neigbhour-index in the buffer 
     269                  IF( nfd_rknei(iin) == -1 .AND. kfillmode == jpfillnothing )   CYCLE 
     270                  ptab(jf)%pt4d(ji,ij1,jk,jl) = psgn(jf) * zbufr(iib,ij2,jk,jl,iin) 
     271               END DO 
     272               !                
     273            END DO   ;   END DO   ! ipl   ; ipk 
     274            !                
     275         END DO   ! ipf 
     276        
     277         ! 
     278         DEALLOCATE( zbufr, ireq_r, ijsnd, ijbuf, ijrcv, ii1st, iiend, ipjfld ) 
     279         ! 
     280         CALL mpi_waitall(nfd_nbnei, ireq_s, MPI_STATUSES_IGNORE, ierr)   ! wait for all Isend 
     281         ! 
     282         DEALLOCATE( zbufs, ireq_s ) 
    220283         ! 
    221284      ELSE                             !==  allgather exchanges  ==! 
     
    265328              ! 
    266329               SELECT CASE ( kfillmode ) 
    267                CASE ( jpfillnothing )               ! no filling  
     330               CASE ( jpfillnothing )               ! no filling 
     331                  CALL ctl_stop( 'STOP', 'mpp_nfd_generic : cannot use jpfillnothing with ln_nnogather = F') 
    268332               CASE ( jpfillcopy    )               ! filling with inner domain values 
    269333                  DO jf = 1, ipf   ;   DO jl = 1, ipl   ;   DO jk = 1, ipk 
     
    329393         DEALLOCATE( ztabglo ) 
    330394         ! 
    331       ENDIF   ! l_north_nogather 
     395      ENDIF   ! ln_nnogather 
    332396      ! 
    333397   END SUBROUTINE mpp_nfd_/**/PRECISION 
  • NEMO/trunk/src/OCE/LBC/mppini.F90

    r15033 r15267  
    2323   USE bdy_oce        ! open BounDarY 
    2424   ! 
    25    USE lbcnfd  , ONLY : isendto, nsndto ! Setup of north fold exchanges 
     25   USE lbcnfd        ! Setup of north fold exchanges 
    2626   USE lib_mpp        ! distribued memory computing library 
    2727   USE iom            ! nemo I/O library 
     
    331331      njmpp = ijmppt(ii,ij) 
    332332      ! 
    333       CALL init_doloop                          ! set start/end indices of do-loop, depending on the halo width value (nn_hls) 
     333      CALL init_doloop    ! set start/end indices of do-loop, depending on the halo width value (nn_hls) 
     334      CALL init_locglo    ! define now functions needed to convert indices from/to global to/from local domains 
    334335      ! 
    335336      IF(lwp) THEN 
     
    503504         WRITE(numout,*) '      mpi nei no-we = ', mpinei(jpnw)  , '   mpi nei no-ea = ', mpinei(jpne) 
    504505      ENDIF 
    505       !                          ! Prepare mpp north fold 
    506       ! 
    507       llmpiNFold =          jpni  > 1 .AND. l_NFold   ! is the North fold done with an MPI communication? 
    508       l_IdoNFold = ijn(narea) == jpnj .AND. l_NFold   ! is this process doing North fold? 
    509       ! 
    510       IF( llmpiNFold ) THEN 
    511          CALL mpp_ini_north 
    512          IF (lwp) THEN 
    513             WRITE(numout,*) 
    514             WRITE(numout,*) '   ==>>>   North fold boundary prepared for jpni >1' 
    515          ENDIF 
    516          IF (llwrtlay) THEN      ! additional prints in layout.dat 
    517             WRITE(inum,*) 
    518             WRITE(inum,*) 
    519             WRITE(inum,*) 'Number of subdomains located along the north fold : ', ndim_rank_north 
    520             WRITE(inum,*) 'Rank of the subdomains located along the north fold : ', ndim_rank_north 
    521             DO jp = 1, ndim_rank_north, 5 
    522                WRITE(inum,*) nrank_north( jp:MINVAL( (/jp+4,ndim_rank_north/) ) ) 
    523             END DO 
    524          ENDIF 
    525          IF ( l_IdoNFold .AND. ln_nnogather ) THEN 
    526             CALL init_nfdcom     ! northfold neighbour lists 
    527             IF (llwrtlay) THEN 
    528                WRITE(inum,*) 
    529                WRITE(inum,*) 'north fold exchanges with explicit point-to-point messaging :' 
    530                WRITE(inum,*) '   nsndto  : ', nsndto 
    531                WRITE(inum,*) '   isendto : ', isendto(1:nsndto) 
    532             ENDIF 
    533          ENDIF 
    534       ENDIF 
    535506      ! 
    536507      CALL mpp_ini_nc(nn_hls)    ! Initialize communicator for neighbourhood collective communications 
     
    539510         mpi_nc_com8(jh) = mpi_nc_com8(nn_hls) 
    540511      END DO 
    541       ! 
    542       IF( jpnij > 1) CALL init_excl_landpt      ! exclude exchanges which contain only land points 
    543       ! 
    544       ! Save processor layout changes in ascii file 
     512      !                          ! Exclude exchanges which contain only land points 
     513      ! 
     514      IF( jpnij > 1 ) CALL init_excl_landpt 
     515      ! 
     516      !                          ! Prepare mpp north fold 
     517      ! 
     518      llmpiNFold =          jpni  > 1 .AND. l_NFold           ! is the North fold done with an MPI communication? 
     519      l_IdoNFold = ijn(narea) == jpnj .AND. l_NFold           ! is this process doing North fold? 
     520      ! 
     521      IF( llmpiNFold )   CALL init_nfdcom( llwrtlay, inum )   ! init northfold communication, must be done after init_excl_landpt 
     522      ! 
     523      !                          !  Save processor layout changes in ascii file 
     524      ! 
    545525      DO jh = 1, n_hlsmax    ! different halo size 
    546526         DO ji = 1, 8 
     
    11651145      INTEGER ::   iiwe, iiea, iist, iisz  
    11661146      INTEGER ::   ijso, ijno, ijst, ijsz  
    1167       LOGICAL ::   llsave 
    11681147      REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   zmsk 
    11691148      LOGICAL , DIMENSION(Ni_0,Nj_0,1)      ::   lloce 
     
    11741153      ! 
    11751154      ! Here we look only at communications excluding the NP folding. 
    1176       ! As lbcnfd not validated if halo size /= nn_hls, we switch if off temporary... 
    1177       llsave = l_IdoNFold 
     1155      !   --> we switch off lbcnfd at this stage (init_nfdcom called after init_excl_landpt)... 
    11781156      l_IdoNFold = .FALSE. 
    11791157      ! 
     
    12331211         ! 
    12341212      END DO 
    1235       l_IdoNFold = llsave 
    12361213 
    12371214   END SUBROUTINE init_excl_landpt 
     
    12781255 
    12791256 
    1280    SUBROUTINE init_nfdcom 
     1257   SUBROUTINE init_nfdcom( ldwrtlay, knum ) 
    12811258      !!---------------------------------------------------------------------- 
    12821259      !!                     ***  ROUTINE  init_nfdcom  *** 
     
    12881265      !!    1.0  ! 2011-10  (A. C. Coward, NOCS & J. Donners, PRACE) 
    12891266      !!    2.0  ! 2013-06 Setup avoiding MPI communication (I. Epicoco, S. Mocavero, CMCC) 
    1290       !!---------------------------------------------------------------------- 
    1291       INTEGER  ::   sxM, dxM, sxT, dxT, jn 
    1292       !!---------------------------------------------------------------------- 
    1293       ! 
    1294       !sxM is the first point (in the global domain) needed to compute the north-fold for the current process 
    1295       sxM = jpiglo - nimpp - jpi + 1 
    1296       !dxM is the last point (in the global domain) needed to compute the north-fold for the current process 
    1297       dxM = jpiglo - nimpp + 2 
    1298       ! 
    1299       ! loop over the other north-fold processes to find the processes 
    1300       ! managing the points belonging to the sxT-dxT range 
    1301       ! 
    1302       nsndto = 0 
    1303       DO jn = 1, jpni 
     1267      !!    3.0  ! 2021-09 complete rewrite using informations from gather north fold 
     1268      !!---------------------------------------------------------------------- 
     1269      LOGICAL, INTENT(in   ) ::   ldwrtlay   ! true if additional prints in layout.dat 
     1270      INTEGER, INTENT(in   ) ::   knum       ! layout.dat unit 
     1271      ! 
     1272      REAL(wp), DIMENSION(jpi,jpj,2,4) ::   zinfo 
     1273      INTEGER , DIMENSION(10) ::   irknei ! too many elements but safe... 
     1274      INTEGER                 ::   ji, jj, jg, jn   ! dummy loop indices 
     1275      LOGICAL                 ::   lnew 
     1276      !!---------------------------------------------------------------------- 
     1277      ! 
     1278      IF (lwp) THEN 
     1279         WRITE(numout,*) 
     1280         WRITE(numout,*) '   ==>>>   North fold boundary prepared for jpni >1' 
     1281      ENDIF 
     1282      ! 
     1283      CALL mpp_ini_northgather   ! we need to init the nfd with gathering in all cases as it is used to define the no-gather case 
     1284      ! 
     1285      IF(ldwrtlay) THEN      ! additional prints in layout.dat 
     1286         WRITE(knum,*) 
     1287         WRITE(knum,*) 
     1288         WRITE(knum,*) 'Number of subdomains located along the north fold : ', ndim_rank_north 
     1289         WRITE(knum,*) 'Rank of the subdomains located along the north fold : ', ndim_rank_north 
     1290         DO jn = 1, ndim_rank_north, 5 
     1291            WRITE(knum,*) nrank_north( jn:MINVAL( (/jn+4,ndim_rank_north/) ) ) 
     1292         END DO 
     1293      ENDIF 
     1294       
     1295      nfd_nbnei = 0   ! defaul def (useless?) 
     1296      IF( ln_nnogather ) THEN 
    13041297         ! 
    1305          sxT = nfimpp(jn)                    ! sxT = 1st  point (in the global domain) of the jn process 
    1306          dxT = nfimpp(jn) + nfjpi(jn) - 1    ! dxT = last point (in the global domain) of the jn process 
     1298         ! Use the "gather nfd" to know how to do the nfd: for ji point, which process send data from which of its ji-index? 
     1299         ! Note that nfd is perfectly symetric: I receive data from X <=> I send data to X  (-> no deadlock) 
    13071300         ! 
    1308          IF    ( sxT < sxM  .AND.  sxM < dxT ) THEN 
    1309             nsndto          = nsndto + 1 
    1310             isendto(nsndto) = jn 
    1311          ELSEIF( sxM <= sxT  .AND.  dxM >= dxT ) THEN 
    1312             nsndto          = nsndto + 1 
    1313             isendto(nsndto) = jn 
    1314          ELSEIF( dxM <  dxT  .AND.  sxT <  dxM ) THEN 
    1315             nsndto          = nsndto + 1 
    1316             isendto(nsndto) = jn 
    1317          ENDIF 
     1301         zinfo(:,:,:,:) = HUGE(0._wp)   ! default def to make sure we don't use the halos 
     1302         DO jg = 1, 4   ! grid type: T, U, V, F 
     1303            DO jj = nn_hls+1, jpj-nn_hls                ! inner domain (warning do_loop_substitute not yet defined) 
     1304               DO ji = nn_hls+1, jpi-nn_hls             ! inner domain (warning do_loop_substitute not yet defined) 
     1305                  zinfo(ji,jj,1,jg) = REAL(narea, wp)   ! mpi_rank + 1 (as default lbc_lnk fill is 0 
     1306                  zinfo(ji,jj,2,jg) = REAL(ji, wp)      ! ji of this proc 
     1307               END DO 
     1308            END DO 
     1309         END DO 
    13181310         ! 
    1319       END DO 
     1311         ln_nnogather = .FALSE.   ! force "classical" North pole folding to fill all halos -> should be no more HUGE values... 
     1312         CALL lbc_lnk( 'mppini', zinfo(:,:,:,1), 'T', 1._wp )   ! Do 4 calls instead of 1 to save memory as the nogather version 
     1313         CALL lbc_lnk( 'mppini', zinfo(:,:,:,2), 'U', 1._wp )   ! creates buffer arrays with jpiglo as the first dimension 
     1314         CALL lbc_lnk( 'mppini', zinfo(:,:,:,3), 'V', 1._wp )   !  
     1315         CALL lbc_lnk( 'mppini', zinfo(:,:,:,4), 'F', 1._wp )   !  
     1316         ln_nnogather = .TRUE. 
     1317          
     1318         IF( l_IdoNFold ) THEN   ! only the procs involed in the NFD must take care of this 
     1319             
     1320            ALLOCATE( nfd_rksnd(jpi,4), nfd_jisnd(jpi,4) )          ! neighbour rand and remote ji-index for each grid (T, U, V, F) 
     1321            nfd_rksnd(:,:) = NINT( zinfo(:, jpj, 1, :) ) - 1        ! neighbour MPI rank 
     1322            nfd_jisnd(:,:) = NINT( zinfo(:, jpj, 2, :) ) - nn_hls   ! neighbour ji index (shifted as we don't send the halos) 
     1323            WHERE( nfd_rksnd == -1 )   nfd_jisnd = 1                ! use ji=1 if no neighbour, see mpp_nfd_generic.h90 
     1324             
     1325            nfd_nbnei = 1                ! Number of neighbour sending data for the nfd. We have at least 1 neighbour! 
     1326            irknei(1) = nfd_rksnd(1,1)   ! which is the 1st one (I can be neighbour of myself, exclude land-proc are also ok) 
     1327            DO jg = 1, 4 
     1328               DO ji = 1, jpi     ! we must be able to fill the full line including halos 
     1329                  lnew = .TRUE.   ! new neighbour? 
     1330                  DO jn = 1, nfd_nbnei 
     1331                     IF( irknei(jn) == nfd_rksnd(ji,jg) )   lnew = .FALSE.   ! already found 
     1332                  END DO 
     1333                  IF( lnew ) THEN 
     1334                     jn = nfd_nbnei + 1 
     1335                     nfd_nbnei = jn 
     1336                     irknei(jn) = nfd_rksnd(ji,jg) 
     1337                  ENDIF 
     1338               END DO 
     1339            END DO 
     1340             
     1341            ALLOCATE( nfd_rknei(nfd_nbnei) ) 
     1342            nfd_rknei(:) = irknei(1:nfd_nbnei) 
     1343            ! re-number nfd_rksnd according to the indexes of nfd_rknei 
     1344            DO jn = 1, nfd_nbnei 
     1345               WHERE( nfd_rksnd == nfd_rknei(jn) )   nfd_rksnd = jn 
     1346            END DO 
     1347             
     1348            IF( ldwrtlay ) THEN 
     1349               WRITE(knum,*) 
     1350               WRITE(knum,*) 'north fold exchanges with explicit point-to-point messaging :' 
     1351               WRITE(knum,*) '   number of neighbours for the NF: nfd_nbnei : ', nfd_nbnei 
     1352               IF( nfd_nbnei > 0 )   WRITE(knum,*) '   neighbours MPI ranks                       : ', nfd_rknei 
     1353            ENDIF 
     1354             
     1355         ENDIF   ! l_IdoNFold 
     1356         ! 
     1357      ENDIF   ! ln_nnogather 
    13201358      ! 
    13211359   END SUBROUTINE init_nfdcom 
     
    13421380   END SUBROUTINE init_doloop 
    13431381 
     1382    
     1383   SUBROUTINE init_locglo 
     1384      !!---------------------------------------------------------------------- 
     1385      !!                     ***  ROUTINE init_locglo  *** 
     1386      !! 
     1387      !! ** Purpose :   initialization of global domain <--> local domain indices 
     1388      !! 
     1389      !! ** Method  : 
     1390      !! 
     1391      !! ** Action  : - mig , mjg : local  domain indices ==> global domain, including halos, indices 
     1392      !!              - mig0, mjg0: local  domain indices ==> global domain, excluding halos, indices 
     1393      !!              - mi0 , mi1 : global domain indices ==> local  domain indices 
     1394      !!              - mj0 , mj1   (if global point not in the local domain ==> mi0>mi1 and/or mj0>mj1) 
     1395      !!---------------------------------------------------------------------- 
     1396      INTEGER ::   ji, jj   ! dummy loop argument 
     1397      !!---------------------------------------------------------------------- 
     1398      ! 
     1399      ALLOCATE( mig(jpi), mjg(jpj), mig0(jpi), mjg0(jpj) ) 
     1400      ALLOCATE( mi0(jpiglo), mi1(jpiglo), mj0(jpjglo), mj1(jpjglo) ) 
     1401      ! 
     1402      DO ji = 1, jpi                 ! local domain indices ==> global domain indices, including halos 
     1403        mig(ji) = ji + nimpp - 1 
     1404      END DO 
     1405      DO jj = 1, jpj 
     1406        mjg(jj) = jj + njmpp - 1 
     1407      END DO 
     1408      !                              ! local domain indices ==> global domain indices, excluding halos 
     1409      ! 
     1410      mig0(:) = mig(:) - nn_hls 
     1411      mjg0(:) = mjg(:) - nn_hls 
     1412      !                              ! global domain, including halos, indices ==> local domain indices 
     1413      !                                   ! (return (m.0,m.1)=(1,0) if data domain gridpoint is to the west/south of the 
     1414      !                                   ! local domain, or (m.0,m.1)=(jp.+1,jp.) to the east/north of local domain. 
     1415      DO ji = 1, jpiglo 
     1416        mi0(ji) = MAX( 1 , MIN( ji - nimpp + 1, jpi+1 ) ) 
     1417        mi1(ji) = MAX( 0 , MIN( ji - nimpp + 1, jpi   ) ) 
     1418      END DO 
     1419      DO jj = 1, jpjglo 
     1420        mj0(jj) = MAX( 1 , MIN( jj - njmpp + 1, jpj+1 ) ) 
     1421        mj1(jj) = MAX( 0 , MIN( jj - njmpp + 1, jpj   ) ) 
     1422      END DO 
     1423      ! 
     1424   END SUBROUTINE init_locglo 
     1425    
    13441426   !!====================================================================== 
    13451427END MODULE mppini 
  • NEMO/trunk/src/OCE/nemogcm.F90

    r15023 r15267  
    7373   USE lib_mpp        ! distributed memory computing 
    7474   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
    75    USE lbcnfd  , ONLY : isendto, nsndto  ! Setup of north fold exchanges 
    7675   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    7776   USE halo_mng       ! halo manager 
  • NEMO/trunk/src/OFF/nemogcm.F90

    r15023 r15267  
    6363   USE timing         ! Timing 
    6464   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    65    USE lbcnfd  , ONLY : isendto, nsndto   ! Setup of north fold exchanges 
    6665#if defined key_qco 
    6766   USE stpmlf , ONLY : Nbb, Nnn, Naa, Nrhs   ! time level indices 
  • NEMO/trunk/src/SAO/nemogcm.F90

    r14433 r15267  
    3333   USE lib_mpp        ! distributed memory computing 
    3434   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
    35    USE lbcnfd  , ONLY : isendto, nsndto ! Setup of north fold exchanges  
    3635   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    3736#if defined key_xios 
  • NEMO/trunk/src/SAS/nemogcm.F90

    r14834 r15267  
    4040   USE lib_mpp        ! distributed memory computing 
    4141   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
    42    USE lbcnfd  , ONLY : isendto, nsndto ! Setup of north fold exchanges 
    4342   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    4443#if defined key_xios 
  • NEMO/trunk/src/SWE/nemogcm.F90

    r15249 r15267  
    3434   USE lib_mpp        ! distributed memory computing 
    3535   USE mppini         ! shared/distributed memory setting (mpp_init routine) 
    36    USE lbcnfd  , ONLY : isendto, nsndto  ! Setup of north fold exchanges  
    3736   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    3837   USE halo_mng       ! halo manager 
  • NEMO/trunk/tests/BENCH/MY_SRC/usrdef_zgr.F90

    r14433 r15267  
    198198       
    199199!!$      IF( c_NFtype == 'T' ) THEN   ! add a small island in the upper corners to avoid model instabilities... 
    200 !!$         z2d(mi0(       nn_hls):mi1(                  nn_hls+2 ),mj0(jpjglo-nn_hls-1):mj1(jpjglo-nn_hls+1)) = 0. 
    201 !!$         z2d(mi0(jpiglo-nn_hls):mi1(MIN(jpiglo,jpiglo-nn_hls+2)),mj0(jpjglo-nn_hls-1):mj1(jpjglo-nn_hls+1)) = 0. 
    202 !!$         z2d(mi0(jpiglo/2     ):mi1(           jpiglo/2     +2 ),mj0(jpjglo-nn_hls-1):mj1(jpjglo-nn_hls+1)) = 0. 
     200!!$         z2d(mi0(       nn_hls):mi1(                  nn_hls+2 ),mj0(jpjglo-nn_hls-1):mj1(jpjglo-nn_hls+1)) = 0._wp 
     201!!$         z2d(mi0(jpiglo-nn_hls):mi1(MIN(jpiglo,jpiglo-nn_hls+2)),mj0(jpjglo-nn_hls-1):mj1(jpjglo-nn_hls+1)) = 0._wp 
     202!!$         z2d(mi0(jpiglo/2     ):mi1(           jpiglo/2     +2 ),mj0(jpjglo-nn_hls-1):mj1(jpjglo-nn_hls+1)) = 0._wp 
    203203!!$      ENDIF 
    204204!!$      ! 
    205 !!$      IF( c_NFtype == 'F' ) THEN   ! add a small island in the upper corners to avoid model instabilities... 
    206 !!$         z2d(mi0(       nn_hls):mi1(       nn_hls+1),mj0(jpjglo-nn_hls):mj1(jpjglo-nn_hls+1)) = 0. 
    207 !!$         z2d(mi0(jpiglo-nn_hls):mi1(jpiglo-nn_hls+1),mj0(jpjglo-nn_hls):mj1(jpjglo-nn_hls+1)) = 0. 
    208 !!$         z2d(mi0(jpiglo/2     ):mi1(jpiglo/2     +1),mj0(jpjglo-nn_hls):mj1(jpjglo-nn_hls+1)) = 0. 
    209 !!$      ENDIF 
    210  
    211       ! 
    212       CALL lbc_lnk( 'usrdef_zgr', z2d, 'T', 1. )           ! set surrounding land to zero (closed boundaries) 
     205      IF( c_NFtype == 'F' ) THEN   ! Must mask the 2 pivot-points  
     206         z2d(mi0(nn_hls+1):mi1(nn_hls+1),mj0(jpjglo-nn_hls):mj1(jpjglo-nn_hls)) = 0._wp 
     207         z2d(mi0(jpiglo/2):mi1(jpiglo/2),mj0(jpjglo-nn_hls):mj1(jpjglo-nn_hls)) = 0._wp 
     208      ENDIF 
     209      ! 
     210      CALL lbc_lnk( 'usrdef_zgr', z2d, 'T', 1._wp )           ! set surrounding land to zero (closed boundaries) 
    213211      ! 
    214212      k_bot(:,:) = INT( z2d(:,:) )           ! =jpkm1 over the ocean point, =0 elsewhere 
Note: See TracChangeset for help on using the changeset viewer.