MODULE domtile !!====================================================================== !! *** MODULE domtile *** !! Tiling utilities !!====================================================================== !! History : 4.2 ! 2020-12 (D. Calvert) Original code !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! dom_tile : Set/initialise the current tile and domain indices !!---------------------------------------------------------------------- USE dom_oce ! ocean space and time domain ! USE prtctl ! Print control (prt_ctl_info routine) USE in_out_manager ! I/O manager IMPLICIT NONE PRIVATE PUBLIC dom_tile ! called by step.F90 !!---------------------------------------------------------------------- !! NEMO/OCE 4.2 , NEMO Consortium (2020) !! $Id: domtile.F90 13982 2020-12-04 10:57:05Z hadcv $ !! Software governed by the CeCILL license (see ./LICENSE) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE dom_tile( ktsi, ktsj, ktei, ktej, ktile ) !!---------------------------------------------------------------------- !! *** ROUTINE dom_tile *** !! !! ** Purpose : Set tile domain variables !! !! ** Action : - ktsi, ktsj : start of internal part of domain !! - ktei, ktej : end of internal part of domain !! - ntile : current tile number !! - nijtile : total number of tiles !!---------------------------------------------------------------------- INTEGER, INTENT(out) :: ktsi, ktsj, ktei, ktej ! Tile domain indices INTEGER, INTENT(in), OPTIONAL :: ktile ! Tile number INTEGER :: jt ! dummy loop argument INTEGER :: iitile, ijtile ! Local integers CHARACTER (len=11) :: charout !!---------------------------------------------------------------------- IF( PRESENT(ktile) .AND. ln_tile ) THEN ntile = ktile ! Set domain indices for tile ktsi = ntsi_a(ktile) ktsj = ntsj_a(ktile) ktei = ntei_a(ktile) ktej = ntej_a(ktile) IF(sn_cfctl%l_prtctl) THEN WRITE(charout, FMT="('ntile =', I4)") ktile CALL prt_ctl_info( charout ) ENDIF ELSE ntile = 0 ! Initialise to full domain nijtile = 1 ktsi = Nis0 ktsj = Njs0 ktei = Nie0 ktej = Nje0 IF( ln_tile ) THEN ! Calculate tile domain indices iitile = Ni_0 / nn_ltile_i ! Number of tiles ijtile = Nj_0 / nn_ltile_j IF( MOD( Ni_0, nn_ltile_i ) /= 0 ) iitile = iitile + 1 IF( MOD( Nj_0, nn_ltile_j ) /= 0 ) ijtile = ijtile + 1 nijtile = iitile * ijtile ALLOCATE( ntsi_a(0:nijtile), ntsj_a(0:nijtile), ntei_a(0:nijtile), ntej_a(0:nijtile) ) ntsi_a(0) = ktsi ! Full domain ntsj_a(0) = ktsj ntei_a(0) = ktei ntej_a(0) = ktej DO jt = 1, nijtile ! Tile domains ntsi_a(jt) = Nis0 + nn_ltile_i * MOD(jt - 1, iitile) ntsj_a(jt) = Njs0 + nn_ltile_j * ((jt - 1) / iitile) ntei_a(jt) = MIN(ntsi_a(jt) + nn_ltile_i - 1, Nie0) ntej_a(jt) = MIN(ntsj_a(jt) + nn_ltile_j - 1, Nje0) ENDDO ENDIF IF(lwp) THEN ! control print WRITE(numout,*) WRITE(numout,*) 'dom_tile : Domain tiling decomposition' WRITE(numout,*) '~~~~~~~~' IF( ln_tile ) THEN WRITE(numout,*) iitile, 'tiles in i' WRITE(numout,*) ' Starting indices' WRITE(numout,*) ' ', (ntsi_a(jt), jt=1, iitile) WRITE(numout,*) ' Ending indices' WRITE(numout,*) ' ', (ntei_a(jt), jt=1, iitile) WRITE(numout,*) ijtile, 'tiles in j' WRITE(numout,*) ' Starting indices' WRITE(numout,*) ' ', (ntsj_a(jt), jt=1, nijtile, iitile) WRITE(numout,*) ' Ending indices' WRITE(numout,*) ' ', (ntej_a(jt), jt=1, nijtile, iitile) ELSE WRITE(numout,*) 'No domain tiling' WRITE(numout,*) ' i indices =', ktsi, ':', ktei WRITE(numout,*) ' j indices =', ktsj, ':', ktej ENDIF ENDIF ENDIF END SUBROUTINE dom_tile !!====================================================================== END MODULE domtile