Go to the documentation of this file.
14 use fruit,
only : assert_equals, assert_true, assert_false, assert_not_equals
24 integer :: i, requests
25 allocate(halo_swap_neigh(6))
28 halo_swap_neigh(i)%recv_size = i
29 halo_swap_neigh(i)%recv_corner_size = i
33 call assert_equals(10, requests,
"Test number of requests")
39 integer :: halo_corner_size
47 call assert_equals(10*2*2, halo_corner_size,
"Test halo_size")
54 allocate( local_grid%corner_neighbours(4,2) )
55 local_grid%corner_neighbours(1,1) = 5
56 local_grid%corner_neighbours(2,2) = 5
61 call assert_equals(10*4, elemnumber,
"Test number elements halo_size")
64 call assert_equals(0, elemnumber,
"Test number elements halo_size")
73 allocate(halo_swap_neigh(6))
75 halo_swap_neigh(i)%pid = i
77 halo_swap_neigh(2)%pid=1
80 call assert_equals(1, pids,
"Test unique pid location")
82 call assert_equals(-1, pids,
"Test pid location")
87 integer :: neigh_pids(8)
94 call assert_false( seen,
"Test not seen if pid=-1" )
97 call assert_true( seen,
"Test seen when there is not -1 if pid=2" )
99 call assert_true( seen,
"Test seen when there is not - 1 if pid=0" )
106 logical,
dimension(3) :: retrieve_same_neigh_info
108 allocate( local_grid%neighbours(3,4) )
109 local_grid%halo_size(
y_index) = 2
110 local_grid%halo_size(
x_index) = 2
112 local_grid%neighbours(
y_index,1) = 2
113 local_grid%neighbours(
y_index,2) = 2
114 local_grid%neighbours(
y_index,3) = 8
115 local_grid%neighbours(
y_index,4) = 8
117 local_grid%neighbours(
x_index,1) = 4
118 local_grid%neighbours(
x_index,2) = 4
119 local_grid%neighbours(
x_index,3) = 6
120 local_grid%neighbours(
x_index,4) = 6
124 call assert_true(retrieve_same_neigh_info(1),
"Test Z")
125 call assert_false(retrieve_same_neigh_info(2),
"Test Y")
126 call assert_false(retrieve_same_neigh_info(3),
"Test Z")
129 local_grid%neighbours(
y_index,3) = 2
130 local_grid%neighbours(
y_index,4) = 2
134 call assert_true(retrieve_same_neigh_info(1),
"Test Z")
135 call assert_true(retrieve_same_neigh_info(2),
"Test Y")
136 call assert_false(retrieve_same_neigh_info(3),
"Test Z")
146 allocate( local_grid%neighbours(3,4) )
148 local_grid%local_domain_start_index(
z_index) = 1
149 local_grid%local_domain_start_index(
y_index) = 11
150 local_grid%local_domain_start_index(
x_index) = 21
152 local_grid%local_domain_end_index(
z_index) = 10
153 local_grid%local_domain_end_index(
y_index) = 20
154 local_grid%local_domain_end_index(
x_index) = 30
158 local_grid%neighbours(j,i) = i*j
164 field_data(i,j,k) = i*j*k
169 call assert_not_equals(field_data(1,2,21), field_data(1,10,21)&
170 ,
"Test fields are not equal before calling")
179 use fruit,
only : init_fruit, run_test_case, fruit_summary
Conversion between common inbuilt FORTRAN data types.
Map data structure that holds string (length 20 maximum) key value pairs.
Describes the neighbours of a process in a specific dimension and contains the communication buffers ...
subroutine test_retrieve_same_neighbour_information
Provides the mechanism for halo swapping. This module contains the functionality required to determin...
subroutine test_get_pid_neighbour_location
subroutine test_determine_halo_corner
Converts data types to integers.
Collection data structures.
subroutine test_determine_halo_corner_elements
integer, parameter, public x_index
integer, parameter, public y_index
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...
Converts data types to logical.
Contains the types used for communication, holding the state of communications and supporting activit...
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, 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 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.
Maintains the state of a halo swap and contains buffers, neighbours etc.
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...
program test_halo_communication_driver
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.
subroutine test_perform_local_data_copy_for_dimension
integer, parameter, public z_index
Grid index parameters.
Converts data types to strings.
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 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...
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.
Defined the local grid, i.e. the grid held on this process after decomposition.
subroutine test_pid_been_seen
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...
Contains common definitions for the data and datatypes used by MONC.
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...
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.
Converts data types to real.
subroutine test_get_number_communication_requests
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...
Functionality to support the different types of grid and abstraction between global grids and local o...
subroutine, public finalise_halo_communication(halo_swap_state)
Finalises the halo swap represented by the state by freeing up all the allocated memory.
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...
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.
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.