MONC
Data Types | Functions/Subroutines
halo_communication_mod Module Reference

Provides the mechanism for halo swapping. This module contains the functionality required to determine what messages get sent where, pack data up en mass, send and receive from neighbouring processes and unpack into the appropriate data locations. There is also some functionality to help local copying of data to corresponding locations. The idea is that the caller will determine the policy (i.e. exactly what fields are to be communicated) through procedure arguments and this mechanism can be used again and again. It implements bulk sending of all field data in one large message to reduce communication overhead. More...

Data Types

interface  copy_corners_to_halo_buffer_proc_interface
 
interface  copy_fields_to_halo_buffer_proc_interface
 
interface  copy_halo_buffer_to_corner_proc_interface
 
interface  copy_halo_buffer_to_field_proc_interface
 
interface  get_fields_per_halo_cell_proc_interface
 Procedure interfaces used to determine the policy (i.e. the fields) of halo swapping and. More...
 
interface  perform_local_data_copy_proc_interface
 

Functions/Subroutines

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 complete procedures and saves the programmer from calling these directly if they do not wish to interleave any computation. More...
 
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 themselves are modified. This will perform local data copying, wait for all outstanding receives to complete, copy the received buffer data into the data and then wait for all outstanding sends to complete. More...
 
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 buffer and then registers send requests for these. As far as this call is concerned, the data is immutable - no modifications will take place to it until the matching completion call is made. More...
 
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 allocating the required memory - which are communication buffers and request handles. All this information is placed into the returned state. More...
 
subroutine, public finalise_halo_communication (halo_swap_state)
 Finalises the halo swap represented by the state by freeing up all the allocated memory. More...
 
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. More...
 
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 field. More...
 
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. More...
 
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. More...
 
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. More...
 
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 potentially corner swaps too if that is enabled. More...
 
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 in a halo swap. More...
 
integer function determine_halo_corner_size (local_grid)
 Determine the halo corner size in elements. More...
 
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 communication. More...
 
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 used to then allocate the appropriate amount of memory to store the neighbour halo swapping data structure. More...
 
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 numbers. More...
 
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 appropriate data structures. We call a "page" of data the contiguous data of a field that we are going to send, such as halo 1 of w, halo 1 of zw and halo 2 of w (assuming these go to the same neighbour as 1) More...
 
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 the number of halo pages per neighbour, only it is for corners. More...
 
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. More...
 
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 initialisation stage. More...
 
subroutine recv_all_halos (current_state, halo_swap_state)
 Registers receive requests for all prognostic fields from the appropriate neighbouring processes (that we have already deduced in the initialisation stage.) More...
 
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 asynchronous sends of this data. More...
 
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 prognostic fields. More...
 
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) More...
 
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. More...
 
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. More...
 
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 seen process ids. More...
 
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 corresponding data for that. Note that this is an O(n) operation so ideally call once and reuse the results. If no process id is found then returns -1. More...
 
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 one field. More...
 

Detailed Description

Provides the mechanism for halo swapping. This module contains the functionality required to determine what messages get sent where, pack data up en mass, send and receive from neighbouring processes and unpack into the appropriate data locations. There is also some functionality to help local copying of data to corresponding locations. The idea is that the caller will determine the policy (i.e. exactly what fields are to be communicated) through procedure arguments and this mechanism can be used again and again. It implements bulk sending of all field data in one large message to reduce communication overhead.

Function/Subroutine Documentation

◆ allocate_halo_buffers_for_each_neighbour()

subroutine halo_communication_mod::allocate_halo_buffers_for_each_neighbour ( type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  number_distinct_neighbours,
type(neighbour_description_type), dimension(:), intent(inout), allocatable  halo_swap_neighbours 
)
private

Allocates the locally stored halo buffers (send and receive) for each neighbouring process.

Parameters
local_gridDescription of the local grid
halo_swap_neighboursMy neighbouring PIDs
number_distinct_neighboursThe number of distinct neighbours that I have

Definition at line 757 of file halocommunication.F90.

757  type(local_grid_type), intent(inout) :: local_grid
758  integer, intent(in) :: number_distinct_neighbours
759  type(neighbour_description_type), dimension(:), allocatable, intent(inout) :: &
760  halo_swap_neighbours
761 
762  integer :: i
763 
764  do i=1,number_distinct_neighbours
765  if (halo_swap_neighbours(i)%halo_pages .gt. 0) then
766  !depending on the direction of the swapping, the send and recv buffer size would change
767  allocate(halo_swap_neighbours(i)%send_halo_buffer(local_grid%size(z_index), &
768  merge(local_grid%size(y_index), local_grid%size(x_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), &
771  merge(local_grid%size(y_index), local_grid%size(x_index), &
772  halo_swap_neighbours(i)%dimension==x_index), halo_swap_neighbours(i)%halo_pages))
773  end if
774  if (halo_swap_neighbours(i)%halo_corners .gt. 0) then
775  ! is 4 because of the 4 cells to swap since the halo_depth is 2(2 in x and 2 in y)??
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))
780  end if
781  end do
Here is the caller graph for this function:

◆ blocking_halo_swap()

