17 use mpi
, only : mpi_request_null, mpi_statuses_ignore
34 pid_location, current_page, source_data)
37 integer,
intent(in) :: dim, pid_location, source_index
38 integer,
intent(inout) :: current_page(:)
44 dim, x_source_index, &
45 y_source_index, pid_location, current_page, source_data)
48 integer,
intent(in) :: dim, pid_location, x_source_index, y_source_index
49 integer,
intent(inout) :: current_page(:)
56 neighbour_location, current_page, source_data)
59 integer,
intent(in) :: dim, target_index, neighbour_location
60 integer,
intent(inout) :: current_page(:)
66 corner_loc, x_target_index, &
67 y_target_index, neighbour_location, current_page, source_data)
70 integer,
intent(in) :: corner_loc, x_target_index, y_target_index, neighbour_location
71 integer,
intent(inout) :: current_page(:)
77 involve_corners, source_data)
80 integer,
intent(in) :: halo_depth
81 logical,
intent(in) :: involve_corners
110 perform_local_data_copy, copy_from_halo_buffer, copy_corners_to_halo_buffer, &
111 copy_from_halo_buffer_to_corner, source_data)
119 copy_corners_to_halo_buffer
121 copy_from_halo_buffer_to_corner
124 if ((
present(copy_corners_to_halo_buffer) .and. .not. &
125 present(copy_from_halo_buffer_to_corner)) .or. &
126 (.not.
present(copy_corners_to_halo_buffer) .and. &
127 present(copy_from_halo_buffer_to_corner)))
then 128 call log_log(
log_error,
"You must either provie no or both corner optional arguments to the halo swap function")
131 if (
present(source_data))
then 132 if (
present(copy_corners_to_halo_buffer))
then 134 copy_to_halo_buffer, copy_corners_to_halo_buffer, source_data)
137 copy_to_halo_buffer, source_data=source_data)
139 if (
present(copy_from_halo_buffer_to_corner))
then 141 perform_local_data_copy, copy_from_halo_buffer, copy_from_halo_buffer_to_corner, &
145 perform_local_data_copy, copy_from_halo_buffer, source_data=source_data)
148 if (
present(copy_corners_to_halo_buffer))
then 150 copy_corners_to_halo_buffer)
154 if (
present(copy_from_halo_buffer_to_corner))
then 156 perform_local_data_copy, &
157 copy_from_halo_buffer, copy_from_halo_buffer_to_corner)
160 perform_local_data_copy, copy_from_halo_buffer)
180 perform_local_data_copy, copy_from_halo_buffer, copy_from_halo_buffer_to_corner, &
187 copy_from_halo_buffer_to_corner
192 if ((
present(copy_from_halo_buffer_to_corner) .and. .not. halo_swap_state%involve_corners) &
193 .or. (.not.
present(copy_from_halo_buffer_to_corner) .and. &
194 halo_swap_state%involve_corners))
then 195 call log_log(
log_warn,
"Inconsistent halo swap corner state and corner subroutine call arguments")
198 if (
present(source_data))
then 199 call perform_local_data_copy(current_state, halo_swap_state%halo_depth, &
200 halo_swap_state%involve_corners, source_data)
202 call perform_local_data_copy(current_state, halo_swap_state%halo_depth, &
203 halo_swap_state%involve_corners)
205 if (halo_swap_state%number_distinct_neighbours .gt. 0)
then 206 call mpi_waitall(
size(halo_swap_state%recv_requests), halo_swap_state%recv_requests, &
207 mpi_statuses_ignore, ierr)
208 if (
present(source_data))
then 209 if (halo_swap_state%involve_corners .and.
present(copy_from_halo_buffer_to_corner))
then 211 copy_from_halo_buffer, copy_from_halo_buffer_to_corner, source_data)
214 copy_from_halo_buffer, source_data=source_data)
217 if (halo_swap_state%involve_corners .and.
present(copy_from_halo_buffer_to_corner))
then 219 copy_from_halo_buffer, copy_from_halo_buffer_to_corner)
222 copy_from_halo_buffer)
225 call mpi_waitall(
size(halo_swap_state%send_requests), halo_swap_state%send_requests, &
226 mpi_statuses_ignore, ierr)
228 halo_swap_state%swap_in_progress=.false.
244 copy_corners_to_halo_buffer, source_data)
252 halo_swap_state%swap_in_progress = .true.
253 if (halo_swap_state%number_distinct_neighbours .gt. 0)
then 254 halo_swap_state%send_requests(:) = mpi_request_null
255 halo_swap_state%recv_requests(:) = mpi_request_null
257 if ((
present(copy_corners_to_halo_buffer) .and. .not. halo_swap_state%involve_corners) &
258 .or. (.not.
present(copy_corners_to_halo_buffer) .and. &
259 halo_swap_state%involve_corners))
then 260 call log_log(
log_warn,
"Inconsistent halo swap corner state and corner subroutine call arguments")
266 if (
present(source_data))
then 267 if (halo_swap_state%involve_corners .and.
present(copy_corners_to_halo_buffer))
then 268 call send_all_halos(current_state, halo_swap_state, copy_to_halo_buffer, &
269 copy_corners_to_halo_buffer, source_data)
271 call send_all_halos(current_state, halo_swap_state, copy_to_halo_buffer, &
272 source_data=source_data)
275 if (halo_swap_state%involve_corners .and.
present(copy_corners_to_halo_buffer))
then 276 call send_all_halos(current_state, halo_swap_state, copy_to_halo_buffer, &
277 copy_corners_to_halo_buffer)
279 call send_all_halos(current_state, halo_swap_state, copy_to_halo_buffer)
294 halo_depth, involve_corners)
297 logical,
intent(in) :: involve_corners
298 integer,
intent(in) :: halo_depth
301 integer :: number_comm_requests
303 halo_state%involve_corners = involve_corners
304 halo_state%halo_depth = halo_depth
306 current_state%local_grid, current_state%parallel%my_rank, involve_corners)
307 if (halo_state%number_distinct_neighbours .gt. 0)
then 308 allocate(halo_state%halo_swap_neighbours(halo_state%number_distinct_neighbours))
310 current_state%parallel%my_rank, halo_state%number_distinct_neighbours, involve_corners)
312 halo_state%number_distinct_neighbours, get_fields_per_halo_cell, &
313 halo_state%fields_per_cell, halo_depth)
315 halo_state%halo_swap_neighbours, &
316 halo_state%number_distinct_neighbours, halo_state%fields_per_cell)
318 halo_state%number_distinct_neighbours, &
319 halo_state%halo_swap_neighbours)
321 halo_state%halo_swap_neighbours, &
322 halo_state%number_distinct_neighbours, involve_corners)
324 halo_state%cell_match)
328 halo_state%number_distinct_neighbours)
329 allocate(halo_state%send_requests(number_comm_requests), &
330 halo_state%recv_requests(number_comm_requests))
331 halo_state%send_requests(:) = mpi_request_null
332 halo_state%recv_requests(:) = mpi_request_null
334 halo_state%initialised = .true.
345 if (
allocated(halo_swap_state%send_requests))
deallocate(halo_swap_state%send_requests)
346 if (
allocated(halo_swap_state%recv_requests))
deallocate(halo_swap_state%recv_requests)
348 do i=1,halo_swap_state%number_distinct_neighbours
349 if (
allocated(halo_swap_state%halo_swap_neighbours(i)%send_halo_buffer)) &
350 deallocate(halo_swap_state%halo_swap_neighbours(i)%send_halo_buffer)
351 if (
allocated(halo_swap_state%halo_swap_neighbours(i)%recv_halo_buffer)) &
352 deallocate(halo_swap_state%halo_swap_neighbours(i)%recv_halo_buffer)
354 if (
allocated(halo_swap_state%halo_swap_neighbours))
deallocate(halo_swap_state%halo_swap_neighbours)
355 halo_swap_state%initialised=.false.
367 x_target_index, y_target_index, halo_page)
369 integer,
intent(in) :: corner_loc, x_target_index, y_target_index, halo_page
370 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: halo_buffer
371 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
373 field_data(local_grid%local_domain_start_index(
z_index):&
374 local_grid%local_domain_end_index(
z_index),&
375 y_target_index, x_target_index)=halo_buffer(:,4,halo_page)
376 field_data(local_grid%local_domain_start_index(
z_index):&
377 local_grid%local_domain_end_index(
z_index),&
378 merge(y_target_index-1, y_target_index+1, corner_loc .lt. 3), x_target_index)=&
379 halo_buffer(:,3,halo_page)
380 field_data(local_grid%local_domain_start_index(
z_index):&
381 local_grid%local_domain_end_index(
z_index),&
382 y_target_index, merge(x_target_index-1, x_target_index+1, corner_loc == 1 .or.&
383 corner_loc == 3))= halo_buffer(:,2,halo_page)
384 field_data(local_grid%local_domain_start_index(
z_index):&
385 local_grid%local_domain_end_index(
z_index),&
386 merge(y_target_index-1, y_target_index+1, corner_loc .lt. 3), &
387 merge(x_target_index-1, x_target_index+1, corner_loc == 1 .or.&
388 corner_loc == 3))= halo_buffer(:,1,halo_page)
402 integer,
intent(in) :: dim, target_index, halo_page
403 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: halo_buffer
404 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
410 field_data(local_grid%local_domain_start_index(
z_index):&
411 local_grid%local_domain_end_index(
z_index),&
412 local_grid%local_domain_start_index(
y_index):&
413 local_grid%local_domain_end_index(
y_index), &
414 target_index) = halo_buffer(:,:,halo_page)
416 field_data(local_grid%local_domain_start_index(
z_index):&
417 local_grid%local_domain_end_index(
z_index), target_index, &
418 local_grid%local_domain_start_index(
x_index):&
419 local_grid%local_domain_end_index(
x_index)) = &
420 halo_buffer(:,:,halo_page)
434 integer,
intent(in) :: dim, source_index, halo_page
435 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: halo_buffer
436 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
439 halo_buffer(:,:,halo_page) = field_data(local_grid%local_domain_start_index(
z_index):&
440 local_grid%local_domain_end_index(
z_index), &
441 local_grid%local_domain_start_index(
y_index):&
442 local_grid%local_domain_end_index(
y_index), source_index)
444 halo_buffer(:,:,halo_page)=field_data(local_grid%local_domain_start_index(
z_index):&
445 local_grid%local_domain_end_index(
z_index), source_index, &
446 local_grid%local_domain_start_index(
x_index):&
447 local_grid%local_domain_end_index(
x_index))
460 x_source_index, y_source_index, halo_page)
462 integer,
intent(in) :: corner_loc, x_source_index, y_source_index, halo_page
463 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: halo_buffer
464 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
467 halo_buffer(:,1,halo_page) = field_data(:,y_source_index,x_source_index)
468 halo_buffer(:,2,halo_page) = field_data(:, merge(y_source_index-1, y_source_index+1, &
469 corner_loc .lt. 3), x_source_index)
470 halo_buffer(:,3,halo_page) = field_data(:, y_source_index, &
471 merge(x_source_index-1,x_source_index+1, corner_loc == 1 .or. corner_loc == 3))
472 halo_buffer(:,4,halo_page) = field_data(:, merge(y_source_index-1, y_source_index+1, &
473 corner_loc .lt. 3), merge(x_source_index-1,x_source_index+1, corner_loc == 1 .or. &
483 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
485 integer,
intent(in) :: my_rank, halo_depth
486 logical,
intent(in) :: involve_corners
506 number_distinct_neighbours)
507 integer,
intent(in) :: number_distinct_neighbours
513 do i=1, number_distinct_neighbours
514 if (halo_swap_neighbours(i)%recv_size .gt. 0) &
516 if (halo_swap_neighbours(i)%recv_corner_size .gt. 0)&
528 number_distinct_neighbours, involve_corners)
530 integer,
intent(in) :: number_distinct_neighbours
531 logical,
intent(in) :: involve_corners
534 integer :: i, normal_size, corner_size
536 do i=1, number_distinct_neighbours
537 if (halo_swap_neighbours(i)%halo_pages .gt. 0)
then 538 if (halo_swap_neighbours(i)%dimension == 0)
then 539 call log_log(
log_error,
"Halo swapping with neighbour needed but dimension is 0 which suggests corner only")
542 halo_swap_neighbours(i)%dimension==
x_index)*halo_swap_neighbours(i)%halo_pages
546 halo_swap_neighbours(i)%recv_size=normal_size
547 halo_swap_neighbours(i)%send_size=normal_size
548 if (involve_corners .and. halo_swap_neighbours(i)%halo_corners .gt. 0)
then 550 corner_size=local_grid%size(
z_index)*4*halo_swap_neighbours(i)%halo_corners
554 halo_swap_neighbours(i)%recv_corner_size=corner_size
555 halo_swap_neighbours(i)%send_corner_size=corner_size
577 integer,
intent(in) :: pid
582 do i=1,
size(local_grid%corner_neighbours, 2)
583 do j=1,
size(local_grid%corner_neighbours, 1)
584 if (local_grid%corner_neighbours(j,i) == pid)
then 604 integer,
intent(in) :: my_rank
605 logical,
intent(in) :: include_corners
607 integer :: i, j, temp_neighbour_pids(merge(16, 8, include_corners))
609 temp_neighbour_pids(:)=-1
613 if (local_grid%neighbours(i,j) .ne. my_rank .and. .not. &
615 local_grid%neighbours(i,j)))
then 619 local_grid%neighbours(i,j)
624 if (include_corners)
then 625 do i=1,
size(local_grid%corner_neighbours, 2)
626 do j=1,
size(local_grid%corner_neighbours, 1)
627 if (local_grid%corner_neighbours(j,i) .ne. my_rank .and. .not. &
629 local_grid%corner_neighbours(j,i)))
then 633 local_grid%corner_neighbours(j,i)
649 integer,
intent(in) :: my_rank, number_distinct_neighbours
650 logical,
intent(in) :: involve_corners
653 populate_halo_swap_neighbours
654 integer :: i, j, current_pid_location, temp_neighbour_pids(merge(16, 8, involve_corners))
656 current_pid_location=0
657 temp_neighbour_pids(:)=-1
660 if (local_grid%neighbours(i,j) .ne. my_rank .and. .not. &
662 local_grid%neighbours(i,j)))
then 663 current_pid_location=current_pid_location+1
664 populate_halo_swap_neighbours(current_pid_location)%pid=local_grid%neighbours(i,j)
665 temp_neighbour_pids(current_pid_location)=local_grid%neighbours(i,j)
666 populate_halo_swap_neighbours(current_pid_location)%dimension=i
671 if (involve_corners)
then 672 do i=1,
size(local_grid%corner_neighbours, 2)
673 do j=1,
size(local_grid%corner_neighbours, 1)
674 if (local_grid%corner_neighbours(j,i) .ne. my_rank .and. .not. &
676 local_grid%corner_neighbours(j,i)))
then 677 current_pid_location=current_pid_location+1
678 populate_halo_swap_neighbours(current_pid_location)%pid = &
679 local_grid%corner_neighbours(j,i)
680 temp_neighbour_pids(current_pid_location)=local_grid%corner_neighbours(j,i)
681 populate_halo_swap_neighbours(current_pid_location)%dimension=0
697 number_distinct_neighbours, get_fields_per_halo_cell, fields_per_cell, halo_depth)
700 integer,
intent(in) :: number_distinct_neighbours, halo_depth
702 integer,
intent(out) :: fields_per_cell
704 integer :: i, j, pid_location, halo_start, halo_end
706 fields_per_cell = get_fields_per_halo_cell(current_state)
707 halo_start = merge(2, 1, halo_depth==1)
708 halo_end = merge(3, 4, halo_depth==1)
712 do j = halo_start, halo_end
713 if (current_state%parallel%my_rank .ne. current_state%local_grid%neighbours(i,j))
then 715 current_state%local_grid%neighbours(i,j), number_distinct_neighbours)
716 halo_swap_neighbours(pid_location)%halo_pages = &
717 halo_swap_neighbours(pid_location)%halo_pages + fields_per_cell
731 number_distinct_neighbours, fields_per_cell)
734 integer,
intent(in) :: number_distinct_neighbours, fields_per_cell
736 integer :: i, j, pid_location
738 do i=1,
size(current_state%local_grid%corner_neighbours, 2)
739 do j=1,
size(current_state%local_grid%corner_neighbours, 1)
740 if (current_state%parallel%my_rank .ne. &
741 current_state%local_grid%corner_neighbours(j,i))
then 743 current_state%local_grid%corner_neighbours(j,i), number_distinct_neighbours)
744 halo_swap_neighbours(pid_location)%halo_corners = &
745 halo_swap_neighbours(pid_location)%halo_corners + fields_per_cell
756 halo_swap_neighbours)
758 integer,
intent(in) :: number_distinct_neighbours
764 do i=1,number_distinct_neighbours
765 if (halo_swap_neighbours(i)%halo_pages .gt. 0)
then 767 allocate(halo_swap_neighbours(i)%send_halo_buffer(local_grid%size(
z_index), &
769 halo_swap_neighbours(i)%dimension==
x_index), halo_swap_neighbours(i)%halo_pages))
770 allocate(halo_swap_neighbours(i)%recv_halo_buffer(local_grid%size(
z_index), &
772 halo_swap_neighbours(i)%dimension==
x_index), halo_swap_neighbours(i)%halo_pages))
774 if (halo_swap_neighbours(i)%halo_corners .gt. 0)
then 776 allocate(halo_swap_neighbours(i)%send_corner_buffer(local_grid%size(
z_index), 4, &
777 halo_swap_neighbours(i)%halo_corners))
778 allocate(halo_swap_neighbours(i)%recv_corner_buffer(local_grid%size(
z_index), 4, &
779 halo_swap_neighbours(i)%halo_corners))
791 integer,
intent(in) :: halo_depth
792 integer,
intent(out) :: cell_match(:,:)
794 logical,
dimension(3) :: same_neighbours
800 if (halo_depth == 1)
then 802 cell_match(i,2)=merge(2, 3, .not. same_neighbours(i))
803 cell_match(i,3)=merge(3, 2, .not. same_neighbours(i))
807 cell_match(i,j) = merge(j, j+halo_depth, .not. same_neighbours(i))
808 cell_match(i,j+halo_depth) = merge(j+halo_depth, j, .not. same_neighbours(i))
825 integer :: i, request_counter, ierr
829 do i = 1, halo_swap_state%number_distinct_neighbours
830 if (halo_swap_state%halo_swap_neighbours(i)%recv_size .gt. 0)
then 831 call mpi_irecv(halo_swap_state%halo_swap_neighbours(i)%recv_halo_buffer, &
832 halo_swap_state%halo_swap_neighbours(i)%recv_size,
precision_type, &
833 halo_swap_state%halo_swap_neighbours(i)%pid, 0, &
834 current_state%parallel%neighbour_comm, &
835 halo_swap_state%recv_requests(request_counter), ierr)
836 request_counter = request_counter + 1
838 if (halo_swap_state%halo_swap_neighbours(i)%recv_corner_size .gt. 0)
then 839 call mpi_irecv(halo_swap_state%halo_swap_neighbours(i)%recv_corner_buffer, &
840 halo_swap_state%halo_swap_neighbours(i)%recv_corner_size,
precision_type, &
841 halo_swap_state%halo_swap_neighbours(i)%pid, 0, &
842 current_state%parallel%neighbour_comm, &
843 halo_swap_state%recv_requests(request_counter), ierr)
844 request_counter = request_counter + 1
858 subroutine send_all_halos(current_state, halo_swap_state, copy_fields_to_halo_buffer, &
859 copy_corner_fields_to_halo_buffer, source_data)
865 copy_corner_fields_to_halo_buffer
868 integer :: i, j, ierr,hstart, hend, pid_location, source_index, request_number, &
869 x_source_index, y_source_index, current_page(halo_swap_state%number_distinct_neighbours),&
872 halo_depth = current_state%local_grid%halo_size(
y_index)
876 hstart = merge(2, 1, halo_swap_state%halo_depth==1)
877 hend = merge(3, halo_depth*2, halo_swap_state%halo_depth==1)
881 if (current_state%parallel%my_rank .ne. current_state%local_grid%neighbours(i,j))
then 884 source_index = current_state%local_grid%local_domain_start_index(i)
886 source_index = current_state%local_grid%local_domain_start_index(i) + &
887 merge(1, 0, halo_swap_state%halo_depth .ne. 1)
889 source_index = current_state%local_grid%local_domain_end_index(i) - &
890 merge(1, 0, halo_swap_state%halo_depth .ne. 1)
892 source_index = current_state%local_grid%local_domain_end_index(i)
896 current_state%local_grid%neighbours(i,j), &
897 halo_swap_state%number_distinct_neighbours)
899 if (
present(source_data))
then 900 call copy_fields_to_halo_buffer(current_state, &
901 halo_swap_state%halo_swap_neighbours(pid_location), i, source_index, &
902 pid_location, current_page, source_data)
904 call copy_fields_to_halo_buffer(current_state, &
905 halo_swap_state%halo_swap_neighbours(pid_location), i, source_index, &
906 pid_location, current_page)
917 if (
present(copy_corner_fields_to_halo_buffer))
then 919 do j = 1,
size(current_state%local_grid%corner_neighbours, 1)
920 if (current_state%parallel%my_rank .ne. &
921 current_state%local_grid%corner_neighbours(j,1))
then 922 x_source_index = merge(current_state%local_grid%local_domain_start_index(
x_index)+1,&
923 current_state%local_grid%local_domain_end_index(
x_index)-1, j==1 .or. j==3)
924 y_source_index =merge(current_state%local_grid%local_domain_start_index(
y_index)+1,&
925 current_state%local_grid%local_domain_end_index(
y_index)-1, j==1 .or. j==2)
927 current_state%local_grid%corner_neighbours(j,1), &
928 halo_swap_state%number_distinct_neighbours)
929 if (
present(source_data))
then 930 call copy_corner_fields_to_halo_buffer(current_state, &
931 halo_swap_state%halo_swap_neighbours(pid_location), j, &
932 x_source_index, y_source_index, pid_location, current_page, source_data)
934 call copy_corner_fields_to_halo_buffer(current_state, &
935 halo_swap_state%halo_swap_neighbours(pid_location), j, &
936 x_source_index, y_source_index, pid_location, current_page)
942 do i=1,halo_swap_state%number_distinct_neighbours
943 if (halo_swap_state%halo_swap_neighbours(i)%send_size .gt. 0)
then 944 call mpi_isend(halo_swap_state%halo_swap_neighbours(i)%send_halo_buffer, &
945 halo_swap_state%halo_swap_neighbours(i)%send_size,
precision_type, &
946 halo_swap_state%halo_swap_neighbours(i)%pid, 0, &
947 current_state%parallel%neighbour_comm, &
948 halo_swap_state%send_requests(request_number), ierr)
949 request_number = request_number+1
951 if (halo_swap_state%halo_swap_neighbours(i)%send_corner_size .gt. 0)
then 952 call mpi_isend(halo_swap_state%halo_swap_neighbours(i)%send_corner_buffer, &
953 halo_swap_state%halo_swap_neighbours(i)%send_corner_size,
precision_type, &
954 halo_swap_state%halo_swap_neighbours(i)%pid, 0, &
955 current_state%parallel%neighbour_comm, &
956 halo_swap_state%send_requests(request_number), ierr)
957 request_number = request_number+1
970 copy_halo_buffer_to_corner, source_data)
977 integer :: i, j, hstart, hend, pid_location, target_index, x_target_index, &
978 y_target_index, current_page(halo_swap_state%number_distinct_neighbours)
980 hstart=merge(2, 1, halo_swap_state%halo_depth==1)
981 hend=merge(3, 4, halo_swap_state%halo_depth==1)
986 if (current_state%parallel%my_rank .ne. current_state%local_grid%neighbours(i,j))
then 987 if (j==halo_swap_state%cell_match(i, 1))
then 989 else if (j==halo_swap_state%cell_match(i, 2))
then 991 else if (j==halo_swap_state%cell_match(i, 3))
then 992 target_index=current_state%local_grid%local_domain_end_index(i)+1
993 else if (j==halo_swap_state%cell_match(i, 4))
then 994 target_index=current_state%local_grid%local_domain_end_index(i)+2
997 current_state%local_grid%neighbours(i,j),&
998 halo_swap_state%number_distinct_neighbours)
999 if (
present(source_data))
then 1000 call copy_halo_buffer_to_field(current_state, &
1001 halo_swap_state%halo_swap_neighbours(pid_location), i, target_index,&
1002 pid_location, current_page, source_data)
1004 call copy_halo_buffer_to_field(current_state, &
1005 halo_swap_state%halo_swap_neighbours(pid_location), i, target_index,&
1006 pid_location, current_page)
1011 if (
present(copy_halo_buffer_to_corner))
then 1013 do j=
size(current_state%local_grid%corner_neighbours, 1),1,-1
1014 if (current_state%parallel%my_rank .ne. &
1015 current_state%local_grid%corner_neighbours(j,1))
then 1016 x_target_index=merge(2, current_state%local_grid%local_domain_end_index(
x_index)+1,&
1018 y_target_index=merge(2, current_state%local_grid%local_domain_end_index(
y_index)+1, &
1021 current_state%local_grid%corner_neighbours(j,1), &
1022 halo_swap_state%number_distinct_neighbours)
1023 if (
present(source_data))
then 1024 call copy_halo_buffer_to_corner(current_state, &
1025 halo_swap_state%halo_swap_neighbours(pid_location), j, &
1026 x_target_index, y_target_index, pid_location, current_page, source_data)
1028 call copy_halo_buffer_to_corner(current_state, &
1029 halo_swap_state%halo_swap_neighbours(pid_location), j, &
1030 x_target_index, y_target_index, pid_location, current_page)
1043 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
1044 integer,
intent(in) :: my_rank
1046 integer :: i, y_source_index, x_source_index, y_target_index, x_target_index
1048 do i=1,
size(local_grid%corner_neighbours, 1)
1049 if (my_rank .eq. local_grid%corner_neighbours(i,1))
then 1051 y_source_index=local_grid%local_domain_end_index(
y_index)-1
1052 x_source_index=local_grid%local_domain_end_index(
x_index)-1
1056 y_source_index=local_grid%local_domain_end_index(
y_index)-1
1057 x_source_index=local_grid%local_domain_start_index(
x_index)
1059 x_target_index=local_grid%local_domain_end_index(
x_index)+1
1061 y_source_index=local_grid%local_domain_start_index(
y_index)
1062 x_source_index=local_grid%local_domain_end_index(
x_index)-1
1063 y_target_index=local_grid%local_domain_end_index(
y_index)+1
1066 y_source_index=local_grid%local_domain_start_index(
y_index)
1067 x_source_index=local_grid%local_domain_start_index(
x_index)
1068 y_target_index=local_grid%local_domain_end_index(
y_index)+1
1069 x_target_index=local_grid%local_domain_end_index(
x_index)+1
1072 field_data(local_grid%local_domain_start_index(
z_index):&
1073 local_grid%local_domain_end_index(
z_index),&
1074 y_target_index, x_target_index)=&
1075 field_data(local_grid%local_domain_start_index(
z_index):&
1076 local_grid%local_domain_end_index(
z_index),y_source_index, x_source_index)
1078 field_data(local_grid%local_domain_start_index(
z_index):&
1079 local_grid%local_domain_end_index(
z_index),&
1080 y_target_index+1, x_target_index)=&
1081 field_data(local_grid%local_domain_start_index(
z_index):&
1082 local_grid%local_domain_end_index(
z_index),y_source_index+1, x_source_index)
1084 field_data(local_grid%local_domain_start_index(
z_index):&
1085 local_grid%local_domain_end_index(
z_index),&
1086 y_target_index, x_target_index+1)=&
1087 field_data(local_grid%local_domain_start_index(
z_index):&
1088 local_grid%local_domain_end_index(
z_index),y_source_index, x_source_index+1)
1090 field_data(local_grid%local_domain_start_index(
z_index):&
1091 local_grid%local_domain_end_index(
z_index),&
1092 y_target_index+1, x_target_index+1)=&
1093 field_data(local_grid%local_domain_start_index(
z_index):&
1094 local_grid%local_domain_end_index(
z_index),y_source_index+1, x_source_index+1)
1107 real(kind=DEFAULT_PRECISION),
dimension(:,:,:),
intent(inout) :: field_data
1108 integer,
intent(in) :: dim, my_rank, halo_depth
1110 integer i, target_index, source_index, hstart, hend
1112 hstart=merge(2,1, halo_depth==1)
1113 hend=merge(3,4, halo_depth==1)
1116 if (local_grid%neighbours(dim,i) .eq. my_rank)
then 1119 source_index=local_grid%local_domain_end_index(dim)-1
1122 source_index=local_grid%local_domain_end_index(dim)
1124 target_index=local_grid%local_domain_end_index(dim)+1
1125 source_index=local_grid%local_domain_start_index(dim)
1127 target_index=local_grid%local_domain_end_index(dim)+2
1128 source_index=local_grid%local_domain_start_index(dim)+1
1131 field_data(local_grid%local_domain_start_index(
z_index):&
1132 local_grid%local_domain_end_index(
z_index),&
1133 local_grid%local_domain_start_index(
y_index):&
1134 local_grid%local_domain_end_index(
y_index), target_index) = &
1135 field_data(local_grid%local_domain_start_index(
z_index):&
1136 local_grid%local_domain_end_index(
z_index),&
1137 local_grid%local_domain_start_index(
y_index):&
1138 local_grid%local_domain_end_index(
y_index), source_index)
1140 field_data(local_grid%local_domain_start_index(
z_index):&
1141 local_grid%local_domain_end_index(
z_index),&
1142 target_index, local_grid%local_domain_start_index(
x_index):&
1143 local_grid%local_domain_end_index(
x_index)) = &
1144 field_data(local_grid%local_domain_start_index(
z_index):&
1145 local_grid%local_domain_end_index(
z_index),&
1146 source_index, local_grid%local_domain_start_index(
x_index):&
1147 local_grid%local_domain_end_index(
x_index))
1162 logical,
dimension(3) :: retrieve_same_neighbour_information
1164 integer :: i, nd1, nd2
1166 retrieve_same_neighbour_information=(/.true., .true., .true./)
1170 do i = 1,local_grid%halo_size(
y_index)*2
1172 nd1=local_grid%neighbours(
y_index,i)
1173 nd2=local_grid%neighbours(
x_index,i)
1175 if (local_grid%neighbours(
y_index,i) .ne. nd1) &
1176 retrieve_same_neighbour_information(
y_index) = .false.
1177 if (local_grid%neighbours(
x_index,i) .ne. nd2) &
1178 retrieve_same_neighbour_information(
x_index) = .false.
1188 integer,
intent(in) :: pid, temp_neighbour_pids(8)
1194 if (temp_neighbour_pids(i) == pid)
return 1195 if (temp_neighbour_pids(i) == -1)
then 1210 number_distinct_neighbours)
1212 integer,
intent(in) :: pid, number_distinct_neighbours
1216 do i=1, number_distinct_neighbours
1217 if (halo_swap_neighbours(i)%pid == pid)
then subroutine perform_local_data_copy_for_corners(my_rank, local_grid, field_data)
Performs a local data copy for corners when the neighbour is local (me)
integer, public precision_type
integer function get_pid_neighbour_location(halo_swap_neighbours, pid, number_distinct_neighbours)
Given the process id of a neighbour this determines the location in the data structure of correspondi...
subroutine, public init_halo_communication(current_state, get_fields_per_halo_cell, halo_state, halo_depth, involve_corners)
Initialises a halo swapping state, by determining the neighbours, size of data in each swap and alloc...
subroutine, public blocking_halo_swap(current_state, halo_swap_state, copy_to_halo_buffer, perform_local_data_copy, copy_from_halo_buffer, copy_corners_to_halo_buffer, copy_from_halo_buffer_to_corner, source_data)
Performs the entire halo swap operation, this is simply a wrapper around the nonblocking initiate and...
subroutine, public perform_local_data_copy_for_field(field_data, local_grid, my_rank, halo_depth, involve_corners)
Will perform a a local copy for the halo data of a field.
logical function, dimension(3) retrieve_same_neighbour_information(local_grid)
Retrieves whether we have the same neighbours for L and R halo swaps in each dimension.
integer, parameter, public log_error
Only log ERROR messages.
Contains the types used for communication, holding the state of communications and supporting activit...
Contains prognostic field definitions and functions.
A prognostic field which is assumed to be 3D.
Procedure interfaces used to determine the policy (i.e. the fields) of halo swapping and...
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
integer, parameter, public z_index
Grid index parameters.
subroutine, public complete_nonblocking_halo_swap(current_state, halo_swap_state, perform_local_data_copy, copy_from_halo_buffer, copy_from_halo_buffer_to_corner, source_data)
This completes a nonblocking halo swap and it is only during this call that the data fields themselve...
Contains common definitions for the data and datatypes used by MONC.
The ModelState which represents the current state of a run.
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
subroutine, public copy_buffer_to_corner(local_grid, halo_buffer, field_data, corner_loc, x_target_index, y_target_index, halo_page)
Copies the received buffer for a specific field to the corresponding corner of that field...
Conversion between common inbuilt FORTRAN data types.
subroutine allocate_halo_buffers_for_each_neighbour(local_grid, number_distinct_neighbours, halo_swap_neighbours)
Allocates the locally stored halo buffers (send and receive) for each neighbouring process...
integer function determine_halo_corner_element_sizes(local_grid, pid)
For a specific process id this determines the number of halo swap corner elements to involve in a com...
Converts data types to strings.
subroutine, public copy_buffer_to_field(local_grid, halo_buffer, field_data, dim, target_index, halo_page)
Copies the received buffer for a specific field to the corresponding halo data of that prognostic fie...
subroutine determine_recv_and_send_sizes(local_grid, halo_swap_neighbours, number_distinct_neighbours, involve_corners)
Determines the amount (in elements) of data that each neighbour will be sent and I will receive from ...
subroutine, public initiate_nonblocking_halo_swap(current_state, halo_swap_state, copy_to_halo_buffer, copy_corners_to_halo_buffer, source_data)
Initiates a non blocking halo swap, this registers the receive requests, copies data into the send bu...
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Maintains the state of a halo swap and contains buffers, neighbours etc.
subroutine, public copy_corner_to_buffer(local_grid, halo_buffer, field_data, corner_loc, x_source_index, y_source_index, halo_page)
Copies prognostic field corner data to send buffer for specific field.
integer function, public get_single_field_per_halo_cell(current_state)
A very common function, which returns a single field per halo cell which is used to halo swap just on...
type(neighbour_description_type) function, dimension(number_distinct_neighbours) populate_halo_swap_neighbours(local_grid, my_rank, number_distinct_neighbours, involve_corners)
Will populate the halo swap neighbour data strutures with appropriate neighbour pid and dimension num...
Defined the local grid, i.e. the grid held on this process after decomposition.
integer function get_number_of_processes_involved_in_communication(local_grid, my_rank, include_corners)
Deduces the number of distinct neighbours that will be involved in a halo swap. This information is u...
subroutine copy_buffer_data_for_prognostics(current_state, halo_swap_state, copy_halo_buffer_to_field, copy_halo_buffer_to_corner, source_data)
Copies the received data (held in buffers) from neighbours into the correct halo location in the prog...
integer function get_number_communication_requests(halo_swap_neighbours, number_distinct_neighbours)
Determines the overall number of communication requests, which is made up of normal halo swaps and po...
Provides the mechanism for halo swapping. This module contains the functionality required to determin...
integer, parameter, public log_warn
Log WARNING and ERROR messages.
Describes the neighbours of a process in a specific dimension and contains the communication buffers ...
subroutine, public copy_field_to_buffer(local_grid, halo_buffer, field_data, dim, source_index, halo_page)
Copies prognostic field data to send buffer for specific field, dimension, halo cell.
subroutine recv_all_halos(current_state, halo_swap_state)
Registers receive requests for all prognostic fields from the appropriate neighbouring processes (tha...
Functionality to support the different types of grid and abstraction between global grids and local o...
logical function has_pid_already_been_seen(temp_neighbour_pids, pid)
Returns whether or not a specific process id has already been "seen" by searching a list of already s...
integer function determine_halo_corner_size(local_grid)
Determine the halo corner size in elements.
subroutine deduce_halo_pages_per_neighbour(current_state, halo_swap_neighbours, number_distinct_neighbours, get_fields_per_halo_cell, fields_per_cell, halo_depth)
Deduces the number of halo pages per neighbour halo swap and places this information in the appropria...
subroutine deduce_halo_corners_per_neighbour(current_state, halo_swap_neighbours, number_distinct_neighbours, fields_per_cell)
Determines the number of halo corners to swap between specific neighours, this is similar to deducing...
subroutine, public finalise_halo_communication(halo_swap_state)
Finalises the halo swap represented by the state by freeing up all the allocated memory.
subroutine perform_local_data_copy_for_dimension(dim, my_rank, halo_depth, local_grid, field_data)
Performs a local data copy for a specific dimension of a prognostic field.
subroutine generate_recv_field_buffer_matches(current_state, halo_depth, cell_match)
Precalculates the received buffer to field halo cell matches for each dimension and called from the i...
The model state which represents the current state of a run.
integer, parameter, public y_index
integer, parameter, public x_index
subroutine send_all_halos(current_state, halo_swap_state, copy_fields_to_halo_buffer, copy_corner_fields_to_halo_buffer, source_data)
Copies all applicable bits of the prognostics into a send buffer for each neighbour and then issues a...