1 !WRF:DRIVER_LAYER:TILING
9 MODULE PROCEDURE set_tiles1 , set_tiles2, set_tiles3, set_tiles_once
14 ! CPP macro for error checking
15 #define ERROR_TEST(A,O,B) IF( A O B )THEN;WRITE(mess,'(3A4)')'A','O','B';CALL WRF_ERROR_FATAL(mess);ENDIF
17 # define MIN_TILE_SIZE 1
20 ! this version is used to compute only on a boundary of some width
21 ! The ids, ide, jds, and jde arguments specify the edge of the boundary (a way of
22 ! accounting for staggering, and the bdyw gives the number of cells
23 ! (idea: if bdyw is negative, have it do the reverse and specify the
24 ! interior, less the boundary.
26 SUBROUTINE set_tiles1 ( grid , ids , ide , jds , jde , bdyw )
28 USE module_domain, ONLY : domain
29 USE module_driver_constants
37 TYPE(domain) , INTENT(INOUT) :: grid
38 INTEGER , INTENT(IN) :: ids , ide , jds , jde , bdyw
42 INTEGER :: spx, epx, spy, epy, t, tt, ts, te
43 INTEGER :: smx, emx, smy, emy
44 INTEGER :: ntiles , num_tiles
48 data_ordering : SELECT CASE ( model_data_order )
49 CASE ( DATA_ORDER_XYZ )
50 spx = grid%sp31 ; epx = grid%ep31 ; spy = grid%sp32 ; epy = grid%ep32
51 CASE ( DATA_ORDER_YXZ )
52 spx = grid%sp32 ; epx = grid%ep32 ; spy = grid%sp31 ; epy = grid%ep31
53 CASE ( DATA_ORDER_ZXY )
54 spx = grid%sp32 ; epx = grid%ep32 ; spy = grid%sp33 ; epy = grid%ep33
55 CASE ( DATA_ORDER_ZYX )
56 spx = grid%sp33 ; epx = grid%ep33 ; spy = grid%sp32 ; epy = grid%ep32
57 CASE ( DATA_ORDER_XZY )
58 spx = grid%sp31 ; epx = grid%ep31 ; spy = grid%sp33 ; epy = grid%ep33
59 CASE ( DATA_ORDER_YZX )
60 spx = grid%sp33 ; epx = grid%ep33 ; spy = grid%sp31 ; epy = grid%ep31
61 END SELECT data_ordering
65 IF ( num_tiles > grid%max_tiles ) THEN
66 IF ( ASSOCIATED(grid%i_start) ) THEN ; DEALLOCATE( grid%i_start ) ; NULLIFY( grid%i_start ) ; ENDIF
67 IF ( ASSOCIATED(grid%i_end) ) THEN ; DEALLOCATE( grid%i_end ) ; NULLIFY( grid%i_end ) ; ENDIF
68 IF ( ASSOCIATED(grid%j_start) ) THEN ; DEALLOCATE( grid%j_start ) ; NULLIFY( grid%j_start ) ; ENDIF
69 IF ( ASSOCIATED(grid%j_end) ) THEN ; DEALLOCATE( grid%j_end ) ; NULLIFY( grid%j_end ) ; ENDIF
70 ALLOCATE(grid%i_start(num_tiles))
71 ALLOCATE(grid%i_end(num_tiles))
72 ALLOCATE(grid%j_start(num_tiles))
73 ALLOCATE(grid%j_end(num_tiles))
74 grid%max_tiles = num_tiles
78 IF ( ids .ge. spx .and. ids .le. epx ) THEN
80 grid%i_end(1) = min( ids+bdyw-1 , epx )
81 grid%j_start(1) = max( spy , jds )
82 grid%j_end(1) = min( epy , jde )
83 ELSEIF ( (ids+bdyw-1) .ge. spx .and. (ids+bdyw-1) .le. epx ) THEN
84 grid%i_start(1) = max( ids , spx )
85 grid%i_end(1) = ids+bdyw-1
86 grid%j_start(1) = max( spy , jds )
87 grid%j_end(1) = min( epy , jde )
96 IF ( ide .ge. spx .and. ide .le. epx ) THEN
97 grid%i_start(2) = max( ide-bdyw+1 , spx )
99 grid%j_start(2) = max( spy , jds )
100 grid%j_end(2) = min( epy , jde )
101 ELSEIF ( (ide-bdyw+1) .ge. spx .and. (ide-bdyw+1) .le. epx ) THEN
102 grid%i_start(2) = ide-bdyw+1
103 grid%i_end(2) = min( ide , epx )
104 grid%j_start(2) = max( spy , jds )
105 grid%j_end(2) = min( epy , jde )
113 ! YS boundary (note that the corners may already be done by XS and XE)
114 IF ( jds .ge. spy .and. jds .le. epy ) THEN
115 grid%j_start(3) = jds
116 grid%j_end(3) = min( jds+bdyw-1 , epy )
117 grid%i_start(3) = max( spx , ids+bdyw )
118 grid%i_end(3) = min( epx , ide-bdyw )
119 ELSEIF ( (jds+bdyw-1) .ge. spy .and. (jds+bdyw-1) .le. epy ) THEN
120 grid%j_start(3) = max( jds , spy )
121 grid%j_end(3) = jds+bdyw-1
122 grid%i_start(3) = max( spx , ids+bdyw )
123 grid%i_end(3) = min( epx , ide-bdyw )
131 ! YE boundary (note that the corners may already be done by XS and XE)
132 IF ( jde .ge. spy .and. jde .le. epy ) THEN
133 grid%j_start(4) = max( jde-bdyw+1 , spy )
135 grid%i_start(4) = max( spx , ids+bdyw )
136 grid%i_end(4) = min( epx , ide-bdyw )
137 ELSEIF ( (jde-bdyw+1) .ge. spy .and. (jde-bdyw+1) .le. epy ) THEN
138 grid%j_start(4) = jde-bdyw+1
139 grid%j_end(4) = min( jde , epy )
140 grid%i_start(4) = max( spx , ids+bdyw )
141 grid%i_end(4) = min( epx , ide-bdyw )
149 grid%num_tiles = num_tiles
152 END SUBROUTINE set_tiles1
154 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
155 ! this version callset set_tiles2 but only if the zone hasn't been cached already
156 ! Up to MAX_TILING_ZONES allowed
158 SUBROUTINE set_tiles_once ( zone, grid , ids , ide , jds , jde , ips , ipe , jps , jpe )
159 USE module_domain, ONLY : domain, MAX_TILING_ZONES
161 TYPE(domain) , INTENT(INOUT) :: grid
162 INTEGER , INTENT(IN) :: zone
163 INTEGER , INTENT(IN) :: ids , ide , jds , jde
164 INTEGER , INTENT(IN) :: ips , ipe , jps , jpe
166 INTEGER num_tiles, num_tiles_x, num_tiles_y
167 IF ( zone .LT. 1 .OR. zone .GT. MAX_TILING_ZONES ) THEN
168 CALL wrf_error_fatal('set_tiles_once: zone out of range, increase MAX_TILE_ZONES in module_domain_type')
170 IF ( .NOT. grid%tiling_latch(zone) ) THEN
171 grid%tiling_latch(zone) = .TRUE.
172 CALL set_tiles2 ( grid, ids, ide, jds, jde, ips, ipe, jps, jpe )
173 num_tiles = grid%num_tiles
174 num_tiles_x = grid%num_tiles_x
175 num_tiles_y = grid%num_tiles_y
176 ALLOCATE(grid%tile_zones(zone)%i_start(num_tiles))
177 ALLOCATE(grid%tile_zones(zone)%i_end(num_tiles))
178 ALLOCATE(grid%tile_zones(zone)%j_start(num_tiles))
179 ALLOCATE(grid%tile_zones(zone)%j_end(num_tiles))
180 grid%tile_zones(zone)%i_start = grid%i_start
181 grid%tile_zones(zone)%i_end = grid%i_end
182 grid%tile_zones(zone)%j_start = grid%j_start
183 grid%tile_zones(zone)%j_end = grid%j_end
184 grid%tile_zones(zone)%num_tiles = num_tiles
185 grid%tile_zones(zone)%num_tiles_x = num_tiles_x
186 grid%tile_zones(zone)%num_tiles_y = num_tiles_y
188 grid%i_start = grid%tile_zones(zone)%i_start
189 grid%i_end = grid%tile_zones(zone)%i_end
190 grid%j_start = grid%tile_zones(zone)%j_start
191 grid%j_end = grid%tile_zones(zone)%j_end
192 grid%num_tiles = grid%tile_zones(zone)%num_tiles
193 grid%num_tiles_x = grid%tile_zones(zone)%num_tiles_x
194 grid%num_tiles_y = grid%tile_zones(zone)%num_tiles_y
196 END SUBROUTINE set_tiles_once
198 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
199 ! this version is used to limit the domain or compute onto halos
200 SUBROUTINE set_tiles2 ( grid , ids , ide , jds , jde , ips , ipe , jps , jpe )
201 USE module_domain, ONLY : domain
202 USE module_driver_constants
210 TYPE(domain) , INTENT(INOUT) :: grid
211 INTEGER , INTENT(IN) :: ids , ide , jds , jde
212 INTEGER , INTENT(IN) :: ips , ipe , jps , jpe
218 INTEGER :: num_tiles_x, num_tiles_y, num_tiles_inc,num_tiles
219 INTEGER :: tile_strategy,tile_strategy_spec
220 INTEGER :: tile_sz_x, tile_sz_y
221 INTEGER :: spx, epx, spy, epy, t, tt, ts, te
222 INTEGER :: smx, emx, smy, emy
227 INTEGER , EXTERNAL :: omp_get_max_threads
229 CHARACTER*255 :: mess
230 CHARACTER*255 :: envval
231 INTEGER :: tnum_tiles, istat
232 LOGICAL :: verbose ! whether to output tile info messages
234 data_ordering : SELECT CASE ( model_data_order )
235 CASE ( DATA_ORDER_XYZ )
236 spx = grid%sp31 ; epx = grid%ep31 ; spy = grid%sp32 ; epy = grid%ep32
237 smx = grid%sm31 ; emx = grid%em31 ; smy = grid%sm32 ; emy = grid%em32
238 CASE ( DATA_ORDER_YXZ )
239 spx = grid%sp32 ; epx = grid%ep32 ; spy = grid%sp31 ; epy = grid%ep31
240 smx = grid%sm32 ; emx = grid%em32 ; smy = grid%sm31 ; emy = grid%em31
241 CASE ( DATA_ORDER_ZXY )
242 spx = grid%sp32 ; epx = grid%ep32 ; spy = grid%sp33 ; epy = grid%ep33
243 smx = grid%sm32 ; emx = grid%em32 ; smy = grid%sm33 ; emy = grid%em33
244 CASE ( DATA_ORDER_ZYX )
245 spx = grid%sp33 ; epx = grid%ep33 ; spy = grid%sp32 ; epy = grid%ep32
246 smx = grid%sm33 ; emx = grid%em33 ; smy = grid%sm32 ; emy = grid%em32
247 CASE ( DATA_ORDER_XZY )
248 spx = grid%sp31 ; epx = grid%ep31 ; spy = grid%sp33 ; epy = grid%ep33
249 smx = grid%sm31 ; emx = grid%em31 ; smy = grid%sm33 ; emy = grid%em33
250 CASE ( DATA_ORDER_YZX )
251 spx = grid%sp33 ; epx = grid%ep33 ; spy = grid%sp31 ; epy = grid%ep31
252 smx = grid%sm33 ; emx = grid%em33 ; smy = grid%sm31 ; emy = grid%em31
253 END SELECT data_ordering
255 ERROR_TEST(ips,<,smx)
256 ERROR_TEST(ipe,>,emx)
257 ERROR_TEST(jps,<,smy)
258 ERROR_TEST(jpe,>,emy)
260 ! Here's how the number of tiles is arrived at:
262 ! if tile sizes are specified use those otherwise
263 ! if num_tiles is specified use that otherwise
264 ! if omp provides a value use that otherwise
269 IF ( grid%num_tiles_spec .EQ. 0 ) THEN
271 CALL nl_get_numtiles( 1, num_tiles )
272 IF ( num_tiles .EQ. 1 ) THEN
274 num_tiles = omp_get_max_threads()
275 WRITE(mess,'("WRF NUMBER OF TILES FROM OMP_GET_MAX_THREADS = ",I3)')num_tiles
276 CALL WRF_MESSAGE ( mess )
280 CALL get_environment_variable("WRF_NUM_TILES",envval, status=istat)
281 IF ( envval .NE. "" .and. istat .eq. 0) THEN
282 READ (envval,*) tnum_tiles
283 IF ( tnum_tiles .GT. 0 ) THEN
285 WRITE(mess,'("WRF NUMBER OF TILES FROM ENV WRF_NUM_TILES = ",I3)')num_tiles
286 CALL WRF_MESSAGE ( mess )
290 ! override num_tiles setting (however gotten) if tile sizes are specified
291 CALL nl_get_tile_sz_x( 1, tile_sz_x )
292 CALL nl_get_tile_sz_y( 1, tile_sz_y )
293 CALL nl_get_tile_strategy( 1, tile_strategy_spec )
294 CALL nl_get_numtiles_inc( 1, num_tiles_inc )
295 CALL nl_get_numtiles_x( 1, num_tiles_x )
296 CALL nl_get_numtiles_y( 1, num_tiles_y )
297 IF ( num_tiles_x .EQ. 0 ) THEN
298 CALL get_environment_variable ("WRF_NUM_TILES_X",envval, status=istat)
299 IF ( envval .NE. "" .and. istat .eq. 0) THEN
300 READ (envval,*) tnum_tiles
301 IF ( tnum_tiles .GT. 0 ) THEN
302 num_tiles_x=tnum_tiles
303 WRITE(mess,'("WRF NUMBER OF TILES X FROM ENV WRF_NUM_TILES_X = ",I3)')num_tiles_x
304 CALL WRF_MESSAGE ( mess )
309 IF ( num_tiles_y .EQ. 0 ) THEN
310 CALL get_environment_variable ("WRF_NUM_TILES_Y",envval, status=istat)
311 IF ( envval .NE. "" .and. istat .eq. 0) THEN
312 READ (envval,*) tnum_tiles
313 IF ( tnum_tiles .GT. 0 ) THEN
314 num_tiles_y=tnum_tiles
315 WRITE(mess,'("WRF NUMBER OF TILES Y FROM ENV WRF_NUM_TILES_Y = ",I3)')num_tiles_y
316 CALL WRF_MESSAGE ( mess )
321 IF ( num_tiles_inc .EQ. 0 ) THEN
322 CALL get_environment_variable ("WRF_NUM_TILES_INC",envval, status=istat)
323 IF ( envval .NE. "" .and. istat .eq. 0) THEN
324 READ (envval,*) tnum_tiles
325 IF ( tnum_tiles .GT. 0 ) THEN
326 num_tiles_inc=tnum_tiles
327 WRITE(mess,'("WRF NUMBER OF TILES INCREMENT FROM ENV WRF_NUM_TILES_INC = ",I3)')num_tiles_inc
328 CALL WRF_MESSAGE ( mess )
332 num_tiles_inc=max( num_tiles_inc , 1)
333 IF ( tile_strategy_spec == TILE_NONE) then
335 WRITE(mess,*)'Tile Strategy is not specified. Assuming 1D-Y'
336 CALL WRF_MESSAGE ( mess )
337 IF ( num_tiles > (epy-spy+1)/MIN_TILE_SIZE .and. num_tiles_x == 0 .and. num_tiles_y == 0) THEN ! number of tiles is too high. Trying to adjust
339 num_tiles_y=(epy-spy+1)/MIN_TILE_SIZE
340 DO WHILE (num_tiles_x*num_tiles_inc*num_tiles_y < num_tiles)
341 num_tiles_x=num_tiles_x+1
343 num_tiles_x=num_tiles_x*num_tiles_inc
344 WRITE(mess,'("Total number of tiles is too big for 1D-Y tiling. Going 2D. New tiling is ",I3,"x",I3)') &
345 num_tiles_x,num_tiles_y
346 CALL WRF_MESSAGE ( mess )
347 tile_strategy=TILE_XY
350 tile_strategy = tile_strategy_spec
353 IF ( tile_sz_x >= 1 .and. tile_sz_y >= 1 ) THEN
354 ! figure number of whole tiles and add 1 for any partials in each dim
355 num_tiles_x = (epx-spx+1) / tile_sz_x
356 if ( tile_sz_x*num_tiles_x < epx-spx+1 ) num_tiles_x = num_tiles_x + 1
357 num_tiles_y = (epy-spy+1) / tile_sz_y
358 if ( tile_sz_y*num_tiles_y < epy-spy+1 ) num_tiles_y = num_tiles_y + 1
359 num_tiles = num_tiles_x * num_tiles_y
361 IF ( num_tiles_x >= 1 .or. num_tiles_y >= 1 ) THEN
362 ! adjust num_tiles_? if several ones are omited
363 IF ( num_tiles_x >= 1 .and. num_tiles_y >= 1 ) THEN
365 num_tiles=num_tiles_x*num_tiles_y;
366 ELSE IF ( num_tiles_x >= 1 ) THEN
367 num_tiles_y=num_tiles/num_tiles_x;
368 ELSE ! IF ( num_tiles_y >= 1 ) THEN
369 num_tiles_x=num_tiles/num_tiles_y;
372 IF ( tile_strategy == TILE_X ) THEN
373 num_tiles_x = num_tiles
375 ELSE IF ( tile_strategy == TILE_Y ) THEN
377 num_tiles_y = num_tiles
378 ELSE ! ( tile_strategy == TILE_XY ) THEN
380 call least_aspect( num_tiles, one, one, num_tiles_y, num_tiles_x )
385 num_tiles=max( num_tiles , 1)
386 num_tiles_x=max( num_tiles_x , 1)
387 num_tiles_y=max( num_tiles_y , 1)
388 grid%num_tiles_spec = num_tiles
389 grid%num_tiles_x = num_tiles_x
390 grid%num_tiles_y = num_tiles_y
393 num_tiles = grid%num_tiles_spec
394 num_tiles_x = grid%num_tiles_x
395 num_tiles_y = grid%num_tiles_y
397 IF ( num_tiles > grid%max_tiles ) THEN
398 IF ( ASSOCIATED(grid%i_start) ) THEN ; DEALLOCATE( grid%i_start ) ; NULLIFY( grid%i_start ) ; ENDIF
399 IF ( ASSOCIATED(grid%i_end) ) THEN ; DEALLOCATE( grid%i_end ) ; NULLIFY( grid%i_end ) ; ENDIF
400 IF ( ASSOCIATED(grid%j_start) ) THEN ; DEALLOCATE( grid%j_start ) ; NULLIFY( grid%j_start ) ; ENDIF
401 IF ( ASSOCIATED(grid%j_end) ) THEN ; DEALLOCATE( grid%j_end ) ; NULLIFY( grid%j_end ) ; ENDIF
402 ALLOCATE(grid%i_start(num_tiles))
403 ALLOCATE(grid%i_end(num_tiles))
404 ALLOCATE(grid%j_start(num_tiles))
405 ALLOCATE(grid%j_end(num_tiles))
406 grid%max_tiles = num_tiles
410 DO t = 0, num_tiles-1
413 ntiles = t / num_tiles_x
414 CALL region_bounds( spy, epy, &
415 num_tiles_y, ntiles, &
417 ! first y (major dimension)
418 IF ( ts .LE. te ) THEN ! converse happens if number of tiles > number of points in dim
420 ! This bit allows the user to specify execution out onto the halo region
421 ! in the call to set_tiles. If the low patch boundary specified by the arguments
422 ! is less than what the model already knows to be the patch boundary and if
423 ! the user hasn't erred by specifying something that would fall off memory
424 ! (safety tests are higher up in this routine, outside the IF) then adjust
425 ! the tile boundary of the low edge tiles accordingly. Likewise for high edges.
426 IF ( jps .lt. spy .and. ts .eq. spy ) ts = jps ;
427 IF ( jpe .gt. epy .and. te .eq. epy ) te = jpe ;
429 grid%j_start(nt) = max ( ts , jds )
430 grid%j_end(nt) = min ( te , jde )
433 ntiles = mod(t,num_tiles_x)
434 CALL region_bounds( spx, epx, &
435 num_tiles_x, ntiles, &
437 IF ( ts .LE. te ) THEN ! converse happens if number of tiles > number of points in dim
438 IF ( ips .lt. spx .and. ts .eq. spx ) ts = ips ;
439 IF ( ipe .gt. epx .and. te .eq. epx ) te = ipe ;
441 grid%i_start(nt) = max ( ts , ids )
442 grid%i_end(nt) = min ( te , ide )
444 WRITE(mess,'("WRF TILE ",I3," IS ",I6," IE ",I6," JS ",I6," JE ",I6)') &
445 nt,grid%i_start(nt),grid%i_end(nt),grid%j_start(nt),grid%j_end(nt)
446 CALL WRF_MESSAGE ( mess )
454 WRITE(mess,'("WRF NUMBER OF TILES = ",I3)')num_tiles
455 CALL WRF_MESSAGE ( mess )
457 grid%num_tiles = num_tiles
460 END SUBROUTINE set_tiles2
462 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
463 ! this version sets the tiles based on a passed in integer mask
464 ! the assumption here is that the mask is relatively straigthforward
465 ! and coverable with 2 or three rectangles. No weird stuff...
467 SUBROUTINE set_tiles3 ( grid , imask, ims, ime, jms, jme, ips, ipe, jps, jpe )
468 USE module_domain, ONLY : domain
469 USE module_driver_constants
477 TYPE(domain) , INTENT(INOUT) :: grid
478 INTEGER , INTENT(IN) :: ims , ime , jms , jme
479 INTEGER , INTENT(IN) :: ips , ipe , jps , jpe
480 INTEGER, DIMENSION(ims:ime,jms:jme), INTENT(IN) :: imask
482 INTEGER, DIMENSION(50) :: i_start, i_end, j_start, j_end
490 CALL set_tiles_masked ( imask, ims, ime, jms, jme, ips, ipe, jps, jpe, &
491 num_tiles, i_start, i_end, j_start, j_end )
493 IF ( num_tiles > grid%max_tiles ) THEN
494 IF ( ASSOCIATED(grid%i_start) ) THEN ; DEALLOCATE( grid%i_start ) ; NULLIFY( grid%i_start ) ; ENDIF
495 IF ( ASSOCIATED(grid%i_end) ) THEN ; DEALLOCATE( grid%i_end ) ; NULLIFY( grid%i_end ) ; ENDIF
496 IF ( ASSOCIATED(grid%j_start) ) THEN ; DEALLOCATE( grid%j_start ) ; NULLIFY( grid%j_start ) ; ENDIF
497 IF ( ASSOCIATED(grid%j_end) ) THEN ; DEALLOCATE( grid%j_end ) ; NULLIFY( grid%j_end ) ; ENDIF
498 ALLOCATE(grid%i_start(num_tiles))
499 ALLOCATE(grid%i_end(num_tiles))
500 ALLOCATE(grid%j_start(num_tiles))
501 ALLOCATE(grid%j_end(num_tiles))
502 grid%max_tiles = num_tiles
504 grid%num_tiles = num_tiles
505 grid%i_start(1:num_tiles) = i_start(1:num_tiles)
506 grid%i_end(1:num_tiles) = i_end(1:num_tiles)
507 grid%j_start(1:num_tiles) = j_start(1:num_tiles)
508 grid%j_end(1:num_tiles) = j_end(1:num_tiles)
510 WRITE(mess,'("WRF TILE ",I3," IS ",I6," IE ",I6," JS ",I6," JE ",I6)') &
511 nt,grid%i_start(nt),grid%i_end(nt),grid%j_start(nt),grid%j_end(nt)
512 CALL wrf_debug ( 1, mess )
514 WRITE(mess,'("set_tiles3: NUMBER OF TILES = ",I3)')num_tiles
515 CALL wrf_debug ( 1, mess )
518 END SUBROUTINE set_tiles3
520 SUBROUTINE set_tiles_masked ( imask, ims, ime, jms, jme, ips, ipe, jps, jpe, &
521 num_tiles, istarts, iends, jstarts, jends )
527 INTEGER , INTENT(IN) :: ims , ime , jms , jme
528 INTEGER, DIMENSION(ims:ime,jms:jme), INTENT(IN) :: imask
529 INTEGER , INTENT(IN) :: ips , ipe , jps , jpe
530 INTEGER , INTENT(OUT) :: num_tiles
531 INTEGER, DIMENSION(*) , INTENT(OUT) :: istarts, iends
532 INTEGER, DIMENSION(*) , INTENT(OUT) :: jstarts, jends
538 INTEGER :: i, j, ir, jr
539 INTEGER :: imaskcopy(ips:ipe,jps:jpe) ! copy of imask to write on
541 imaskcopy = imask(ips:ipe,jps:jpe)
543 ! simple multi-pass scheme, optimize later...
544 DO WHILE (ANY(imaskcopy == 1))
547 ! find first "1" and build a rectangle from it
548 IF ( imaskcopy(i,j) == 1 ) THEN
549 num_tiles = num_tiles + 1
550 istarts(num_tiles) = i
552 jstarts(num_tiles) = j
554 ! don't check this point again
556 ! find length of first row
557 DO ir = istarts(num_tiles)+1,ipe
558 IF ( imaskcopy(ir,j) == 1 ) THEN
559 iends(num_tiles) = ir
560 ! don't check this point again
566 ! find number of rows
567 DO jr = jstarts(num_tiles)+1,jpe
568 IF (ALL(imaskcopy(istarts(num_tiles):iends(num_tiles),jr) == 1)) THEN
569 jends(num_tiles) = jr
570 ! don't check these points again
571 imaskcopy(istarts(num_tiles):iends(num_tiles),jr) = 0
576 ENDIF ! if ( imaskcopy(i,j) == 1 )
581 END SUBROUTINE set_tiles_masked
584 SUBROUTINE init_module_tiles
585 END SUBROUTINE init_module_tiles
587 END MODULE module_tiles