subroutine, public halo_communication_mod::blocking_halo_swap ( type(model_state_type), intent(inout)  current_state,
type(halo_communication_type), intent(inout)  halo_swap_state,
procedure(copy_fields_to_halo_buffer_proc_interface copy_to_halo_buffer,
procedure(perform_local_data_copy_proc_interface perform_local_data_copy,
procedure(copy_halo_buffer_to_field_proc_interface copy_from_halo_buffer,
procedure(copy_corners_to_halo_buffer_proc_interface), optional  copy_corners_to_halo_buffer,
procedure(copy_halo_buffer_to_corner_proc_interface), optional  copy_from_halo_buffer_to_corner,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)

Performs the entire halo swap operation, this is simply a wrapper around the nonblocking initiate and complete procedures and saves the programmer from calling these directly if they do not wish to interleave any computation.

Parameters
current_stateThe current model state
halo_swap_stateThe halo swapping state
copy_to_halo_bufferProcedure pointer which is called to copy the data into the halo send buffer
perform_local_data_copyProcedure pointer which performs local data copying (where the neighbour is the local process)
copy_from_halo_bufferProcedure pointer which copies received data from the halo buffer into the field
copy_corners_to_halo_bufferOptional procedure pointer which copies data in corners to the halo send buffer
copy_from_halo_buffer_to_cornerOptional procedure pointer which copies received data into halo corners
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 112 of file halocommunication.F90.

112 
113  type(model_state_type), intent(inout) :: current_state
114  type(halo_communication_type), intent(inout) :: halo_swap_state
115  procedure(copy_fields_to_halo_buffer_proc_interface) :: copy_to_halo_buffer
116  procedure(copy_halo_buffer_to_field_proc_interface) :: copy_from_halo_buffer
117  procedure(perform_local_data_copy_proc_interface) :: perform_local_data_copy
118  procedure(copy_corners_to_halo_buffer_proc_interface), optional :: &
119  copy_corners_to_halo_buffer
120  procedure(copy_halo_buffer_to_corner_proc_interface), optional :: &
121  copy_from_halo_buffer_to_corner
122  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
123 
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")
129  end if
130 
131  if (present(source_data)) then
132  if (present(copy_corners_to_halo_buffer)) then
133  call initiate_nonblocking_halo_swap(current_state, halo_swap_state, &
134  copy_to_halo_buffer, copy_corners_to_halo_buffer, source_data)
135  else
136  call initiate_nonblocking_halo_swap(current_state, halo_swap_state, &
137  copy_to_halo_buffer, source_data=source_data)
138  end if
139  if (present(copy_from_halo_buffer_to_corner)) then
140  call complete_nonblocking_halo_swap(current_state, halo_swap_state, &
141  perform_local_data_copy, copy_from_halo_buffer, copy_from_halo_buffer_to_corner, &
142  source_data)
143  else
144  call complete_nonblocking_halo_swap(current_state, halo_swap_state, &
145  perform_local_data_copy, copy_from_halo_buffer, source_data=source_data)
146  end if
147  else
148  if (present(copy_corners_to_halo_buffer)) then
149  call initiate_nonblocking_halo_swap(current_state, halo_swap_state, copy_to_halo_buffer,&
150  copy_corners_to_halo_buffer)
151  else
152  call initiate_nonblocking_halo_swap(current_state, halo_swap_state, copy_to_halo_buffer)
153  end if
154  if (present(copy_from_halo_buffer_to_corner)) then
155  call complete_nonblocking_halo_swap(current_state, halo_swap_state, &
156  perform_local_data_copy, &
157  copy_from_halo_buffer, copy_from_halo_buffer_to_corner)
158  else
159  call complete_nonblocking_halo_swap(current_state, halo_swap_state, &
160  perform_local_data_copy, copy_from_halo_buffer)
161  end if
162  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ complete_nonblocking_halo_swap()

subroutine, public halo_communication_mod::complete_nonblocking_halo_swap ( type(model_state_type), intent(inout)  current_state,
type(halo_communication_type), intent(inout)  halo_swap_state,
procedure(perform_local_data_copy_proc_interface perform_local_data_copy,
procedure(copy_halo_buffer_to_field_proc_interface copy_from_halo_buffer,
procedure(copy_halo_buffer_to_corner_proc_interface), optional  copy_from_halo_buffer_to_corner,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)

This completes a nonblocking halo swap and it is only during this call that the data fields themselves are modified. This will perform local data copying, wait for all outstanding receives to complete, copy the received buffer data into the data and then wait for all outstanding sends to complete.

Parameters
current_stateThe current model state
halo_swap_stateThe halo swapping state
perform_local_data_copyProcedure pointer which performs local data copying (where the neighbour is the local process)
copy_from_halo_bufferProcedure pointer which copies received data from the halo buffer into the field
copy_from_halo_buffer_to_cornerOptional procedure pointer which copies halo data into field corners
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 182 of file halocommunication.F90.

182  type(model_state_type), intent(inout) :: current_state
183  type(halo_communication_type), intent(inout) :: halo_swap_state
184  procedure(copy_halo_buffer_to_field_proc_interface) :: copy_from_halo_buffer
185  procedure(perform_local_data_copy_proc_interface) :: perform_local_data_copy
186  procedure(copy_halo_buffer_to_corner_proc_interface), optional :: &
187  copy_from_halo_buffer_to_corner
188  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
189 
190  integer :: ierr
191 
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")
196  end if
197 
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)
201  else
202  call perform_local_data_copy(current_state, halo_swap_state%halo_depth, &
203  halo_swap_state%involve_corners)
204  end if
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
210  call copy_buffer_data_for_prognostics(current_state, halo_swap_state, &
211  copy_from_halo_buffer, copy_from_halo_buffer_to_corner, source_data)
212  else
213  call copy_buffer_data_for_prognostics(current_state, halo_swap_state, &
214  copy_from_halo_buffer, source_data=source_data)
215  end if
216  else
217  if (halo_swap_state%involve_corners .and. present(copy_from_halo_buffer_to_corner)) then
218  call copy_buffer_data_for_prognostics(current_state, halo_swap_state, &
219  copy_from_halo_buffer, copy_from_halo_buffer_to_corner)
220  else
221  call copy_buffer_data_for_prognostics(current_state, halo_swap_state, &
222  copy_from_halo_buffer)
223  end if
224  end if
225  call mpi_waitall(size(halo_swap_state%send_requests), halo_swap_state%send_requests, &
226  mpi_statuses_ignore, ierr)
227  end if
228  halo_swap_state%swap_in_progress=.false.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ copy_buffer_data_for_prognostics()

subroutine halo_communication_mod::copy_buffer_data_for_prognostics ( type(model_state_type), intent(inout)  current_state,
type(halo_communication_type), intent(inout)  halo_swap_state,
procedure(copy_halo_buffer_to_field_proc_interface copy_halo_buffer_to_field,
procedure(copy_halo_buffer_to_corner_proc_interface), optional  copy_halo_buffer_to_corner,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the received data (held in buffers) from neighbours into the correct halo location in the prognostic fields.

Parameters
current_stateThe current model state
halo_swap_stateThe halo swapping state
copy_halo_buffer_to_fieldProcedure pointer which copies the halo buffer into the data field
copy_halo_buffer_to_cornerOptional procedure pointer to copy halo buffer into field corners
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 971 of file halocommunication.F90.

971  type(model_state_type), intent(inout) :: current_state
972  type(halo_communication_type), intent(inout) :: halo_swap_state
973  procedure(copy_halo_buffer_to_field_proc_interface) :: copy_halo_buffer_to_field
974  procedure(copy_halo_buffer_to_corner_proc_interface), optional :: copy_halo_buffer_to_corner
975  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
976 
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)
979 
980  hstart=merge(2, 1, halo_swap_state%halo_depth==1)
981  hend=merge(3, 4, halo_swap_state%halo_depth==1)
982 
983  current_page(:)=1
984  do i=2, 3
985  do j=hstart, hend
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
988  target_index=1
989  else if (j==halo_swap_state%cell_match(i, 2)) then
990  target_index=2
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
995  end if
996  pid_location=get_pid_neighbour_location(halo_swap_state%halo_swap_neighbours, &
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)
1003  else
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)
1007  end if
1008  end if
1009  end do
1010  end do
1011  if (present(copy_halo_buffer_to_corner)) then
1012  current_page(:)=1
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,&
1017  j==1 .or. j==3)
1018  y_target_index=merge(2, current_state%local_grid%local_domain_end_index(y_index)+1, &
1019  j==1 .or. j==2)
1020  pid_location=get_pid_neighbour_location(halo_swap_state%halo_swap_neighbours, &
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)
1027  else
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)
1031  end if
1032  end if
1033  end do
1034  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ copy_buffer_to_corner()

subroutine, public halo_communication_mod::copy_buffer_to_corner ( type(local_grid_type), intent(inout)  local_grid,
real(kind=default_precision), dimension(:,:,:), intent(inout)  halo_buffer,
real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data,
integer, intent(in)  corner_loc,
integer, intent(in)  x_target_index,
integer, intent(in)  y_target_index,
integer, intent(in)  halo_page 
)

Copies the received buffer for a specific field to the corresponding corner of that field.

Parameters
local_gridDescription of the local grid
halo_bufferRaw halo buffer data
field_dataRaw prognostic field data
corner_locLocation of the corner
x_target_indexThe X target index that we copy into
y_target_indexThe Y target index that we copy into
halo_pageThe halo page to read from

Definition at line 368 of file halocommunication.F90.

368  type(local_grid_type), intent(inout) :: local_grid
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
372 
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)
Here is the caller graph for this function:

◆ copy_buffer_to_field()

subroutine, public halo_communication_mod::copy_buffer_to_field ( type(local_grid_type), intent(inout)  local_grid,
real(kind=default_precision), dimension(:,:,:), intent(inout)  halo_buffer,
real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data,
integer, intent(in)  dim,
integer, intent(in)  target_index,
integer, intent(in)  halo_page 
)

Copies the received buffer for a specific field to the corresponding halo data of that prognostic field.

Parameters
local_gridDescription of the local grid
halo_bufferRaw halo buffer data
field_dataRaw prognostic field data
dimThe dimension to copy into
target_indexThe target index in the dimension that we copy into
halo_pageThe halo page to read from

Definition at line 401 of file halocommunication.F90.

401  type(local_grid_type), intent(inout) :: local_grid
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
405 
406  ! If the neighbours are the same then reverse our placement of the data due to wrapping
407  ! around and order of
408  ! messages being sent. This is not an issue if the neighbours are different
409  if (dim == x_index) then
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)
415  else
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)
421  end if
Here is the caller graph for this function:

