Line data Source code
1 : !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
2 :
3 : module ice_gather_scatter
4 :
5 : ! This module contains routines for gathering data to a single
6 : ! processor from a distributed array, and scattering data from a
7 : ! single processor to a distributed array.
8 : !
9 : ! NOTE: The arrays gathered and scattered are assumed to have
10 : ! horizontal dimensions (nx_block, ny_block).
11 : !
12 : ! author: Phil Jones, LANL
13 : ! Oct. 2004: Adapted from POP version by William H. Lipscomb, LANL
14 : ! Jan. 2008: Elizabeth Hunke replaced old routines with new POP
15 : ! infrastructure, added specialized routine scatter_global_stress
16 :
17 : use mpi ! MPI Fortran module
18 : use ice_kinds_mod
19 : use ice_communicate, only: my_task, mpiR8, mpiR4, mpitag_gs, MPI_COMM_ICE, &
20 : ice_barrier, add_mpi_barriers
21 : use ice_constants, only: spval_dbl, c0, &
22 : field_loc_center, field_loc_NEcorner, field_loc_Nface, field_loc_Eface, & ! LCOV_EXCL_LINE
23 : field_loc_noupdate, & ! LCOV_EXCL_LINE
24 : field_type_scalar, field_type_vector, field_type_angle, & ! LCOV_EXCL_LINE
25 : field_type_noupdate
26 : use ice_blocks, only: block, nx_block, ny_block, nblocks_tot, get_block, &
27 : nblocks_x, nblocks_y, nghost
28 : use ice_distribution, only: distrb
29 : use ice_domain_size, only: nx_global, ny_global
30 : use ice_exit, only: abort_ice
31 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
32 :
33 : implicit none
34 : private
35 :
36 : public :: gather_global, &
37 : gather_global_ext, & ! LCOV_EXCL_LINE
38 : scatter_global, & ! LCOV_EXCL_LINE
39 : scatter_global_ext, & ! LCOV_EXCL_LINE
40 : scatter_global_stress
41 :
42 : !-----------------------------------------------------------------------
43 : !
44 : ! overload module functions
45 : !
46 : !-----------------------------------------------------------------------
47 :
48 : interface gather_global_ext
49 : module procedure gather_global_ext_dbl, &
50 : gather_global_ext_int, & ! LCOV_EXCL_LINE
51 : gather_global_ext_log
52 : end interface
53 :
54 : interface gather_global
55 : module procedure gather_global_dbl, &
56 : gather_global_real, & ! LCOV_EXCL_LINE
57 : gather_global_int
58 : end interface
59 :
60 : interface scatter_global
61 : module procedure scatter_global_dbl, &
62 : scatter_global_real, & ! LCOV_EXCL_LINE
63 : scatter_global_int
64 : end interface
65 :
66 : !-----------------------------------------------------------------------
67 : !
68 : ! module variables
69 : !
70 : !-----------------------------------------------------------------------
71 :
72 : !***********************************************************************
73 :
74 : contains
75 :
76 : !***********************************************************************
77 :
78 11856 : subroutine gather_global_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
79 :
80 : ! This subroutine gathers a distributed array to a global-sized
81 : ! array on the processor dst_task.
82 : !
83 : ! This is the specific inteface for double precision arrays
84 : ! corresponding to the generic interface gather_global. It is shown
85 : ! to provide information on the generic interface (the generic
86 : ! interface is identical, but chooses a specific inteface based
87 : ! on the data type of the input argument).
88 :
89 : integer (int_kind), intent(in) :: &
90 : dst_task ! task to which array should be gathered
91 :
92 : type (distrb), intent(in) :: &
93 : src_dist ! distribution of blocks in the source array
94 :
95 : real (dbl_kind), dimension(:,:,:), intent(in) :: &
96 : ARRAY ! array containing horizontal slab of distributed field
97 :
98 : real (dbl_kind), intent(in), optional :: &
99 : spc_val
100 :
101 : real (dbl_kind), dimension(:,:), intent(inout) :: &
102 : ARRAY_G ! array containing global horizontal field on dst_task
103 :
104 : !-----------------------------------------------------------------------
105 : !
106 : ! local variables
107 : !
108 : !-----------------------------------------------------------------------
109 :
110 : integer (int_kind) :: &
111 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
112 : nsends ,&! number of actual sends ! LCOV_EXCL_LINE
113 : src_block ,&! block locator for send ! LCOV_EXCL_LINE
114 : ierr ! MPI error flag
115 :
116 : integer (int_kind), dimension(MPI_STATUS_SIZE) :: &
117 : status
118 :
119 : integer (int_kind), dimension(:), allocatable :: &
120 11856 : snd_request
121 :
122 : integer (int_kind), dimension(:,:), allocatable :: &
123 11856 : snd_status
124 :
125 : real (dbl_kind), dimension(:,:), allocatable :: &
126 11856 : msg_buffer
127 :
128 : real (dbl_kind) :: &
129 2856 : special_value
130 :
131 : type (block) :: &
132 : this_block ! block info for current block
133 :
134 : character(len=*), parameter :: subname = '(gather_global_dbl)'
135 :
136 11856 : if (present(spc_val)) then
137 7968 : special_value = spc_val
138 : else
139 3888 : special_value = spval_dbl
140 : endif
141 :
142 : !-----------------------------------------------------------------------
143 : !
144 : ! if this task is the dst_task, copy local blocks into the global
145 : ! array and post receives for non-local blocks.
146 : !
147 : !-----------------------------------------------------------------------
148 :
149 11856 : if (my_task == dst_task) then
150 :
151 51128 : do n=1,nblocks_tot
152 :
153 : !*** copy local blocks
154 :
155 51128 : if (src_dist%blockLocation(n) == my_task+1) then
156 :
157 8120 : this_block = get_block(n,n)
158 :
159 247884 : do j=this_block%jlo,this_block%jhi
160 4972656 : do i=this_block%ilo,this_block%ihi
161 4141200 : ARRAY_G(this_block%i_glob(i), &
162 : this_block%j_glob(j)) = & ! LCOV_EXCL_LINE
163 9105736 : ARRAY(i,j,src_dist%blockLocalID(n))
164 : end do
165 : end do
166 :
167 : !*** fill land blocks with special values
168 :
169 41032 : else if (src_dist%blockLocation(n) == 0) then
170 :
171 274 : this_block = get_block(n,n)
172 :
173 8220 : do j=this_block%jlo,this_block%jhi
174 47950 : do i=this_block%ilo,this_block%ihi
175 0 : ARRAY_G(this_block%i_glob(i), &
176 47676 : this_block%j_glob(j)) = special_value
177 : end do
178 : end do
179 : endif
180 :
181 : end do
182 :
183 : !*** receive blocks to fill up the rest
184 :
185 1976 : allocate (msg_buffer(nx_block,ny_block))
186 :
187 51128 : do n=1,nblocks_tot
188 49152 : if (src_dist%blockLocation(n) > 0 .and. &
189 1976 : src_dist%blockLocation(n) /= my_task+1) then
190 :
191 40758 : this_block = get_block(n,n)
192 :
193 : call MPI_RECV(msg_buffer, size(msg_buffer), &
194 : mpiR8, src_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
195 40758 : MPI_COMM_ICE, status, ierr)
196 :
197 1252728 : do j=this_block%jlo,this_block%jhi
198 22825602 : do i=this_block%ilo,this_block%ihi
199 12423600 : ARRAY_G(this_block%i_glob(i), &
200 35208444 : this_block%j_glob(j)) = msg_buffer(i,j)
201 : end do
202 : end do
203 : endif
204 : end do
205 :
206 1976 : deallocate(msg_buffer)
207 :
208 : !-----------------------------------------------------------------------
209 : !
210 : ! otherwise send data to dst_task
211 : !
212 : !-----------------------------------------------------------------------
213 :
214 : else
215 :
216 0 : allocate(snd_request(nblocks_tot), &
217 9880 : snd_status (MPI_STATUS_SIZE, nblocks_tot))
218 :
219 9880 : nsends = 0
220 290712 : do n=1,nblocks_tot
221 290712 : if (src_dist%blockLocation(n) == my_task+1) then
222 :
223 40758 : nsends = nsends + 1
224 40758 : src_block = src_dist%blockLocalID(n)
225 8568 : call MPI_ISEND(ARRAY(1,1,src_block), nx_block*ny_block, &
226 : mpiR8, dst_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
227 40758 : MPI_COMM_ICE, snd_request(nsends), ierr)
228 : endif
229 : end do
230 :
231 9880 : if (nsends > 0) &
232 9880 : call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
233 9880 : deallocate(snd_request, snd_status)
234 :
235 : endif
236 :
237 11856 : if (add_mpi_barriers) then
238 0 : call ice_barrier()
239 : endif
240 :
241 : !-----------------------------------------------------------------------
242 :
243 11856 : end subroutine gather_global_dbl
244 :
245 : !***********************************************************************
246 :
247 0 : subroutine gather_global_real(ARRAY_G, ARRAY, dst_task, src_dist)
248 :
249 : !-----------------------------------------------------------------------
250 : !
251 : ! This subroutine gathers a distributed array to a global-sized
252 : ! array on the processor dst_task.
253 : !
254 : !-----------------------------------------------------------------------
255 :
256 : !-----------------------------------------------------------------------
257 : !
258 : ! input variables
259 : !
260 : !-----------------------------------------------------------------------
261 :
262 : integer (int_kind), intent(in) :: &
263 : dst_task ! task to which array should be gathered
264 :
265 : type (distrb), intent(in) :: &
266 : src_dist ! distribution of blocks in the source array
267 :
268 : real (real_kind), dimension(:,:,:), intent(in) :: &
269 : ARRAY ! array containing distributed field
270 :
271 : !-----------------------------------------------------------------------
272 : !
273 : ! output variables
274 : !
275 : !-----------------------------------------------------------------------
276 :
277 : real (real_kind), dimension(:,:), intent(inout) :: &
278 : ARRAY_G ! array containing global field on dst_task
279 :
280 : !-----------------------------------------------------------------------
281 : !
282 : ! local variables
283 : !
284 : !-----------------------------------------------------------------------
285 :
286 : integer (int_kind) :: &
287 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
288 : nsends ,&! number of actual sends ! LCOV_EXCL_LINE
289 : src_block ,&! block locator for send ! LCOV_EXCL_LINE
290 : ierr ! MPI error flag
291 :
292 : integer (int_kind), dimension(MPI_STATUS_SIZE) :: &
293 : status
294 :
295 : integer (int_kind), dimension(:), allocatable :: &
296 0 : snd_request
297 :
298 : integer (int_kind), dimension(:,:), allocatable :: &
299 0 : snd_status
300 :
301 : real (real_kind), dimension(:,:), allocatable :: &
302 0 : msg_buffer
303 :
304 : type (block) :: &
305 : this_block ! block info for current block
306 :
307 : character(len=*), parameter :: subname = '(gather_global_real)'
308 :
309 : !-----------------------------------------------------------------------
310 : !
311 : ! if this task is the dst_task, copy local blocks into the global
312 : ! array and post receives for non-local blocks.
313 : !
314 : !-----------------------------------------------------------------------
315 :
316 0 : if (my_task == dst_task) then
317 :
318 0 : do n=1,nblocks_tot
319 :
320 : !*** copy local blocks
321 :
322 0 : if (src_dist%blockLocation(n) == my_task+1) then
323 :
324 0 : this_block = get_block(n,n)
325 :
326 0 : do j=this_block%jlo,this_block%jhi
327 0 : do i=this_block%ilo,this_block%ihi
328 0 : ARRAY_G(this_block%i_glob(i), &
329 : this_block%j_glob(j)) = & ! LCOV_EXCL_LINE
330 0 : ARRAY(i,j,src_dist%blockLocalID(n))
331 : end do
332 : end do
333 :
334 : !*** fill land blocks with zeroes
335 :
336 0 : else if (src_dist%blockLocation(n) == 0) then
337 :
338 0 : this_block = get_block(n,n)
339 :
340 0 : do j=this_block%jlo,this_block%jhi
341 0 : do i=this_block%ilo,this_block%ihi
342 0 : ARRAY_G(this_block%i_glob(i), &
343 0 : this_block%j_glob(j)) = 0._real_kind
344 : end do
345 : end do
346 : endif
347 :
348 : end do
349 :
350 : !*** receive blocks to fill up the rest
351 :
352 0 : allocate (msg_buffer(nx_block,ny_block))
353 :
354 0 : do n=1,nblocks_tot
355 0 : if (src_dist%blockLocation(n) > 0 .and. &
356 0 : src_dist%blockLocation(n) /= my_task+1) then
357 :
358 0 : this_block = get_block(n,n)
359 :
360 : call MPI_RECV(msg_buffer, size(msg_buffer), &
361 : mpiR4, src_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
362 0 : MPI_COMM_ICE, status, ierr)
363 :
364 0 : do j=this_block%jlo,this_block%jhi
365 0 : do i=this_block%ilo,this_block%ihi
366 0 : ARRAY_G(this_block%i_glob(i), &
367 0 : this_block%j_glob(j)) = msg_buffer(i,j)
368 : end do
369 : end do
370 : endif
371 : end do
372 :
373 0 : deallocate(msg_buffer)
374 :
375 : !-----------------------------------------------------------------------
376 : !
377 : ! otherwise send data to dst_task
378 : !
379 : !-----------------------------------------------------------------------
380 :
381 : else
382 :
383 0 : allocate(snd_request(nblocks_tot), &
384 0 : snd_status (MPI_STATUS_SIZE, nblocks_tot))
385 :
386 0 : nsends = 0
387 0 : do n=1,nblocks_tot
388 0 : if (src_dist%blockLocation(n) == my_task+1) then
389 :
390 0 : nsends = nsends + 1
391 0 : src_block = src_dist%blockLocalID(n)
392 0 : call MPI_ISEND(ARRAY(1,1,src_block), nx_block*ny_block, &
393 : mpiR4, dst_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
394 0 : MPI_COMM_ICE, snd_request(nsends), ierr)
395 : endif
396 : end do
397 :
398 0 : if (nsends > 0) &
399 0 : call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
400 0 : deallocate(snd_request, snd_status)
401 :
402 : endif
403 :
404 0 : if (add_mpi_barriers) then
405 0 : call ice_barrier()
406 : endif
407 :
408 : !-----------------------------------------------------------------------
409 :
410 0 : end subroutine gather_global_real
411 :
412 : !***********************************************************************
413 :
414 0 : subroutine gather_global_int(ARRAY_G, ARRAY, dst_task, src_dist)
415 :
416 : !-----------------------------------------------------------------------
417 : !
418 : ! This subroutine gathers a distributed array to a global-sized
419 : ! array on the processor dst_task.
420 : !
421 : !-----------------------------------------------------------------------
422 :
423 : !-----------------------------------------------------------------------
424 : !
425 : ! input variables
426 : !
427 : !-----------------------------------------------------------------------
428 :
429 : integer (int_kind), intent(in) :: &
430 : dst_task ! task to which array should be gathered
431 :
432 : type (distrb), intent(in) :: &
433 : src_dist ! distribution of blocks in the source array
434 :
435 : integer (int_kind), dimension(:,:,:), intent(in) :: &
436 : ARRAY ! array containing distributed field
437 :
438 : !-----------------------------------------------------------------------
439 : !
440 : ! output variables
441 : !
442 : !-----------------------------------------------------------------------
443 :
444 : integer (int_kind), dimension(:,:), intent(inout) :: &
445 : ARRAY_G ! array containing global field on dst_task
446 :
447 : !-----------------------------------------------------------------------
448 : !
449 : ! local variables
450 : !
451 : !-----------------------------------------------------------------------
452 :
453 : integer (int_kind) :: &
454 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
455 : nsends ,&! number of actual sends ! LCOV_EXCL_LINE
456 : src_block ,&! block locator for send ! LCOV_EXCL_LINE
457 : ierr ! MPI error flag
458 :
459 : integer (int_kind), dimension(MPI_STATUS_SIZE) :: &
460 : status
461 :
462 : integer (int_kind), dimension(:), allocatable :: &
463 0 : snd_request
464 :
465 : integer (int_kind), dimension(:,:), allocatable :: &
466 0 : snd_status
467 :
468 : integer (int_kind), dimension(:,:), allocatable :: &
469 0 : msg_buffer
470 :
471 : type (block) :: &
472 : this_block ! block info for current block
473 :
474 : character(len=*), parameter :: subname = '(gather_global_int)'
475 :
476 : !-----------------------------------------------------------------------
477 : !
478 : ! if this task is the dst_task, copy local blocks into the global
479 : ! array and post receives for non-local blocks.
480 : !
481 : !-----------------------------------------------------------------------
482 :
483 0 : if (my_task == dst_task) then
484 :
485 0 : do n=1,nblocks_tot
486 :
487 : !*** copy local blocks
488 :
489 0 : if (src_dist%blockLocation(n) == my_task+1) then
490 :
491 0 : this_block = get_block(n,n)
492 :
493 0 : do j=this_block%jlo,this_block%jhi
494 0 : do i=this_block%ilo,this_block%ihi
495 0 : ARRAY_G(this_block%i_glob(i), &
496 : this_block%j_glob(j)) = & ! LCOV_EXCL_LINE
497 0 : ARRAY(i,j,src_dist%blockLocalID(n))
498 : end do
499 : end do
500 :
501 : !*** fill land blocks with zeroes
502 :
503 0 : else if (src_dist%blockLocation(n) == 0) then
504 :
505 0 : this_block = get_block(n,n)
506 :
507 0 : do j=this_block%jlo,this_block%jhi
508 0 : do i=this_block%ilo,this_block%ihi
509 0 : ARRAY_G(this_block%i_glob(i), &
510 0 : this_block%j_glob(j)) = 0
511 : end do
512 : end do
513 : endif
514 :
515 : end do
516 :
517 : !*** receive blocks to fill up the rest
518 :
519 0 : allocate (msg_buffer(nx_block,ny_block))
520 :
521 0 : do n=1,nblocks_tot
522 0 : if (src_dist%blockLocation(n) > 0 .and. &
523 0 : src_dist%blockLocation(n) /= my_task+1) then
524 :
525 0 : this_block = get_block(n,n)
526 :
527 : call MPI_RECV(msg_buffer, size(msg_buffer), &
528 : mpi_integer, src_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
529 0 : MPI_COMM_ICE, status, ierr)
530 :
531 0 : do j=this_block%jlo,this_block%jhi
532 0 : do i=this_block%ilo,this_block%ihi
533 0 : ARRAY_G(this_block%i_glob(i), &
534 0 : this_block%j_glob(j)) = msg_buffer(i,j)
535 : end do
536 : end do
537 : endif
538 : end do
539 :
540 0 : deallocate(msg_buffer)
541 :
542 : !-----------------------------------------------------------------------
543 : !
544 : ! otherwise send data to dst_task
545 : !
546 : !-----------------------------------------------------------------------
547 :
548 : else
549 :
550 0 : allocate(snd_request(nblocks_tot), &
551 0 : snd_status (MPI_STATUS_SIZE, nblocks_tot))
552 :
553 0 : nsends = 0
554 0 : do n=1,nblocks_tot
555 0 : if (src_dist%blockLocation(n) == my_task+1) then
556 :
557 0 : nsends = nsends + 1
558 0 : src_block = src_dist%blockLocalID(n)
559 0 : call MPI_ISEND(ARRAY(1,1,src_block), nx_block*ny_block, &
560 : mpi_integer, dst_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
561 0 : MPI_COMM_ICE, snd_request(nsends), ierr)
562 : endif
563 : end do
564 :
565 0 : if (nsends > 0) &
566 0 : call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
567 0 : deallocate(snd_request, snd_status)
568 :
569 : endif
570 :
571 0 : if (add_mpi_barriers) then
572 0 : call ice_barrier()
573 : endif
574 :
575 : !-----------------------------------------------------------------------
576 :
577 0 : end subroutine gather_global_int
578 :
579 : !***********************************************************************
580 :
581 0 : subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
582 :
583 : ! This subroutine gathers a distributed array to a global-sized
584 : ! array on the processor dst_task, including ghost cells.
585 :
586 : integer (int_kind), intent(in) :: &
587 : dst_task ! task to which array should be gathered
588 :
589 : type (distrb), intent(in) :: &
590 : src_dist ! distribution of blocks in the source array
591 :
592 : real (dbl_kind), dimension(:,:,:), intent(in) :: &
593 : ARRAY ! array containing horizontal slab of distributed field
594 :
595 : real (dbl_kind), dimension(:,:), intent(inout) :: &
596 : ARRAY_G ! array containing global horizontal field on dst_task
597 :
598 : real (dbl_kind), intent(in), optional :: &
599 : spc_val
600 :
601 : !-----------------------------------------------------------------------
602 : !
603 : ! local variables
604 : !
605 : !-----------------------------------------------------------------------
606 :
607 : integer (int_kind) :: &
608 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
609 : nx, ny ,&! global dimensions ! LCOV_EXCL_LINE
610 : nsends ,&! number of actual sends ! LCOV_EXCL_LINE
611 : src_block ,&! block locator for send ! LCOV_EXCL_LINE
612 : ierr ! MPI error flag
613 :
614 : integer (int_kind), dimension(MPI_STATUS_SIZE) :: &
615 : status
616 :
617 : integer (int_kind), dimension(:), allocatable :: &
618 0 : snd_request
619 :
620 : integer (int_kind), dimension(:,:), allocatable :: &
621 0 : snd_status
622 :
623 : real (dbl_kind), dimension(:,:), allocatable :: &
624 0 : msg_buffer
625 :
626 : real (dbl_kind) :: &
627 0 : special_value
628 :
629 : type (block) :: &
630 : this_block ! block info for current block
631 :
632 : character(len=*), parameter :: subname = '(gather_global_ext_dbl)'
633 :
634 0 : if (present(spc_val)) then
635 0 : special_value = spc_val
636 : else
637 0 : special_value = spval_dbl
638 : endif
639 0 : ARRAY_G = special_value
640 :
641 0 : nx = nx_global + 2*nghost
642 0 : ny = ny_global + 2*nghost
643 :
644 : !-----------------------------------------------------------------------
645 : !
646 : ! if this task is the dst_task, copy local blocks into the global
647 : ! array and post receives for non-local blocks.
648 : !
649 : !-----------------------------------------------------------------------
650 :
651 0 : if (my_task == dst_task) then
652 :
653 0 : do n=1,nblocks_tot
654 :
655 : !*** copy local blocks
656 :
657 0 : if (src_dist%blockLocation(n) == my_task+1) then
658 :
659 0 : this_block = get_block(n,n)
660 :
661 : ! interior
662 0 : do j=this_block%jlo,this_block%jhi
663 0 : do i=this_block%ilo,this_block%ihi
664 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
665 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
666 0 : ARRAY (i,j,src_dist%blockLocalID(n))
667 : end do
668 : end do
669 :
670 : ! fill ghost cells
671 0 : if (this_block%jblock == 1) then
672 : ! south block
673 0 : do j=1, nghost
674 0 : do i=this_block%ilo,this_block%ihi
675 0 : ARRAY_G(this_block%i_glob(i)+nghost,j) = &
676 0 : ARRAY (i,j,src_dist%blockLocalID(n))
677 : end do
678 : end do
679 0 : if (this_block%iblock == 1) then
680 : ! southwest corner
681 0 : do j=1, nghost
682 0 : do i=1, nghost
683 0 : ARRAY_G(i,j) = &
684 0 : ARRAY (i,j,src_dist%blockLocalID(n))
685 : end do
686 : end do
687 : endif
688 : endif
689 0 : if (this_block%jblock == nblocks_y) then
690 : ! north block
691 0 : do j=1, nghost
692 0 : do i=this_block%ilo,this_block%ihi
693 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
694 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
695 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
696 : end do
697 : end do
698 0 : if (this_block%iblock == nblocks_x) then
699 : ! northeast corner
700 0 : do j=1, nghost
701 0 : do i=1, nghost
702 0 : ARRAY_G(nx-i+1, ny-j+1) = &
703 : ARRAY (this_block%ihi+nghost-i+1, & ! LCOV_EXCL_LINE
704 : this_block%jhi+nghost-j+1, & ! LCOV_EXCL_LINE
705 0 : src_dist%blockLocalID(n))
706 : end do
707 : end do
708 : endif
709 : endif
710 0 : if (this_block%iblock == 1) then
711 : ! west block
712 0 : do j=this_block%jlo,this_block%jhi
713 0 : do i=1, nghost
714 0 : ARRAY_G(i,this_block%j_glob(j)+nghost) = &
715 0 : ARRAY (i,j,src_dist%blockLocalID(n))
716 : end do
717 : end do
718 0 : if (this_block%jblock == nblocks_y) then
719 : ! northwest corner
720 0 : do j=1, nghost
721 0 : do i=1, nghost
722 0 : ARRAY_G(i, ny-j+1) = &
723 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
724 : end do
725 : end do
726 : endif
727 : endif
728 0 : if (this_block%iblock == nblocks_x) then
729 : ! east block
730 0 : do j=this_block%jlo,this_block%jhi
731 0 : do i=1, nghost
732 0 : ARRAY_G(nx_global + nghost + i, &
733 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
734 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
735 : end do
736 : end do
737 0 : if (this_block%jblock == 1) then
738 : ! southeast corner
739 0 : do j=1, nghost
740 0 : do i=1, nghost
741 0 : ARRAY_G( nx-i+1,j) = &
742 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
743 : end do
744 : end do
745 : endif
746 : endif
747 :
748 : endif ! src_dist%blockLocation
749 :
750 : end do
751 :
752 : !*** receive blocks to fill up the rest
753 :
754 0 : allocate (msg_buffer(nx_block,ny_block))
755 :
756 0 : do n=1,nblocks_tot
757 0 : if (src_dist%blockLocation(n) > 0 .and. &
758 0 : src_dist%blockLocation(n) /= my_task+1) then
759 :
760 0 : this_block = get_block(n,n)
761 :
762 : call MPI_RECV(msg_buffer, size(msg_buffer), &
763 : mpiR8, src_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
764 0 : MPI_COMM_ICE, status, ierr)
765 :
766 : ! block interior
767 0 : do j=this_block%jlo,this_block%jhi
768 0 : do i=this_block%ilo,this_block%ihi
769 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
770 0 : this_block%j_glob(j)+nghost) = msg_buffer(i,j)
771 : end do
772 : end do
773 0 : if (this_block%jblock == 1) then
774 : ! south block
775 0 : do j=1, nghost
776 0 : do i=this_block%ilo,this_block%ihi
777 0 : ARRAY_G (this_block%i_glob(i)+nghost,j) = &
778 0 : msg_buffer(i,j)
779 : end do
780 : end do
781 0 : if (this_block%iblock == 1) then
782 : ! southwest corner
783 0 : do j=1, nghost
784 0 : do i=1, nghost
785 0 : ARRAY_G(i,j) = msg_buffer(i,j)
786 : end do
787 : end do
788 : endif
789 : endif
790 0 : if (this_block%jblock == nblocks_y) then
791 : ! north block
792 0 : do j=1, nghost
793 0 : do i=this_block%ilo,this_block%ihi
794 0 : ARRAY_G (this_block%i_glob(i)+nghost, &
795 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
796 0 : msg_buffer(i, this_block%jhi+j)
797 : end do
798 : end do
799 0 : if (this_block%iblock == nblocks_x) then
800 : ! northeast corner
801 0 : do j=1, nghost
802 0 : do i=1, nghost
803 0 : ARRAY_G (nx-i+1, ny-j+1) = &
804 : msg_buffer(this_block%ihi+nghost-i+1,& ! LCOV_EXCL_LINE
805 0 : this_block%jhi+nghost-j+1)
806 : end do
807 : end do
808 : endif
809 : endif
810 0 : if (this_block%iblock == 1) then
811 : ! west block
812 0 : do j=this_block%jlo,this_block%jhi
813 0 : do i=1, nghost
814 0 : ARRAY_G (i, this_block%j_glob(j)+nghost) = &
815 0 : msg_buffer(i, j)
816 : end do
817 : end do
818 0 : if (this_block%jblock == nblocks_y) then
819 : ! northwest corner
820 0 : do j=1, nghost
821 0 : do i=1, nghost
822 0 : ARRAY_G (i, ny-j+1) = &
823 0 : msg_buffer(i, this_block%jhi+nghost-j+1)
824 : end do
825 : end do
826 : endif
827 : endif
828 0 : if (this_block%iblock == nblocks_x) then
829 : ! east block
830 0 : do j=this_block%jlo,this_block%jhi
831 0 : do i=1, nghost
832 0 : ARRAY_G (nx_global+nghost+i, &
833 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
834 0 : msg_buffer(this_block%ihi+i, j)
835 : end do
836 : end do
837 0 : if (this_block%jblock == 1) then
838 : ! southeast corner
839 0 : do j=1, nghost
840 0 : do i=1, nghost
841 0 : ARRAY_G (nx-i+1, j) = &
842 0 : msg_buffer(this_block%ihi+nghost-i+1, j)
843 : end do
844 : end do
845 : endif
846 : endif
847 : endif
848 : end do
849 :
850 0 : deallocate(msg_buffer)
851 :
852 : !-----------------------------------------------------------------------
853 : !
854 : ! otherwise send data to dst_task
855 : !
856 : !-----------------------------------------------------------------------
857 :
858 : else ! master task
859 :
860 0 : allocate(snd_request(nblocks_tot), &
861 0 : snd_status (MPI_STATUS_SIZE, nblocks_tot))
862 :
863 0 : nsends = 0
864 0 : do n=1,nblocks_tot
865 0 : if (src_dist%blockLocation(n) == my_task+1) then
866 :
867 0 : nsends = nsends + 1
868 0 : src_block = src_dist%blockLocalID(n)
869 0 : call MPI_ISEND(ARRAY(1,1,src_block), nx_block*ny_block, &
870 : mpiR8, dst_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
871 0 : MPI_COMM_ICE, snd_request(nsends), ierr)
872 : endif
873 : end do
874 :
875 0 : if (nsends > 0) &
876 0 : call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
877 0 : deallocate(snd_request, snd_status)
878 :
879 : endif ! master task
880 :
881 0 : if (add_mpi_barriers) then
882 0 : call ice_barrier()
883 : endif
884 :
885 : !-----------------------------------------------------------------------
886 :
887 0 : end subroutine gather_global_ext_dbl
888 :
889 : !***********************************************************************
890 :
891 0 : subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
892 :
893 : ! This subroutine gathers a distributed array to a global-sized
894 : ! array on the processor dst_task, including ghost cells.
895 :
896 : integer (int_kind), intent(in) :: &
897 : dst_task ! task to which array should be gathered
898 :
899 : type (distrb), intent(in) :: &
900 : src_dist ! distribution of blocks in the source array
901 :
902 : integer (int_kind), dimension(:,:,:), intent(in) :: &
903 : ARRAY ! array containing horizontal slab of distributed field
904 :
905 : integer (int_kind), dimension(:,:), intent(inout) :: &
906 : ARRAY_G ! array containing global horizontal field on dst_task
907 :
908 : integer (int_kind), intent(in), optional :: &
909 : spc_val
910 :
911 : !-----------------------------------------------------------------------
912 : !
913 : ! local variables
914 : !
915 : !-----------------------------------------------------------------------
916 :
917 : integer (int_kind) :: &
918 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
919 : nx, ny ,&! global dimensions ! LCOV_EXCL_LINE
920 : nsends ,&! number of actual sends ! LCOV_EXCL_LINE
921 : src_block ,&! block locator for send ! LCOV_EXCL_LINE
922 : ierr ! MPI error flag
923 :
924 : integer (int_kind), dimension(MPI_STATUS_SIZE) :: &
925 : status
926 :
927 : integer (int_kind), dimension(:), allocatable :: &
928 0 : snd_request
929 :
930 : integer (int_kind), dimension(:,:), allocatable :: &
931 0 : snd_status
932 :
933 : integer (int_kind), dimension(:,:), allocatable :: &
934 0 : msg_buffer
935 :
936 : integer (int_kind) :: &
937 : special_value
938 :
939 : type (block) :: &
940 : this_block ! block info for current block
941 :
942 : character(len=*), parameter :: subname = '(gather_global_ext_int)'
943 :
944 0 : if (present(spc_val)) then
945 0 : special_value = spc_val
946 : else
947 0 : special_value = -9999
948 : endif
949 0 : ARRAY_G = special_value
950 :
951 0 : nx = nx_global + 2*nghost
952 0 : ny = ny_global + 2*nghost
953 :
954 : !-----------------------------------------------------------------------
955 : !
956 : ! if this task is the dst_task, copy local blocks into the global
957 : ! array and post receives for non-local blocks.
958 : !
959 : !-----------------------------------------------------------------------
960 :
961 0 : if (my_task == dst_task) then
962 :
963 0 : do n=1,nblocks_tot
964 :
965 : !*** copy local blocks
966 :
967 0 : if (src_dist%blockLocation(n) == my_task+1) then
968 :
969 0 : this_block = get_block(n,n)
970 :
971 : ! interior
972 0 : do j=this_block%jlo,this_block%jhi
973 0 : do i=this_block%ilo,this_block%ihi
974 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
975 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
976 0 : ARRAY (i,j,src_dist%blockLocalID(n))
977 : end do
978 : end do
979 :
980 : ! fill ghost cells
981 0 : if (this_block%jblock == 1) then
982 : ! south block
983 0 : do j=1, nghost
984 0 : do i=this_block%ilo,this_block%ihi
985 0 : ARRAY_G(this_block%i_glob(i)+nghost,j) = &
986 0 : ARRAY (i,j,src_dist%blockLocalID(n))
987 : end do
988 : end do
989 0 : if (this_block%iblock == 1) then
990 : ! southwest corner
991 0 : do j=1, nghost
992 0 : do i=1, nghost
993 0 : ARRAY_G(i,j) = &
994 0 : ARRAY (i,j,src_dist%blockLocalID(n))
995 : end do
996 : end do
997 : endif
998 : endif
999 0 : if (this_block%jblock == nblocks_y) then
1000 : ! north block
1001 0 : do j=1, nghost
1002 0 : do i=this_block%ilo,this_block%ihi
1003 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
1004 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
1005 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
1006 : end do
1007 : end do
1008 0 : if (this_block%iblock == nblocks_x) then
1009 : ! northeast corner
1010 0 : do j=1, nghost
1011 0 : do i=1, nghost
1012 0 : ARRAY_G(nx-i+1, ny-j+1) = &
1013 : ARRAY (this_block%ihi+nghost-i+1, & ! LCOV_EXCL_LINE
1014 : this_block%jhi+nghost-j+1, & ! LCOV_EXCL_LINE
1015 0 : src_dist%blockLocalID(n))
1016 : end do
1017 : end do
1018 : endif
1019 : endif
1020 0 : if (this_block%iblock == 1) then
1021 : ! west block
1022 0 : do j=this_block%jlo,this_block%jhi
1023 0 : do i=1, nghost
1024 0 : ARRAY_G(i,this_block%j_glob(j)+nghost) = &
1025 0 : ARRAY (i,j,src_dist%blockLocalID(n))
1026 : end do
1027 : end do
1028 0 : if (this_block%jblock == nblocks_y) then
1029 : ! northwest corner
1030 0 : do j=1, nghost
1031 0 : do i=1, nghost
1032 0 : ARRAY_G(i, ny-j+1) = &
1033 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
1034 : end do
1035 : end do
1036 : endif
1037 : endif
1038 0 : if (this_block%iblock == nblocks_x) then
1039 : ! east block
1040 0 : do j=this_block%jlo,this_block%jhi
1041 0 : do i=1, nghost
1042 0 : ARRAY_G(nx_global + nghost + i, &
1043 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
1044 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
1045 : end do
1046 : end do
1047 0 : if (this_block%jblock == 1) then
1048 : ! southeast corner
1049 0 : do j=1, nghost
1050 0 : do i=1, nghost
1051 0 : ARRAY_G( nx-i+1,j) = &
1052 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
1053 : end do
1054 : end do
1055 : endif
1056 : endif
1057 :
1058 : endif ! src_dist%blockLocation
1059 :
1060 : end do
1061 :
1062 : !*** receive blocks to fill up the rest
1063 :
1064 0 : allocate (msg_buffer(nx_block,ny_block))
1065 :
1066 0 : do n=1,nblocks_tot
1067 0 : if (src_dist%blockLocation(n) > 0 .and. &
1068 0 : src_dist%blockLocation(n) /= my_task+1) then
1069 :
1070 0 : this_block = get_block(n,n)
1071 :
1072 : call MPI_RECV(msg_buffer, size(msg_buffer), &
1073 : MPI_INTEGER, src_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
1074 0 : MPI_COMM_ICE, status, ierr)
1075 :
1076 : ! block interior
1077 0 : do j=this_block%jlo,this_block%jhi
1078 0 : do i=this_block%ilo,this_block%ihi
1079 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
1080 0 : this_block%j_glob(j)+nghost) = msg_buffer(i,j)
1081 : end do
1082 : end do
1083 0 : if (this_block%jblock == 1) then
1084 : ! south block
1085 0 : do j=1, nghost
1086 0 : do i=this_block%ilo,this_block%ihi
1087 0 : ARRAY_G (this_block%i_glob(i)+nghost,j) = &
1088 0 : msg_buffer(i,j)
1089 : end do
1090 : end do
1091 0 : if (this_block%iblock == 1) then
1092 : ! southwest corner
1093 0 : do j=1, nghost
1094 0 : do i=1, nghost
1095 0 : ARRAY_G(i,j) = msg_buffer(i,j)
1096 : end do
1097 : end do
1098 : endif
1099 : endif
1100 0 : if (this_block%jblock == nblocks_y) then
1101 : ! north block
1102 0 : do j=1, nghost
1103 0 : do i=this_block%ilo,this_block%ihi
1104 0 : ARRAY_G (this_block%i_glob(i)+nghost, &
1105 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
1106 0 : msg_buffer(i, this_block%jhi+j)
1107 : end do
1108 : end do
1109 0 : if (this_block%iblock == nblocks_x) then
1110 : ! northeast corner
1111 0 : do j=1, nghost
1112 0 : do i=1, nghost
1113 0 : ARRAY_G (nx-i+1, ny-j+1) = &
1114 : msg_buffer(this_block%ihi+nghost-i+1,& ! LCOV_EXCL_LINE
1115 0 : this_block%jhi+nghost-j+1)
1116 : end do
1117 : end do
1118 : endif
1119 : endif
1120 0 : if (this_block%iblock == 1) then
1121 : ! west block
1122 0 : do j=this_block%jlo,this_block%jhi
1123 0 : do i=1, nghost
1124 0 : ARRAY_G (i, this_block%j_glob(j)+nghost) = &
1125 0 : msg_buffer(i, j)
1126 : end do
1127 : end do
1128 0 : if (this_block%jblock == nblocks_y) then
1129 : ! northwest corner
1130 0 : do j=1, nghost
1131 0 : do i=1, nghost
1132 0 : ARRAY_G (i, ny-j+1) = &
1133 0 : msg_buffer(i, this_block%jhi+nghost-j+1)
1134 : end do
1135 : end do
1136 : endif
1137 : endif
1138 0 : if (this_block%iblock == nblocks_x) then
1139 : ! east block
1140 0 : do j=this_block%jlo,this_block%jhi
1141 0 : do i=1, nghost
1142 0 : ARRAY_G (nx_global+nghost+i, &
1143 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
1144 0 : msg_buffer(this_block%ihi+i, j)
1145 : end do
1146 : end do
1147 0 : if (this_block%jblock == 1) then
1148 : ! southeast corner
1149 0 : do j=1, nghost
1150 0 : do i=1, nghost
1151 0 : ARRAY_G (nx-i+1, j) = &
1152 0 : msg_buffer(this_block%ihi+nghost-i+1, j)
1153 : end do
1154 : end do
1155 : endif
1156 : endif
1157 : endif
1158 : end do
1159 :
1160 0 : deallocate(msg_buffer)
1161 :
1162 : !-----------------------------------------------------------------------
1163 : !
1164 : ! otherwise send data to dst_task
1165 : !
1166 : !-----------------------------------------------------------------------
1167 :
1168 : else ! master task
1169 :
1170 0 : allocate(snd_request(nblocks_tot), &
1171 0 : snd_status (MPI_STATUS_SIZE, nblocks_tot))
1172 :
1173 0 : nsends = 0
1174 0 : do n=1,nblocks_tot
1175 0 : if (src_dist%blockLocation(n) == my_task+1) then
1176 :
1177 0 : nsends = nsends + 1
1178 0 : src_block = src_dist%blockLocalID(n)
1179 0 : call MPI_ISEND(ARRAY(1,1,src_block), nx_block*ny_block, &
1180 : MPI_INTEGER, dst_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
1181 0 : MPI_COMM_ICE, snd_request(nsends), ierr)
1182 : endif
1183 : end do
1184 :
1185 0 : if (nsends > 0) &
1186 0 : call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
1187 0 : deallocate(snd_request, snd_status)
1188 :
1189 : endif ! master task
1190 :
1191 0 : if (add_mpi_barriers) then
1192 0 : call ice_barrier()
1193 : endif
1194 :
1195 : !-----------------------------------------------------------------------
1196 :
1197 0 : end subroutine gather_global_ext_int
1198 :
1199 : !***********************************************************************
1200 :
1201 0 : subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
1202 :
1203 : ! This subroutine gathers a distributed array to a global-sized
1204 : ! array on the processor dst_task, including ghost cells.
1205 :
1206 : integer (int_kind), intent(in) :: &
1207 : dst_task ! task to which array should be gathered
1208 :
1209 : type (distrb), intent(in) :: &
1210 : src_dist ! distribution of blocks in the source array
1211 :
1212 : logical (log_kind), dimension(:,:,:), intent(in) :: &
1213 : ARRAY ! array containing horizontal slab of distributed field
1214 :
1215 : logical (log_kind), dimension(:,:), intent(inout) :: &
1216 : ARRAY_G ! array containing global horizontal field on dst_task
1217 :
1218 : logical (log_kind), intent(in), optional :: &
1219 : spc_val
1220 :
1221 : !-----------------------------------------------------------------------
1222 : !
1223 : ! local variables
1224 : !
1225 : !-----------------------------------------------------------------------
1226 :
1227 : integer (int_kind) :: &
1228 : i,j,n ,&! dummy loop counters ! LCOV_EXCL_LINE
1229 : nx, ny ,&! global dimensions ! LCOV_EXCL_LINE
1230 : nsends ,&! number of actual sends ! LCOV_EXCL_LINE
1231 : src_block ,&! block locator for send ! LCOV_EXCL_LINE
1232 : ierr ! MPI error flag
1233 :
1234 : integer (int_kind), dimension(MPI_STATUS_SIZE) :: &
1235 : status
1236 :
1237 : integer (int_kind), dimension(:), allocatable :: &
1238 0 : snd_request
1239 :
1240 : integer (int_kind), dimension(:,:), allocatable :: &
1241 0 : snd_status
1242 :
1243 : logical (log_kind), dimension(:,:), allocatable :: &
1244 0 : msg_buffer
1245 :
1246 : logical (log_kind) :: &
1247 : special_value
1248 :
1249 : type (block) :: &
1250 : this_block ! block info for current block
1251 :
1252 : character(len=*), parameter :: subname = '(gather_global_ext_log)'
1253 :
1254 0 : if (present(spc_val)) then
1255 0 : special_value = spc_val
1256 : else
1257 0 : special_value = .false.
1258 : endif
1259 0 : ARRAY_G = special_value
1260 :
1261 0 : nx = nx_global + 2*nghost
1262 0 : ny = ny_global + 2*nghost
1263 :
1264 : !-----------------------------------------------------------------------
1265 : !
1266 : ! if this task is the dst_task, copy local blocks into the global
1267 : ! array and post receives for non-local blocks.
1268 : !
1269 : !-----------------------------------------------------------------------
1270 :
1271 0 : if (my_task == dst_task) then
1272 :
1273 0 : do n=1,nblocks_tot
1274 :
1275 : !*** copy local blocks
1276 :
1277 0 : if (src_dist%blockLocation(n) == my_task+1) then
1278 :
1279 0 : this_block = get_block(n,n)
1280 :
1281 : ! interior
1282 0 : do j=this_block%jlo,this_block%jhi
1283 0 : do i=this_block%ilo,this_block%ihi
1284 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
1285 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
1286 0 : ARRAY (i,j,src_dist%blockLocalID(n))
1287 : end do
1288 : end do
1289 :
1290 : ! fill ghost cells
1291 0 : if (this_block%jblock == 1) then
1292 : ! south block
1293 0 : do j=1, nghost
1294 0 : do i=this_block%ilo,this_block%ihi
1295 0 : ARRAY_G(this_block%i_glob(i)+nghost,j) = &
1296 0 : ARRAY (i,j,src_dist%blockLocalID(n))
1297 : end do
1298 : end do
1299 0 : if (this_block%iblock == 1) then
1300 : ! southwest corner
1301 0 : do j=1, nghost
1302 0 : do i=1, nghost
1303 0 : ARRAY_G(i,j) = &
1304 0 : ARRAY (i,j,src_dist%blockLocalID(n))
1305 : end do
1306 : end do
1307 : endif
1308 : endif
1309 0 : if (this_block%jblock == nblocks_y) then
1310 : ! north block
1311 0 : do j=1, nghost
1312 0 : do i=this_block%ilo,this_block%ihi
1313 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
1314 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
1315 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
1316 : end do
1317 : end do
1318 0 : if (this_block%iblock == nblocks_x) then
1319 : ! northeast corner
1320 0 : do j=1, nghost
1321 0 : do i=1, nghost
1322 0 : ARRAY_G(nx-i+1, ny-j+1) = &
1323 : ARRAY (this_block%ihi+nghost-i+1, & ! LCOV_EXCL_LINE
1324 : this_block%jhi+nghost-j+1, & ! LCOV_EXCL_LINE
1325 0 : src_dist%blockLocalID(n))
1326 : end do
1327 : end do
1328 : endif
1329 : endif
1330 0 : if (this_block%iblock == 1) then
1331 : ! west block
1332 0 : do j=this_block%jlo,this_block%jhi
1333 0 : do i=1, nghost
1334 0 : ARRAY_G(i,this_block%j_glob(j)+nghost) = &
1335 0 : ARRAY (i,j,src_dist%blockLocalID(n))
1336 : end do
1337 : end do
1338 0 : if (this_block%jblock == nblocks_y) then
1339 : ! northwest corner
1340 0 : do j=1, nghost
1341 0 : do i=1, nghost
1342 0 : ARRAY_G(i, ny-j+1) = &
1343 0 : ARRAY (i,this_block%jhi+nghost-j+1,src_dist%blockLocalID(n))
1344 : end do
1345 : end do
1346 : endif
1347 : endif
1348 0 : if (this_block%iblock == nblocks_x) then
1349 : ! east block
1350 0 : do j=this_block%jlo,this_block%jhi
1351 0 : do i=1, nghost
1352 0 : ARRAY_G(nx_global + nghost + i, &
1353 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
1354 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
1355 : end do
1356 : end do
1357 0 : if (this_block%jblock == 1) then
1358 : ! southeast corner
1359 0 : do j=1, nghost
1360 0 : do i=1, nghost
1361 0 : ARRAY_G( nx-i+1,j) = &
1362 0 : ARRAY (this_block%ihi+nghost-i+1,j,src_dist%blockLocalID(n))
1363 : end do
1364 : end do
1365 : endif
1366 : endif
1367 :
1368 : endif ! src_dist%blockLocation
1369 :
1370 : end do
1371 :
1372 : !*** receive blocks to fill up the rest
1373 :
1374 0 : allocate (msg_buffer(nx_block,ny_block))
1375 :
1376 0 : do n=1,nblocks_tot
1377 0 : if (src_dist%blockLocation(n) > 0 .and. &
1378 0 : src_dist%blockLocation(n) /= my_task+1) then
1379 :
1380 0 : this_block = get_block(n,n)
1381 :
1382 : call MPI_RECV(msg_buffer, size(msg_buffer), &
1383 : MPI_LOGICAL, src_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
1384 0 : MPI_COMM_ICE, status, ierr)
1385 :
1386 : ! block interior
1387 0 : do j=this_block%jlo,this_block%jhi
1388 0 : do i=this_block%ilo,this_block%ihi
1389 0 : ARRAY_G(this_block%i_glob(i)+nghost, &
1390 0 : this_block%j_glob(j)+nghost) = msg_buffer(i,j)
1391 : end do
1392 : end do
1393 0 : if (this_block%jblock == 1) then
1394 : ! south block
1395 0 : do j=1, nghost
1396 0 : do i=this_block%ilo,this_block%ihi
1397 0 : ARRAY_G (this_block%i_glob(i)+nghost,j) = &
1398 0 : msg_buffer(i,j)
1399 : end do
1400 : end do
1401 0 : if (this_block%iblock == 1) then
1402 : ! southwest corner
1403 0 : do j=1, nghost
1404 0 : do i=1, nghost
1405 0 : ARRAY_G(i,j) = msg_buffer(i,j)
1406 : end do
1407 : end do
1408 : endif
1409 : endif
1410 0 : if (this_block%jblock == nblocks_y) then
1411 : ! north block
1412 0 : do j=1, nghost
1413 0 : do i=this_block%ilo,this_block%ihi
1414 0 : ARRAY_G (this_block%i_glob(i)+nghost, &
1415 : ny_global + nghost + j) = & ! LCOV_EXCL_LINE
1416 0 : msg_buffer(i, this_block%jhi+j)
1417 : end do
1418 : end do
1419 0 : if (this_block%iblock == nblocks_x) then
1420 : ! northeast corner
1421 0 : do j=1, nghost
1422 0 : do i=1, nghost
1423 0 : ARRAY_G (nx-i+1, ny-j+1) = &
1424 : msg_buffer(this_block%ihi+nghost-i+1,& ! LCOV_EXCL_LINE
1425 0 : this_block%jhi+nghost-j+1)
1426 : end do
1427 : end do
1428 : endif
1429 : endif
1430 0 : if (this_block%iblock == 1) then
1431 : ! west block
1432 0 : do j=this_block%jlo,this_block%jhi
1433 0 : do i=1, nghost
1434 0 : ARRAY_G (i, this_block%j_glob(j)+nghost) = &
1435 0 : msg_buffer(i, j)
1436 : end do
1437 : end do
1438 0 : if (this_block%jblock == nblocks_y) then
1439 : ! northwest corner
1440 0 : do j=1, nghost
1441 0 : do i=1, nghost
1442 0 : ARRAY_G (i, ny-j+1) = &
1443 0 : msg_buffer(i, this_block%jhi+nghost-j+1)
1444 : end do
1445 : end do
1446 : endif
1447 : endif
1448 0 : if (this_block%iblock == nblocks_x) then
1449 : ! east block
1450 0 : do j=this_block%jlo,this_block%jhi
1451 0 : do i=1, nghost
1452 0 : ARRAY_G (nx_global+nghost+i, &
1453 : this_block%j_glob(j)+nghost) = & ! LCOV_EXCL_LINE
1454 0 : msg_buffer(this_block%ihi+i, j)
1455 : end do
1456 : end do
1457 0 : if (this_block%jblock == 1) then
1458 : ! southeast corner
1459 0 : do j=1, nghost
1460 0 : do i=1, nghost
1461 0 : ARRAY_G (nx-i+1, j) = &
1462 0 : msg_buffer(this_block%ihi+nghost-i+1, j)
1463 : end do
1464 : end do
1465 : endif
1466 : endif
1467 : endif
1468 : end do
1469 :
1470 0 : deallocate(msg_buffer)
1471 :
1472 : !-----------------------------------------------------------------------
1473 : !
1474 : ! otherwise send data to dst_task
1475 : !
1476 : !-----------------------------------------------------------------------
1477 :
1478 : else ! master task
1479 :
1480 0 : allocate(snd_request(nblocks_tot), &
1481 0 : snd_status (MPI_STATUS_SIZE, nblocks_tot))
1482 :
1483 0 : nsends = 0
1484 0 : do n=1,nblocks_tot
1485 0 : if (src_dist%blockLocation(n) == my_task+1) then
1486 :
1487 0 : nsends = nsends + 1
1488 0 : src_block = src_dist%blockLocalID(n)
1489 0 : call MPI_ISEND(ARRAY(1,1,src_block), nx_block*ny_block, &
1490 : MPI_LOGICAL, dst_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
1491 0 : MPI_COMM_ICE, snd_request(nsends), ierr)
1492 : endif
1493 : end do
1494 :
1495 0 : if (nsends > 0) &
1496 0 : call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
1497 0 : deallocate(snd_request, snd_status)
1498 :
1499 : endif ! master task
1500 :
1501 0 : if (add_mpi_barriers) then
1502 0 : call ice_barrier()
1503 : endif
1504 :
1505 : !-----------------------------------------------------------------------
1506 :
1507 0 : end subroutine gather_global_ext_log
1508 :
1509 : !***********************************************************************
1510 :
1511 12156 : subroutine scatter_global_dbl(ARRAY, ARRAY_G, src_task, dst_dist, &
1512 : field_loc, field_type)
1513 :
1514 : ! This subroutine scatters a global-sized array to a distributed array.
1515 : !
1516 : ! This is the specific interface for double precision arrays
1517 : ! corresponding to the generic interface scatter_global.
1518 :
1519 : integer (int_kind), intent(in) :: &
1520 : src_task ! task from which array should be scattered
1521 :
1522 : type (distrb), intent(in) :: &
1523 : dst_dist ! distribution of resulting blocks
1524 :
1525 : real (dbl_kind), dimension(:,:), intent(in) :: &
1526 : ARRAY_G ! array containing global field on src_task
1527 :
1528 : integer (int_kind), intent(in) :: &
1529 : field_type, &! id for type of field (scalar, vector, angle) ! LCOV_EXCL_LINE
1530 : field_loc ! id for location on horizontal grid
1531 : ! (center, NEcorner, Nface, Eface)
1532 :
1533 : real (dbl_kind), dimension(:,:,:), intent(inout) :: &
1534 : ARRAY ! array containing distributed field
1535 :
1536 : !-----------------------------------------------------------------------
1537 : !
1538 : ! local variables
1539 : !
1540 : !-----------------------------------------------------------------------
1541 :
1542 : integer (int_kind) :: &
1543 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
1544 : nrecvs, &! actual number of messages received ! LCOV_EXCL_LINE
1545 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
1546 : dst_block, &! location of block in dst array ! LCOV_EXCL_LINE
1547 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
1548 : yoffset2, &! ! LCOV_EXCL_LINE
1549 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
1550 : ierr ! MPI error flag
1551 :
1552 : type (block) :: &
1553 : this_block ! block info for current block
1554 :
1555 : integer (int_kind), dimension(:), allocatable :: &
1556 12156 : rcv_request ! request array for receives
1557 :
1558 : integer (int_kind), dimension(:,:), allocatable :: &
1559 12156 : rcv_status ! status array for receives
1560 :
1561 : real (dbl_kind), dimension(:,:), allocatable :: &
1562 12156 : msg_buffer ! buffer for sending blocks
1563 :
1564 : character(len=*), parameter :: subname = '(scatter_global_dbl)'
1565 :
1566 : !-----------------------------------------------------------------------
1567 : !
1568 : ! initialize return array to zero and set up tripole quantities
1569 : !
1570 : !-----------------------------------------------------------------------
1571 :
1572 36723084 : ARRAY = c0
1573 :
1574 12156 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
1575 12156 : if (this_block%tripoleTFlag) then
1576 0 : select case (field_loc)
1577 : case (field_loc_center) ! cell center location
1578 0 : xoffset = 2
1579 0 : yoffset = 0
1580 : case (field_loc_NEcorner) ! cell corner (velocity) location
1581 0 : xoffset = 1
1582 0 : yoffset = -1
1583 : case (field_loc_Eface) ! cell face location
1584 0 : xoffset = 1
1585 0 : yoffset = 0
1586 : case (field_loc_Nface) ! cell face location
1587 0 : xoffset = 2
1588 0 : yoffset = -1
1589 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1590 0 : xoffset = 1
1591 0 : yoffset = 1
1592 : end select
1593 : else
1594 23212 : select case (field_loc)
1595 : case (field_loc_center) ! cell center location
1596 11056 : xoffset = 1
1597 11056 : yoffset = 1
1598 : case (field_loc_NEcorner) ! cell corner (velocity) location
1599 1028 : xoffset = 0
1600 1028 : yoffset = 0
1601 : case (field_loc_Eface) ! cell face location
1602 36 : xoffset = 0
1603 36 : yoffset = 1
1604 : case (field_loc_Nface) ! cell face location
1605 36 : xoffset = 1
1606 36 : yoffset = 0
1607 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1608 0 : xoffset = 1
1609 12156 : yoffset = 1
1610 : end select
1611 : endif
1612 :
1613 24180 : select case (field_type)
1614 : case (field_type_scalar)
1615 12024 : isign = 1
1616 : case (field_type_vector)
1617 112 : isign = -1
1618 : case (field_type_angle)
1619 20 : isign = -1
1620 : case (field_type_noupdate) ! ghost cells never used - use cell center
1621 0 : isign = 1
1622 : case default
1623 12156 : call abort_ice(subname//'ERROR: Unknown field type in scatter')
1624 : end select
1625 :
1626 : !-----------------------------------------------------------------------
1627 : !
1628 : ! if this task is the src_task, copy blocks of global array into
1629 : ! message buffer and send to other processors. also copy local blocks
1630 : !
1631 : !-----------------------------------------------------------------------
1632 :
1633 12156 : if (my_task == src_task) then
1634 :
1635 : !*** send non-local blocks away
1636 :
1637 2361 : allocate (msg_buffer(nx_block,ny_block))
1638 :
1639 69065 : do n=1,nblocks_tot
1640 66704 : if (dst_dist%blockLocation(n) > 0 .and. &
1641 2361 : dst_dist%blockLocation(n)-1 /= my_task) then
1642 :
1643 29084460 : msg_buffer = c0
1644 54548 : this_block = get_block(n,n)
1645 :
1646 : !*** if this is an interior block, then there is no
1647 : !*** padding or update checking required
1648 :
1649 : if (this_block%iblock > 1 .and. &
1650 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
1651 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
1652 : this_block%jblock < nblocks_y) then
1653 :
1654 627128 : do j=1,ny_block
1655 8840685 : do i=1,nx_block
1656 6182082 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
1657 15003254 : this_block%j_glob(j))
1658 : end do
1659 : end do
1660 :
1661 : !*** if this is an edge block but not a northern edge
1662 : !*** we only need to check for closed boundaries and
1663 : !*** padding (global index = 0)
1664 :
1665 35035 : else if (this_block%jblock /= nblocks_y) then
1666 :
1667 620484 : do j=1,ny_block
1668 620484 : if (this_block%j_glob(j) /= 0) then
1669 11271876 : do i=1,nx_block
1670 11271876 : if (this_block%i_glob(i) /= 0) then
1671 10303470 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
1672 20974125 : this_block%j_glob(j))
1673 : endif
1674 : end do
1675 : endif
1676 : end do
1677 :
1678 : !*** if this is a northern edge block, we need to check
1679 : !*** for and properly deal with tripole boundaries
1680 :
1681 : else
1682 :
1683 507416 : do j=1,ny_block
1684 507416 : if (this_block%j_glob(j) > 0) then ! normal boundary
1685 :
1686 8936864 : do i=1,nx_block
1687 8936864 : if (this_block%i_glob(i) /= 0) then
1688 8242776 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
1689 16687996 : this_block%j_glob(j))
1690 : endif
1691 : end do
1692 :
1693 0 : else if (this_block%j_glob(j) < 0) then ! tripole
1694 :
1695 : ! for yoffset=0 or 1, yoffset2=0,0
1696 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
1697 0 : do yoffset2=0,max(yoffset,0)-yoffset
1698 : jsrc = ny_global + yoffset + yoffset2 + &
1699 0 : (this_block%j_glob(j) + ny_global)
1700 0 : do i=1,nx_block
1701 0 : if (this_block%i_glob(i) /= 0) then
1702 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
1703 0 : if (isrc < 1) isrc = isrc + nx_global
1704 0 : if (isrc > nx_global) isrc = isrc - nx_global
1705 0 : msg_buffer(i,j-yoffset2) = isign * ARRAY_G(isrc,jsrc)
1706 : endif
1707 : end do
1708 : end do
1709 :
1710 : endif
1711 : end do
1712 :
1713 : endif
1714 :
1715 : call MPI_SEND(msg_buffer, nx_block*ny_block, &
1716 : mpiR8, dst_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
1717 54548 : MPI_COMM_ICE, ierr)
1718 :
1719 : endif
1720 : end do
1721 :
1722 2361 : deallocate(msg_buffer)
1723 :
1724 : !*** copy any local blocks
1725 :
1726 69065 : do n=1,nblocks_tot
1727 69065 : if (dst_dist%blockLocation(n) == my_task+1) then
1728 11704 : dst_block = dst_dist%blockLocalID(n)
1729 11704 : this_block = get_block(n,n)
1730 :
1731 : !*** if this is an interior block, then there is no
1732 : !*** padding or update checking required
1733 :
1734 : if (this_block%iblock > 1 .and. &
1735 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
1736 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
1737 : this_block%jblock < nblocks_y) then
1738 :
1739 126176 : do j=1,ny_block
1740 2025267 : do i=1,nx_block
1741 3091041 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1742 4082018 : this_block%j_glob(j))
1743 : end do
1744 : end do
1745 :
1746 : !*** if this is an edge block but not a northern edge
1747 : !*** we only need to check for closed boundaries and
1748 : !*** padding (global index = 0)
1749 :
1750 7761 : else if (this_block%jblock /= nblocks_y) then
1751 :
1752 220780 : do j=1,ny_block
1753 220780 : if (this_block%j_glob(j) /= 0) then
1754 5256700 : do i=1,nx_block
1755 5256700 : if (this_block%i_glob(i) /= 0) then
1756 9273123 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1757 11224859 : this_block%j_glob(j))
1758 : endif
1759 : end do
1760 : endif
1761 : end do
1762 :
1763 : !*** if this is a northern edge block, we need to check
1764 : !*** for and properly deal with tripole boundaries
1765 :
1766 : else
1767 :
1768 28928 : do j=1,ny_block
1769 28928 : if (this_block%j_glob(j) > 0) then ! normal boundary
1770 :
1771 224192 : do i=1,nx_block
1772 224192 : if (this_block%i_glob(i) /= 0) then
1773 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
1774 196168 : this_block%j_glob(j))
1775 : endif
1776 : end do
1777 :
1778 0 : else if (this_block%j_glob(j) < 0) then ! tripole
1779 :
1780 : ! for yoffset=0 or 1, yoffset2=0,0
1781 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
1782 0 : do yoffset2=0,max(yoffset,0)-yoffset
1783 : jsrc = ny_global + yoffset + yoffset2 + &
1784 0 : (this_block%j_glob(j) + ny_global)
1785 0 : do i=1,nx_block
1786 0 : if (this_block%i_glob(i) /= 0) then
1787 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
1788 0 : if (isrc < 1) isrc = isrc + nx_global
1789 0 : if (isrc > nx_global) isrc = isrc - nx_global
1790 0 : ARRAY(i,j-yoffset2,dst_block) &
1791 0 : = isign * ARRAY_G(isrc,jsrc)
1792 : endif
1793 : end do
1794 : end do
1795 :
1796 : endif
1797 : end do
1798 :
1799 : endif
1800 : endif
1801 : end do
1802 :
1803 : !-----------------------------------------------------------------------
1804 : !
1805 : ! otherwise receive data from src_task
1806 : !
1807 : !-----------------------------------------------------------------------
1808 :
1809 : else
1810 :
1811 0 : allocate (rcv_request(nblocks_tot), &
1812 9795 : rcv_status(MPI_STATUS_SIZE, nblocks_tot))
1813 :
1814 369011 : rcv_request = 0
1815 2524307 : rcv_status = 0
1816 :
1817 9795 : nrecvs = 0
1818 369011 : do n=1,nblocks_tot
1819 369011 : if (dst_dist%blockLocation(n) == my_task+1) then
1820 54548 : nrecvs = nrecvs + 1
1821 54548 : dst_block = dst_dist%blockLocalID(n)
1822 14772 : call MPI_IRECV(ARRAY(1,1,dst_block), nx_block*ny_block, &
1823 : mpiR8, src_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
1824 54548 : MPI_COMM_ICE, rcv_request(nrecvs), ierr)
1825 : endif
1826 : end do
1827 :
1828 9795 : if (nrecvs > 0) &
1829 9795 : call MPI_WAITALL(nrecvs, rcv_request, rcv_status, ierr)
1830 :
1831 9795 : deallocate(rcv_request, rcv_status)
1832 : endif
1833 :
1834 : !-----------------------------------------------------------------
1835 : ! Ensure unused ghost cell values are 0
1836 : !-----------------------------------------------------------------
1837 :
1838 12156 : if (field_loc == field_loc_noupdate) then
1839 0 : do n=1,nblocks_tot
1840 0 : dst_block = dst_dist%blockLocalID(n)
1841 0 : this_block = get_block(n,n)
1842 :
1843 0 : if (dst_block > 0) then
1844 :
1845 : ! north edge
1846 0 : do j = this_block%jhi+1,ny_block
1847 0 : do i = 1, nx_block
1848 0 : ARRAY (i,j,dst_block) = c0
1849 : enddo
1850 : enddo
1851 : ! east edge
1852 0 : do j = 1, ny_block
1853 0 : do i = this_block%ihi+1,nx_block
1854 0 : ARRAY (i,j,dst_block) = c0
1855 : enddo
1856 : enddo
1857 : ! south edge
1858 0 : do j = 1, this_block%jlo-1
1859 0 : do i = 1, nx_block
1860 0 : ARRAY (i,j,dst_block) = c0
1861 : enddo
1862 : enddo
1863 : ! west edge
1864 0 : do j = 1, ny_block
1865 0 : do i = 1, this_block%ilo-1
1866 0 : ARRAY (i,j,dst_block) = c0
1867 : enddo
1868 : enddo
1869 :
1870 : endif
1871 : enddo
1872 : endif
1873 :
1874 12156 : if (add_mpi_barriers) then
1875 0 : call ice_barrier()
1876 : endif
1877 :
1878 : !-----------------------------------------------------------------------
1879 :
1880 12156 : end subroutine scatter_global_dbl
1881 :
1882 : !***********************************************************************
1883 :
1884 0 : subroutine scatter_global_real(ARRAY, ARRAY_G, src_task, dst_dist, &
1885 : field_loc, field_type)
1886 :
1887 : !-----------------------------------------------------------------------
1888 : !
1889 : ! This subroutine scatters a global-sized array to a distributed array.
1890 : !
1891 : !-----------------------------------------------------------------------
1892 :
1893 : !-----------------------------------------------------------------------
1894 : !
1895 : ! input variables
1896 : !
1897 : !-----------------------------------------------------------------------
1898 :
1899 : integer (int_kind), intent(in) :: &
1900 : src_task ! task from which array should be scattered
1901 :
1902 : type (distrb), intent(in) :: &
1903 : dst_dist ! distribution of resulting blocks
1904 :
1905 : real (real_kind), dimension(:,:), intent(in) :: &
1906 : ARRAY_G ! array containing global field on src_task
1907 :
1908 : integer (int_kind), intent(in) :: &
1909 : field_type, &! id for type of field (scalar, vector, angle) ! LCOV_EXCL_LINE
1910 : field_loc ! id for location on horizontal grid
1911 : ! (center, NEcorner, Nface, Eface)
1912 :
1913 : !-----------------------------------------------------------------------
1914 : !
1915 : ! output variables
1916 : !
1917 : !-----------------------------------------------------------------------
1918 :
1919 : real (real_kind), dimension(:,:,:), intent(inout) :: &
1920 : ARRAY ! array containing distributed field
1921 :
1922 : !-----------------------------------------------------------------------
1923 : !
1924 : ! local variables
1925 : !
1926 : !-----------------------------------------------------------------------
1927 :
1928 : integer (int_kind) :: &
1929 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
1930 : nrecvs, &! actual number of messages received ! LCOV_EXCL_LINE
1931 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
1932 : dst_block, &! location of block in dst array ! LCOV_EXCL_LINE
1933 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
1934 : yoffset2, &! ! LCOV_EXCL_LINE
1935 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
1936 : ierr ! MPI error flag
1937 :
1938 : type (block) :: &
1939 : this_block ! block info for current block
1940 :
1941 : integer (int_kind), dimension(:), allocatable :: &
1942 0 : rcv_request ! request array for receives
1943 :
1944 : integer (int_kind), dimension(:,:), allocatable :: &
1945 0 : rcv_status ! status array for receives
1946 :
1947 : real (real_kind), dimension(:,:), allocatable :: &
1948 0 : msg_buffer ! buffer for sending blocks
1949 :
1950 : character(len=*), parameter :: subname = '(scatter_global_real)'
1951 :
1952 : !-----------------------------------------------------------------------
1953 : !
1954 : ! initialize return array to zero and set up tripole quantities
1955 : !
1956 : !-----------------------------------------------------------------------
1957 :
1958 0 : ARRAY = 0._real_kind
1959 :
1960 0 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
1961 0 : if (this_block%tripoleTFlag) then
1962 0 : select case (field_loc)
1963 : case (field_loc_center) ! cell center location
1964 0 : xoffset = 2
1965 0 : yoffset = 0
1966 : case (field_loc_NEcorner) ! cell corner (velocity) location
1967 0 : xoffset = 1
1968 0 : yoffset = 1
1969 : case (field_loc_Eface) ! cell face location
1970 0 : xoffset = 1
1971 0 : yoffset = 0
1972 : case (field_loc_Nface) ! cell face location
1973 0 : xoffset = 2
1974 0 : yoffset = 1
1975 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1976 0 : xoffset = 1
1977 0 : yoffset = 1
1978 : end select
1979 : else
1980 0 : select case (field_loc)
1981 : case (field_loc_center) ! cell center location
1982 0 : xoffset = 1
1983 0 : yoffset = 1
1984 : case (field_loc_NEcorner) ! cell corner (velocity) location
1985 0 : xoffset = 0
1986 0 : yoffset = 0
1987 : case (field_loc_Eface) ! cell face location
1988 0 : xoffset = 0
1989 0 : yoffset = 1
1990 : case (field_loc_Nface) ! cell face location
1991 0 : xoffset = 1
1992 0 : yoffset = 0
1993 : case (field_loc_noupdate) ! ghost cells never used - use cell center
1994 0 : xoffset = 1
1995 0 : yoffset = 1
1996 : end select
1997 : endif
1998 :
1999 0 : select case (field_type)
2000 : case (field_type_scalar)
2001 0 : isign = 1
2002 : case (field_type_vector)
2003 0 : isign = -1
2004 : case (field_type_angle)
2005 0 : isign = -1
2006 : case (field_type_noupdate) ! ghost cells never used - use cell center
2007 0 : isign = 1
2008 : case default
2009 0 : call abort_ice(subname//'ERROR: Unknown field type in scatter')
2010 : end select
2011 :
2012 : !-----------------------------------------------------------------------
2013 : !
2014 : ! if this task is the src_task, copy blocks of global array into
2015 : ! message buffer and send to other processors. also copy local blocks
2016 : !
2017 : !-----------------------------------------------------------------------
2018 :
2019 0 : if (my_task == src_task) then
2020 :
2021 : !*** send non-local blocks away
2022 :
2023 0 : allocate (msg_buffer(nx_block,ny_block))
2024 :
2025 0 : do n=1,nblocks_tot
2026 0 : if (dst_dist%blockLocation(n) > 0 .and. &
2027 0 : dst_dist%blockLocation(n)-1 /= my_task) then
2028 :
2029 0 : msg_buffer = 0._real_kind
2030 0 : this_block = get_block(n,n)
2031 :
2032 : !*** if this is an interior block, then there is no
2033 : !*** padding or update checking required
2034 :
2035 : if (this_block%iblock > 1 .and. &
2036 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
2037 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
2038 : this_block%jblock < nblocks_y) then
2039 :
2040 0 : do j=1,ny_block
2041 0 : do i=1,nx_block
2042 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
2043 0 : this_block%j_glob(j))
2044 : end do
2045 : end do
2046 :
2047 : !*** if this is an edge block but not a northern edge
2048 : !*** we only need to check for closed boundaries and
2049 : !*** padding (global index = 0)
2050 :
2051 0 : else if (this_block%jblock /= nblocks_y) then
2052 :
2053 0 : do j=1,ny_block
2054 0 : if (this_block%j_glob(j) /= 0) then
2055 0 : do i=1,nx_block
2056 0 : if (this_block%i_glob(i) /= 0) then
2057 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
2058 0 : this_block%j_glob(j))
2059 : endif
2060 : end do
2061 : endif
2062 : end do
2063 :
2064 : !*** if this is a northern edge block, we need to check
2065 : !*** for and properly deal with tripole boundaries
2066 :
2067 : else
2068 :
2069 0 : do j=1,ny_block
2070 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
2071 :
2072 0 : do i=1,nx_block
2073 0 : if (this_block%i_glob(i) /= 0) then
2074 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
2075 0 : this_block%j_glob(j))
2076 : endif
2077 : end do
2078 :
2079 0 : else if (this_block%j_glob(j) < 0) then ! tripole
2080 :
2081 : ! for yoffset=0 or 1, yoffset2=0,0
2082 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
2083 0 : do yoffset2=0,max(yoffset,0)-yoffset
2084 : jsrc = ny_global + yoffset + yoffset2 + &
2085 0 : (this_block%j_glob(j) + ny_global)
2086 0 : do i=1,nx_block
2087 0 : if (this_block%i_glob(i) /= 0) then
2088 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
2089 0 : if (isrc < 1) isrc = isrc + nx_global
2090 0 : if (isrc > nx_global) isrc = isrc - nx_global
2091 0 : msg_buffer(i,j-yoffset2) = isign * ARRAY_G(isrc,jsrc)
2092 : endif
2093 : end do
2094 : end do
2095 :
2096 : endif
2097 : end do
2098 :
2099 : endif
2100 :
2101 : call MPI_SEND(msg_buffer, nx_block*ny_block, &
2102 : mpiR4, dst_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
2103 0 : MPI_COMM_ICE, ierr)
2104 :
2105 : endif
2106 : end do
2107 :
2108 0 : deallocate(msg_buffer)
2109 :
2110 : !*** copy any local blocks
2111 :
2112 0 : do n=1,nblocks_tot
2113 0 : if (dst_dist%blockLocation(n) == my_task+1) then
2114 0 : dst_block = dst_dist%blockLocalID(n)
2115 0 : this_block = get_block(n,n)
2116 :
2117 : !*** if this is an interior block, then there is no
2118 : !*** padding or update checking required
2119 :
2120 : if (this_block%iblock > 1 .and. &
2121 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
2122 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
2123 : this_block%jblock < nblocks_y) then
2124 :
2125 0 : do j=1,ny_block
2126 0 : do i=1,nx_block
2127 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
2128 0 : this_block%j_glob(j))
2129 : end do
2130 : end do
2131 :
2132 : !*** if this is an edge block but not a northern edge
2133 : !*** we only need to check for closed boundaries and
2134 : !*** padding (global index = 0)
2135 :
2136 0 : else if (this_block%jblock /= nblocks_y) then
2137 :
2138 0 : do j=1,ny_block
2139 0 : if (this_block%j_glob(j) /= 0) then
2140 0 : do i=1,nx_block
2141 0 : if (this_block%i_glob(i) /= 0) then
2142 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
2143 0 : this_block%j_glob(j))
2144 : endif
2145 : end do
2146 : endif
2147 : end do
2148 :
2149 : !*** if this is a northern edge block, we need to check
2150 : !*** for and properly deal with tripole boundaries
2151 :
2152 : else
2153 :
2154 0 : do j=1,ny_block
2155 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
2156 :
2157 0 : do i=1,nx_block
2158 0 : if (this_block%i_glob(i) /= 0) then
2159 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
2160 0 : this_block%j_glob(j))
2161 : endif
2162 : end do
2163 :
2164 0 : else if (this_block%j_glob(j) < 0) then ! tripole
2165 :
2166 : ! for yoffset=0 or 1, yoffset2=0,0
2167 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
2168 0 : do yoffset2=0,max(yoffset,0)-yoffset
2169 : jsrc = ny_global + yoffset + yoffset2 + &
2170 0 : (this_block%j_glob(j) + ny_global)
2171 0 : do i=1,nx_block
2172 0 : if (this_block%i_glob(i) /= 0) then
2173 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
2174 0 : if (isrc < 1) isrc = isrc + nx_global
2175 0 : if (isrc > nx_global) isrc = isrc - nx_global
2176 0 : ARRAY(i,j-yoffset2,dst_block) &
2177 0 : = isign * ARRAY_G(isrc,jsrc)
2178 : endif
2179 : end do
2180 : end do
2181 :
2182 : endif
2183 : end do
2184 :
2185 : endif
2186 : endif
2187 : end do
2188 :
2189 : !-----------------------------------------------------------------------
2190 : !
2191 : ! otherwise receive data from src_task
2192 : !
2193 : !-----------------------------------------------------------------------
2194 :
2195 : else
2196 :
2197 0 : allocate (rcv_request(nblocks_tot), &
2198 0 : rcv_status(MPI_STATUS_SIZE, nblocks_tot))
2199 :
2200 0 : rcv_request = 0
2201 0 : rcv_status = 0
2202 :
2203 0 : nrecvs = 0
2204 0 : do n=1,nblocks_tot
2205 0 : if (dst_dist%blockLocation(n) == my_task+1) then
2206 0 : nrecvs = nrecvs + 1
2207 0 : dst_block = dst_dist%blockLocalID(n)
2208 0 : call MPI_IRECV(ARRAY(1,1,dst_block), nx_block*ny_block, &
2209 : mpiR4, src_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
2210 0 : MPI_COMM_ICE, rcv_request(nrecvs), ierr)
2211 : endif
2212 : end do
2213 :
2214 0 : if (nrecvs > 0) &
2215 0 : call MPI_WAITALL(nrecvs, rcv_request, rcv_status, ierr)
2216 :
2217 0 : deallocate(rcv_request, rcv_status)
2218 : endif
2219 :
2220 : !-----------------------------------------------------------------
2221 : ! Ensure unused ghost cell values are 0
2222 : !-----------------------------------------------------------------
2223 :
2224 0 : if (field_loc == field_loc_noupdate) then
2225 0 : do n=1,nblocks_tot
2226 0 : dst_block = dst_dist%blockLocalID(n)
2227 0 : this_block = get_block(n,n)
2228 :
2229 0 : if (dst_block > 0) then
2230 :
2231 : ! north edge
2232 0 : do j = this_block%jhi+1,ny_block
2233 0 : do i = 1, nx_block
2234 0 : ARRAY (i,j,dst_block) = 0._real_kind
2235 : enddo
2236 : enddo
2237 : ! east edge
2238 0 : do j = 1, ny_block
2239 0 : do i = this_block%ihi+1,nx_block
2240 0 : ARRAY (i,j,dst_block) = 0._real_kind
2241 : enddo
2242 : enddo
2243 : ! south edge
2244 0 : do j = 1, this_block%jlo-1
2245 0 : do i = 1, nx_block
2246 0 : ARRAY (i,j,dst_block) = 0._real_kind
2247 : enddo
2248 : enddo
2249 : ! west edge
2250 0 : do j = 1, ny_block
2251 0 : do i = 1, this_block%ilo-1
2252 0 : ARRAY (i,j,dst_block) = 0._real_kind
2253 : enddo
2254 : enddo
2255 :
2256 : endif
2257 : enddo
2258 : endif
2259 :
2260 0 : if (add_mpi_barriers) then
2261 0 : call ice_barrier()
2262 : endif
2263 :
2264 : !-----------------------------------------------------------------------
2265 :
2266 0 : end subroutine scatter_global_real
2267 :
2268 : !***********************************************************************
2269 :
2270 0 : subroutine scatter_global_int(ARRAY, ARRAY_G, src_task, dst_dist, &
2271 : field_loc, field_type)
2272 :
2273 : !-----------------------------------------------------------------------
2274 : !
2275 : ! This subroutine scatters a global-sized array to a distributed array.
2276 : !
2277 : !-----------------------------------------------------------------------
2278 :
2279 : !-----------------------------------------------------------------------
2280 : !
2281 : ! input variables
2282 : !
2283 : !-----------------------------------------------------------------------
2284 :
2285 : integer (int_kind), intent(in) :: &
2286 : src_task ! task from which array should be scattered
2287 :
2288 : integer (int_kind), intent(in) :: &
2289 : field_type, &! id for type of field (scalar, vector, angle) ! LCOV_EXCL_LINE
2290 : field_loc ! id for location on horizontal grid
2291 : ! (center, NEcorner, Nface, Eface)
2292 :
2293 : type (distrb), intent(in) :: &
2294 : dst_dist ! distribution of resulting blocks
2295 :
2296 : integer (int_kind), dimension(:,:), intent(in) :: &
2297 : ARRAY_G ! array containing global field on src_task
2298 :
2299 : !-----------------------------------------------------------------------
2300 : !
2301 : ! output variables
2302 : !
2303 : !-----------------------------------------------------------------------
2304 :
2305 : integer (int_kind), dimension(:,:,:), intent(inout) :: &
2306 : ARRAY ! array containing distributed field
2307 :
2308 : !-----------------------------------------------------------------------
2309 : !
2310 : ! local variables
2311 : !
2312 : !-----------------------------------------------------------------------
2313 :
2314 : integer (int_kind) :: &
2315 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
2316 : nrecvs, &! actual number of messages received ! LCOV_EXCL_LINE
2317 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
2318 : dst_block, &! location of block in dst array ! LCOV_EXCL_LINE
2319 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
2320 : yoffset2, &! ! LCOV_EXCL_LINE
2321 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
2322 : ierr ! MPI error flag
2323 :
2324 : type (block) :: &
2325 : this_block ! block info for current block
2326 :
2327 : integer (int_kind), dimension(:), allocatable :: &
2328 0 : rcv_request ! request array for receives
2329 :
2330 : integer (int_kind), dimension(:,:), allocatable :: &
2331 0 : rcv_status ! status array for receives
2332 :
2333 : integer (int_kind), dimension(:,:), allocatable :: &
2334 0 : msg_buffer ! buffer for sending blocks
2335 :
2336 : character(len=*), parameter :: subname = '(scatter_global_int)'
2337 :
2338 : !-----------------------------------------------------------------------
2339 : !
2340 : ! initialize return array to zero and set up tripole quantities
2341 : !
2342 : !-----------------------------------------------------------------------
2343 :
2344 0 : ARRAY = 0
2345 :
2346 0 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
2347 0 : if (this_block%tripoleTFlag) then
2348 0 : select case (field_loc)
2349 : case (field_loc_center) ! cell center location
2350 0 : xoffset = 2
2351 0 : yoffset = 0
2352 : case (field_loc_NEcorner) ! cell corner (velocity) location
2353 0 : xoffset = 1
2354 0 : yoffset = 1
2355 : case (field_loc_Eface) ! cell face location
2356 0 : xoffset = 1
2357 0 : yoffset = 0
2358 : case (field_loc_Nface) ! cell face location
2359 0 : xoffset = 2
2360 0 : yoffset = 1
2361 : case (field_loc_noupdate) ! ghost cells never used - use cell center
2362 0 : xoffset = 1
2363 0 : yoffset = 1
2364 : end select
2365 : else
2366 0 : select case (field_loc)
2367 : case (field_loc_center) ! cell center location
2368 0 : xoffset = 1
2369 0 : yoffset = 1
2370 : case (field_loc_NEcorner) ! cell corner (velocity) location
2371 0 : xoffset = 0
2372 0 : yoffset = 0
2373 : case (field_loc_Eface) ! cell face location
2374 0 : xoffset = 0
2375 0 : yoffset = 1
2376 : case (field_loc_Nface) ! cell face location
2377 0 : xoffset = 1
2378 0 : yoffset = 0
2379 : case (field_loc_noupdate) ! ghost cells never used - use cell center
2380 0 : xoffset = 1
2381 0 : yoffset = 1
2382 : end select
2383 : endif
2384 :
2385 0 : select case (field_type)
2386 : case (field_type_scalar)
2387 0 : isign = 1
2388 : case (field_type_vector)
2389 0 : isign = -1
2390 : case (field_type_angle)
2391 0 : isign = -1
2392 : case (field_type_noupdate) ! ghost cells never used - use cell center
2393 0 : isign = 1
2394 : case default
2395 0 : call abort_ice(subname//'ERROR: Unknown field type in scatter')
2396 : end select
2397 :
2398 : !-----------------------------------------------------------------------
2399 : !
2400 : ! if this task is the src_task, copy blocks of global array into
2401 : ! message buffer and send to other processors. also copy local blocks
2402 : !
2403 : !-----------------------------------------------------------------------
2404 :
2405 0 : if (my_task == src_task) then
2406 :
2407 : !*** send non-local blocks away
2408 :
2409 0 : allocate (msg_buffer(nx_block,ny_block))
2410 :
2411 0 : do n=1,nblocks_tot
2412 0 : if (dst_dist%blockLocation(n) > 0 .and. &
2413 0 : dst_dist%blockLocation(n)-1 /= my_task) then
2414 :
2415 0 : msg_buffer = 0
2416 0 : this_block = get_block(n,n)
2417 :
2418 : !*** if this is an interior block, then there is no
2419 : !*** padding or update checking required
2420 :
2421 : if (this_block%iblock > 1 .and. &
2422 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
2423 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
2424 : this_block%jblock < nblocks_y) then
2425 :
2426 0 : do j=1,ny_block
2427 0 : do i=1,nx_block
2428 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
2429 0 : this_block%j_glob(j))
2430 : end do
2431 : end do
2432 :
2433 : !*** if this is an edge block but not a northern edge
2434 : !*** we only need to check for closed boundaries and
2435 : !*** padding (global index = 0)
2436 :
2437 0 : else if (this_block%jblock /= nblocks_y) then
2438 :
2439 0 : do j=1,ny_block
2440 0 : if (this_block%j_glob(j) /= 0) then
2441 0 : do i=1,nx_block
2442 0 : if (this_block%i_glob(i) /= 0) then
2443 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
2444 0 : this_block%j_glob(j))
2445 : endif
2446 : end do
2447 : endif
2448 : end do
2449 :
2450 : !*** if this is a northern edge block, we need to check
2451 : !*** for and properly deal with tripole boundaries
2452 :
2453 : else
2454 :
2455 0 : do j=1,ny_block
2456 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
2457 :
2458 0 : do i=1,nx_block
2459 0 : if (this_block%i_glob(i) /= 0) then
2460 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i),&
2461 0 : this_block%j_glob(j))
2462 : endif
2463 : end do
2464 :
2465 0 : else if (this_block%j_glob(j) < 0) then ! tripole
2466 :
2467 : ! for yoffset=0 or 1, yoffset2=0,0
2468 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
2469 0 : do yoffset2=0,max(yoffset,0)-yoffset
2470 : jsrc = ny_global + yoffset + yoffset2 + &
2471 0 : (this_block%j_glob(j) + ny_global)
2472 0 : do i=1,nx_block
2473 0 : if (this_block%i_glob(i) /= 0) then
2474 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
2475 0 : if (isrc < 1) isrc = isrc + nx_global
2476 0 : if (isrc > nx_global) isrc = isrc - nx_global
2477 0 : msg_buffer(i,j-yoffset2) = isign * ARRAY_G(isrc,jsrc)
2478 : endif
2479 : end do
2480 : end do
2481 :
2482 : endif
2483 : end do
2484 :
2485 : endif
2486 :
2487 : call MPI_SEND(msg_buffer, nx_block*ny_block, &
2488 : mpi_integer, dst_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
2489 0 : MPI_COMM_ICE, ierr)
2490 :
2491 : endif
2492 : end do
2493 :
2494 0 : deallocate(msg_buffer)
2495 :
2496 : !*** copy any local blocks
2497 :
2498 0 : do n=1,nblocks_tot
2499 0 : if (dst_dist%blockLocation(n) == my_task+1) then
2500 0 : dst_block = dst_dist%blockLocalID(n)
2501 0 : this_block = get_block(n,n)
2502 :
2503 : !*** if this is an interior block, then there is no
2504 : !*** padding or update checking required
2505 :
2506 : if (this_block%iblock > 1 .and. &
2507 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
2508 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
2509 : this_block%jblock < nblocks_y) then
2510 :
2511 0 : do j=1,ny_block
2512 0 : do i=1,nx_block
2513 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
2514 0 : this_block%j_glob(j))
2515 : end do
2516 : end do
2517 :
2518 : !*** if this is an edge block but not a northern edge
2519 : !*** we only need to check for closed boundaries and
2520 : !*** padding (global index = 0)
2521 :
2522 0 : else if (this_block%jblock /= nblocks_y) then
2523 :
2524 0 : do j=1,ny_block
2525 0 : if (this_block%j_glob(j) /= 0) then
2526 0 : do i=1,nx_block
2527 0 : if (this_block%i_glob(i) /= 0) then
2528 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
2529 0 : this_block%j_glob(j))
2530 : endif
2531 : end do
2532 : endif
2533 : end do
2534 :
2535 : !*** if this is a northern edge block, we need to check
2536 : !*** for and properly deal with tripole boundaries
2537 :
2538 : else
2539 :
2540 0 : do j=1,ny_block
2541 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
2542 :
2543 0 : do i=1,nx_block
2544 0 : if (this_block%i_glob(i) /= 0) then
2545 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i),&
2546 0 : this_block%j_glob(j))
2547 : endif
2548 : end do
2549 :
2550 0 : else if (this_block%j_glob(j) < 0) then ! tripole
2551 :
2552 : ! for yoffset=0 or 1, yoffset2=0,0
2553 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
2554 0 : do yoffset2=0,max(yoffset,0)-yoffset
2555 : jsrc = ny_global + yoffset + yoffset2 + &
2556 0 : (this_block%j_glob(j) + ny_global)
2557 0 : do i=1,nx_block
2558 0 : if (this_block%i_glob(i) /= 0) then
2559 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
2560 0 : if (isrc < 1) isrc = isrc + nx_global
2561 0 : if (isrc > nx_global) isrc = isrc - nx_global
2562 0 : ARRAY(i,j-yoffset2,dst_block) &
2563 0 : = isign * ARRAY_G(isrc,jsrc)
2564 : endif
2565 : end do
2566 : end do
2567 :
2568 : endif
2569 : end do
2570 :
2571 : endif
2572 : endif
2573 : end do
2574 :
2575 : !-----------------------------------------------------------------------
2576 : !
2577 : ! otherwise receive data from src_task
2578 : !
2579 : !-----------------------------------------------------------------------
2580 :
2581 : else
2582 :
2583 0 : allocate (rcv_request(nblocks_tot), &
2584 0 : rcv_status(MPI_STATUS_SIZE, nblocks_tot))
2585 :
2586 0 : rcv_request = 0
2587 0 : rcv_status = 0
2588 :
2589 0 : nrecvs = 0
2590 0 : do n=1,nblocks_tot
2591 0 : if (dst_dist%blockLocation(n) == my_task+1) then
2592 0 : nrecvs = nrecvs + 1
2593 0 : dst_block = dst_dist%blockLocalID(n)
2594 0 : call MPI_IRECV(ARRAY(1,1,dst_block), nx_block*ny_block, &
2595 : mpi_integer, src_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
2596 0 : MPI_COMM_ICE, rcv_request(nrecvs), ierr)
2597 : endif
2598 : end do
2599 :
2600 0 : if (nrecvs > 0) &
2601 0 : call MPI_WAITALL(nrecvs, rcv_request, rcv_status, ierr)
2602 :
2603 0 : deallocate(rcv_request, rcv_status)
2604 : endif
2605 :
2606 : !-----------------------------------------------------------------
2607 : ! Ensure unused ghost cell values are 0
2608 : !-----------------------------------------------------------------
2609 :
2610 0 : if (field_loc == field_loc_noupdate) then
2611 0 : do n=1,nblocks_tot
2612 0 : dst_block = dst_dist%blockLocalID(n)
2613 0 : this_block = get_block(n,n)
2614 :
2615 0 : if (dst_block > 0) then
2616 :
2617 : ! north edge
2618 0 : do j = this_block%jhi+1,ny_block
2619 0 : do i = 1, nx_block
2620 0 : ARRAY (i,j,dst_block) = 0
2621 : enddo
2622 : enddo
2623 : ! east edge
2624 0 : do j = 1, ny_block
2625 0 : do i = this_block%ihi+1,nx_block
2626 0 : ARRAY (i,j,dst_block) = 0
2627 : enddo
2628 : enddo
2629 : ! south edge
2630 0 : do j = 1, this_block%jlo-1
2631 0 : do i = 1, nx_block
2632 0 : ARRAY (i,j,dst_block) = 0
2633 : enddo
2634 : enddo
2635 : ! west edge
2636 0 : do j = 1, ny_block
2637 0 : do i = 1, this_block%ilo-1
2638 0 : ARRAY (i,j,dst_block) = 0
2639 : enddo
2640 : enddo
2641 :
2642 : endif
2643 : enddo
2644 : endif
2645 :
2646 0 : if (add_mpi_barriers) then
2647 0 : call ice_barrier()
2648 : endif
2649 :
2650 : !-----------------------------------------------------------------------
2651 :
2652 0 : end subroutine scatter_global_int
2653 :
2654 : !***********************************************************************
2655 :
2656 0 : subroutine scatter_global_ext(ARRAY, ARRAY_G, src_task, dst_dist)
2657 :
2658 : ! This subroutine scatters a global-sized array to a distributed array.
2659 : !
2660 : ! This is the specific interface for double precision arrays
2661 : ! corresponding to the generic interface scatter_global.
2662 :
2663 : integer (int_kind), intent(in) :: &
2664 : src_task ! task from which array should be scattered
2665 :
2666 : type (distrb), intent(in) :: &
2667 : dst_dist ! distribution of resulting blocks
2668 :
2669 : real (dbl_kind), dimension(:,:), intent(in) :: &
2670 : ARRAY_G ! array containing global field on src_task
2671 :
2672 : real (dbl_kind), dimension(:,:,:), intent(inout) :: &
2673 : ARRAY ! array containing distributed field
2674 :
2675 : !-----------------------------------------------------------------------
2676 : !
2677 : ! local variables
2678 : !
2679 : !-----------------------------------------------------------------------
2680 :
2681 : integer (int_kind) :: &
2682 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
2683 : iblk, jblk, &! block indices ! LCOV_EXCL_LINE
2684 : iglb, jglb, &! global indices ! LCOV_EXCL_LINE
2685 : nrecvs, &! actual number of messages received ! LCOV_EXCL_LINE
2686 : dst_block, &! location of block in dst array ! LCOV_EXCL_LINE
2687 : ierr ! MPI error flag
2688 :
2689 : type (block) :: &
2690 : this_block ! block info for current block
2691 :
2692 : integer (int_kind), dimension(:), allocatable :: &
2693 0 : rcv_request ! request array for receives
2694 :
2695 : integer (int_kind), dimension(:,:), allocatable :: &
2696 0 : rcv_status ! status array for receives
2697 :
2698 : real (dbl_kind), dimension(:,:), allocatable :: &
2699 0 : msg_buffer ! buffer for sending blocks
2700 :
2701 : character(len=*), parameter :: subname = '(scatter_global_ext)'
2702 :
2703 : !-----------------------------------------------------------------------
2704 : !
2705 : ! initialize return array to zero
2706 : !
2707 : !-----------------------------------------------------------------------
2708 :
2709 0 : ARRAY = c0
2710 :
2711 : !-----------------------------------------------------------------------
2712 : !
2713 : ! if this task is the src_task, copy blocks of global array into
2714 : ! message buffer and send to other processors. also copy local blocks
2715 : !
2716 : !-----------------------------------------------------------------------
2717 :
2718 0 : if (my_task == src_task) then
2719 :
2720 : !*** send non-local blocks away
2721 :
2722 0 : allocate (msg_buffer(nx_block,ny_block))
2723 :
2724 0 : do n=1,nblocks_tot
2725 0 : if (dst_dist%blockLocation(n) > 0 .and. &
2726 0 : dst_dist%blockLocation(n)-1 /= my_task) then
2727 :
2728 0 : msg_buffer = c0
2729 0 : this_block = get_block(n,n)
2730 :
2731 : ! interior
2732 0 : do j = 1, ny_block
2733 0 : do i = 1, nx_block
2734 0 : if (this_block%j_glob(j) > 0) then
2735 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i)+nghost,&
2736 0 : this_block%j_glob(j)+nghost)
2737 : endif
2738 : end do
2739 : end do
2740 :
2741 0 : if (this_block%jblock == 1) then
2742 : ! south edge
2743 0 : do j = 1, nghost
2744 0 : do i = this_block%ilo,this_block%ihi
2745 0 : msg_buffer(i,j) = ARRAY_G(this_block%i_glob(i)+nghost,j)
2746 : enddo
2747 0 : do i = 1, nghost
2748 : ! southwest corner
2749 0 : iblk = i
2750 0 : jblk = j
2751 0 : iglb = this_block%i_glob(this_block%ilo)+i-1
2752 0 : jglb = j
2753 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2754 : ! southeast corner
2755 0 : iblk = this_block%ihi+i
2756 0 : iglb = this_block%i_glob(this_block%ihi)+nghost+i
2757 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2758 : enddo
2759 : enddo
2760 : endif
2761 0 : if (this_block%jblock == nblocks_y) then
2762 : ! north edge
2763 0 : do j = 1, nghost
2764 0 : do i = this_block%ilo,this_block%ihi
2765 0 : msg_buffer(i,this_block%jhi+j) = ARRAY_G(this_block%i_glob(i)+nghost,&
2766 0 : ny_global+nghost+j)
2767 : enddo
2768 0 : do i = 1, nghost
2769 : ! northwest corner
2770 0 : iblk = i
2771 0 : jblk = this_block%jhi+j
2772 0 : iglb = this_block%i_glob(this_block%ilo)+i-1
2773 0 : jglb = ny_global+nghost+j
2774 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2775 : ! northeast corner
2776 0 : iblk = this_block%ihi+i
2777 0 : iglb = this_block%i_glob(this_block%ihi)+nghost+i
2778 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2779 : enddo
2780 : enddo
2781 : endif
2782 0 : if (this_block%iblock == 1) then
2783 : ! west edge
2784 0 : do j = this_block%jlo,this_block%jhi
2785 0 : do i = 1, nghost
2786 0 : msg_buffer(i,j) = ARRAY_G(i,this_block%j_glob(j)+nghost)
2787 : enddo
2788 : enddo
2789 0 : do j = 1, nghost
2790 0 : do i = 1, nghost
2791 : ! northwest corner
2792 0 : iblk = i
2793 0 : jblk = this_block%jhi+j
2794 0 : iglb = i
2795 0 : jglb = this_block%j_glob(this_block%jhi)+nghost+j
2796 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2797 : ! southwest corner
2798 0 : jblk = j
2799 0 : jglb = this_block%j_glob(this_block%jlo)+j-1
2800 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2801 : enddo
2802 : enddo
2803 : endif
2804 0 : if (this_block%iblock == nblocks_x) then
2805 : ! east edge
2806 0 : do j = this_block%jlo,this_block%jhi
2807 0 : do i = 1, nghost
2808 0 : msg_buffer(this_block%ihi+i,j) = ARRAY_G(nx_global+nghost+i, &
2809 0 : this_block%j_glob(j)+nghost)
2810 : enddo
2811 : enddo
2812 0 : do j = 1, nghost
2813 0 : do i = 1, nghost
2814 : ! northeast corner
2815 0 : iblk = this_block%ihi+i
2816 0 : jblk = this_block%jhi+j
2817 0 : iglb = nx_global+nghost+i
2818 0 : jglb = this_block%j_glob(this_block%jhi)+nghost+j
2819 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2820 : ! southeast corner
2821 0 : jblk = j
2822 0 : jglb = this_block%j_glob(this_block%jlo)+j-1
2823 0 : msg_buffer(iblk,jblk) = ARRAY_G(iglb,jglb)
2824 : enddo
2825 : enddo
2826 : endif
2827 :
2828 : call MPI_SEND(msg_buffer, nx_block*ny_block, &
2829 : mpiR8, dst_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
2830 0 : MPI_COMM_ICE, ierr)
2831 :
2832 : endif
2833 : end do
2834 :
2835 0 : deallocate(msg_buffer)
2836 :
2837 : !*** copy any local blocks
2838 :
2839 0 : do n=1,nblocks_tot
2840 0 : if (dst_dist%blockLocation(n) == my_task+1) then
2841 0 : dst_block = dst_dist%blockLocalID(n)
2842 0 : this_block = get_block(n,n)
2843 :
2844 : ! interior
2845 0 : do j = 1, ny_block
2846 0 : do i = 1, nx_block
2847 0 : if (this_block%j_glob(j) > 0) then
2848 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i)+nghost,&
2849 0 : this_block%j_glob(j)+nghost)
2850 : endif
2851 : end do
2852 : end do
2853 :
2854 0 : if (this_block%jblock == 1) then
2855 : ! south edge
2856 0 : do j = 1, nghost
2857 0 : do i = this_block%ilo,this_block%ihi
2858 0 : ARRAY(i,j,dst_block) = ARRAY_G(this_block%i_glob(i)+nghost,j)
2859 : enddo
2860 0 : do i = 1, nghost
2861 : ! southwest corner
2862 0 : iblk = i
2863 0 : jblk = j
2864 0 : iglb = this_block%i_glob(this_block%ilo)+i-1
2865 0 : jglb = j
2866 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2867 : ! southeast corner
2868 0 : iblk = this_block%ihi+i
2869 0 : iglb = this_block%i_glob(this_block%ihi)+nghost+i
2870 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2871 : enddo
2872 : enddo
2873 : endif
2874 0 : if (this_block%jblock == nblocks_y) then
2875 : ! north edge
2876 0 : do j = 1, nghost
2877 0 : do i = this_block%ilo,this_block%ihi
2878 0 : ARRAY(i,this_block%jhi+j,dst_block) = ARRAY_G(this_block%i_glob(i)+nghost,&
2879 0 : ny_global+nghost+j)
2880 : enddo
2881 0 : do i = 1, nghost
2882 : ! northwest corner
2883 0 : iblk = i
2884 0 : jblk = this_block%jhi+j
2885 0 : iglb = this_block%i_glob(this_block%ilo)+i-1
2886 0 : jglb = ny_global+nghost+j
2887 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2888 : ! northeast corner
2889 0 : iblk = this_block%ihi+i
2890 0 : iglb = this_block%i_glob(this_block%ihi)+nghost+i
2891 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2892 : enddo
2893 : enddo
2894 : endif
2895 0 : if (this_block%iblock == 1) then
2896 : ! west edge
2897 0 : do j = this_block%jlo,this_block%jhi
2898 0 : do i = 1, nghost
2899 0 : ARRAY(i,j,dst_block) = ARRAY_G(i,this_block%j_glob(j)+nghost)
2900 : enddo
2901 : enddo
2902 0 : do j = 1, nghost
2903 0 : do i = 1, nghost
2904 : ! northwest corner
2905 0 : iblk = i
2906 0 : jblk = this_block%jhi+j
2907 0 : iglb = i
2908 0 : jglb = this_block%j_glob(this_block%jhi)+nghost+j
2909 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2910 : ! southwest corner
2911 0 : jblk = j
2912 0 : jglb = this_block%j_glob(this_block%jlo)+j-1
2913 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2914 : enddo
2915 : enddo
2916 : endif
2917 0 : if (this_block%iblock == nblocks_x) then
2918 : ! east edge
2919 0 : do j = this_block%jlo,this_block%jhi
2920 0 : do i = 1, nghost
2921 0 : ARRAY(this_block%ihi+i,j,dst_block) = ARRAY_G(nx_global+nghost+i, &
2922 0 : this_block%j_glob(j)+nghost)
2923 : enddo
2924 : enddo
2925 0 : do j = 1, nghost
2926 0 : do i = 1, nghost
2927 : ! northeast corner
2928 0 : iblk = this_block%ihi+i
2929 0 : jblk = this_block%jhi+j
2930 0 : iglb = nx_global+nghost+i
2931 0 : jglb = this_block%j_glob(this_block%jhi)+nghost+j
2932 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2933 : ! southeast corner
2934 0 : jblk = j
2935 0 : jglb = this_block%j_glob(this_block%jlo)+j-1
2936 0 : ARRAY(iblk,jblk,dst_block) = ARRAY_G(iglb,jglb)
2937 : enddo
2938 : enddo
2939 : endif
2940 :
2941 : endif
2942 : end do
2943 :
2944 : !-----------------------------------------------------------------------
2945 : !
2946 : ! otherwise receive data from src_task
2947 : !
2948 : !-----------------------------------------------------------------------
2949 :
2950 : else
2951 :
2952 0 : allocate (rcv_request(nblocks_tot), &
2953 0 : rcv_status(MPI_STATUS_SIZE, nblocks_tot))
2954 :
2955 0 : rcv_request = 0
2956 0 : rcv_status = 0
2957 :
2958 0 : nrecvs = 0
2959 0 : do n=1,nblocks_tot
2960 0 : if (dst_dist%blockLocation(n) == my_task+1) then
2961 0 : nrecvs = nrecvs + 1
2962 0 : dst_block = dst_dist%blockLocalID(n)
2963 0 : call MPI_IRECV(ARRAY(1,1,dst_block), nx_block*ny_block, &
2964 : mpiR8, src_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
2965 0 : MPI_COMM_ICE, rcv_request(nrecvs), ierr)
2966 : endif
2967 : end do
2968 :
2969 0 : if (nrecvs > 0) &
2970 0 : call MPI_WAITALL(nrecvs, rcv_request, rcv_status, ierr)
2971 :
2972 0 : deallocate(rcv_request, rcv_status)
2973 : endif
2974 :
2975 0 : if (add_mpi_barriers) then
2976 0 : call ice_barrier()
2977 : endif
2978 :
2979 : !-----------------------------------------------------------------------
2980 :
2981 0 : end subroutine scatter_global_ext
2982 :
2983 : !***********************************************************************
2984 :
2985 0 : subroutine scatter_global_stress(ARRAY, ARRAY_G1, ARRAY_G2, &
2986 : src_task, dst_dist)
2987 :
2988 : ! This subroutine scatters global stresses to a distributed array.
2989 : !
2990 : ! Ghost cells in the stress tensor must be handled separately on tripole
2991 : ! grids, because matching the corner values requires 2 different arrays.
2992 :
2993 : integer (int_kind), intent(in) :: &
2994 : src_task ! task from which array should be scattered
2995 :
2996 : type (distrb), intent(in) :: &
2997 : dst_dist ! distribution of resulting blocks
2998 :
2999 : real (dbl_kind), dimension(:,:), intent(in) :: &
3000 : ARRAY_G1, &! array containing global field on src_task ! LCOV_EXCL_LINE
3001 : ARRAY_G2 ! array containing global field on src_task
3002 :
3003 : real (dbl_kind), dimension(:,:,:), intent(inout) :: &
3004 : ARRAY ! array containing distributed field
3005 :
3006 : !-----------------------------------------------------------------------
3007 : !
3008 : ! local variables
3009 : !
3010 : !-----------------------------------------------------------------------
3011 :
3012 : integer (int_kind) :: &
3013 : i,j,n, &! dummy loop indices ! LCOV_EXCL_LINE
3014 : nrecvs, &! actual number of messages received ! LCOV_EXCL_LINE
3015 : isrc, jsrc, &! source addresses ! LCOV_EXCL_LINE
3016 : dst_block, &! location of block in dst array ! LCOV_EXCL_LINE
3017 : xoffset, yoffset, &! offsets for tripole boundary conditions ! LCOV_EXCL_LINE
3018 : yoffset2, &! ! LCOV_EXCL_LINE
3019 : isign, &! sign factor for tripole boundary conditions ! LCOV_EXCL_LINE
3020 : ierr ! MPI error flag
3021 :
3022 : type (block) :: &
3023 : this_block ! block info for current block
3024 :
3025 : integer (int_kind), dimension(:), allocatable :: &
3026 0 : rcv_request ! request array for receives
3027 :
3028 : integer (int_kind), dimension(:,:), allocatable :: &
3029 0 : rcv_status ! status array for receives
3030 :
3031 : real (dbl_kind), dimension(:,:), allocatable :: &
3032 0 : msg_buffer ! buffer for sending blocks
3033 :
3034 : character(len=*), parameter :: subname = '(scatter_global_stress)'
3035 :
3036 : !-----------------------------------------------------------------------
3037 : !
3038 : ! initialize return array to zero and set up tripole quantities
3039 : !
3040 : !-----------------------------------------------------------------------
3041 :
3042 0 : ARRAY = c0
3043 :
3044 0 : this_block = get_block(1,1) ! for the tripoleTflag - all blocks have it
3045 0 : if (this_block%tripoleTFlag) then
3046 0 : xoffset = 2 ! treat stresses as cell-centered scalars (they are not
3047 0 : yoffset = 0 ! shared with neighboring grid cells)
3048 : else
3049 0 : xoffset = 1 ! treat stresses as cell-centered scalars (they are not
3050 0 : yoffset = 1 ! shared with neighboring grid cells)
3051 : endif
3052 0 : isign = 1
3053 :
3054 : !-----------------------------------------------------------------------
3055 : !
3056 : ! if this task is the src_task, copy blocks of global array into
3057 : ! message buffer and send to other processors. also copy local blocks
3058 : !
3059 : !-----------------------------------------------------------------------
3060 :
3061 0 : if (my_task == src_task) then
3062 :
3063 : !*** send non-local blocks away
3064 :
3065 0 : allocate (msg_buffer(nx_block,ny_block))
3066 :
3067 0 : do n=1,nblocks_tot
3068 0 : if (dst_dist%blockLocation(n) > 0 .and. &
3069 0 : dst_dist%blockLocation(n)-1 /= my_task) then
3070 :
3071 0 : msg_buffer = c0
3072 0 : this_block = get_block(n,n)
3073 :
3074 : !*** if this is an interior block, then there is no
3075 : !*** padding or update checking required
3076 :
3077 : if (this_block%iblock > 1 .and. &
3078 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
3079 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
3080 : this_block%jblock < nblocks_y) then
3081 :
3082 0 : do j=1,ny_block
3083 0 : do i=1,nx_block
3084 0 : msg_buffer(i,j) = ARRAY_G1(this_block%i_glob(i),&
3085 0 : this_block%j_glob(j))
3086 : end do
3087 : end do
3088 :
3089 : !*** if this is an edge block but not a northern edge
3090 : !*** we only need to check for closed boundaries and
3091 : !*** padding (global index = 0)
3092 :
3093 0 : else if (this_block%jblock /= nblocks_y) then
3094 :
3095 0 : do j=1,ny_block
3096 0 : if (this_block%j_glob(j) /= 0) then
3097 0 : do i=1,nx_block
3098 0 : if (this_block%i_glob(i) /= 0) then
3099 0 : msg_buffer(i,j) = ARRAY_G1(this_block%i_glob(i),&
3100 0 : this_block%j_glob(j))
3101 : endif
3102 : end do
3103 : endif
3104 : end do
3105 :
3106 : !*** if this is a northern edge block, we need to check
3107 : !*** for and properly deal with tripole boundaries
3108 :
3109 : else
3110 :
3111 0 : do j=1,ny_block
3112 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
3113 :
3114 0 : do i=1,nx_block
3115 0 : if (this_block%i_glob(i) /= 0) then
3116 0 : msg_buffer(i,j) = ARRAY_G1(this_block%i_glob(i),&
3117 0 : this_block%j_glob(j))
3118 : endif
3119 : end do
3120 :
3121 0 : else if (this_block%j_glob(j) < 0) then ! tripole
3122 :
3123 : jsrc = ny_global + yoffset + &
3124 0 : (this_block%j_glob(j) + ny_global)
3125 0 : do i=1,nx_block
3126 0 : if (this_block%i_glob(i) /= 0) then
3127 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
3128 0 : if (isrc < 1) isrc = isrc + nx_global
3129 0 : if (isrc > nx_global) isrc = isrc - nx_global
3130 0 : msg_buffer(i,j) = isign * ARRAY_G2(isrc,jsrc)
3131 : endif
3132 : end do
3133 :
3134 : endif
3135 : end do
3136 :
3137 : endif
3138 :
3139 : call MPI_SEND(msg_buffer, nx_block*ny_block, &
3140 : mpiR8, dst_dist%blockLocation(n)-1, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
3141 0 : MPI_COMM_ICE, ierr)
3142 :
3143 : endif
3144 : end do
3145 :
3146 0 : deallocate(msg_buffer)
3147 :
3148 : !*** copy any local blocks
3149 :
3150 0 : do n=1,nblocks_tot
3151 0 : if (dst_dist%blockLocation(n) == my_task+1) then
3152 0 : dst_block = dst_dist%blockLocalID(n)
3153 0 : this_block = get_block(n,n)
3154 :
3155 : !*** if this is an interior block, then there is no
3156 : !*** padding or update checking required
3157 :
3158 : if (this_block%iblock > 1 .and. &
3159 : this_block%iblock < nblocks_x .and. & ! LCOV_EXCL_LINE
3160 : this_block%jblock > 1 .and. & ! LCOV_EXCL_LINE
3161 : this_block%jblock < nblocks_y) then
3162 :
3163 0 : do j=1,ny_block
3164 0 : do i=1,nx_block
3165 0 : ARRAY(i,j,dst_block) = ARRAY_G1(this_block%i_glob(i),&
3166 0 : this_block%j_glob(j))
3167 : end do
3168 : end do
3169 :
3170 : !*** if this is an edge block but not a northern edge
3171 : !*** we only need to check for closed boundaries and
3172 : !*** padding (global index = 0)
3173 :
3174 0 : else if (this_block%jblock /= nblocks_y) then
3175 :
3176 0 : do j=1,ny_block
3177 0 : if (this_block%j_glob(j) /= 0) then
3178 0 : do i=1,nx_block
3179 0 : if (this_block%i_glob(i) /= 0) then
3180 0 : ARRAY(i,j,dst_block) = ARRAY_G1(this_block%i_glob(i),&
3181 0 : this_block%j_glob(j))
3182 : endif
3183 : end do
3184 : endif
3185 : end do
3186 :
3187 : !*** if this is a northern edge block, we need to check
3188 : !*** for and properly deal with tripole boundaries
3189 :
3190 : else
3191 :
3192 0 : do j=1,ny_block
3193 0 : if (this_block%j_glob(j) > 0) then ! normal boundary
3194 :
3195 0 : do i=1,nx_block
3196 0 : if (this_block%i_glob(i) /= 0) then
3197 0 : ARRAY(i,j,dst_block) = ARRAY_G1(this_block%i_glob(i),&
3198 0 : this_block%j_glob(j))
3199 : endif
3200 : end do
3201 :
3202 0 : else if (this_block%j_glob(j) < 0) then ! tripole
3203 :
3204 : ! for yoffset=0 or 1, yoffset2=0,0
3205 : ! for yoffset=-1, yoffset2=0,1, for u-rows on T-fold grid
3206 0 : do yoffset2=0,max(yoffset,0)-yoffset
3207 : jsrc = ny_global + yoffset + yoffset2 + &
3208 0 : (this_block%j_glob(j) + ny_global)
3209 0 : do i=1,nx_block
3210 0 : if (this_block%i_glob(i) /= 0) then
3211 0 : isrc = nx_global + xoffset - this_block%i_glob(i)
3212 0 : if (isrc < 1) isrc = isrc + nx_global
3213 0 : if (isrc > nx_global) isrc = isrc - nx_global
3214 0 : ARRAY(i,j-yoffset2,dst_block) &
3215 0 : = isign * ARRAY_G2(isrc,jsrc)
3216 : endif
3217 : end do
3218 : end do
3219 :
3220 : endif
3221 : end do
3222 :
3223 : endif
3224 : endif
3225 : end do
3226 :
3227 : !-----------------------------------------------------------------------
3228 : !
3229 : ! otherwise receive data from src_task
3230 : !
3231 : !-----------------------------------------------------------------------
3232 :
3233 : else
3234 :
3235 0 : allocate (rcv_request(nblocks_tot), &
3236 0 : rcv_status(MPI_STATUS_SIZE, nblocks_tot))
3237 :
3238 0 : rcv_request = 0
3239 0 : rcv_status = 0
3240 :
3241 0 : nrecvs = 0
3242 0 : do n=1,nblocks_tot
3243 0 : if (dst_dist%blockLocation(n) == my_task+1) then
3244 0 : nrecvs = nrecvs + 1
3245 0 : dst_block = dst_dist%blockLocalID(n)
3246 0 : call MPI_IRECV(ARRAY(1,1,dst_block), nx_block*ny_block, &
3247 : mpiR8, src_task, 3*mpitag_gs+n, & ! LCOV_EXCL_LINE
3248 0 : MPI_COMM_ICE, rcv_request(nrecvs), ierr)
3249 : endif
3250 : end do
3251 :
3252 0 : if (nrecvs > 0) &
3253 0 : call MPI_WAITALL(nrecvs, rcv_request, rcv_status, ierr)
3254 :
3255 0 : deallocate(rcv_request, rcv_status)
3256 : endif
3257 :
3258 0 : if (add_mpi_barriers) then
3259 0 : call ice_barrier()
3260 : endif
3261 :
3262 : !-----------------------------------------------------------------------
3263 :
3264 0 : end subroutine scatter_global_stress
3265 :
3266 : !***********************************************************************
3267 :
3268 : end module ice_gather_scatter
3269 :
3270 : !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
|