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 12979 for NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src – NEMO

Ignore:
Timestamp:
2020-05-27T14:38:32+02:00 (4 years ago)
Author:
hadcv
Message:

Replace references to tile index arrays with scalars

Location:
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/DOM/dom_oce.F90

    r12958 r12979  
    7777   LOGICAL, PUBLIC ::   ln_tile 
    7878   INTEGER         ::   nn_ltile_i, nn_ltile_j 
     79 
     80   ! Domain tiling (all tiles) 
     81   INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntsi_a       !: start of internal part of tile domain 
     82   INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntsj_a       ! 
     83   INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntei_a       !: end of internal part of tile domain 
     84   INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntej_a       ! 
    7985 
    8086   !                                 !  domain MPP decomposition parameters 
  • NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/DOM/domain.F90

    r12879 r12979  
    5050   PUBLIC   dom_init     ! called by nemogcm.F90 
    5151   PUBLIC   domain_cfg   ! called by nemogcm.F90 
     52   PUBLIC   dom_tile     ! called by step.F90 
    5253 
    5354   !!------------------------------------------------------------------------- 
     
    120121      !           !==  Reference coordinate system  ==! 
    121122      ! 
    122       CALL dom_glo                     ! global domain versus local domain 
    123       CALL dom_nam                     ! read namelist ( namrun, namdom ) 
    124       CALL dom_tile                    ! Tile domains 
     123      CALL dom_glo                            ! global domain versus local domain 
     124      CALL dom_nam                            ! read namelist ( namrun, namdom ) 
     125      CALL dom_tile( ntsi, ntsj, ntei, ntej ) ! Tile domain 
    125126 
    126127      ! 
     
    272273 
    273274 
    274    SUBROUTINE dom_tile 
     275   SUBROUTINE dom_tile( ktsi, ktsj, ktei, ktej, ktile ) 
    275276      !!---------------------------------------------------------------------- 
    276277      !!                     ***  ROUTINE dom_tile  *** 
     
    278279      !! ** Purpose :   Set tile domain variables 
    279280      !! 
    280       !! ** Action  : - ntsi, ntsj     : start of internal part of domain 
    281       !!              - ntei, ntej     : end of internal part of domain 
     281      !! ** Action  : - ktsi, ktsj     : start of internal part of domain 
     282      !!              - ktei, ktej     : end of internal part of domain 
     283      !!              - ntile          : current tile number 
    282284      !!              - nijtile        : total number of tiles 
    283285      !!---------------------------------------------------------------------- 
    284       INTEGER ::   jt               ! dummy loop argument 
    285       INTEGER ::   iitile, ijtile   ! Local integers 
    286       !!---------------------------------------------------------------------- 
    287       ntile = 0                     ! Initialise to full domain indices 
    288  
    289       IF( ln_tile ) THEN            ! Set tile decomposition 
    290          iitile = (jpi - 2 * nn_hls) / nn_ltile_i 
    291          ijtile = (jpj - 2 * nn_hls) / nn_ltile_j 
    292          IF( MOD( jpi - 2 * nn_hls, nn_ltile_i ) /= 0 ) iitile = iitile + 1 
    293          IF( MOD( jpj - 2 * nn_hls, nn_ltile_j ) /= 0 ) ijtile = ijtile + 1 
    294  
    295          nijtile = iitile * ijtile 
    296          ALLOCATE( ntsi(0:nijtile), ntsj(0:nijtile), ntei(0:nijtile), ntej(0:nijtile) ) 
     286      INTEGER, INTENT(out) :: ktsi, ktsj, ktei, ktej      ! Tile domain indices 
     287      INTEGER, INTENT(in), OPTIONAL :: ktile              ! Tile number 
     288      INTEGER ::   jt                                     ! dummy loop argument 
     289      INTEGER ::   iitile, ijtile                         ! Local integers 
     290      !!---------------------------------------------------------------------- 
     291      IF( PRESENT(ktile) .AND. ln_tile ) THEN 
     292         ntile = ktile                 ! Set domain indices for tile 
     293         ktsi = ntsi_a(ktile) 
     294         ktsj = ntsj_a(ktile) 
     295         ktei = ntei_a(ktile) 
     296         ktej = ntej_a(ktile) 
    297297      ELSE 
     298         ntile = 0                     ! Initialise to full domain 
    298299         nijtile = 1 
    299          ALLOCATE( ntsi(0:0), ntsj(0:0), ntei(0:0), ntej(0:0) ) 
    300       ENDIF 
    301  
    302       ntsi(0) = 1 + nn_hls          ! Full domain 
    303       ntsj(0) = 1 + nn_hls 
    304       ntei(0) = jpi - nn_hls 
    305       ntej(0) = jpj - nn_hls 
    306  
    307       IF( ln_tile ) THEN            ! Tile domains 
    308          DO jt = 1, nijtile 
    309             ntsi(jt) = ntsi(0) + nn_ltile_i * MOD(jt - 1, iitile) 
    310             ntsj(jt) = ntsj(0) + nn_ltile_j * ((jt - 1) / iitile) 
    311             ntei(jt) = MIN(ntsi(jt) + nn_ltile_i - 1, ntei(0)) 
    312             ntej(jt) = MIN(ntsj(jt) + nn_ltile_j - 1, ntej(0)) 
    313          ENDDO 
    314       ENDIF 
    315  
    316       IF(lwp) THEN                  ! control print 
    317          WRITE(numout,*) 
    318          WRITE(numout,*) 'dom_tile : Domain tiling decomposition' 
    319          WRITE(numout,*) '~~~~~~~~' 
    320          IF( ln_tile ) THEN 
    321             WRITE(numout,*) iitile, 'tiles in i' 
    322             WRITE(numout,*) '    Starting indices' 
    323             WRITE(numout,*) '        ', (ntsi(jt), jt=1, iitile) 
    324             WRITE(numout,*) '    Ending indices' 
    325             WRITE(numout,*) '        ', (ntei(jt), jt=1, iitile) 
    326             WRITE(numout,*) ijtile, 'tiles in j' 
    327             WRITE(numout,*) '    Starting indices' 
    328             WRITE(numout,*) '        ', (ntsj(jt), jt=1, nijtile, iitile) 
    329             WRITE(numout,*) '    Ending indices' 
    330             WRITE(numout,*) '        ', (ntej(jt), jt=1, nijtile, iitile) 
    331          ELSE 
    332             WRITE(numout,*) 'No domain tiling' 
    333             WRITE(numout,*) '    i indices =', ntsi(0), ':', ntei(0) 
    334             WRITE(numout,*) '    j indices =', ntsj(0), ':', ntej(0) 
     300         ktsi = 1 + nn_hls 
     301         ktsj = 1 + nn_hls 
     302         ktei = jpi - nn_hls 
     303         ktej = jpj - nn_hls 
     304 
     305         IF( ln_tile ) THEN            ! Calculate tile domain indices 
     306            iitile = (jpi - 2 * nn_hls) / nn_ltile_i       ! Number of tiles 
     307            ijtile = (jpj - 2 * nn_hls) / nn_ltile_j 
     308            IF( MOD( jpi - 2 * nn_hls, nn_ltile_i ) /= 0 ) iitile = iitile + 1 
     309            IF( MOD( jpj - 2 * nn_hls, nn_ltile_j ) /= 0 ) ijtile = ijtile + 1 
     310 
     311            nijtile = iitile * ijtile 
     312            ALLOCATE( ntsi_a(0:nijtile), ntsj_a(0:nijtile), ntei_a(0:nijtile), ntej_a(0:nijtile) ) 
     313 
     314            ntsi_a(0) = ktsi                               ! Full domain 
     315            ntsj_a(0) = ktsj 
     316            ntei_a(0) = ktei 
     317            ntej_a(0) = ktej 
     318 
     319            DO jt = 1, nijtile                             ! Tile domains 
     320               ntsi_a(jt) = ktsi + nn_ltile_i * MOD(jt - 1, iitile) 
     321               ntsj_a(jt) = ktsj + nn_ltile_j * ((jt - 1) / iitile) 
     322               ntei_a(jt) = MIN(ntsi_a(jt) + nn_ltile_i - 1, ktei) 
     323               ntej_a(jt) = MIN(ntsj_a(jt) + nn_ltile_j - 1, ktej) 
     324            ENDDO 
     325         ENDIF 
     326 
     327         IF(lwp) THEN                  ! control print 
     328            WRITE(numout,*) 
     329            WRITE(numout,*) 'dom_tile : Domain tiling decomposition' 
     330            WRITE(numout,*) '~~~~~~~~' 
     331            IF( ln_tile ) THEN 
     332               WRITE(numout,*) iitile, 'tiles in i' 
     333               WRITE(numout,*) '    Starting indices' 
     334               WRITE(numout,*) '        ', (ntsi_a(jt), jt=1, iitile) 
     335               WRITE(numout,*) '    Ending indices' 
     336               WRITE(numout,*) '        ', (ntei_a(jt), jt=1, iitile) 
     337               WRITE(numout,*) ijtile, 'tiles in j' 
     338               WRITE(numout,*) '    Starting indices' 
     339               WRITE(numout,*) '        ', (ntsj_a(jt), jt=1, nijtile, iitile) 
     340               WRITE(numout,*) '    Ending indices' 
     341               WRITE(numout,*) '        ', (ntej_a(jt), jt=1, nijtile, iitile) 
     342            ELSE 
     343               WRITE(numout,*) 'No domain tiling' 
     344               WRITE(numout,*) '    i indices =', ktsi, ':', ktei 
     345               WRITE(numout,*) '    j indices =', ktsj, ':', ktej 
     346            ENDIF 
    335347         ENDIF 
    336348      ENDIF 
  • NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/do_loop_substitute.h90

    r12879 r12979  
    5555! 
    5656#endif 
    57 #define __kIs_     ntsi(ntile) 
    58 #define __kJs_     ntsj(ntile) 
     57#define __kIs_     ntsi 
     58#define __kJs_     ntsj 
    5959#define __kIsm1_   __kIs_ - nn_hls 
    6060#define __kJsm1_   __kJs_ - nn_hls 
    6161 
    62 #define __kIe_     ntei(ntile) 
    63 #define __kJe_     ntej(ntile) 
     62#define __kIe_     ntei 
     63#define __kJe_     ntej 
    6464#define __kIep1_   __kIe_ + nn_hls 
    6565#define __kJep1_   __kJe_ + nn_hls 
  • NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/par_oce.F90

    r12879 r12979  
    6363 
    6464   ! Domain tiling 
    65    INTEGER, PUBLIC                                  ::   nijtile    !: number of tiles in total 
    66    INTEGER, PUBLIC                                  ::   ntile      !: current tile number 
    67    INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntsi       !: start of internal part of tile domain 
    68    INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntsj       ! 
    69    INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntei       !: end of internal part of tile domain 
    70    INTEGER, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) ::   ntej       ! 
     65   INTEGER, PUBLIC ::   nijtile    !: number of tiles in total 
     66   INTEGER, PUBLIC ::   ntile      !: current tile number 
     67   INTEGER, PUBLIC ::   ntsi       !: start of internal part of tile domain 
     68   INTEGER, PUBLIC ::   ntsj       ! 
     69   INTEGER, PUBLIC ::   ntei       !: end of internal part of tile domain 
     70   INTEGER, PUBLIC ::   ntej       ! 
    7171 
    7272   !!--------------------------------------------------------------------- 
  • NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/step.F90

    r12958 r12979  
    264264      ! Loop over tile domains 
    265265      DO jtile = 1, nijtile 
    266          IF( ln_tile ) ntile = jtile 
     266         IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = jtile ) 
    267267                         CALL tra_ldf    ( kstp, Nbb, Nnn, ts, Nrhs )  ! lateral mixing 
    268268      END DO 
    269       IF( ln_tile ) ntile = 0                                          ! Revert to tile over full domain 
     269      IF( ln_tile ) CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 ) ! Revert to tile over full domain 
    270270                         CALL tra_zdf    ( kstp, Nbb, Nnn, Nrhs, ts, Naa  )  ! vertical mixing and after tracer fields 
    271271      IF( ln_zdfnpc  )   CALL tra_npc    ( kstp,      Nnn, Nrhs, ts, Naa  )  ! update after fields by non-penetrative convection 
  • NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/step_oce.F90

    r12879 r12979  
    99   USE oce             ! ocean dynamics and tracers variables 
    1010   USE dom_oce         ! ocean space and time domain variables 
     11   USE domain, ONLY : dom_tile 
    1112   USE zdf_oce         ! ocean vertical physics variables 
    1213   USE zdfdrg  ,  ONLY : ln_drgimp   ! implicit top/bottom friction 
Note: See TracChangeset for help on using the changeset viewer.