◆ copy_corner_to_buffer()

subroutine, public halo_communication_mod::copy_corner_to_buffer ( type(local_grid_type), intent(inout)  local_grid,
real(kind=default_precision), dimension(:,:,:), intent(inout)  halo_buffer,
real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data,
integer, intent(in)  corner_loc,
integer, intent(in)  x_source_index,
integer, intent(in)  y_source_index,
integer, intent(in)  halo_page 
)

Copies prognostic field corner data to send buffer for specific field.

Parameters
local_gridDescription of the local grid
halo_bufferRaw halo_buffer data that is written to
field_dataRaw prognostic field data that is read from
corner_locLocation of the corner
x_source_indexThe X index in the read dimension
y_source_indexThe Y index in the read dimension
halo_pageThe halo buffer page to write to

Definition at line 461 of file halocommunication.F90.

461  type(local_grid_type), intent(inout) :: local_grid
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
465 
466  !! TODO: hardcoded size of the corners
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. &
474  corner_loc == 3))
Here is the caller graph for this function:

◆ copy_field_to_buffer()

subroutine, public halo_communication_mod::copy_field_to_buffer ( type(local_grid_type), intent(inout)  local_grid,
real(kind=default_precision), dimension(:,:,:), intent(inout)  halo_buffer,
real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data,
integer, intent(in)  dim,
integer, intent(in)  source_index,
integer, intent(in)  halo_page 
)

Copies prognostic field data to send buffer for specific field, dimension, halo cell.

Parameters
local_gridDescription of the local grid
halo_bufferRaw halo_buffer data that is written to
field_dataRaw prognostic field data that is read from
dimDimension that we are reading from
source_indexThe index in the read dimension
halo_pageThe halo buffer page to write to

Definition at line 433 of file halocommunication.F90.

