Line data Source code
1 : !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
2 :
3 : module ice_blocks
4 :
5 : ! This module contains data types and tools for decomposing a global
6 : ! horizontal domain into a set of blocks. It contains a data type
7 : ! for describing each block and contains routines for creating and
8 : ! querying the block decomposition for a global domain.
9 : !
10 : ! author: Phil Jones, LANL
11 : ! Oct. 2004: Adapted from POP by William H. Lipscomb, LANL
12 :
13 : use ice_kinds_mod
14 : use ice_domain_size, only: block_size_x, block_size_y
15 : use ice_fileunits, only: nu_diag
16 : use ice_exit, only: abort_ice
17 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
18 :
19 : implicit none
20 : private
21 :
22 : type, public :: block ! block data type
23 : integer (int_kind) :: &
24 : block_id ,&! global block number ! LCOV_EXCL_LINE
25 : local_id ,&! local address of block in current distrib ! LCOV_EXCL_LINE
26 : ilo, ihi, jlo, jhi ,&! begin, end indices for physical domain ! LCOV_EXCL_LINE
27 : iblock, jblock ! cartesian i,j position for block
28 :
29 : logical (log_kind) :: &
30 : tripole, & ! flag is true if block is at tripole bndy ! LCOV_EXCL_LINE
31 : tripoleTFlag ! tripole boundary is a T-fold
32 :
33 : integer (int_kind), dimension(:), pointer :: &
34 : i_glob, j_glob ! global domain location for each point
35 : end type
36 :
37 : public :: create_blocks ,&
38 : get_block ,& ! LCOV_EXCL_LINE
39 : get_block_parameter ,& ! LCOV_EXCL_LINE
40 : ice_blocksGetNbrID
41 :
42 : integer (int_kind), parameter, public :: &
43 : nghost = 1 ! number of ghost cells around each block
44 :
45 : integer (int_kind), public :: &! size of block domain in
46 : nx_block, ny_block ! x,y dir including ghost
47 :
48 : ! predefined directions for neighbor id routine
49 : ! Note: the directions that are commented out are implemented in
50 : ! POP but not in CICE. If the tripole cut were in the south
51 : ! instead of the north, these would need to be used (and also
52 : ! implemented in ice_boundary.F90).
53 : integer (int_kind), parameter, public :: &
54 : ice_blocksNorth = 1, & ! (i ,j+1) ! LCOV_EXCL_LINE
55 : ice_blocksSouth = 2, & ! (i ,j-1) ! LCOV_EXCL_LINE
56 : ice_blocksEast = 3, & ! (i+1,j ) ! LCOV_EXCL_LINE
57 : ice_blocksWest = 4, & ! (i-1,j ) ! LCOV_EXCL_LINE
58 : ice_blocksNorthEast = 5, & ! (i+1,j+1) ! LCOV_EXCL_LINE
59 : ice_blocksNorthWest = 6, & ! (i-1,j+1) ! LCOV_EXCL_LINE
60 : ice_blocksSouthEast = 7, & ! (i+1,j-1) ! LCOV_EXCL_LINE
61 : ice_blocksSouthWest = 8 ! (i-1,j-1)
62 : integer (int_kind), parameter, public :: &
63 : ! ice_blocksNorth2 = 9, & ! (i ,j+2) ! LCOV_EXCL_LINE
64 : ! ice_blocksSouth2 = 10, & ! (i ,j-2) ! LCOV_EXCL_LINE
65 : ice_blocksEast2 = 11, & ! (i+2,j ) ! LCOV_EXCL_LINE
66 : ice_blocksWest2 = 12 ! (i-2,j )
67 : ! ice_blocksNorthEast2 = 13, & ! (i+2,j+2)
68 : ! ice_blocksNorthWest2 = 14, & ! (i-2,j+2) ! LCOV_EXCL_LINE
69 : ! ice_blocksSouthEast2 = 15, & ! (i+2,j-2) ! LCOV_EXCL_LINE
70 : ! ice_blocksSouthWest2 = 16 ! (i-2,j-2)
71 : integer (int_kind), parameter, public :: &
72 : ice_blocksEastNorthEast = 17, & ! (i+2,j+1) ! LCOV_EXCL_LINE
73 : ! ice_blocksEastSouthEast = 18, & ! (i+2,j-1) ! LCOV_EXCL_LINE
74 : ice_blocksWestNorthWest = 19 ! (i-2,j+1)
75 : ! ice_blocksWestSouthWest = 20, & ! (i-2,j-1)
76 : ! ice_blocksNorthNorthEast = 21, & ! (i+1,j-2) ! LCOV_EXCL_LINE
77 : ! ice_blocksSouthSouthEast = 22, & ! (i+1,j-2) ! LCOV_EXCL_LINE
78 : ! ice_blocksNorthNorthWest = 23, & ! (i-1,j+2) ! LCOV_EXCL_LINE
79 : ! ice_blocksSouthSouthWest = 24 ! (i-1,j-2)
80 :
81 : integer (int_kind), public :: &
82 : nblocks_tot ,&! total number of blocks in decomposition ! LCOV_EXCL_LINE
83 : nblocks_x ,&! tot num blocks in i direction ! LCOV_EXCL_LINE
84 : nblocks_y ! tot num blocks in j direction
85 :
86 : logical (kind=log_kind), public :: &
87 : debug_blocks ! print verbose block information
88 :
89 : !-----------------------------------------------------------------------
90 : !
91 : ! module private data
92 : !
93 : !-----------------------------------------------------------------------
94 :
95 : type (block), dimension(:), allocatable, public :: &
96 : all_blocks ! block information for all blocks in domain
97 :
98 : integer (int_kind), dimension(:,:),allocatable, public :: &
99 : all_blocks_ij ! block index stored in Cartesian order
100 : ! useful for determining block index
101 : ! of neighbor blocks
102 :
103 : integer (int_kind), dimension(:,:), allocatable, target, public :: &
104 : i_global, &! global i index for each point in each block ! LCOV_EXCL_LINE
105 : j_global ! global j index for each point in each block
106 :
107 : !***********************************************************************
108 :
109 : contains
110 :
111 : !***********************************************************************
112 :
113 37 : subroutine create_blocks(nx_global, ny_global, ew_boundary_type, &
114 : ns_boundary_type)
115 :
116 : ! This subroutine decomposes the global domain into blocks and
117 : ! fills the data structures with all the necessary block information.
118 :
119 : use ice_communicate, only: my_task, master_task
120 :
121 : integer (int_kind), intent(in) :: &
122 : nx_global, ny_global ! global domain size in x,y
123 :
124 : character (*), intent(in) :: &
125 : ew_boundary_type, &! type of boundary in logical east-west dir ! LCOV_EXCL_LINE
126 : ns_boundary_type ! type of boundary in logical north-south dir
127 :
128 : !----------------------------------------------------------------------
129 : !
130 : ! local variables
131 : !
132 : !----------------------------------------------------------------------
133 :
134 : integer (int_kind) :: &
135 : i, j, n ,&! loop indices ! LCOV_EXCL_LINE
136 : iblock, jblock ,&! block loop indices ! LCOV_EXCL_LINE
137 : is, ie, js, je ! temp start, end indices
138 :
139 : character(len=*), parameter :: subname = '(create_blocks)'
140 :
141 : !----------------------------------------------------------------------
142 : !
143 : ! compute number of blocks and cartesian decomposition
144 : ! if the requested block size does not divide the global domain
145 : ! size evenly, add additional block space to accomodate padding
146 : !
147 : !----------------------------------------------------------------------
148 :
149 37 : nx_block = block_size_x + 2*nghost ! size of block domain in x,y dir
150 37 : ny_block = block_size_y + 2*nghost ! including ghost cells
151 37 : nblocks_x = (nx_global-1)/block_size_x + 1
152 37 : nblocks_y = (ny_global-1)/block_size_y + 1
153 37 : nblocks_tot = nblocks_x*nblocks_y
154 :
155 : !----------------------------------------------------------------------
156 : !
157 : ! allocate block arrays
158 : !
159 : !----------------------------------------------------------------------
160 :
161 37 : if (.not.allocated(all_blocks)) allocate(all_blocks(nblocks_tot))
162 37 : if (.not.allocated(i_global)) allocate(i_global(nx_block,nblocks_tot))
163 37 : if (.not.allocated(j_global)) allocate(j_global(ny_block,nblocks_tot))
164 37 : if (.not.allocated(all_blocks_ij)) allocate(all_blocks_ij(nblocks_x,nblocks_y))
165 :
166 : !----------------------------------------------------------------------
167 : !
168 : ! fill block data structures for all blocks in domain
169 : !
170 : !----------------------------------------------------------------------
171 :
172 37 : n = 0
173 182 : do jblock=1,nblocks_y
174 145 : js = (jblock-1)*block_size_y + 1
175 145 : if (js > ny_global) call abort_ice(subname// &
176 0 : 'ERROR: Bad block decomp: ny_block too large?')
177 145 : je = js + block_size_y - 1
178 145 : if (je > ny_global) je = ny_global ! pad array
179 :
180 1271 : do iblock=1,nblocks_x
181 1089 : n = n + 1 ! global block id
182 :
183 1089 : is = (iblock-1)*block_size_x + 1
184 1089 : if (is > nx_global) call abort_ice(subname// &
185 0 : 'ERROR: Bad block decomp: nx_block too large?')
186 1089 : ie = is + block_size_x - 1
187 1089 : if (ie > nx_global) ie = nx_global
188 :
189 1089 : all_blocks(n)%block_id = n
190 1089 : all_blocks(n)%iblock = iblock
191 1089 : all_blocks(n)%jblock = jblock
192 1089 : all_blocks(n)%ilo = nghost + 1
193 1089 : all_blocks(n)%jlo = nghost + 1
194 1089 : all_blocks(n)%ihi = nx_block - nghost ! default value
195 1089 : all_blocks(n)%jhi = ny_block - nghost ! default value
196 :
197 1089 : if (jblock == nblocks_y .and. &
198 : (ns_boundary_type == 'tripole' .or. & ! LCOV_EXCL_LINE
199 : ns_boundary_type == 'tripoleT')) then
200 0 : all_blocks(n)%tripole = .true.
201 : else
202 1089 : all_blocks(n)%tripole = .false.
203 : endif
204 1089 : all_blocks(n)%tripoleTFlag = (ns_boundary_type == 'tripoleT')
205 :
206 1089 : all_blocks_ij(iblock,jblock) = n
207 :
208 35703 : do j=1,ny_block
209 34614 : j_global(j,n) = js - nghost + j - 1
210 :
211 : !*** southern ghost cells
212 :
213 34614 : if (j_global(j,n) < 1) then
214 0 : select case (ns_boundary_type)
215 : case ('cyclic')
216 0 : j_global(j,n) = j_global(j,n) + ny_global
217 : case ('open')
218 273 : j_global(j,n) = nghost - j + 1
219 : case ('closed')
220 0 : j_global(j,n) = 0
221 : case ('tripole')
222 0 : j_global(j,n) = nghost - j + 1 ! open
223 : case ('tripoleT')
224 0 : j_global(j,n) = -j_global(j,n) + 1 ! open
225 : case default
226 273 : call abort_ice(subname//'ERROR: unknown n-s bndy type')
227 : end select
228 : endif
229 :
230 : !*** padding required
231 :
232 35703 : if (j_global(j,n) > ny_global + nghost) then
233 0 : j_global(j,n) = 0 ! padding
234 :
235 : !*** northern ghost cells
236 :
237 34614 : else if (j_global(j,n) > ny_global) then
238 0 : select case (ns_boundary_type)
239 : case ('cyclic')
240 0 : j_global(j,n) = j_global(j,n) - ny_global
241 : case ('open')
242 273 : j_global(j,n) = 2*ny_global - j_global(j,n) + 1
243 : case ('closed')
244 0 : j_global(j,n) = 0
245 : case ('tripole')
246 0 : j_global(j,n) = -j_global(j,n)
247 : case ('tripoleT')
248 0 : j_global(j,n) = -j_global(j,n)
249 : case default
250 273 : call abort_ice(subname//'ERROR: unknown n-s bndy type')
251 : end select
252 :
253 : !*** set last physical point if padded domain
254 :
255 0 : else if (j_global(j,n) == ny_global .and. &
256 : j >= all_blocks(n)%jlo .and. & ! LCOV_EXCL_LINE
257 3936 : j < all_blocks(n)%jhi) then
258 0 : all_blocks(n)%jhi = j ! last physical point in padded domain
259 : endif
260 : end do
261 :
262 1089 : all_blocks(n)%j_glob => j_global(:,n)
263 :
264 19559 : do i=1,nx_block
265 18470 : i_global(i,n) = is - nghost + i - 1
266 :
267 : !*** western ghost cells
268 :
269 18470 : if (i_global(i,n) < 1) then
270 145 : select case (ew_boundary_type)
271 : case ('cyclic')
272 145 : i_global(i,n) = i_global(i,n) + nx_global
273 : case ('open')
274 0 : i_global(i,n) = nghost - i + 1
275 : case ('closed')
276 0 : i_global(i,n) = 0
277 : case default
278 145 : call abort_ice(subname//'ERROR: unknown e-w bndy type')
279 : end select
280 : endif
281 :
282 : !*** padded domain - fill padded region with zero
283 :
284 19559 : if (i_global(i,n) > nx_global + nghost) then
285 0 : i_global(i,n) = 0
286 :
287 : !*** eastern ghost cells
288 :
289 18470 : else if (i_global(i,n) > nx_global) then
290 145 : select case (ew_boundary_type)
291 : case ('cyclic')
292 145 : i_global(i,n) = i_global(i,n) - nx_global
293 : case ('open')
294 0 : i_global(i,n) = 2*nx_global - i_global(i,n) + 1
295 : case ('closed')
296 0 : i_global(i,n) = 0
297 : case default
298 145 : call abort_ice(subname//'ERROR: unknown e-w bndy type')
299 : end select
300 :
301 : !*** last physical point in padded domain
302 :
303 0 : else if (i_global(i,n) == nx_global .and. &
304 : i >= all_blocks(n)%ilo .and. & ! LCOV_EXCL_LINE
305 3424 : i < all_blocks(n)%ihi) then
306 0 : all_blocks(n)%ihi = i
307 : endif
308 : end do
309 :
310 1234 : all_blocks(n)%i_glob => i_global(:,n)
311 :
312 : end do
313 : end do
314 :
315 37 : if (debug_blocks) then
316 0 : if (my_task == master_task) then
317 0 : write(nu_diag,*) ' '
318 0 : write(nu_diag,'(2a)') subname,' block ID, iblock, jblock Locations:'
319 0 : do n = 1, nblocks_tot
320 0 : write(nu_diag,'(2a,3i8,l4)') subname,' global block ID, iblock, jblock, tripole:', &
321 : all_blocks(n)%block_id, & ! LCOV_EXCL_LINE
322 : all_blocks(n)%iblock, & ! LCOV_EXCL_LINE
323 : all_blocks(n)%jblock, & ! LCOV_EXCL_LINE
324 0 : all_blocks(n)%tripole
325 : enddo
326 : endif
327 : endif
328 :
329 : !----------------------------------------------------------------------
330 :
331 37 : end subroutine create_blocks
332 :
333 : !***********************************************************************
334 :
335 17424 : function ice_blocksGetNbrID(blockID, direction, iBoundary, jBoundary) &
336 : result (nbrID)
337 :
338 : ! This function returns the block id of a neighboring block in a
339 : ! requested direction. Directions:
340 : ! ice\_blocksNorth (i ,j+1)
341 : ! ice\_blocksSouth (i ,j-1)
342 : ! ice\_blocksEast (i+1,j )
343 : ! ice\_blocksWest (i-1,j )
344 : ! ice\_blocksNorthEast (i+1,j+1)
345 : ! ice\_blocksNorthWest (i-1,j+1)
346 : ! ice\_blocksSouthEast (i ,j-1)
347 : ! ice\_blocksSouthWest (i-1,j-1)
348 : ! ice\_blocksNorth2 (i ,j+2)
349 : ! ice\_blocksSouth2 (i ,j-2)
350 : ! ice\_blocksEast2 (i+2,j )
351 : ! ice\_blocksWest2 (i-2,j )
352 : ! ice\_blocksNorthEast2 (i+2,j+2)
353 : ! ice\_blocksNorthWest2 (i-2,j+2)
354 : ! ice\_blocksSouthEast2 (i+2,j-2)
355 : ! ice\_blocksSouthWest2 (i-2,j-2)
356 : ! ice\_blocksEastNorthEast (i+2,j+1)
357 : ! ice\_blocksEastSouthEast (i+2,j-1)
358 : ! ice\_blocksWestNorthWest (i-2,j+1)
359 : ! ice\_blocksWestSouthWest (i-2,j-1)
360 : ! ice\_blocksNorthNorthEast (i+1,j-2)
361 : ! ice\_blocksSouthSouthEast (i+1,j-2)
362 : ! ice\_blocksNorthNorthWest (i-1,j+2)
363 : ! ice\_blocksSouthSouthWest (i-1,j-2)
364 : !
365 :
366 : integer (int_kind), intent(in) :: &
367 : blockID, &! id of block for which neighbor id requested ! LCOV_EXCL_LINE
368 : direction ! direction for which to look for neighbor -
369 : ! must be one of the predefined module
370 : ! variables for block direction
371 :
372 : character (*), intent(in) :: &
373 : iBoundary, &! determines what to do at edges of domain ! LCOV_EXCL_LINE
374 : jBoundary ! options are - open, closed, cyclic, tripole, tripoleT
375 :
376 : integer (int_kind) :: &
377 : nbrID ! block ID of neighbor in requested dir
378 :
379 : !----------------------------------------------------------------------
380 : !
381 : ! local variables
382 : !
383 : !----------------------------------------------------------------------
384 :
385 : integer (int_kind) :: &
386 : iBlock, jBlock, &! i,j block location of current block ! LCOV_EXCL_LINE
387 : inbr, jnbr ! i,j block location of neighboring block
388 :
389 : character(len=*), parameter :: subname = '(ice_blocksGetNbrID)'
390 :
391 : !----------------------------------------------------------------------
392 : !
393 : ! retrieve info for current block
394 : !
395 : !----------------------------------------------------------------------
396 :
397 17424 : call get_block_parameter(blockID, iblock=iBlock, jblock=jBlock)
398 17424 : nbrID = 0 ! initial default
399 :
400 : !----------------------------------------------------------------------
401 : !
402 : ! compute i,j block location of neighbor
403 : !
404 : !----------------------------------------------------------------------
405 :
406 19602 : select case(direction)
407 :
408 : case (ice_blocksNorth)
409 :
410 2178 : inbr = iBlock
411 2178 : jnbr = jBlock + 1
412 2178 : if (jnbr > nblocks_y) then
413 546 : select case(jBoundary)
414 : case ('open')
415 546 : jnbr = 0 ! do not write into the neighbor's ghost cells
416 : case ('closed')
417 0 : jnbr = 0
418 : case ('cyclic')
419 0 : jnbr = 1
420 : case ('tripole':'tripoleT')
421 : !*** return negative j value to flag tripole
422 : !*** i index of main northern neighbor across the
423 : !*** tripole cut - may also need i+1,i-1 to get
424 : !*** other points if there has been padding or
425 : !*** if the block size does not divide the domain
426 : !*** evenly
427 0 : inbr = nblocks_x - iBlock + 1
428 0 : jnbr = -jBlock
429 : case default
430 546 : call abort_ice(subname//'ERROR: unknown north boundary')
431 : end select
432 : endif
433 :
434 : case (ice_blocksSouth)
435 :
436 2178 : inbr = iBlock
437 2178 : jnbr = jBlock - 1
438 2724 : if (jnbr < 1) then
439 546 : select case(jBoundary)
440 : case ('open')
441 546 : jnbr = 0 ! do not write into the neighbor's ghost cells
442 : case ('closed')
443 0 : jnbr = 0
444 : case ('cyclic')
445 0 : jnbr = nblocks_y
446 : case ('tripole')
447 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
448 : case ('tripoleT')
449 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
450 : case default
451 546 : call abort_ice(subname//'ERROR: unknown south boundary')
452 : end select
453 : endif
454 :
455 : case (ice_blocksEast )
456 :
457 2178 : inbr = iBlock + 1
458 2178 : jnbr = jBlock
459 2178 : if (inbr > nblocks_x) then
460 0 : select case(iBoundary)
461 : case ('open')
462 0 : inbr = 0 ! do not write into the neighbor's ghost cells
463 : case ('closed')
464 0 : inbr = 0
465 : case ('cyclic')
466 290 : inbr = 1
467 : case default
468 290 : call abort_ice(subname//'ERROR: unknown east boundary')
469 : end select
470 : endif
471 :
472 : case (ice_blocksWest )
473 :
474 2178 : jnbr = jBlock
475 2178 : inbr = iBlock - 1
476 2468 : if (inbr < 1) then
477 0 : select case(iBoundary)
478 : case ('open')
479 0 : inbr = 0 ! do not write into the neighbor's ghost cells
480 : case ('closed')
481 0 : inbr = 0
482 : case ('cyclic')
483 290 : inbr = nblocks_x
484 : case default
485 290 : call abort_ice(subname//'ERROR: unknown west boundary')
486 : end select
487 : endif
488 :
489 : case (ice_blocksNorthEast)
490 :
491 2178 : inbr = iBlock + 1
492 2178 : jnbr = jBlock + 1
493 2178 : if (inbr > nblocks_x) then
494 0 : select case(iBoundary)
495 : case ('open')
496 0 : inbr = 0 ! do not write into the neighbor's ghost cells
497 : case ('closed')
498 0 : inbr = 0
499 : case ('cyclic')
500 290 : inbr = 1
501 : case default
502 290 : call abort_ice(subname//'ERROR: unknown east boundary')
503 : end select
504 : endif
505 2178 : if (jnbr > nblocks_y) then
506 546 : select case(jBoundary)
507 : case ('open')
508 546 : jnbr = 0 ! do not write into the neighbor's ghost cells
509 : case ('closed')
510 0 : jnbr = 0
511 : case ('cyclic')
512 0 : jnbr = 1
513 : case ('tripole':'tripoleT')
514 : !*** return negative j value to flag tripole
515 : !*** i index of main northern neighbor across the
516 : !*** tripole cut - may also need i+1,i-1 to get
517 : !*** other points if there has been padding or
518 : !*** if the block size does not divide the domain
519 : !*** evenly
520 0 : inbr = nblocks_x - iBlock
521 0 : if (inbr == 0) inbr = nblocks_x
522 0 : jnbr = -jBlock
523 : case default
524 546 : call abort_ice(subname//'ERROR: unknown north boundary')
525 : end select
526 : endif
527 :
528 : case (ice_blocksNorthWest)
529 :
530 2178 : inbr = iBlock - 1
531 2178 : jnbr = jBlock + 1
532 2178 : if (inbr < 1) then
533 0 : select case(iBoundary)
534 : case ('open')
535 0 : inbr = 0 ! do not write into the neighbor's ghost cells
536 : case ('closed')
537 0 : inbr = 0
538 : case ('cyclic')
539 290 : inbr = nblocks_x
540 : case default
541 290 : call abort_ice(subname//'ERROR: unknown west boundary')
542 : end select
543 : endif
544 2178 : if (jnbr > nblocks_y) then
545 546 : select case(jBoundary)
546 : case ('open')
547 546 : jnbr = 0 ! do not write into the neighbor's ghost cells
548 : case ('closed')
549 0 : jnbr = 0
550 : case ('cyclic')
551 0 : jnbr = 1
552 : case ('tripole':'tripoleT')
553 : !*** return negative j value to flag tripole
554 : !*** i index of main northern neighbor across the
555 : !*** tripole cut - may also need i+1,i-1 to get
556 : !*** other points if there has been padding or
557 : !*** if the block size does not divide the domain
558 : !*** evenly
559 0 : inbr = nblocks_x - iBlock + 2
560 0 : if (inbr > nblocks_x) inbr = 1
561 0 : jnbr = -jBlock
562 : case default
563 546 : call abort_ice(subname//'ERROR: unknown north boundary')
564 : end select
565 : endif
566 :
567 : case (ice_blocksSouthEast )
568 :
569 2178 : inbr = iBlock + 1
570 2178 : jnbr = jBlock - 1
571 2178 : if (inbr > nblocks_x) then
572 0 : select case(iBoundary)
573 : case ('open')
574 0 : inbr = 0 ! do not write into the neighbor's ghost cells
575 : case ('closed')
576 0 : inbr = 0
577 : case ('cyclic')
578 290 : inbr = 1
579 : case default
580 290 : call abort_ice(subname//'ERROR: unknown east boundary')
581 : end select
582 : endif
583 2724 : if (jnbr < 1) then
584 546 : select case(jBoundary)
585 : case ('open')
586 546 : inbr = 0 ! do not write into the neighbor's ghost cells
587 : case ('closed')
588 0 : jnbr = 0
589 : case ('cyclic')
590 0 : jnbr = nblocks_y
591 : case ('tripole')
592 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
593 : case ('tripoleT')
594 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
595 : case default
596 546 : call abort_ice(subname//'ERROR: unknown south boundary')
597 : end select
598 : endif
599 :
600 : case (ice_blocksSouthWest )
601 2178 : inbr = iBlock - 1
602 2178 : jnbr = jBlock - 1
603 2178 : if (inbr < 1) then
604 0 : select case(iBoundary)
605 : case ('open')
606 0 : inbr = 0 ! do not write into the neighbor's ghost cells
607 : case ('closed')
608 0 : inbr = 0
609 : case ('cyclic')
610 290 : inbr = nblocks_x
611 : case default
612 290 : call abort_ice(subname//'ERROR: unknown west boundary')
613 : end select
614 : endif
615 2724 : if (jnbr < 1) then
616 546 : select case(jBoundary)
617 : case ('open')
618 546 : jnbr = 0 ! do not write into the neighbor's ghost cells
619 : case ('closed')
620 0 : jnbr = 0
621 : case ('cyclic')
622 0 : jnbr = nblocks_y
623 : case ('tripole')
624 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
625 : case ('tripoleT')
626 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
627 : case default
628 546 : call abort_ice(subname//'ERROR: unknown south boundary')
629 : end select
630 : endif
631 :
632 : case (ice_blocksEast2)
633 :
634 0 : inbr = iBlock + 2
635 0 : jnbr = jBlock
636 0 : if (inbr > nblocks_x) then
637 0 : select case(iBoundary)
638 : case ('open')
639 0 : inbr = 0 ! do not write into the neighbor's ghost cells
640 : case ('closed')
641 0 : inbr = 0
642 : case ('cyclic')
643 0 : inbr = inbr - nblocks_x
644 : case default
645 0 : call abort_ice(subname//'ERROR: unknown east boundary')
646 : end select
647 : endif
648 :
649 : case (ice_blocksWest2)
650 0 : jnbr = jBlock
651 0 : inbr = iBlock - 2
652 0 : if (inbr < 1) then
653 0 : select case(iBoundary)
654 : case ('open')
655 0 : inbr = 0 ! do not write into the neighbor's ghost cells
656 : case ('closed')
657 0 : inbr = 0
658 : case ('cyclic')
659 0 : inbr = nblocks_x + inbr
660 : case default
661 0 : call abort_ice(subname//'ERROR: unknown west boundary')
662 : end select
663 : endif
664 :
665 : case (ice_blocksEastNorthEast)
666 :
667 0 : inbr = iBlock + 2
668 0 : jnbr = jBlock + 1
669 0 : if (inbr > nblocks_x) then
670 0 : select case(iBoundary)
671 : case ('open')
672 0 : inbr = 0 ! do not write into the neighbor's ghost cells
673 : case ('closed')
674 0 : inbr = 0
675 : case ('cyclic')
676 0 : inbr = inbr - nblocks_x
677 : case default
678 0 : call abort_ice(subname//'ERROR: unknown east boundary')
679 : end select
680 : endif
681 0 : if (jnbr > nblocks_y) then
682 0 : select case(jBoundary)
683 : case ('open')
684 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
685 : case ('closed')
686 0 : jnbr = 0
687 : case ('cyclic')
688 0 : jnbr = jnbr - nblocks_y
689 : case ('tripole':'tripoleT')
690 : !*** return negative j value to flag tripole
691 : !*** i index of main northern neighbor across the
692 : !*** tripole cut - may also need i+1,i-1 to get
693 : !*** other points if there has been padding or
694 : !*** if the block size does not divide the domain
695 : !*** evenly
696 0 : inbr = nblocks_x - iBlock - 1
697 0 : if (inbr <= 0) inbr = inbr + nblocks_x
698 0 : jnbr = -jBlock
699 : case default
700 0 : call abort_ice(subname//'ERROR: unknown north boundary')
701 : end select
702 : endif
703 :
704 : case (ice_blocksWestNorthWest)
705 :
706 0 : inbr = iBlock - 2
707 0 : jnbr = jBlock + 1
708 0 : if (inbr < 1) then
709 0 : select case(iBoundary)
710 : case ('open')
711 0 : inbr = 0 ! do not write into the neighbor's ghost cells
712 : case ('closed')
713 0 : inbr = 0
714 : case ('cyclic')
715 0 : inbr = nblocks_x + inbr
716 : case default
717 0 : call abort_ice(subname//'ERROR: unknown west boundary')
718 : end select
719 : endif
720 0 : if (jnbr > nblocks_y) then
721 0 : select case(jBoundary)
722 : case ('open')
723 0 : jnbr = 0 ! do not write into the neighbor's ghost cells
724 : case ('closed')
725 0 : jnbr = 0
726 : case ('cyclic')
727 0 : jnbr = jnbr + nblocks_y
728 : case ('tripole':'tripoleT')
729 : !*** return negative j value to flag tripole
730 : !*** i index of main northern neighbor across the
731 : !*** tripole cut - may also need i+1,i-1 to get
732 : !*** other points if there has been padding or
733 : !*** if the block size does not divide the domain
734 : !*** evenly
735 0 : inbr = nblocks_x - iBlock + 3
736 0 : if (inbr > nblocks_x) inbr = inbr - nblocks_x
737 0 : jnbr = -jBlock
738 : case default
739 0 : call abort_ice(subname//'ERROR: unknown north boundary')
740 : end select
741 : endif
742 :
743 : case default
744 :
745 0 : call abort_ice(subname//'ERROR: unknown direction')
746 17424 : return
747 :
748 : end select
749 :
750 : !----------------------------------------------------------------------
751 : !
752 : ! now get block id for this neighbor block
753 : !
754 : !----------------------------------------------------------------------
755 :
756 17424 : if (inbr > 0 .and. jnbr > 0) then
757 14148 : nbrID = all_blocks_ij(inbr,jnbr)
758 3276 : else if (inbr > 0 .and. jnbr < 0) then ! tripole upper boundary
759 : !*** return negative value to flag tripole
760 0 : nbrID = -all_blocks_ij(inbr,abs(jnbr))
761 : else
762 3276 : nbrID = 0 ! neighbor outside domain
763 : endif
764 :
765 : !----------------------------------------------------------------------
766 :
767 34848 : end function ice_blocksGetNbrID
768 :
769 : !**********************************************************************
770 :
771 16196707 : function get_block(block_id,local_id)
772 :
773 : ! This function returns the block data structure for the block
774 : ! associated with the input block id.
775 :
776 : integer (int_kind), intent(in) :: &
777 : block_id, &! global block id for requested block info ! LCOV_EXCL_LINE
778 : local_id ! local block id to assign to this block
779 :
780 : type (block) :: &
781 : get_block ! block information returned for requested block
782 :
783 : character(len=*), parameter :: subname = '(get_block)'
784 :
785 : !----------------------------------------------------------------------
786 : !
787 : ! check for valid id. if valid, return block info for requested block
788 : !
789 : !----------------------------------------------------------------------
790 :
791 16196707 : if (block_id < 1 .or. block_id > nblocks_tot) then
792 0 : call abort_ice(subname//'ERROR: invalid block_id')
793 : endif
794 :
795 16196707 : get_block = all_blocks(block_id)
796 16196707 : get_block%local_id = local_id
797 :
798 : !----------------------------------------------------------------------
799 :
800 16196707 : end function get_block
801 :
802 : !**********************************************************************
803 :
804 9542404 : subroutine get_block_parameter(block_id, local_id, &
805 : ilo, ihi, jlo, jhi, & ! LCOV_EXCL_LINE
806 : iblock, jblock, tripole, & ! LCOV_EXCL_LINE
807 : i_glob, j_glob)
808 :
809 : ! This routine returns requested parts of the block data type
810 : ! for the block associated with the input block id
811 :
812 : integer (int_kind), intent(in) :: &
813 : block_id ! global block id for which parameters are requested
814 :
815 : !(optional) parts of block data type to extract if requested
816 :
817 : integer (int_kind), intent(out), optional :: &
818 : local_id ,&! local id assigned to block in current distrb ! LCOV_EXCL_LINE
819 : ilo, ihi, jlo, jhi ,&! begin,end indices for physical domain ! LCOV_EXCL_LINE
820 : iblock, jblock ! cartesian i,j position for bloc
821 :
822 : logical (log_kind), intent(out), optional :: &
823 : tripole ! flag is true if block on tripole bndy
824 :
825 : integer (int_kind), dimension(:), pointer, optional :: &
826 : i_glob, j_glob ! global domain location for each point
827 :
828 : character(len=*), parameter :: subname = '(get_block_parameter)'
829 :
830 : !----------------------------------------------------------------------
831 : !
832 : ! extract each component of data type if requested
833 : !
834 : !----------------------------------------------------------------------
835 :
836 9542404 : if (block_id < 1 .or. block_id > nblocks_tot) then
837 0 : call abort_ice(subname//'ERROR: invalid block_id')
838 : endif
839 :
840 9542404 : if (present(local_id)) local_id = all_blocks(block_id)%local_id
841 9542404 : if (present(ilo )) ilo = all_blocks(block_id)%ilo
842 9542404 : if (present(ihi )) ihi = all_blocks(block_id)%ihi
843 9542404 : if (present(jlo )) jlo = all_blocks(block_id)%jlo
844 9542404 : if (present(jhi )) jhi = all_blocks(block_id)%jhi
845 9542404 : if (present(iblock )) iblock = all_blocks(block_id)%iblock
846 9542404 : if (present(jblock )) jblock = all_blocks(block_id)%jblock
847 9542404 : if (present(i_glob )) i_glob => all_blocks(block_id)%i_glob
848 9542404 : if (present(j_glob )) j_glob => all_blocks(block_id)%j_glob
849 9542404 : if (present(tripole )) tripole = all_blocks(block_id)%tripole
850 :
851 : !----------------------------------------------------------------------
852 :
853 9542404 : end subroutine get_block_parameter
854 :
855 : !**********************************************************************
856 :
857 0 : end module ice_blocks
858 :
859 : !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
|