Line data Source code
1 : !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
2 :
3 : module ice_gather_scatter
4 :
5 : ! This module contains routines that mimic the behavior of the mpi
6 : ! version in the case of a single processor: gathering data to a single
7 : ! processor from a distributed array, and scattering data from a
8 : ! single processor to a distributed array.
9 : !
10 : ! author: Phil Jones, LANL
11 : ! Oct. 2004: Adapted from POP version by William H. Lipscomb, LANL
12 : ! Jan. 2008: Elizabeth Hunke replaced old routines with new POP
13 : ! infrastructure, added specialized routine scatter_global_stress
14 :
15 : use ice_kinds_mod
16 : use ice_constants
17 : use ice_blocks, only: block, get_block, nghost, nx_block, ny_block, &
18 : nblocks_x, nblocks_y, nblocks_tot
19 : use ice_distribution, only: distrb
20 : use ice_domain_size, only: nx_global, ny_global
21 : use ice_exit, only: abort_ice
22 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
23 :
24 : implicit none
25 : private
26 :
27 : public :: gather_global, &
28 : gather_global_ext, & ! LCOV_EXCL_LINE
29 : scatter_global, & ! LCOV_EXCL_LINE
30 : scatter_global_ext, & ! LCOV_EXCL_LINE
31 : scatter_global_stress
32 :
33 : !-----------------------------------------------------------------------
34 : !
35 : ! overload module functions
36 : !
37 : !-----------------------------------------------------------------------
38 :
39 : interface gather_global_ext
40 : module procedure gather_global_ext_dbl, &
41 : gather_global_ext_int, & ! LCOV_EXCL_LINE
42 : gather_global_ext_log
43 : end interface
44 :
45 : interface gather_global
46 : module procedure gather_global_dbl, &
47 : gather_global_real, & ! LCOV_EXCL_LINE
48 : gather_global_int
49 : end interface
50 :
51 : interface scatter_global
52 : module procedure scatter_global_dbl, &
53 : scatter_global_real, & ! LCOV_EXCL_LINE
54 : scatter_global_int
55 : end interface
56 :
57 : !-----------------------------------------------------------------------
58 : !
59 : ! module variables
60 : !
61 : !-----------------------------------------------------------------------
62 :
63 : !***********************************************************************
64 :
65 : contains
66 :
67 : !***********************************************************************
68 :
69 214 : subroutine gather_global_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
70 :
71 : ! This subroutine gathers a distributed array to a global-sized
72 : ! array on the processor dst_task.
73 : !
74 : ! This is the specific inteface for double precision arrays
75 : ! corresponding to the generic interface gather_global. It is shown
76 : ! to provide information on the generic interface (the generic
77 : ! interface is identical, but chooses a specific inteface based
78 : ! on the data type of the input argument).
79 :
80 : integer (int_kind), intent(in) :: &
81 : dst_task ! task to which array should be gathered
82 :
83 : type (distrb), intent(in) :: &
84 : src_dist ! distribution of blocks in the source array
85 :
86 : real (dbl_kind), dimension(:,:,:), intent(in) :: &
87 : ARRAY ! array containing horizontal slab of distributed field
88 :
89 : real (dbl_kind), intent(in), optional :: &
90 : spc_val
91 :
92 : real (dbl_kind), dimension(:,:), intent(inout) :: &
93 : ARRAY_G ! array containing global horizontal field on dst_task
94 :
95 : !-----------------------------------------------------------------------
96 : !
97 : ! local variables
98 : !
99 : !-----------------------------------------------------------------------
100 :
101 : integer (int_kind) :: &
102 : i,j,n ! dummy loop counters
103 :
104 : real (dbl_kind) :: &
105 : special_value
106 :
107 : type (block) :: &
108 : this_block ! block info for current block
109 :
110 : character(len=*), parameter :: subname = '(gather_global_dbl)'
111 :
112 214 : if (present(spc_val)) then
113 166 : special_value = spc_val
114 : else
115 48 : special_value = spval_dbl
116 : endif
117 :
118 : !-----------------------------------------------------------------------
119 : !
120 : ! copy local array into block decomposition
121 : !
122 : !-----------------------------------------------------------------------
123 :
124 428 : do n=1,nblocks_tot
125 :
126 214 : this_block = get_block(n,n)
127 :
128 : !*** copy local blocks
129 :
130 428 : if (src_dist%blockLocation(n) /= 0) then
131 :
132 25038 : do j=this_block%jlo,this_block%jhi
133 2507438 : do i=this_block%ilo,this_block%ihi
134 : ARRAY_G(this_block%i_glob(i), &
135 : this_block%j_glob(j)) = & ! LCOV_EXCL_LINE
136 2507224 : ARRAY(i,j,src_dist%blockLocalID(n))
137 : end do
138 : end do
139 :
140 : else !*** fill land blocks with special values
141 :
142 0 : do j=this_block%jlo,this_block%jhi
143 0 : do i=this_block%ilo,this_block%ihi
144 : ARRAY_G(this_block%i_glob(i), &
145 0 : this_block%j_glob(j)) = special_value
146 : end do
147 : end do
148 :
149 : endif
150 :
151 : end do
152 :
153 : !-----------------------------------------------------------------------
154 :
155 214 : end subroutine gather_global_dbl
156 :
157 : !***********************************************************************
158 :
159 0 : subroutine gather_global_real(ARRAY_G, ARRAY, dst_task, src_dist)
160 :
161 : !-----------------------------------------------------------------------
162 : !
163 : ! This subroutine gathers a distributed array to a global-sized
164 : ! array on the processor dst_task.
165 : !
166 : !-----------------------------------------------------------------------
167 : !-----------------------------------------------------------------------
168 : !
169 : ! input variables
170 : !
171 : !-----------------------------------------------------------------------
172 :
173 : integer (int_kind), intent(in) :: &
174 : dst_task ! task to which array should be gathered
175 :
176 : type (distrb), intent(in) :: &
177 : src_dist ! distribution of blocks in the source array
178 :
179 : real (real_kind), dimension(:,:,:), intent(in) :: &
180 : ARRAY ! array containing distributed field
181 :
182 : !-----------------------------------------------------------------------
183 : !
184 : ! output variables
185 : !
186 : !-----------------------------------------------------------------------
187 :
188 : real (real_kind), dimension(:,:), intent(inout) :: &
189 : ARRAY_G ! array containing global field on dst_task
190 :
191 : !-----------------------------------------------------------------------
192 : !
193 : ! local variables
194 : !
195 : !-----------------------------------------------------------------------
196 :
197 : integer (int_kind) :: &
198 : i,j,n ! dummy loop counters
199 :
200 : type (block) :: &
201 : this_block ! block info for current block
202 :
203 : character(len=*), parameter :: subname = '(gather_global_real)'
204 :
205 : !-----------------------------------------------------------------------
206 : !
207 : ! copy local array into block decomposition
208 : !
209 : !-----------------------------------------------------------------------
210 :
211 0 : do n=1,nblocks_tot
212 :
213 0 : this_block = get_block(n,n)
214 :
215 : !*** copy local blocks
216 :
217 0 : if (src_dist%blockLocation(n) /= 0) then
218 :
219 0 : do j=this_block%jlo,this_block%jhi
220 0 : do i=this_block%ilo,this_block%ihi
221 : ARRAY_G(this_block%i_glob(i), &
222 : this_block%j_glob(j)) = & ! LCOV_EXCL_LINE
223 0 : ARRAY(i,j,src_dist%blockLocalID(n))
224 : end do
225 : end do
226 :
227 : else !*** fill land blocks with zeroes
228 :
229 0 : do j=this_block%jlo,this_block%jhi
230 0 : do i=this_block%ilo,this_block%ihi
231 : ARRAY_G(this_block%i_glob(i), &
232 0 : this_block%j_glob(j)) = 0._real_kind
233 : end do
234 : end do
235 :
236 : endif
237 :
238 : end do
239 :
240 : !-----------------------------------------------------------------------
241 :
242 0 : end subroutine gather_global_real
243 :
244 : !***********************************************************************
245 :
246 0 : subroutine gather_global_int(ARRAY_G, ARRAY, dst_task, src_dist)
247 :
248 : !-----------------------------------------------------------------------
249 : !
250 : ! This subroutine gathers a distributed array to a global-sized
251 : ! array on the processor dst_task.
252 : !
253 : !-----------------------------------------------------------------------
254 : !-----------------------------------------------------------------------
255 : !
256 : ! input variables
257 : !
258 : !-----------------------------------------------------------------------
259 :
260 : integer (int_kind), intent(in) :: &
261 : dst_task ! task to which array should be gathered
262 :
263 : type (distrb), intent(in) :: &
264 : src_dist ! distribution of blocks in the source array
265 :
266 : integer (int_kind), dimension(:,:,:), intent(in) :: &
267 : ARRAY ! array containing distributed field
268 :
269 : !-----------------------------------------------------------------------
270 : !
271 : ! output variables
272 : !
273 : !-----------------------------------------------------------------------
274 :
275 : integer (int_kind), dimension(:,:), intent(inout) :: &
276 : ARRAY_G ! array containing global field on dst_task
277 :
278 : !-----------------------------------------------------------------------
279 : !
280 : ! local variables
281 : !
282 : !-----------------------------------------------------------------------
283 :
284 : integer (int_kind) :: &
285 : i,j,n ! dummy loop counters
286 :
287 : type (block) :: &
288 : this_block ! block info for current block
289 :
290 : character(len=*), parameter :: subname = '(gather_global_int)'
291 :
292 : !-----------------------------------------------------------------------
293 : !
294 : ! copy local array into block decomposition
295 : !
296 : !-----------------------------------------------------------------------
297 :
298 0 : do n=1,nblocks_tot
299 :
300 0 : this_block = get_block(n,n)
301 :
302 : !*** copy local blocks
303 :
304 0 : if (src_dist%blockLocation(n) /= 0) then
305 :
306 0 : do j=this_block%jlo,this_block%jhi
307 0 : do i=this_block%ilo,this_block%ihi
308 : ARRAY_G(this_block%i_glob(i), &
309 : this_block%j_glob(j)) = & ! LCOV_EXCL_LINE
310 0 : ARRAY(i,j,src_dist%blockLocalID(n))
311 : end do
312 : end do
313 :
314 : else !*** fill land blocks with zeroes
315 :
316 0 : do j=this_block%jlo,this_block%jhi
317 0 : do i=this_block%ilo,this_block%ihi
318 : ARRAY_G(this_block%i_glob(i), &
319 0 : this_block%j_glob(j)) = 0
320 : end do
321 : end do
322 :
323 : endif
324 :
325 : end do
326 :
327 : !-----------------------------------------------------------------------
328 :
329 0 : end subroutine gather_global_int
330 :
331 : !***********************************************************************
332 :
333 0 : subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
334 :
335 : ! This subroutine gathers a distributed array to a global-sized
336 : ! array on the processor dst_task, including ghost cells.
337 :
338 : integer (int_kind), intent(in) :: &
339 : dst_task ! task to which array should be gathered
340 :
341 : type (distrb), intent(in) :: &
342 : src_dist ! distribution of blocks in the source array
343 :
344 : real (dbl_kind), dimension(:,:,:), intent(in) :: &
345 : ARRAY ! array containing horizontal slab of distributed field
346 :
347 : real (dbl_kind), optional :: &
348 : spc_val
349 :
350 : real (dbl_kind), dimension(:,:), intent(inout) :: &
351 : ARRAY_G ! array containing global horizontal field on dst_task
352 :
353 : !-----------------------------------------------------------------------
354 : !
355 : ! local variables
356 : !
357 : !-----------------------------------------------------------------------
358 :
359 : integer (int_kind) :: &
360 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
361 : nx, ny ! global dimensions
362 :
363 : real (dbl_kind) :: &
364 : special_value
365 :
366 : type (block) :: &
367 : this_block ! block info for current block
368 :
369 : character(len=*), parameter :: subname = '(gather_global_ext_dbl)'
370 :
371 0 : if (present(spc_val)) then
372 0 : special_value = spc_val
373 : else
374 0 : special_value = spval_dbl
375 : endif
376 0 : ARRAY_G = special_value
377 :
378 0 : nx = nx_global + 2*nghost
379 0 : ny = ny_global + 2*nghost
380 :
381 : !-----------------------------------------------------------------------
382 : !
383 : ! copy local array into block decomposition
384 : !
385 : !-----------------------------------------------------------------------
386 :
387 0 : do n=1,nblocks_tot
388 :
389 0 : this_block = get_block(n,n)
390 :
391 : !*** copy local blocks
392 :
393 0 : if (src_dist%blockLocation(n) /= 0) then
394 :
395 0 : do j=this_block%jlo,this_block%jhi
396 0 : do i=this_block%ilo,this_block%ihi
397 : ARRAY_G(this_block%i_glob(i)+nghost, &
398 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
399 0 : ARRAY(i,j,src_dist%blockLocalID(n))
400 : end do
401 : end do
402 :
403 : ! fill ghost cells
404 0 : if (this_block%jblock == 1) then
405 : ! south block
406 0 : do j=1, nghost
407 0 : do i=this_block%ilo,this_block%ihi
408 : ARRAY_G(this_block%i_glob(i)+nghost,j) = &
409 0 : ARRAY (i,j,src_dist%blockLocalID(n))
410 : end do
411 : end do
412 0 : if (this_block%iblock == 1) then
413 : ! southwest corner
414 0 : do j=1, nghost
415 0 : do i=1, nghost
416 : ARRAY_G(i,j) = &
417 0 : ARRAY (i,j,src_dist%blockLocalID(n))
418 : end do
419 : end do
420 : endif
421 : endif
422 0 : if (this_block%jblock == nblocks_y) then
423 : ! north block
424 0 : do j=1, nghost
425 0 : do i=this_block%ilo,this_block%ihi
426 : ARRAY_G(this_block%i_glob(i)+nghost, &
427 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
428 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
429 : end do
430 : end do
431 0 : if (this_block%iblock == nblocks_x) then
432 : ! northeast corner
433 0 : do j=1, nghost
434 0 : do i=1, nghost
435 : ARRAY_G(nx-i+1, ny-j+1) = &
436 : ARRAY (this_block%ihi+nghost-i+1, & ! LCOV_EXCL_LINE
437 : this_block%jhi+nghost-j+1, & ! LCOV_EXCL_LINE
438 0 : src_dist%blockLocalID(n))
439 : end do
440 : end do
441 : endif
442 : endif
443 0 : if (this_block%iblock == 1) then
444 : ! west block
445 0 : do j=this_block%jlo,this_block%jhi
446 0 : do i=1, nghost
447 : ARRAY_G(i,this_block%j_glob(j)+nghost) = &
448 0 : ARRAY (i,j,src_dist%blockLocalID(n))
449 : end do
450 : end do
451 0 : if (this_block%jblock == nblocks_y) then
452 : ! northwest corner
453 0 : do j=1, nghost
454 0 : do i=1, nghost
455 : ARRAY_G(i, ny-j+1) = &
456 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
457 : end do
458 : end do
459 : endif
460 : endif
461 0 : if (this_block%iblock == nblocks_x) then
462 : ! east block
463 0 : do j=this_block%jlo,this_block%jhi
464 0 : do i=1, nghost
465 : ARRAY_G(nx_global + nghost + i, &
466 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
467 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
468 : end do
469 : end do
470 0 : if (this_block%jblock == 1) then
471 : ! southeast corner
472 0 : do j=1, nghost
473 0 : do i=1, nghost
474 : ARRAY_G( nx-i+1,j) = &
475 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
476 : end do
477 : end do
478 : endif
479 : endif
480 :
481 : endif ! src_dist%blockLocation
482 :
483 : end do
484 :
485 : !-----------------------------------------------------------------------
486 :
487 0 : end subroutine gather_global_ext_dbl
488 :
489 : !***********************************************************************
490 :
491 0 : subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
492 :
493 : ! This subroutine gathers a distributed array to a global-sized
494 : ! array on the processor dst_task, including ghost cells.
495 :
496 : integer (int_kind), intent(in) :: &
497 : dst_task ! task to which array should be gathered
498 :
499 : type (distrb), intent(in) :: &
500 : src_dist ! distribution of blocks in the source array
501 :
502 : integer (int_kind), dimension(:,:,:), intent(in) :: &
503 : ARRAY ! array containing horizontal slab of distributed field
504 :
505 : integer (int_kind), optional :: &
506 : spc_val
507 :
508 : integer (int_kind), dimension(:,:), intent(inout) :: &
509 : ARRAY_G ! array containing global horizontal field on dst_task
510 :
511 : !-----------------------------------------------------------------------
512 : !
513 : ! local variables
514 : !
515 : !-----------------------------------------------------------------------
516 :
517 : integer (int_kind) :: &
518 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
519 : nx, ny ! global dimensions
520 :
521 : integer (int_kind) :: &
522 : special_value
523 :
524 : type (block) :: &
525 : this_block ! block info for current block
526 :
527 : character(len=*), parameter :: subname = '(gather_global_ext_int)'
528 :
529 0 : if (present(spc_val)) then
530 0 : special_value = spc_val
531 : else
532 0 : special_value = -9999
533 : endif
534 0 : ARRAY_G = special_value
535 :
536 0 : nx = nx_global + 2*nghost
537 0 : ny = ny_global + 2*nghost
538 :
539 : !-----------------------------------------------------------------------
540 : !
541 : ! copy local array into block decomposition
542 : !
543 : !-----------------------------------------------------------------------
544 :
545 0 : do n=1,nblocks_tot
546 :
547 0 : this_block = get_block(n,n)
548 :
549 : !*** copy local blocks
550 :
551 0 : if (src_dist%blockLocation(n) /= 0) then
552 :
553 0 : do j=this_block%jlo,this_block%jhi
554 0 : do i=this_block%ilo,this_block%ihi
555 : ARRAY_G(this_block%i_glob(i)+nghost, &
556 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
557 0 : ARRAY(i,j,src_dist%blockLocalID(n))
558 : end do
559 : end do
560 :
561 : ! fill ghost cells
562 0 : if (this_block%jblock == 1) then
563 : ! south block
564 0 : do j=1, nghost
565 0 : do i=this_block%ilo,this_block%ihi
566 : ARRAY_G(this_block%i_glob(i)+nghost,j) = &
567 0 : ARRAY (i,j,src_dist%blockLocalID(n))
568 : end do
569 : end do
570 0 : if (this_block%iblock == 1) then
571 : ! southwest corner
572 0 : do j=1, nghost
573 0 : do i=1, nghost
574 : ARRAY_G(i,j) = &
575 0 : ARRAY (i,j,src_dist%blockLocalID(n))
576 : end do
577 : end do
578 : endif
579 : endif
580 0 : if (this_block%jblock == nblocks_y) then
581 : ! north block
582 0 : do j=1, nghost
583 0 : do i=this_block%ilo,this_block%ihi
584 : ARRAY_G(this_block%i_glob(i)+nghost, &
585 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
586 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
587 : end do
588 : end do
589 0 : if (this_block%iblock == nblocks_x) then
590 : ! northeast corner
591 0 : do j=1, nghost
592 0 : do i=1, nghost
593 : ARRAY_G(nx-i+1, ny-j+1) = &
594 : ARRAY (this_block%ihi+nghost-i+1, & ! LCOV_EXCL_LINE
595 : this_block%jhi+nghost-j+1, & ! LCOV_EXCL_LINE
596 0 : src_dist%blockLocalID(n))
597 : end do
598 : end do
599 : endif
600 : endif
601 0 : if (this_block%iblock == 1) then
602 : ! west block
603 0 : do j=this_block%jlo,this_block%jhi
604 0 : do i=1, nghost
605 : ARRAY_G(i,this_block%j_glob(j)+nghost) = &
606 0 : ARRAY (i,j,src_dist%blockLocalID(n))
607 : end do
608 : end do
609 0 : if (this_block%jblock == nblocks_y) then
610 : ! northwest corner
611 0 : do j=1, nghost
612 0 : do i=1, nghost
613 : ARRAY_G(i, ny-j+1) = &
614 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
615 : end do
616 : end do
617 : endif
618 : endif
619 0 : if (this_block%iblock == nblocks_x) then
620 : ! east block
621 0 : do j=this_block%jlo,this_block%jhi
622 0 : do i=1, nghost
623 : ARRAY_G(nx_global + nghost + i, &
624 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
625 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
626 : end do
627 : end do
628 0 : if (this_block%jblock == 1) then
629 : ! southeast corner
630 0 : do j=1, nghost
631 0 : do i=1, nghost
632 : ARRAY_G( nx-i+1,j) = &
633 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
634 : end do
635 : end do
636 : endif
637 : endif
638 :
639 : endif ! src_dist%blockLocation
640 :
641 : end do
642 :
643 : !-----------------------------------------------------------------------
644 :
645 0 : end subroutine gather_global_ext_int
646 :
647 : !***********************************************************************
648 :
649 0 : subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
650 :
651 : ! This subroutine gathers a distributed array to a global-sized
652 : ! array on the processor dst_task, including ghost cells.
653 :
654 : integer (int_kind), intent(in) :: &
655 : dst_task ! task to which array should be gathered
656 :
657 : type (distrb), intent(in) :: &
658 : src_dist ! distribution of blocks in the source array
659 :
660 : logical (log_kind), dimension(:,:,:), intent(in) :: &
661 : ARRAY ! array containing horizontal slab of distributed field
662 :
663 : logical (log_kind), optional :: &
664 : spc_val
665 :
666 : logical (log_kind), dimension(:,:), intent(inout) :: &
667 : ARRAY_G ! array containing global horizontal field on dst_task
668 :
669 : !-----------------------------------------------------------------------
670 : !
671 : ! local variables
672 : !
673 : !-----------------------------------------------------------------------
674 :
675 : integer (int_kind) :: &
676 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
677 : nx, ny ! global dimensions
678 :
679 : logical (log_kind) :: &
680 : special_value
681 :
682 : type (block) :: &
683 : this_block ! block info for current block
684 :
685 : character(len=*), parameter :: subname = '(gather_global_ext_log)'
686 :
687 0 : if (present(spc_val)) then
688 0 : special_value = spc_val
689 : else
690 0 : special_value = .false.
691 : endif
692 0 : ARRAY_G = special_value
693 :
694 0 : nx = nx_global + 2*nghost
695 0 : ny = ny_global + 2*nghost
696 :
697 : !-----------------------------------------------------------------------
698 : !
699 : ! copy local array into block decomposition
700 : !
701 : !-----------------------------------------------------------------------
702 :
703 0 : do n=1,nblocks_tot
704 :
705 0 : this_block = get_block(n,n)
706 :
707 : !*** copy local blocks
708 :
709 0 : if (src_dist%blockLocation(n) /= 0) then
710 :
711 0 : do j=this_block%jlo,this_block%jhi
712 0 : do i=this_block%ilo,this_block%ihi
713 : ARRAY_G(this_block%i_glob(i)+nghost, &
714 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
715 0 : ARRAY(i,j,src_dist%blockLocalID(n))
716 : end do
717 : end do
718 :
719 : ! fill ghost cells
720 0 : if (this_block%jblock == 1) then
721 : ! south block
722 0 : do j=1, nghost
723 0 : do i=this_block%ilo,this_block%ihi
724 : ARRAY_G(this_block%i_glob(i)+nghost,j) = &
725 0 : ARRAY (i,j,src_dist%blockLocalID(n))
726 : end do
727 : end do
728 0 : if (this_block%iblock == 1) then
729 : ! southwest corner
730 0 : do j=1, nghost
731 0 : do i=1, nghost
732 : ARRAY_G(i,j) = &
733 0 : ARRAY (i,j,src_dist%blockLocalID(n))
734 : end do
735 : end do
736 : endif
737 : endif
738 0 : if (this_block%jblock == nblocks_y) then
739 : ! north block
740 0 : do j=1, nghost
741 0 : do i=this_block%ilo,this_block%ihi
742 : ARRAY_G(this_block%i_glob(i)+nghost, &
743 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
744 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
745 : end do
746 : end do
747 0 : if (this_block%iblock == nblocks_x) then
748 : ! northeast corner
749 0 : do j=1, nghost
750 0 : do i=1, nghost
751 : ARRAY_G(nx-i+1, ny-j+1) = &
752 : ARRAY (this_block%ihi+nghost-i+1, & ! LCOV_EXCL_LINE
753 : this_block%jhi+nghost-j+1, & ! LCOV_EXCL_LINE
754 0 : src_dist%blockLocalID(n))
755 : end do
756 : end do
757 : endif
758 : endif
759 0 : if (this_block%iblock == 1) then
760 : ! west block
761 0 : do j=this_block%jlo,this_block%jhi
762 0 : do i=1, nghost
763 : ARRAY_G(i,this_block%j_glob(j)+nghost) = &
764 0 : ARRAY (i,j,src_dist%blockLocalID(n))
765 : end do
766 : end do
767 0 : if (this_block%jblock == nblocks_y) then
768 : ! northwest corner
769 0 : do j=1, nghost
770 0 : do i=1, nghost
771 : ARRAY_G(i, ny-j+1) = &
772 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
773 : end do
774 : end do
775 : endif
776 : endif
777 0 : if (this_block%iblock == nblocks_x) then
778 : ! east block
779 0 : do j=this_block%jlo,this_block%jhi
780 0 : do i=1, nghost
781 : ARRAY_G(nx_global + nghost + i, &
782 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
783 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
784 : end do
785 : end do
786 0 : if (this_block%jblock == 1) then
787 : ! southeast corner
788 0 : do j=1, nghost
789 0 : do i=1, nghost
790 : ARRAY_G( nx-i+1,j) = &
791 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
792 : end do
793 : end do
794 : endif
795 : endif
796 :
797 : endif ! src_dist%blockLocation
798 :
799 : end do
800 :
801 : !-----------------------------------------------------------------------
802 :
803 0 : end subroutine gather_global_ext_log
804 :
805 : !***********************************************************************
806 :
807 228 : subroutine scatter_global_dbl(ARRAY, ARRAY_G, src_task, dst_dist, &
808 : field_loc, field_type)
809 :
810 : ! This subroutine scatters a global-sized array on the processor
811 : ! src\_task to a distribution of blocks given by dst\_dist.
812 : !
813 : ! This is the specific interface for double precision arrays
814 : ! corresponding to the generic interface scatter_global. It is shown
815 : ! to provide information on the generic interface (the generic
816 : ! interface is identical, but chooses a specific interface based
817 : ! on the data type of the input argument).
818 :
819 : integer (int_kind), intent(in) :: &
820 : src_task ! task from which array should be scattered
821 :
822 : type (distrb), intent(in) :: &
823 : dst_dist ! distribution of resulting blocks
824 :
825 : real (dbl_kind), dimension(:,:), intent(in) :: &
826 : ARRAY_G ! array containing global field on src_task
827 :
828 : integer (int_kind), intent(in) :: &
829 : field_type, &! id for type of field (scalar, vector, angle) ! LCOV_EXCL_LINE
830 : field_loc ! id for location on horizontal grid
831 : ! (center, NEcorner, Nface, Eface)
832 :
833 : real (dbl_kind), dimension(:,:,:), intent(inout) :: &
834 : ARRAY ! array containing distributed field
835 :
836 : !-----------------------------------------------------------------------
837 : !
838 : ! local variables
839 : !
840 : !-----------------------------------------------------------------------
841 :
842 : integer (int_kind) :: &
843 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
844 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
845 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
846 : yoffset2, &! ! LCOV_EXCL_LINE
847 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
848 : dst_block ! local block index in dest distribution
849 :
850 : type (block) :: &
851 : this_block ! block info for current block
852 :
853 : character(len=*), parameter :: subname = '(scatter_global_dbl)'
854 :
855 : !-----------------------------------------------------------------------
856 : !
857 : ! initialize return array to zero and set up tripole quantities
858 : !
859 : !-----------------------------------------------------------------------
860 :
861 2771568 : ARRAY = c0
862 :
863 228 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
864 228 : if (this_block%tripoleTFlag) then
865 0 : select case (field_loc)
866 : case (field_loc_center) ! cell center location
867 0 : xoffset = 2
868 0 : yoffset = 0
869 : case (field_loc_NEcorner) ! cell corner (velocity) location
870 0 : xoffset = 1
871 0 : yoffset = -1
872 : case (field_loc_Eface) ! cell face location
873 0 : xoffset = 1
874 0 : yoffset = 0
875 : case (field_loc_Nface) ! cell face location
876 0 : xoffset = 2
877 0 : yoffset = -1
878 : case (field_loc_noupdate) ! ghost cells never used - use cell center
879 0 : xoffset = 1
880 0 : yoffset = 1
881 : end select
882 : else
883 421 : select case (field_loc)
884 : case (field_loc_center) ! cell center location
885 193 : xoffset = 1
886 193 : yoffset = 1
887 : case (field_loc_NEcorner) ! cell corner (velocity) location
888 33 : xoffset = 0
889 33 : yoffset = 0
890 : case (field_loc_Eface) ! cell face location
891 1 : xoffset = 0
892 1 : yoffset = 1
893 : case (field_loc_Nface) ! cell face location
894 1 : xoffset = 1
895 1 : yoffset = 0
896 : case (field_loc_noupdate) ! ghost cells never used - use cell center
897 0 : xoffset = 1
898 228 : yoffset = 1
899 : end select
900 : endif
901 :
902 451 : select case (field_type)
903 : case (field_type_scalar)
904 223 : isign = 1
905 : case (field_type_vector)
906 4 : isign = -1
907 : case (field_type_angle)
908 1 : isign = -1
909 : case (field_type_noupdate) ! ghost cells not needed - use cell center
910 0 : isign = 1
911 : case default
912 228 : call abort_ice(subname//'ERROR: Unknown field type in scatter')
913 : end select
914 :
915 : !-----------------------------------------------------------------------
916 : !
917 : ! copy blocks of global array into local block distribution
918 : !
919 : !-----------------------------------------------------------------------
920 :
921 456 : do n=1,nblocks_tot
922 :
923 456 : if (dst_dist%blockLocation(n) /= 0) then
924 :
925 228 : this_block = get_block(n,n)
926 228 : dst_block = dst_dist%blockLocalID(n)
927 :
928 : !*** if this is an interior block, then there is no
929 : !*** padding or update checking required
930 :
931 : if (this_block%iblock > 1 .and. &
932 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
933 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
934 : this_block%jblock < nblocks_y) then
935 :
936 0 : do j=1,ny_block
937 0 : do i=1,nx_block
938 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
939 0 : this_block%j_glob(j))
940 : end do
941 : end do
942 :
943 : !*** if this is an edge block but not a northern edge
944 : !*** we only need to check for closed boundaries and
945 : !*** padding (global index = 0)
946 :
947 228 : else if (this_block%jblock /= nblocks_y) then
948 :
949 0 : do j=1,ny_block
950 0 : if (this_block%j_glob(j) /= 0) then
951 0 : do i=1,nx_block
952 0 : if (this_block%i_glob(i) /= 0) then
953 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
954 0 : this_block%j_glob(j))
955 : endif
956 : end do
957 : endif
958 : end do
959 :
960 : !*** if this is a northern edge block, we need to check
961 : !*** for and properly deal with tripole boundaries
962 :
963 : else
964 :
965 27132 : do j=1,ny_block
966 27132 : if (this_block%j_glob(j) > 0) then ! normal boundary
967 :
968 2771112 : do i=1,nx_block
969 2771112 : if (this_block%i_glob(i) /= 0) then
970 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
971 2744208 : this_block%j_glob(j))
972 : endif
973 : end do
974 :
975 0 : else if (this_block%j_glob(j) < 0) then ! tripole
976 :
977 : ! for yoffset=0 or 1, yoffset2=0,0
978 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
979 0 : do yoffset2=0,max(yoffset,0)-yoffset
980 : jsrc = ny_global + yoffset + yoffset2 + &
981 0 : (this_block%j_glob(j) + ny_global)
982 0 : do i=1,nx_block
983 0 : if (this_block%i_glob(i) /= 0) then
984 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
985 0 : if (isrc < 1) isrc = isrc + nx_global
986 0 : if (isrc > nx_global) isrc = isrc - nx_global
987 : ARRAY(i,j-yoffset2,dst_block) &
988 0 : = isign * ARRAY_G(isrc,jsrc)
989 : endif
990 : end do
991 : end do
992 :
993 : endif
994 : end do
995 :
996 : endif
997 : endif ! dst block not land
998 : end do ! block loop
999 :
1000 : !-----------------------------------------------------------------
1001 : ! Ensure unused ghost cell values are 0
1002 : !-----------------------------------------------------------------
1003 :
1004 228 : if (field_loc == field_loc_noupdate) then
1005 0 : do n=1,nblocks_tot
1006 0 : dst_block = dst_dist%blockLocalID(n)
1007 0 : this_block = get_block(n,n)
1008 :
1009 0 : if (dst_block > 0) then
1010 :
1011 : ! north edge
1012 0 : do j = this_block%jhi+1,ny_block
1013 0 : do i = 1, nx_block
1014 0 : ARRAY (i,j,dst_block) = c0
1015 : enddo
1016 : enddo
1017 : ! east edge
1018 0 : do j = 1, ny_block
1019 0 : do i = this_block%ihi+1,nx_block
1020 0 : ARRAY (i,j,dst_block) = c0
1021 : enddo
1022 : enddo
1023 : ! south edge
1024 0 : do j = 1, this_block%jlo-1
1025 0 : do i = 1, nx_block
1026 0 : ARRAY (i,j,dst_block) = c0
1027 : enddo
1028 : enddo
1029 : ! west edge
1030 0 : do j = 1, ny_block
1031 0 : do i = 1, this_block%ilo-1
1032 0 : ARRAY (i,j,dst_block) = c0
1033 : enddo
1034 : enddo
1035 :
1036 : endif
1037 : enddo
1038 : endif
1039 :
1040 : !-----------------------------------------------------------------------
1041 :
1042 228 : end subroutine scatter_global_dbl
1043 :
1044 : !***********************************************************************
1045 :
1046 0 : subroutine scatter_global_real(ARRAY, ARRAY_G, src_task, dst_dist, &
1047 : field_loc, field_type)
1048 :
1049 : !-----------------------------------------------------------------------
1050 : !
1051 : ! This subroutine scatters a global-sized array on the processor
1052 : ! src\_task to a distribution of blocks given by dst\_dist.
1053 : !
1054 : !-----------------------------------------------------------------------
1055 : !-----------------------------------------------------------------------
1056 : !
1057 : ! input variables
1058 : !
1059 : !-----------------------------------------------------------------------
1060 :
1061 : integer (int_kind), intent(in) :: &
1062 : src_task ! task from which array should be scattered
1063 :
1064 : type (distrb), intent(in) :: &
1065 : dst_dist ! distribution of resulting blocks
1066 :
1067 : real (real_kind), dimension(:,:), intent(in) :: &
1068 : ARRAY_G ! array containing global field on src_task
1069 :
1070 : integer (int_kind), intent(in) :: &
1071 : field_type, &! id for type of field (scalar, vector, angle) ! LCOV_EXCL_LINE
1072 : field_loc ! id for location on horizontal grid
1073 : ! (center, NEcorner, Nface, Eface)
1074 :
1075 : !-----------------------------------------------------------------------
1076 : !
1077 : ! output variables
1078 : !
1079 : !-----------------------------------------------------------------------
1080 :
1081 : real (real_kind), dimension(:,:,:), intent(inout) :: &
1082 : ARRAY ! array containing distributed field
1083 :
1084 : !-----------------------------------------------------------------------
1085 : !
1086 : ! local variables
1087 : !
1088 : !-----------------------------------------------------------------------
1089 :
1090 : integer (int_kind) :: &
1091 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
1092 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
1093 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
1094 : yoffset2, &! ! LCOV_EXCL_LINE
1095 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
1096 : dst_block ! local block index in dest distribution
1097 :
1098 : type (block) :: &
1099 : this_block ! block info for current block
1100 :
1101 : character(len=*), parameter :: subname = '(scatter_global_real)'
1102 :
1103 : !-----------------------------------------------------------------------
1104 : !
1105 : ! initialize return array to zero and set up tripole quantities
1106 : !
1107 : !-----------------------------------------------------------------------
1108 :
1109 0 : ARRAY = 0._real_kind
1110 :
1111 0 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
1112 0 : if (this_block%tripoleTFlag) then
1113 0 : select case (field_loc)
1114 : case (field_loc_center) ! cell center location
1115 0 : xoffset = 2
1116 0 : yoffset = 0
1117 : case (field_loc_NEcorner) ! cell corner (velocity) location
1118 0 : xoffset = 1
1119 0 : yoffset = 1
1120 : case (field_loc_Eface) ! cell face location
1121 0 : xoffset = 1
1122 0 : yoffset = 0
1123 : case (field_loc_Nface) ! cell face location
1124 0 : xoffset = 2
1125 0 : yoffset = 1
1126 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1127 0 : xoffset = 1
1128 0 : yoffset = 1
1129 : end select
1130 : else
1131 0 : select case (field_loc)
1132 : case (field_loc_center) ! cell center location
1133 0 : xoffset = 1
1134 0 : yoffset = 1
1135 : case (field_loc_NEcorner) ! cell corner (velocity) location
1136 0 : xoffset = 0
1137 0 : yoffset = 0
1138 : case (field_loc_Eface) ! cell face location
1139 0 : xoffset = 0
1140 0 : yoffset = 1
1141 : case (field_loc_Nface) ! cell face location
1142 0 : xoffset = 1
1143 0 : yoffset = 0
1144 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1145 0 : xoffset = 1
1146 0 : yoffset = 1
1147 : end select
1148 : endif
1149 :
1150 0 : select case (field_type)
1151 : case (field_type_scalar)
1152 0 : isign = 1
1153 : case (field_type_vector)
1154 0 : isign = -1
1155 : case (field_type_angle)
1156 0 : isign = -1
1157 : case (field_type_noupdate) ! ghost cells not needed - use cell center
1158 0 : isign = 1
1159 : case default
1160 0 : call abort_ice(subname//'ERROR: Unknown field type in scatter')
1161 : end select
1162 :
1163 : !-----------------------------------------------------------------------
1164 : !
1165 : ! copy blocks of global array into local block distribution
1166 : !
1167 : !-----------------------------------------------------------------------
1168 :
1169 0 : do n=1,nblocks_tot
1170 :
1171 0 : if (dst_dist%blockLocation(n) /= 0) then
1172 :
1173 0 : this_block = get_block(n,n)
1174 0 : dst_block = dst_dist%blockLocalID(n)
1175 :
1176 : !*** if this is an interior block, then there is no
1177 : !*** padding or update checking required
1178 :
1179 : if (this_block%iblock > 1 .and. &
1180 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
1181 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
1182 : this_block%jblock < nblocks_y) then
1183 :
1184 0 : do j=1,ny_block
1185 0 : do i=1,nx_block
1186 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1187 0 : this_block%j_glob(j))
1188 : end do
1189 : end do
1190 :
1191 : !*** if this is an edge block but not a northern edge
1192 : !*** we only need to check for closed boundaries and
1193 : !*** padding (global index = 0)
1194 :
1195 0 : else if (this_block%jblock /= nblocks_y) then
1196 :
1197 0 : do j=1,ny_block
1198 0 : if (this_block%j_glob(j) /= 0) then
1199 0 : do i=1,nx_block
1200 0 : if (this_block%i_glob(i) /= 0) then
1201 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1202 0 : this_block%j_glob(j))
1203 : endif
1204 : end do
1205 : endif
1206 : end do
1207 :
1208 : !*** if this is a northern edge block, we need to check
1209 : !*** for and properly deal with tripole boundaries
1210 :
1211 : else
1212 :
1213 0 : do j=1,ny_block
1214 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
1215 :
1216 0 : do i=1,nx_block
1217 0 : if (this_block%i_glob(i) /= 0) then
1218 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1219 0 : this_block%j_glob(j))
1220 : endif
1221 : end do
1222 :
1223 0 : else if (this_block%j_glob(j) < 0) then ! tripole
1224 :
1225 : ! for yoffset=0 or 1, yoffset2=0,0
1226 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
1227 0 : do yoffset2=0,max(yoffset,0)-yoffset
1228 : jsrc = ny_global + yoffset + yoffset2 + &
1229 0 : (this_block%j_glob(j) + ny_global)
1230 0 : do i=1,nx_block
1231 0 : if (this_block%i_glob(i) /= 0) then
1232 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
1233 0 : if (isrc < 1) isrc = isrc + nx_global
1234 0 : if (isrc > nx_global) isrc = isrc - nx_global
1235 : ARRAY(i,j-yoffset2,dst_block) &
1236 0 : = isign * ARRAY_G(isrc,jsrc)
1237 : endif
1238 : end do
1239 : end do
1240 :
1241 : endif
1242 : end do
1243 :
1244 : endif
1245 : endif ! dst block not land
1246 : end do ! block loop
1247 :
1248 : !-----------------------------------------------------------------
1249 : ! Ensure unused ghost cell values are 0
1250 : !-----------------------------------------------------------------
1251 :
1252 0 : if (field_loc == field_loc_noupdate) then
1253 0 : do n=1,nblocks_tot
1254 0 : dst_block = dst_dist%blockLocalID(n)
1255 0 : this_block = get_block(n,n)
1256 :
1257 0 : if (dst_block > 0) then
1258 :
1259 : ! north edge
1260 0 : do j = this_block%jhi+1,ny_block
1261 0 : do i = 1, nx_block
1262 0 : ARRAY (i,j,dst_block) = 0._real_kind
1263 : enddo
1264 : enddo
1265 : ! east edge
1266 0 : do j = 1, ny_block
1267 0 : do i = this_block%ihi+1,nx_block
1268 0 : ARRAY (i,j,dst_block) = 0._real_kind
1269 : enddo
1270 : enddo
1271 : ! south edge
1272 0 : do j = 1, this_block%jlo-1
1273 0 : do i = 1, nx_block
1274 0 : ARRAY (i,j,dst_block) = 0._real_kind
1275 : enddo
1276 : enddo
1277 : ! west edge
1278 0 : do j = 1, ny_block
1279 0 : do i = 1, this_block%ilo-1
1280 0 : ARRAY (i,j,dst_block) = 0._real_kind
1281 : enddo
1282 : enddo
1283 :
1284 : endif
1285 : enddo
1286 : endif
1287 :
1288 : !-----------------------------------------------------------------------
1289 :
1290 0 : end subroutine scatter_global_real
1291 :
1292 : !***********************************************************************
1293 :
1294 0 : subroutine scatter_global_int(ARRAY, ARRAY_G, src_task, dst_dist, &
1295 : field_loc, field_type)
1296 :
1297 : !-----------------------------------------------------------------------
1298 : !
1299 : ! This subroutine scatters a global-sized array on the processor
1300 : ! src\_task to a distribution of blocks given by dst\_dist.
1301 : !
1302 : !-----------------------------------------------------------------------
1303 : !-----------------------------------------------------------------------
1304 : !
1305 : ! input variables
1306 : !
1307 : !-----------------------------------------------------------------------
1308 :
1309 : integer (int_kind), intent(in) :: &
1310 : src_task ! task from which array should be scattered
1311 :
1312 : type (distrb), intent(in) :: &
1313 : dst_dist ! distribution of resulting blocks
1314 :
1315 : integer (int_kind), dimension(:,:), intent(in) :: &
1316 : ARRAY_G ! array containing global field on src_task
1317 :
1318 : integer (int_kind), intent(in) :: &
1319 : field_type, &! id for type of field (scalar, vector, angle) ! LCOV_EXCL_LINE
1320 : field_loc ! id for location on horizontal grid
1321 : ! (center, NEcorner, Nface, Eface)
1322 :
1323 : !-----------------------------------------------------------------------
1324 : !
1325 : ! output variables
1326 : !
1327 : !-----------------------------------------------------------------------
1328 :
1329 : integer (int_kind), dimension(:,:,:), intent(inout) :: &
1330 : ARRAY ! array containing distributed field
1331 :
1332 : !-----------------------------------------------------------------------
1333 : !
1334 : ! local variables
1335 : !
1336 : !-----------------------------------------------------------------------
1337 :
1338 : integer (int_kind) :: &
1339 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
1340 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
1341 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
1342 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
1343 : yoffset2, &! ! LCOV_EXCL_LINE
1344 : dst_block ! local block index in dest distribution
1345 :
1346 : type (block) :: &
1347 : this_block ! block info for current block
1348 :
1349 : character(len=*), parameter :: subname = '(scatter_global_int)'
1350 :
1351 : !-----------------------------------------------------------------------
1352 : !
1353 : ! initialize return array to zero and set up tripole quantities
1354 : !
1355 : !-----------------------------------------------------------------------
1356 :
1357 0 : ARRAY = 0
1358 :
1359 0 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
1360 0 : if (this_block%tripoleTFlag) then
1361 0 : select case (field_loc)
1362 : case (field_loc_center) ! cell center location
1363 0 : xoffset = 2
1364 0 : yoffset = 0
1365 : case (field_loc_NEcorner) ! cell corner (velocity) location
1366 0 : xoffset = 1
1367 0 : yoffset = 1
1368 : case (field_loc_Eface) ! cell face location
1369 0 : xoffset = 1
1370 0 : yoffset = 0
1371 : case (field_loc_Nface) ! cell face location
1372 0 : xoffset = 2
1373 0 : yoffset = 1
1374 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1375 0 : xoffset = 1
1376 0 : yoffset = 1
1377 : end select
1378 : else
1379 0 : select case (field_loc)
1380 : case (field_loc_center) ! cell center location
1381 0 : xoffset = 1
1382 0 : yoffset = 1
1383 : case (field_loc_NEcorner) ! cell corner (velocity) location
1384 0 : xoffset = 0
1385 0 : yoffset = 0
1386 : case (field_loc_Eface) ! cell face location
1387 0 : xoffset = 0
1388 0 : yoffset = 1
1389 : case (field_loc_Nface) ! cell face location
1390 0 : xoffset = 1
1391 0 : yoffset = 0
1392 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1393 0 : xoffset = 1
1394 0 : yoffset = 1
1395 : end select
1396 : endif
1397 :
1398 0 : select case (field_type)
1399 : case (field_type_scalar)
1400 0 : isign = 1
1401 : case (field_type_vector)
1402 0 : isign = -1
1403 : case (field_type_angle)
1404 0 : isign = -1
1405 : case (field_type_noupdate) ! ghost cells not needed - use cell center
1406 0 : isign = 1
1407 : case default
1408 0 : call abort_ice(subname//'ERROR: Unknown field type in scatter')
1409 : end select
1410 :
1411 : !-----------------------------------------------------------------------
1412 : !
1413 : ! copy blocks of global array into local block distribution
1414 : !
1415 : !-----------------------------------------------------------------------
1416 :
1417 0 : do n=1,nblocks_tot
1418 :
1419 0 : if (dst_dist%blockLocation(n) /= 0) then
1420 :
1421 0 : this_block = get_block(n,n)
1422 0 : dst_block = dst_dist%blockLocalID(n)
1423 :
1424 : !*** if this is an interior block, then there is no
1425 : !*** padding or update checking required
1426 :
1427 : if (this_block%iblock > 1 .and. &
1428 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
1429 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
1430 : this_block%jblock < nblocks_y) then
1431 :
1432 0 : do j=1,ny_block
1433 0 : do i=1,nx_block
1434 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1435 0 : this_block%j_glob(j))
1436 : end do
1437 : end do
1438 :
1439 : !*** if this is an edge block but not a northern edge
1440 : !*** we only need to check for closed boundaries and
1441 : !*** padding (global index = 0)
1442 :
1443 0 : else if (this_block%jblock /= nblocks_y) then
1444 :
1445 0 : do j=1,ny_block
1446 0 : if (this_block%j_glob(j) /= 0) then
1447 0 : do i=1,nx_block
1448 0 : if (this_block%i_glob(i) /= 0) then
1449 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1450 0 : this_block%j_glob(j))
1451 : endif
1452 : end do
1453 : endif
1454 : end do
1455 :
1456 : !*** if this is a northern edge block, we need to check
1457 : !*** for and properly deal with tripole boundaries
1458 :
1459 : else
1460 :
1461 0 : do j=1,ny_block
1462 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
1463 :
1464 0 : do i=1,nx_block
1465 0 : if (this_block%i_glob(i) /= 0) then
1466 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1467 0 : this_block%j_glob(j))
1468 : endif
1469 : end do
1470 :
1471 0 : else if (this_block%j_glob(j) < 0) then ! tripole
1472 :
1473 : ! for yoffset=0 or 1, yoffset2=0,0
1474 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
1475 0 : do yoffset2=0,max(yoffset,0)-yoffset
1476 : jsrc = ny_global + yoffset + yoffset2 + &
1477 0 : (this_block%j_glob(j) + ny_global)
1478 0 : do i=1,nx_block
1479 0 : if (this_block%i_glob(i) /= 0) then
1480 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
1481 0 : if (isrc < 1) isrc = isrc + nx_global
1482 0 : if (isrc > nx_global) isrc = isrc - nx_global
1483 : ARRAY(i,j-yoffset2,dst_block) &
1484 0 : = isign * ARRAY_G(isrc,jsrc)
1485 : endif
1486 : end do
1487 : end do
1488 :
1489 : endif
1490 : end do
1491 :
1492 : endif
1493 : endif ! dst block not land
1494 : end do ! block loop
1495 :
1496 : !-----------------------------------------------------------------
1497 : ! Ensure unused ghost cell values are 0
1498 : !-----------------------------------------------------------------
1499 :
1500 0 : if (field_loc == field_loc_noupdate) then
1501 0 : do n=1,nblocks_tot
1502 0 : dst_block = dst_dist%blockLocalID(n)
1503 0 : this_block = get_block(n,n)
1504 :
1505 0 : if (dst_block > 0) then
1506 :
1507 : ! north edge
1508 0 : do j = this_block%jhi+1,ny_block
1509 0 : do i = 1, nx_block
1510 0 : ARRAY (i,j,dst_block) = 0
1511 : enddo
1512 : enddo
1513 : ! east edge
1514 0 : do j = 1, ny_block
1515 0 : do i = this_block%ihi+1,nx_block
1516 0 : ARRAY (i,j,dst_block) = 0
1517 : enddo
1518 : enddo
1519 : ! south edge
1520 0 : do j = 1, this_block%jlo-1
1521 0 : do i = 1, nx_block
1522 0 : ARRAY (i,j,dst_block) = 0
1523 : enddo
1524 : enddo
1525 : ! west edge
1526 0 : do j = 1, ny_block
1527 0 : do i = 1, this_block%ilo-1
1528 0 : ARRAY (i,j,dst_block) = 0
1529 : enddo
1530 : enddo
1531 :
1532 : endif
1533 : enddo
1534 : endif
1535 :
1536 : !-----------------------------------------------------------------------
1537 :
1538 0 : end subroutine scatter_global_int
1539 :
1540 : !***********************************************************************
1541 :
1542 0 : subroutine scatter_global_ext(ARRAY, ARRAY_G, src_task, dst_dist)
1543 :
1544 : ! This subroutine scatters a global-sized array on the processor
1545 : ! src\_task to a distribution of blocks given by dst\_dist.
1546 : !
1547 : ! This is the specific interface for double precision arrays
1548 : ! corresponding to the generic interface scatter_global. It is shown
1549 : ! to provide information on the generic interface (the generic
1550 : ! interface is identical, but chooses a specific interface based
1551 : ! on the data type of the input argument).
1552 :
1553 : integer (int_kind), intent(in) :: &
1554 : src_task ! task from which array should be scattered
1555 :
1556 : type (distrb), intent(in) :: &
1557 : dst_dist ! distribution of resulting blocks
1558 :
1559 : real (dbl_kind), dimension(:,:), intent(in) :: &
1560 : ARRAY_G ! array containing global field on src_task
1561 :
1562 : real (dbl_kind), dimension(:,:,:), intent(inout) :: &
1563 : ARRAY ! array containing distributed field
1564 :
1565 : !-----------------------------------------------------------------------
1566 : !
1567 : ! local variables
1568 : !
1569 : !-----------------------------------------------------------------------
1570 :
1571 : integer (int_kind) :: &
1572 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
1573 : iblk, jblk, &! source addresses ! LCOV_EXCL_LINE
1574 : iglb, jglb, &! global indices ! LCOV_EXCL_LINE
1575 : dst_block ! local block index in dest distribution
1576 :
1577 : type (block) :: &
1578 : this_block ! block info for current block
1579 :
1580 : character(len=*), parameter :: subname = '(scatter_global_ext)'
1581 :
1582 : !-----------------------------------------------------------------------
1583 : !
1584 : ! initialize return array to zero
1585 : !
1586 : !-----------------------------------------------------------------------
1587 :
1588 0 : ARRAY = c0
1589 :
1590 : !-----------------------------------------------------------------------
1591 : !
1592 : ! copy blocks of global array into local block distribution
1593 : !
1594 : !-----------------------------------------------------------------------
1595 :
1596 0 : do n=1,nblocks_tot
1597 :
1598 0 : if (dst_dist%blockLocation(n) /= 0) then
1599 :
1600 0 : this_block = get_block(n,n)
1601 0 : dst_block = dst_dist%blockLocalID(n)
1602 :
1603 : ! interior
1604 0 : do j = 1, ny_block
1605 0 : do i = 1, nx_block
1606 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i)+nghost,&
1607 0 : this_block%j_glob(j)+nghost)
1608 : end do
1609 : end do
1610 :
1611 0 : if (this_block%jblock == 1) then
1612 : ! south edge
1613 0 : do j = 1, nghost
1614 0 : do i = this_block%ilo,this_block%ihi
1615 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i)+nghost,j)
1616 : enddo
1617 0 : do i = 1, nghost
1618 : ! southwest corner
1619 0 : iblk = i
1620 0 : jblk = j
1621 0 : iglb = this_block%i_glob(this_block%ilo)+i-1
1622 0 : jglb = j
1623 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1624 : ! southeast corner
1625 0 : iblk = this_block%ihi+i
1626 0 : iglb = this_block%i_glob(this_block%ihi)+nghost+i
1627 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1628 : enddo
1629 : enddo
1630 : endif
1631 0 : if (this_block%jblock == nblocks_y) then
1632 : ! north edge
1633 0 : do j = 1, nghost
1634 0 : do i = this_block%ilo,this_block%ihi
1635 : ARRAY(i,this_block%jhi+j,dst_block) = ARRAY_G(this_block%i_glob(i)+nghost,&
1636 0 : ny_global+nghost+j)
1637 : enddo
1638 0 : do i = 1, nghost
1639 : ! northwest corner
1640 0 : iblk = i
1641 0 : jblk = this_block%jhi+j
1642 0 : iglb = this_block%i_glob(this_block%ilo)+i-1
1643 0 : jglb = ny_global+nghost+j
1644 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1645 : ! northeast corner
1646 0 : iblk = this_block%ihi+i
1647 0 : iglb = this_block%i_glob(this_block%ihi)+nghost+i
1648 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1649 : enddo
1650 : enddo
1651 : endif
1652 0 : if (this_block%iblock == 1) then
1653 : ! west edge
1654 0 : do j = this_block%jlo,this_block%jhi
1655 0 : do i = 1, nghost
1656 0 : ARRAY(i,j,dst_block) = ARRAY_G(i,this_block%j_glob(j)+nghost)
1657 : enddo
1658 : enddo
1659 0 : do j = 1, nghost
1660 0 : do i = 1, nghost
1661 : ! northwest corner
1662 0 : iblk = i
1663 0 : jblk = this_block%jhi+j
1664 0 : iglb = i
1665 0 : jglb = this_block%j_glob(this_block%jhi)+nghost+j
1666 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1667 : ! southwest corner
1668 0 : jblk = j
1669 0 : jglb = this_block%j_glob(this_block%jlo)+j-1
1670 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1671 : enddo
1672 : enddo
1673 : endif
1674 0 : if (this_block%iblock == nblocks_x) then
1675 : ! east edge
1676 0 : do j = this_block%jlo,this_block%jhi
1677 0 : do i = 1, nghost
1678 : ARRAY(this_block%ihi+i,j,dst_block) = ARRAY_G(nx_global+nghost+i, &
1679 0 : this_block%j_glob(j)+nghost)
1680 : enddo
1681 : enddo
1682 0 : do j = 1, nghost
1683 0 : do i = 1, nghost
1684 : ! northeast corner
1685 0 : iblk = this_block%ihi+i
1686 0 : jblk = this_block%jhi+j
1687 0 : iglb = nx_global+nghost+i
1688 0 : jglb = this_block%j_glob(this_block%jhi)+nghost+j
1689 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1690 : ! southeast corner
1691 0 : jblk = j
1692 0 : jglb = this_block%j_glob(this_block%jlo)+j-1
1693 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
1694 : enddo
1695 : enddo
1696 : endif
1697 :
1698 : endif ! dst block not land
1699 : end do ! block loop
1700 :
1701 : !-----------------------------------------------------------------------
1702 :
1703 0 : end subroutine scatter_global_ext
1704 :
1705 : !***********************************************************************
1706 :
1707 0 : subroutine scatter_global_stress(ARRAY, ARRAY_G1, ARRAY_G2, &
1708 : src_task, dst_dist)
1709 :
1710 : ! This subroutine scatters global stresses to a distributed array.
1711 : !
1712 : ! Ghost cells in the stress tensor must be handled separately on tripole
1713 : ! grids, because matching the corner values requires 2 different arrays.
1714 :
1715 : integer (int_kind), intent(in) :: &
1716 : src_task ! task from which array should be scattered
1717 :
1718 : type (distrb), intent(in) :: &
1719 : dst_dist ! distribution of resulting blocks
1720 :
1721 : real (dbl_kind), dimension(:,:), intent(in) :: &
1722 : ARRAY_G1, &! array containing global field on src_task ! LCOV_EXCL_LINE
1723 : ARRAY_G2 ! array containing global field on src_task
1724 :
1725 : real (dbl_kind), dimension(:,:,:), intent(inout) :: &
1726 : ARRAY ! array containing distributed field
1727 :
1728 : !-----------------------------------------------------------------------
1729 : !
1730 : ! local variables
1731 : !
1732 : !-----------------------------------------------------------------------
1733 :
1734 : integer (int_kind) :: &
1735 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
1736 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
1737 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
1738 : yoffset2, &! ! LCOV_EXCL_LINE
1739 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
1740 : dst_block ! local block index in dest distribution
1741 :
1742 : type (block) :: &
1743 : this_block ! block info for current block
1744 :
1745 : character(len=*), parameter :: subname = '(scatter_global_stress)'
1746 :
1747 : !-----------------------------------------------------------------------
1748 : !
1749 : ! initialize return array to zero and set up tripole quantities
1750 : !
1751 : !-----------------------------------------------------------------------
1752 :
1753 0 : ARRAY = c0
1754 :
1755 0 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
1756 0 : if (this_block%tripoleTFlag) then
1757 0 : xoffset = 2 ! treat stresses as cell-centered scalars (they are not
1758 0 : yoffset = 0 ! shared with neighboring grid cells)
1759 : else
1760 0 : xoffset = 1 ! treat stresses as cell-centered scalars (they are not
1761 0 : yoffset = 1 ! shared with neighboring grid cells)
1762 : endif
1763 0 : isign = 1
1764 :
1765 : !-----------------------------------------------------------------------
1766 : !
1767 : ! copy blocks of global array into local block distribution
1768 : !
1769 : !-----------------------------------------------------------------------
1770 :
1771 0 : do n=1,nblocks_tot
1772 :
1773 0 : if (dst_dist%blockLocation(n) /= 0) then
1774 :
1775 0 : this_block = get_block(n,n)
1776 0 : dst_block = dst_dist%blockLocalID(n)
1777 :
1778 : !*** if this is an interior block, then there is no
1779 : !*** padding or update checking required
1780 :
1781 : if (this_block%iblock > 1 .and. &
1782 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
1783 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
1784 : this_block%jblock < nblocks_y) then
1785 :
1786 0 : do j=1,ny_block
1787 0 : do i=1,nx_block
1788 : ARRAY(i,j,dst_block) = ARRAY_G1(this_block%i_glob(i),&
1789 0 : this_block%j_glob(j))
1790 : end do
1791 : end do
1792 :
1793 : !*** if this is an edge block but not a northern edge
1794 : !*** we only need to check for closed boundaries and
1795 : !*** padding (global index = 0)
1796 :
1797 0 : else if (this_block%jblock /= nblocks_y) then
1798 :
1799 0 : do j=1,ny_block
1800 0 : if (this_block%j_glob(j) /= 0) then
1801 0 : do i=1,nx_block
1802 0 : if (this_block%i_glob(i) /= 0) then
1803 : ARRAY(i,j,dst_block) = ARRAY_G1(this_block%i_glob(i),&
1804 0 : this_block%j_glob(j))
1805 : endif
1806 : end do
1807 : endif
1808 : end do
1809 :
1810 : !*** if this is a northern edge block, we need to check
1811 : !*** for and properly deal with tripole boundaries
1812 :
1813 : else
1814 :
1815 0 : do j=1,ny_block
1816 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
1817 :
1818 0 : do i=1,nx_block
1819 0 : if (this_block%i_glob(i) /= 0) then
1820 : ARRAY(i,j,dst_block) = ARRAY_G1(this_block%i_glob(i),&
1821 0 : this_block%j_glob(j))
1822 : endif
1823 : end do
1824 :
1825 0 : else if (this_block%j_glob(j) < 0) then ! tripole
1826 :
1827 : ! for yoffset=0 or 1, yoffset2=0,0
1828 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
1829 0 : do yoffset2=0,max(yoffset,0)-yoffset
1830 : jsrc = ny_global + yoffset + yoffset2 + &
1831 0 : (this_block%j_glob(j) + ny_global)
1832 0 : do i=1,nx_block
1833 0 : if (this_block%i_glob(i) /= 0) then
1834 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
1835 0 : if (isrc < 1) isrc = isrc + nx_global
1836 0 : if (isrc > nx_global) isrc = isrc - nx_global
1837 : ARRAY(i,j-yoffset2,dst_block) &
1838 0 : = isign * ARRAY_G2(isrc,jsrc)
1839 : endif
1840 : end do
1841 : end do
1842 :
1843 : endif
1844 : end do
1845 :
1846 : endif
1847 : endif ! dst block not land
1848 : end do ! block loop
1849 :
1850 : !-----------------------------------------------------------------------
1851 :
1852 0 : end subroutine scatter_global_stress
1853 :
1854 : !***********************************************************************
1855 :
1856 : end module ice_gather_scatter
1857 :
1858 : !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
|