433  type(local_grid_type), intent(inout) :: local_grid
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
437 
438  if (dim == x_index) then
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)
443  else
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))
448  end if
Here is the caller graph for this function:

◆ deduce_halo_corners_per_neighbour()

subroutine halo_communication_mod::deduce_halo_corners_per_neighbour ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), dimension(:), allocatable  halo_swap_neighbours,
integer, intent(in)  number_distinct_neighbours,
integer, intent(in)  fields_per_cell 
)
private

Determines the number of halo corners to swap between specific neighours, this is similar to deducing the number of halo pages per neighbour, only it is for corners.

Parameters
current_stateThe current model state
halo_swap_neighboursDescription of neighbouring halo swap state
number_distinct_neighboursThe number of distinct neighbours that I swap with
fields_per_cellThe number of fields per cell is written into here

Definition at line 732 of file halocommunication.F90.

732  type(model_state_type), intent(inout) :: current_state
733  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neighbours
734  integer, intent(in) :: number_distinct_neighbours, fields_per_cell
735 
736  integer :: i, j, pid_location
737  ! i moves in x ,j in y
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
742  pid_location = get_pid_neighbour_location(halo_swap_neighbours, &
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
746  end if
747  end do
748  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ deduce_halo_pages_per_neighbour()

subroutine halo_communication_mod::deduce_halo_pages_per_neighbour ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), dimension(:), allocatable  halo_swap_neighbours,
integer, intent(in)  number_distinct_neighbours,
procedure(get_fields_per_halo_cell_proc_interface get_fields_per_halo_cell,
integer, intent(out)  fields_per_cell,
integer, intent(in)  halo_depth 
)
private

Deduces the number of halo pages per neighbour halo swap and places this information in the appropriate data structures. We call a "page" of data the contiguous data of a field that we are going to send, such as halo 1 of w, halo 1 of zw and halo 2 of w (assuming these go to the same neighbour as 1)

Parameters
current_stateThe current model state
halo_swap_neighboursMy neighbouring PIDs
number_distinct_neighboursThe number of distinct neighbours that I have
get_fields_per_halo_cellProcedure pointer to get the number of fields per halo cell
fields_per_cellThe number of fields per cell is written into here

Definition at line 698 of file halocommunication.F90.

698  type(model_state_type), intent(inout) :: current_state
699  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neighbours
700  integer, intent(in) :: number_distinct_neighbours, halo_depth
701  procedure(get_fields_per_halo_cell_proc_interface) :: get_fields_per_halo_cell
702  integer, intent(out) :: fields_per_cell
703 
704  integer :: i, j, pid_location, halo_start, halo_end
705 
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)
709 
710  ! i moves in x and y. z is 1
711  do i = 2, 3
712  do j = halo_start, halo_end
713  if (current_state%parallel%my_rank .ne. current_state%local_grid%neighbours(i,j)) then
714  pid_location = get_pid_neighbour_location(halo_swap_neighbours, &
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
718  end if
719  end do
720  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ determine_halo_corner_element_sizes()

integer function halo_communication_mod::determine_halo_corner_element_sizes ( type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  pid 
)
private

For a specific process id this determines the number of halo swap corner elements to involve in a communication.

Parameters
local_gridThe local grid
pidThe process id that this calculation is based upon The number of elements involved in a specific corner halo swap

Definition at line 576 of file halocommunication.F90.

576  type(local_grid_type), intent(inout) :: local_grid
577  integer, intent(in) :: pid
578 
579  integer :: i,j
580  determine_halo_corner_element_sizes=0
581 
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
585  determine_halo_corner_element_sizes=&
586  determine_halo_corner_element_sizes+local_grid%size(z_index)
587  ! If second halo then there are 3 corner elements, therefore add extra two to this
588  if (i==2) determine_halo_corner_element_sizes=determine_halo_corner_element_sizes+&
589  local_grid%size(z_index)*2
590  end if
591  end do
592  end do
Here is the caller graph for this function:

◆ determine_halo_corner_size()

integer function halo_communication_mod::determine_halo_corner_size ( type(local_grid_type), intent(inout)  local_grid)
private

Determine the halo corner size in elements.

Parameters
local_gridThe local grid
pidThe process id that this calculation is based upon The number of elements involved in a specific corner halo swap

Definition at line 564 of file halocommunication.F90.

564  type(local_grid_type), intent(inout) :: local_grid
565 
566  determine_halo_corner_size = local_grid%halo_size(x_index)*local_grid%halo_size(y_index)*&
567  local_grid%size(z_index)
Here is the caller graph for this function:

◆ determine_recv_and_send_sizes()

subroutine halo_communication_mod::determine_recv_and_send_sizes ( type(local_grid_type), intent(inout)  local_grid,
type(neighbour_description_type), dimension(:), allocatable  halo_swap_neighbours,
integer, intent(in)  number_distinct_neighbours,
logical, intent(in)  involve_corners 
)
private

Determines the amount (in elements) of data that each neighbour will be sent and I will receive from in a halo swap.

Parameters
local_gridThe local grid
halo_swap_neighboursStructure describing state of halo swap neighbours
number_distinct_neighboursThe number of distinct neighbours I swap with
involve_cornersWhether or not to involve corners in a halo swap

Definition at line 529 of file halocommunication.F90.

529  type(local_grid_type), intent(inout) :: local_grid
530  integer, intent(in) :: number_distinct_neighbours
531  logical, intent(in) :: involve_corners
532  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neighbours
533 
534  integer :: i, normal_size, corner_size
535 
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")
540  end if
541  normal_size=local_grid%size(z_index) * merge(local_grid%size(y_index), local_grid%size(x_index), &
542  halo_swap_neighbours(i)%dimension==x_index)*halo_swap_neighbours(i)%halo_pages
543  else
544  normal_size=0
545  end if
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
549  ! For the moment assume both halos are the same neighbour - hence the 4, otherwise should call determine_halo_corner_element_sizes
550  corner_size=local_grid%size(z_index)*4*halo_swap_neighbours(i)%halo_corners
551  else
552  corner_size=0
553  end if
554  halo_swap_neighbours(i)%recv_corner_size=corner_size
555  halo_swap_neighbours(i)%send_corner_size=corner_size
556  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalise_halo_communication()

subroutine, public halo_communication_mod::finalise_halo_communication ( type(halo_communication_type), intent(inout)  halo_swap_state)

Finalises the halo swap represented by the state by freeing up all the allocated memory.

Parameters
halo_swap_stateState of the specific halo swap

Definition at line 340 of file halocommunication.F90.

340  type(halo_communication_type), intent(inout) :: halo_swap_state
341 
342  integer :: i
343 
344  ! TODO - issue cancel of outstanding requests here
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)
347 
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)
353  end do
354  if (allocated(halo_swap_state%halo_swap_neighbours)) deallocate(halo_swap_state%halo_swap_neighbours)
355  halo_swap_state%initialised=.false.
Here is the caller graph for this function:

