Changeset 12879
- Timestamp:
- 2020-05-06T14:38:01+02:00 (3 years ago)
- Location:
- NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public
- Files:
-
- 9 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/cfgs/SHARED/namelist_ref
r12765 r12879 99 99 !----------------------------------------------------------------------- 100 100 ln_tile = .false. ! Use tiling (T) or not (F) 101 nn_ tile_i = 10! Length of tiles in i102 nn_ tile_j = 10! Length of tiles in j101 nn_ltile_i = 10 ! Length of tiles in i 102 nn_ltile_j = 10 ! Length of tiles in j 103 103 / 104 104 !----------------------------------------------------------------------- -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/DOM/dom_oce.F90
r12765 r12879 75 75 ! Tiling namelist 76 76 LOGICAL, PUBLIC :: ln_tile 77 INTEGER :: nn_ tile_i, nn_tile_j77 INTEGER :: nn_ltile_i, nn_ltile_j 78 78 79 79 ! ! domain MPP decomposition parameters -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/DOM/domain.F90
r12765 r12879 49 49 50 50 PUBLIC dom_init ! called by nemogcm.F90 51 PUBLIC dom_tile ! called by step.F9052 51 PUBLIC domain_cfg ! called by nemogcm.F90 53 52 … … 123 122 CALL dom_glo ! global domain versus local domain 124 123 CALL dom_nam ! read namelist ( namrun, namdom ) 125 126 ! Initialise tile to full domain 127 CALL dom_tile(0) 124 CALL dom_tile ! Tile domains 128 125 129 126 ! … … 275 272 276 273 277 SUBROUTINE dom_tile (kntile)274 SUBROUTINE dom_tile 278 275 !!---------------------------------------------------------------------- 279 276 !! *** ROUTINE dom_tile *** 280 277 !! 281 !! ** Purpose : Set domain indices for specified tile 282 !! 283 !! ** Action : - ntile : current tile number 284 !! - ntsi, ntsj : start of internal part of domain 278 !! ** Purpose : Set tile domain variables 279 !! 280 !! ** Action : - ntsi, ntsj : start of internal part of domain 285 281 !! - ntei, ntej : end of internal part of domain 286 !! - ntsim1, ntsjm1 : start of domain 287 !! - nteip1, ntejp1 : end of domain 288 !!---------------------------------------------------------------------- 289 INTEGER , INTENT(in ) :: kntile ! Tile number 290 INTEGER :: iitile, ijtile ! Tile number in i and j 291 !!---------------------------------------------------------------------- 292 293 IF( ln_tile .AND. kntile > 0 ) THEN ! Tile domain 294 iitile = 1 + MOD( kntile - 1, jpnitile ) 295 ijtile = 1 + (kntile - 1) / jpnitile 296 297 ntile = kntile 298 ntsi = 2 + (iitile - 1) * nn_tile_i 299 ntsj = 2 + (ijtile - 1) * nn_tile_j 300 ntei = MIN(ntsi + nn_tile_i - 1, jpim1) ! Size of last tile limited by full domain 301 ntej = MIN(ntsj + nn_tile_j - 1, jpjm1) ! 302 ntsim1 = ntsi - 1 303 ntsjm1 = ntsj - 1 304 nteip1 = ntei + 1 305 ntejp1 = ntej + 1 306 ELSE ! Full domain 307 ntile = 1 308 ntsi = 2 309 ntsj = 2 310 ntei = jpim1 311 ntej = jpjm1 312 ntsim1 = 1 313 ntsjm1 = 1 314 nteip1 = jpi 315 ntejp1 = jpj 282 !! - nijtile : total number of tiles 283 !!---------------------------------------------------------------------- 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) ) 297 ELSE 298 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) 335 ENDIF 316 336 ENDIF 317 337 END SUBROUTINE dom_tile … … 339 359 & ln_cfmeta, ln_xios_read, nn_wxios 340 360 NAMELIST/namdom/ ln_linssh, rn_Dt, rn_atfp, ln_crs, ln_meshmask 341 NAMELIST/namtile/ ln_tile, nn_ tile_i, nn_tile_j361 NAMELIST/namtile/ ln_tile, nn_ltile_i, nn_ltile_j 342 362 #if defined key_netcdf4 343 363 NAMELIST/namnc4/ nn_nchunks_i, nn_nchunks_j, nn_nchunks_k, ln_nc4zip … … 473 493 IF(lwm) WRITE( numond, namtile ) 474 494 475 ! Set tile decomposition476 IF( ln_tile ) THEN477 jpnitile = (jpi - 2) / nn_tile_i478 jpnjtile = (jpj - 2) / nn_tile_j479 IF( MOD( jpi - 2, nn_tile_i ) /= 0 ) jpnitile = jpnitile + 1480 IF( MOD( jpj - 2, nn_tile_j ) /= 0 ) jpnjtile = jpnjtile + 1481 ELSE482 jpnitile = 1483 jpnjtile = 1484 ENDIF485 jpnijtile = jpnitile * jpnjtile486 487 495 IF(lwp) THEN 488 496 WRITE(numout,*) 489 WRITE(numout,*) ' Namelist : namtile --- tiling decomposition'497 WRITE(numout,*) ' Namelist : namtile --- Domain tiling decomposition' 490 498 WRITE(numout,*) ' Tiling (T) or not (F) ln_tile = ', ln_tile 491 WRITE(numout,*) ' Length of tile in i nn_ tile_i = ', nn_tile_i492 WRITE(numout,*) ' Length of tile in j nn_ tile_j = ', nn_tile_j499 WRITE(numout,*) ' Length of tile in i nn_ltile_i = ', nn_ltile_i 500 WRITE(numout,*) ' Length of tile in j nn_ltile_j = ', nn_ltile_j 493 501 WRITE(numout,*) 494 502 IF( ln_tile ) THEN 495 WRITE(numout,*) ' The domain will be decomposed into ', jpnijtile, 'tiles of size', nn_tile_i, 'x', nn_tile_j503 WRITE(numout,*) ' The domain will be decomposed into tiles of size', nn_ltile_i, 'x', nn_ltile_j 496 504 ELSE 497 505 WRITE(numout,*) ' Domain tiling will NOT be used' -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/TRA/traldf.F90
r12765 r12879 63 63 ENDIF 64 64 ! 65 IF( ntile == jpnijtile ) THEN! Do only after all tiles finish65 IF( ntile == nijtile ) THEN ! Do only after all tiles finish 66 66 IF( l_trdtra ) THEN !* Save ta and sa trends 67 67 ! TODO: TO BE TILED … … 83 83 END SELECT 84 84 ! 85 IF( ntile == jpnijtile ) THEN! Do only after all tiles finish85 IF( ntile == nijtile ) THEN ! Do only after all tiles finish 86 86 IF( l_trdtra ) THEN !* save the horizontal diffusive trends for further diagnostics 87 87 ! TODO: TO BE TILED -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/TRA/traldf_iso.F90
r12765 r12879 109 109 REAL(wp) :: zmskv, zahv_w, zabe2, zcof2, zcoef4 ! - - 110 110 REAL(wp) :: zcoef0, ze3w_2, zsign ! - - 111 REAL(wp), DIMENSION( IND_2D) :: zdkt, zdk1t, z2d112 REAL(wp), DIMENSION( IND_2D,jpk) :: zdit, zdjt, zftu, zftv, ztfw111 REAL(wp), DIMENSION(A2D) :: zdkt, zdk1t, z2d 112 REAL(wp), DIMENSION(A2D,jpk) :: zdit, zdjt, zftu, zftv, ztfw 113 113 !!---------------------------------------------------------------------- 114 114 ! … … 321 321 END_3D 322 322 ! 323 IF( ntile == jpnijtile ) THEN! Do only after all tiles finish323 IF( ntile == nijtile ) THEN ! Do only after all tiles finish 324 324 IF( ( kpass == 1 .AND. ln_traldf_lap ) .OR. & !== first pass only ( laplacian) ==! 325 325 ( kpass == 2 .AND. ln_traldf_blp ) ) THEN !== 2nd pass (bilaplacian) ==! -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/do_loop_substitute.h90
r12765 r12879 55 55 ! 56 56 #endif 57 #define __kIs_ ntsi 58 #define __kJs_ ntsj 59 #define __kIsm1_ ntsim160 #define __kJsm1_ ntsjm157 #define __kIs_ ntsi(ntile) 58 #define __kJs_ ntsj(ntile) 59 #define __kIsm1_ __kIs_ - nn_hls 60 #define __kJsm1_ __kJs_ - nn_hls 61 61 62 #define __kIe_ ntei 63 #define __kJe_ ntej 64 #define __kIep1_ nteip165 #define __kJep1_ ntejp162 #define __kIe_ ntei(ntile) 63 #define __kJe_ ntej(ntile) 64 #define __kIep1_ __kIe_ + nn_hls 65 #define __kJep1_ __kJe_ + nn_hls 66 66 67 #define IND_2D__kIsm1_:__kIep1_,__kJsm1_:__kJep1_67 #define A2D __kIsm1_:__kIep1_,__kJsm1_:__kJep1_ 68 68 69 69 #define DO_2D_00_00 DO jj = __kJs_, __kJe_ ; DO ji = __kIs_, __kIe_ -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/par_oce.F90
r12765 r12879 62 62 INTEGER, PUBLIC :: jpjmax! = ( jpjglo-2*nn_hls + (jpnj-1) ) / jpnj + 2*nn_hls !: maximum jpj 63 63 64 ! Tiling decomposition 65 INTEGER, PUBLIC :: jpnitile !: number of tiles following i 66 INTEGER, PUBLIC :: jpnjtile !: number of tiles following j 67 INTEGER, PUBLIC :: jpnijtile !: number of tiles in total (jpnitile x jpnjtile) 68 69 ! Tile indices 70 INTEGER, PUBLIC :: ntsi !: start of internal part of tile domain 71 INTEGER, PUBLIC :: ntsj ! 72 INTEGER, PUBLIC :: ntei !: end of internal part of tile domain 73 INTEGER, PUBLIC :: ntej ! 74 INTEGER, PUBLIC :: ntsim1 !: start of tile domain 75 INTEGER, PUBLIC :: ntsjm1 ! 76 INTEGER, PUBLIC :: nteip1 !: end of tile domain 77 INTEGER, PUBLIC :: ntejp1 ! 78 INTEGER, PUBLIC :: ntile !: current tile number 64 ! 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 ! 79 71 80 72 !!--------------------------------------------------------------------- -
NEMO/branches/UKMO/dev_r12745_HPC-02_Daley_Tiling_trial_public/src/OCE/step.F90
r12765 r12879 266 266 267 267 ! Loop over tile domains 268 DO jtile = 1, jpnijtile269 IF( ln_tile ) CALL dom_tile( jtile )268 DO jtile = 1, nijtile 269 IF( ln_tile ) ntile = jtile 270 270 CALL tra_ldf ( kstp, Nbb, Nnn, ts, Nrhs ) ! lateral mixing 271 271 END DO 272 IF( ln_tile ) CALL dom_tile( 0 ) ! Revert to tile over full domain 273 272 IF( ln_tile ) ntile = 0 ! Revert to tile over full domain 274 273 CALL tra_zdf ( kstp, Nbb, Nnn, Nrhs, ts, Naa ) ! vertical mixing and after tracer fields 275 274 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
r12765 r12879 62 62 USE domvvl ! variable vertical scale factors (dom_vvl_sf_nxt routine) 63 63 ! (dom_vvl_sf_swp routine) 64 USE domain , ONLY : dom_tile65 64 66 65 USE ldfslp ! iso-neutral slopes (ldf_slp routine)
Note: See TracChangeset
for help on using the changeset viewer.