◆ generate_recv_field_buffer_matches()

subroutine halo_communication_mod::generate_recv_field_buffer_matches ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  halo_depth,
integer, dimension(:,:), intent(out)  cell_match 
)
private

Precalculates the received buffer to field halo cell matches for each dimension and called from the initialisation stage.

Parameters
current_stateThe current model state
halo_depthThe halo depth
cell_matchThe matching cells are written into here

Definition at line 790 of file halocommunication.F90.

790  type(model_state_type), intent(inout) :: current_state
791  integer, intent(in) :: halo_depth
792  integer, intent(out) :: cell_match(:,:)
793 
794  logical, dimension(3) :: same_neighbours
795  integer :: i,j
796 
797  same_neighbours = retrieve_same_neighbour_information(current_state%local_grid)
798 
799  do i = 2,3
800  if (halo_depth == 1) then
801  cell_match(i,1)=1
802  cell_match(i,2)=merge(2, 3, .not. same_neighbours(i))
803  cell_match(i,3)=merge(3, 2, .not. same_neighbours(i))
804  cell_match(i,4)=4
805  else
806  do j = 1, halo_depth
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))
809  end do
810  ! cell_match(i,j) = merge(2, 4, .not. same_neighbours(i))
811  ! cell_match(i,j) = merge(3, 1, .not. same_neighbours(i))
812 ! cell_match(i,j) = merge(4, 2, .not. same_neighbours(i))
813  end if
814  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_number_communication_requests()

integer function halo_communication_mod::get_number_communication_requests ( type(neighbour_description_type), dimension(:), allocatable  halo_swap_neighbours,
integer, intent(in)  number_distinct_neighbours 
)
private

Determines the overall number of communication requests, which is made up of normal halo swaps and potentially corner swaps too if that is enabled.

Parameters
halo_swap_neighboursNeighbouring halo swap state
number_distinct_neighboursThe number of distinct neighours that I will swap with
Returns
The number of communication requests

Definition at line 507 of file halocommunication.F90.

507  integer, intent(in) :: number_distinct_neighbours
508  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neighbours
509 
510  integer :: i
511 
512  get_number_communication_requests=0
513  do i=1, number_distinct_neighbours
514  if (halo_swap_neighbours(i)%recv_size .gt. 0) &
515  get_number_communication_requests=get_number_communication_requests+1
516  if (halo_swap_neighbours(i)%recv_corner_size .gt. 0)&
517  get_number_communication_requests=get_number_communication_requests+1
518  end do
Here is the caller graph for this function:

◆ get_number_of_processes_involved_in_communication()

integer function halo_communication_mod::get_number_of_processes_involved_in_communication ( type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  my_rank,
logical, intent(in)  include_corners 
)
private

Deduces the number of distinct neighbours that will be involved in a halo swap. This information is used to then allocate the appropriate amount of memory to store the neighbour halo swapping data structure.

Parameters
local_gridDescription of the local grid
my_rankMy global PID
include_cornersWhether to include corners or not

Definition at line 603 of file halocommunication.F90.

603  type(local_grid_type), intent(inout) :: local_grid
604  integer, intent(in) :: my_rank
605  logical, intent(in) :: include_corners
606 
607  integer :: i, j, temp_neighbour_pids(merge(16, 8, include_corners))
608 
609  temp_neighbour_pids(:)=-1
610  get_number_of_processes_involved_in_communication=0
611  do i=2,3
612  do j=1,4
613  if (local_grid%neighbours(i,j) .ne. my_rank .and. .not. &
614  has_pid_already_been_seen(temp_neighbour_pids, &
615  local_grid%neighbours(i,j))) then
616  get_number_of_processes_involved_in_communication = &
617  get_number_of_processes_involved_in_communication+1
618  temp_neighbour_pids(get_number_of_processes_involved_in_communication) = &
619  local_grid%neighbours(i,j)
620  end if
621  end do
622  end do
623 
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. &
628  has_pid_already_been_seen(temp_neighbour_pids, &
629  local_grid%corner_neighbours(j,i))) then
630  get_number_of_processes_involved_in_communication = &
631  get_number_of_processes_involved_in_communication+1
632  temp_neighbour_pids(get_number_of_processes_involved_in_communication) = &
633  local_grid%corner_neighbours(j,i)
634  end if
635  end do
636  end do
637  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_pid_neighbour_location()

integer function halo_communication_mod::get_pid_neighbour_location ( type(neighbour_description_type), dimension(:), allocatable  halo_swap_neighbours,
integer, intent(in)  pid,
integer, intent(in)  number_distinct_neighbours 
)
private

Given the process id of a neighbour this determines the location in the data structure of corresponding data for that. Note that this is an O(n) operation so ideally call once and reuse the results. If no process id is found then returns -1.

Parameters
halo_swap_neighboursPIDs of the neighbouring processes to me
pidThe process id to find
number_distinct_neighboursThe number of distinct neighbours I have

Definition at line 1211 of file halocommunication.F90.

1211  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neighbours
1212  integer, intent(in) :: pid, number_distinct_neighbours
1213 
1214  integer :: i
1215 
1216  do i=1, number_distinct_neighbours
1217  if (halo_swap_neighbours(i)%pid == pid) then
1218  get_pid_neighbour_location = i
1219  return
1220  end if
1221  end do
1222  ! Not found
1223  get_pid_neighbour_location=-1
Here is the caller graph for this function:

◆ get_single_field_per_halo_cell()

integer function, public halo_communication_mod::get_single_field_per_halo_cell ( type(model_state_type), intent(inout)  current_state)

A very common function, which returns a single field per halo cell which is used to halo swap just one field.

Parameters
current_stateThe current model state

Definition at line 1230 of file halocommunication.F90.

1230  type(model_state_type), intent(inout) :: current_state
1231 
1232  get_single_field_per_halo_cell=1
Here is the caller graph for this function:

◆ has_pid_already_been_seen()

logical function halo_communication_mod::has_pid_already_been_seen ( integer, dimension(8), intent(in)  temp_neighbour_pids,
integer, intent(in)  pid 
)
private

Returns whether or not a specific process id has already been "seen" by searching a list of already seen process ids.

Parameters
temp_neighbour_pidsThe process Ids already seen
pidThe PID to search for

Definition at line 1188 of file halocommunication.F90.

1188  integer, intent(in) :: pid, temp_neighbour_pids(8)
1189 
1190  integer :: i
1191 
1192  has_pid_already_been_seen=.true.
1193  do i=1,8
1194  if (temp_neighbour_pids(i) == pid) return
1195  if (temp_neighbour_pids(i) == -1) then
1196  has_pid_already_been_seen=.false.
1197  return
1198  end if
1199  end do
1200  has_pid_already_been_seen=.false.
Here is the caller graph for this function:

◆ init_halo_communication()

subroutine, public halo_communication_mod::init_halo_communication ( type(model_state_type), intent(inout)  current_state,
procedure(get_fields_per_halo_cell_proc_interface get_fields_per_halo_cell,
type(halo_communication_type), intent(out)  halo_state,
integer, intent(in)  halo_depth,
logical, intent(in)  involve_corners 
)

Initialises a halo swapping state, by determining the neighbours, size of data in each swap and allocating the required memory - which are communication buffers and request handles. All this information is placed into the returned state.

Parameters
current_stateThe current model state
get_fields_per_halo_cellProcedure pointer to the procedure which determines how many fields per halo cell there are
involve_cornersWhether this involves corners or not
Returns
A halo swapping state which can be used for halo swapping operations

Definition at line 295 of file halocommunication.F90.

295  type(model_state_type), intent(inout) :: current_state
296  procedure(get_fields_per_halo_cell_proc_interface) :: get_fields_per_halo_cell
297  logical, intent(in) :: involve_corners
298  integer, intent(in) :: halo_depth
299  type(halo_communication_type), intent(out) :: halo_state
300 
301  integer :: number_comm_requests
302 
303  halo_state%involve_corners = involve_corners
304  halo_state%halo_depth = halo_depth
305  halo_state%number_distinct_neighbours = get_number_of_processes_involved_in_communication(&
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))
309  halo_state%halo_swap_neighbours = populate_halo_swap_neighbours(current_state%local_grid, &
310  current_state%parallel%my_rank, halo_state%number_distinct_neighbours, involve_corners)
311  call deduce_halo_pages_per_neighbour(current_state, halo_state%halo_swap_neighbours, &
312  halo_state%number_distinct_neighbours, get_fields_per_halo_cell, &
313  halo_state%fields_per_cell, halo_depth)
314  if (involve_corners) call deduce_halo_corners_per_neighbour(current_state, &
315  halo_state%halo_swap_neighbours, &
316  halo_state%number_distinct_neighbours, halo_state%fields_per_cell)
317  call allocate_halo_buffers_for_each_neighbour(current_state%local_grid, &
318  halo_state%number_distinct_neighbours, &
319  halo_state%halo_swap_neighbours)
320  call determine_recv_and_send_sizes(current_state%local_grid, &
321  halo_state%halo_swap_neighbours, &
322  halo_state%number_distinct_neighbours, involve_corners)
323  call generate_recv_field_buffer_matches(current_state, halo_state%halo_depth, &
324  halo_state%cell_match)
325 
326  ! required for nonblocking MPI communcations
327  number_comm_requests = get_number_communication_requests(halo_state%halo_swap_neighbours, &
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
333  end if
334  halo_state%initialised = .true.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initiate_nonblocking_halo_swap()

subroutine, public halo_communication_mod::initiate_nonblocking_halo_swap ( type(model_state_type), intent(inout)  current_state,
type(halo_communication_type), intent(inout)  halo_swap_state,
procedure(copy_fields_to_halo_buffer_proc_interface copy_to_halo_buffer,
procedure(copy_corners_to_halo_buffer_proc_interface), optional  copy_corners_to_halo_buffer,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)

Initiates a non blocking halo swap, this registers the receive requests, copies data into the send buffer and then registers send requests for these. As far as this call is concerned, the data is immutable - no modifications will take place to it until the matching completion call is made.

Parameters
current_stateThe current model state
halo_swap_stateThe halo swapping state
copy_to_halo_bufferProcedure pointer which is called to copy the data into the halo send buffer
copy_corners_to_halo_bufferOptional procedure pointer which copies field corners into halo buffer
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 245 of file halocommunication.F90.

245 
246  type(model_state_type), intent(inout) :: current_state
247  type(halo_communication_type), intent(inout) :: halo_swap_state
248  procedure(copy_fields_to_halo_buffer_proc_interface) :: copy_to_halo_buffer
249  procedure(copy_corners_to_halo_buffer_proc_interface), optional :: copy_corners_to_halo_buffer
250  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
251 
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
256 
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")
261  end if
262 
263  ! we call recv before send
264  call recv_all_halos(current_state, halo_swap_state)
265 
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)
270  else
271  call send_all_halos(current_state, halo_swap_state, copy_to_halo_buffer, &
272  source_data=source_data)
273  end if
274  else
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)
278  else
279  call send_all_halos(current_state, halo_swap_state, copy_to_halo_buffer)
280  end if
281  end if
282  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ perform_local_data_copy_for_corners()

subroutine halo_communication_mod::perform_local_data_copy_for_corners ( integer, intent(in)  my_rank,
type(local_grid_type), intent(inout)  local_grid,
real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data 
)
private

Performs a local data copy for corners when the neighbour is local (me)

Parameters
my_rankMy global rank
local_gridThe local grid
field_dataThe field data top copy into

Definition at line 1042 of file halocommunication.F90.

1042  type(local_grid_type), intent(inout) :: local_grid
1043  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout) :: field_data
1044  integer, intent(in) :: my_rank
1045 
1046  integer :: i, y_source_index, x_source_index, y_target_index, x_target_index
1047 
1048  do i=1,size(local_grid%corner_neighbours, 1)
1049  if (my_rank .eq. local_grid%corner_neighbours(i,1)) then
1050  if (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
1053  y_target_index=1
1054  x_target_index=1
1055  else if (i==2) then
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)
1058  y_target_index=1
1059  x_target_index=local_grid%local_domain_end_index(x_index)+1
1060  else if (i==3) then
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
1064  x_target_index=1
1065  else if (i==4) then
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
1070  end if
1071 
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)
1077 
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)
1083 
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)
1089 
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)
1095  end if
1096  end do
Here is the caller graph for this function:

◆ perform_local_data_copy_for_dimension()

subroutine halo_communication_mod::perform_local_data_copy_for_dimension ( integer, intent(in)  dim,
integer, intent(in)  my_rank,
integer, intent(in)  halo_depth,
type(local_grid_type), intent(inout)  local_grid,
real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data 
)
private

Performs a local data copy for a specific dimension of a prognostic field.

Parameters
dimThe dimension
my_rankMy process Id
local_gridThe local grid description
field_dataThe field data that we are going to copy from and to

Definition at line 1106 of file halocommunication.F90.

1106  type(local_grid_type), intent(inout) :: local_grid
1107  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout) :: field_data
1108  integer, intent(in) :: dim, my_rank, halo_depth
1109 
1110  integer i, target_index, source_index, hstart, hend
1111 
1112  hstart=merge(2,1, halo_depth==1)
1113  hend=merge(3,4, halo_depth==1)
1114 
1115  do i=hstart, hend
1116  if (local_grid%neighbours(dim,i) .eq. my_rank) then
1117  if (i==1) then
1118  target_index=1
1119  source_index=local_grid%local_domain_end_index(dim)-1
1120  else if (i==2) then
1121  target_index=2
1122  source_index=local_grid%local_domain_end_index(dim)
1123  else if (i==3) then
1124  target_index=local_grid%local_domain_end_index(dim)+1
1125  source_index=local_grid%local_domain_start_index(dim)
1126  else if (i==4) then
1127  target_index=local_grid%local_domain_end_index(dim)+2
1128  source_index=local_grid%local_domain_start_index(dim)+1
1129  end if
1130  if (dim == x_index) then
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)
1139  else
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))
1148  end if
1149  end if
1150  end do
Here is the caller graph for this function:

◆ perform_local_data_copy_for_field()

subroutine, public halo_communication_mod::perform_local_data_copy_for_field ( real(kind=default_precision), dimension(:,:,:), intent(inout)  field_data,
type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  my_rank,
integer, intent(in)  halo_depth,
logical, intent(in)  involve_corners 
)

Will perform a a local copy for the halo data of a field.

Parameters
field_dataThe field data to perform the copy for
local_gridLocal grid information
my_rankMy process id

Definition at line 483 of file halocommunication.F90.

483  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout) :: field_data
484  type(local_grid_type), intent(inout) :: local_grid
485  integer, intent(in) :: my_rank, halo_depth
486  logical, intent(in) :: involve_corners
487 
488  call perform_local_data_copy_for_dimension(y_index, my_rank, halo_depth, local_grid, &
489  field_data)
490  call perform_local_data_copy_for_dimension(x_index, my_rank, halo_depth, local_grid, &
491  field_data)
492  if (involve_corners) call perform_local_data_copy_for_corners(my_rank, local_grid, field_data)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ populate_halo_swap_neighbours()

type(neighbour_description_type) function, dimension(number_distinct_neighbours) halo_communication_mod::populate_halo_swap_neighbours ( type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  my_rank,
integer, intent(in)  number_distinct_neighbours,
logical, intent(in)  involve_corners 
)
private

Will populate the halo swap neighbour data strutures with appropriate neighbour pid and dimension numbers.

Parameters
local_gridDescription of the local grid
my_rankMy global PID
number_distinct_neighboursThe number of distinct neighbours that I have
include_cornersWhether to include corners or not

Definition at line 648 of file halocommunication.F90.

648  type(local_grid_type), intent(inout) :: local_grid
649  integer, intent(in) :: my_rank, number_distinct_neighbours
650  logical, intent(in) :: involve_corners
651 
652  type(neighbour_description_type), dimension(number_distinct_neighbours) :: &
653  populate_halo_swap_neighbours
654  integer :: i, j, current_pid_location, temp_neighbour_pids(merge(16, 8, involve_corners))
655 
656  current_pid_location=0
657  temp_neighbour_pids(:)=-1
658  do i=2,3
659  do j=1,4
660  if (local_grid%neighbours(i,j) .ne. my_rank .and. .not. &
661  has_pid_already_been_seen(temp_neighbour_pids, &
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
667  end if
668  end do
669  end do
670 
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. &
675  has_pid_already_been_seen(temp_neighbour_pids, &
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
682  end if
683  end do
684  end do
685  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ recv_all_halos()

subroutine halo_communication_mod::recv_all_halos ( type(model_state_type), intent(inout)  current_state,
type(halo_communication_type), intent(inout)  halo_swap_state 
)
private

Registers receive requests for all prognostic fields from the appropriate neighbouring processes (that we have already deduced in the initialisation stage.)

Parameters
current_stateThe current model state
halo_swap_stateThe halo swap state

Definition at line 822 of file halocommunication.F90.

822  type(model_state_type), intent(inout) :: current_state
823  type(halo_communication_type), intent(inout) :: halo_swap_state
824 
825  integer :: i, request_counter, ierr
826 
827  request_counter = 1
828 
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
837  end if
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
845  end if
846  end do
Here is the caller graph for this function:

◆ retrieve_same_neighbour_information()

logical function, dimension(3) halo_communication_mod::retrieve_same_neighbour_information ( type(local_grid_type), intent(inout)  local_grid)
private

Retrieves whether we have the same neighbours for L and R halo swaps in each dimension.

This is because if we do then the data storage needs to be reversed - i.e. the first two messages sent to this process are Left and the second two are Right, but due to wrapping around, the L messages need to be put on the right halo and R messages on the left halo

Parameters
neighboursThe neighbouring processes in each dimension

Definition at line 1160 of file halocommunication.F90.

1160  type(local_grid_type), intent(inout) :: local_grid
1161 ! integer, dimension(:,:), intent(in) :: neighbours
1162  logical, dimension(3) :: retrieve_same_neighbour_information
1163 
1164  integer :: i, nd1, nd2
1165 
1166  retrieve_same_neighbour_information=(/.true., .true., .true./)
1167 
1168  ! halo_size in X and Y are the same, therefore it does not matter which one we take
1169  ! we multiply by 2 since there are 2 sides Up&Down or Left&Right
1170  do i = 1,local_grid%halo_size(y_index)*2
1171  if (i==1) then
1172  nd1=local_grid%neighbours(y_index,i)
1173  nd2=local_grid%neighbours(x_index,i)
1174  else
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.
1179  end if
1180  end do
Here is the caller graph for this function:

◆ send_all_halos()

subroutine halo_communication_mod::send_all_halos ( type(model_state_type), intent(inout)  current_state,
type(halo_communication_type), intent(inout)  halo_swap_state,
procedure(copy_fields_to_halo_buffer_proc_interface copy_fields_to_halo_buffer,
procedure(copy_corners_to_halo_buffer_proc_interface), optional  copy_corner_fields_to_halo_buffer,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies all applicable bits of the prognostics into a send buffer for each neighbour and then issues asynchronous sends of this data.

Parameters
current_stateThe current model state
halo_swap_stateThe halo swap state
copy_fields_to_halo_bufferProcedure pointer to copy the data from the fields into the send buffer
copy_corner_fields_to_halo_bufferOptional procedure pointer to copy corner fields into halo send buffer
source_dataOptional field data that can be copied from (passed to the user procedure)

Definition at line 860 of file halocommunication.F90.

860 
861  type(model_state_type), intent(inout) :: current_state
862  type(halo_communication_type), intent(inout) :: halo_swap_state
863  procedure(copy_fields_to_halo_buffer_proc_interface) :: copy_fields_to_halo_buffer
864  procedure(copy_corners_to_halo_buffer_proc_interface), optional :: &
865  copy_corner_fields_to_halo_buffer
866  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
867 
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),&
870  halo_depth
871  ! halo_size(Y_INDEX) = halo_size(X_INDEX), just pick one
872  halo_depth = current_state%local_grid%halo_size(y_index)
873  current_page(:) = 1
874  request_number = 1
875  ! TODO: hardcoded to halodepth 1 or 2
876  hstart = merge(2, 1, halo_swap_state%halo_depth==1)
877  hend = merge(3, halo_depth*2, halo_swap_state%halo_depth==1)
878 
879  do i=2, 3
880  do j=hstart, hend
881  if (current_state%parallel%my_rank .ne. current_state%local_grid%neighbours(i,j)) then
882 
883  if (j==1) then
884  source_index = current_state%local_grid%local_domain_start_index(i)
885  else if (j==2) then
886  source_index = current_state%local_grid%local_domain_start_index(i) + &
887  merge(1, 0, halo_swap_state%halo_depth .ne. 1)
888  else if (j==3) then
889  source_index = current_state%local_grid%local_domain_end_index(i) - &
890  merge(1, 0, halo_swap_state%halo_depth .ne. 1)
891  else if (j==4) then
892  source_index = current_state%local_grid%local_domain_end_index(i)
893  end if
894 
895  pid_location = get_pid_neighbour_location(halo_swap_state%halo_swap_neighbours, &
896  current_state%local_grid%neighbours(i,j), &
897  halo_swap_state%number_distinct_neighbours)
898 
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)
903  else
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)
907  end if
908  ! call log_log(LOG_DEBUG, "PID ="//trim(conv_to_String(&
909  !current_state%parallel%my_rank))//" source_index = "//&
910  ! trim(conv_to_string(source_index))//" PID location ="//trim(&
911  !conv_to_string(pid_location))//&
912  ! " i= "//trim(conv_to_string(i))//" j="//trim(conv_to_string(j)))
913  end if
914  end do
915  end do
916 
917  if (present(copy_corner_fields_to_halo_buffer)) then
918  current_page(:)=1
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)
926  pid_location = get_pid_neighbour_location(halo_swap_state%halo_swap_neighbours, &
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)
933  else
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)
937  end if
938  end if
939  end do
940  end if
941 
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
950  end if
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
958  end if
959  end do
Here is the call graph for this function:
Here is the caller graph for this function: