MONC
test_halocommunication.F90
Go to the documentation of this file.
1 ! Unit tests for the options database functionality
9  use collections_mod, only : map_type, c_get
14  use fruit, only : assert_equals, assert_true, assert_false, assert_not_equals
16  implicit none
17 
18  contains
19 
20 
21  ! Test the number of requests
23  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neigh
24  integer :: i, requests
25  allocate(halo_swap_neigh(6))
26 
27  do i=1,5
28  halo_swap_neigh(i)%recv_size = i
29  halo_swap_neigh(i)%recv_corner_size = i
30  enddo
31 
32  requests = get_number_communication_requests(halo_swap_neigh, 5)
33  call assert_equals(10, requests, "Test number of requests")
35 
36  ! Test the size of halo corner
38  type(local_grid_type) :: local_grid
39  integer :: halo_corner_size
40 
41  local_grid%halo_size(x_index)=2
42  local_grid%halo_size(y_index)=2
43  local_grid%size(z_index) = 10
44 
45  halo_corner_size = determine_halo_corner_size(local_grid)
46 
47  call assert_equals(10*2*2, halo_corner_size, "Test halo_size")
48  end subroutine test_determine_halo_corner
49 
50  ! Test the size of halo corner elements
52  type(local_grid_type) :: local_grid
53  integer :: elemNumber
54  allocate( local_grid%corner_neighbours(4,2) )
55  local_grid%corner_neighbours(1,1) = 5
56  local_grid%corner_neighbours(2,2) = 5
57  local_grid%size(z_index) = 10
58  elemnumber = determine_halo_corner_element_sizes(local_grid, 5)
59 
60  !4 elements per corner*10 in Z direction
61  call assert_equals(10*4, elemnumber, "Test number elements halo_size")
62  ! if not PID elemNumber should be 0
63  elemnumber = determine_halo_corner_element_sizes(local_grid, 15)
64  call assert_equals(0, elemnumber, "Test number elements halo_size")
66 
67  ! Test get neighbour pid
69  type(local_grid_type) :: local_grid
70  type(neighbour_description_type), dimension(:), allocatable :: halo_swap_neigh
71  integer :: i, pids
72 
73  allocate(halo_swap_neigh(6))
74  do i=1,5
75  halo_swap_neigh(i)%pid = i
76  enddo
77  halo_swap_neigh(2)%pid=1
78  pids = get_pid_neighbour_location(halo_swap_neigh, 1, 5)
79 
80  call assert_equals(1, pids, "Test unique pid location")
81  pids = get_pid_neighbour_location(halo_swap_neigh, 10, 5)
82  call assert_equals(-1, pids, "Test pid location")
83  end subroutine test_get_pid_neighbour_location
84 
85  ! Test has_pid_already_been_seen which actually tells if a list contains a pid
86  subroutine test_pid_been_seen
87  integer :: neigh_pids(8)
88  logical :: seen
89  neigh_pids(1) = -1
90  neigh_pids(2) = 2
91  neigh_pids(3) = 0
92 
93  seen = has_pid_already_been_seen(neigh_pids, 1)
94  call assert_false( seen, "Test not seen if pid=-1" )
95  neigh_pids(1) = 0
96  seen = has_pid_already_been_seen(neigh_pids, 2)
97  call assert_true( seen, "Test seen when there is not -1 if pid=2" )
98  seen = has_pid_already_been_seen(neigh_pids, 0)
99  call assert_true( seen, "Test seen when there is not - 1 if pid=0" )
100 
101  end subroutine test_pid_been_seen
102 
103  ! Test if we have same neighbourd L,R,U,D
105  type(local_grid_type) :: local_grid
106  logical, dimension(3) :: retrieve_same_neigh_info
107  integer :: i
108  allocate( local_grid%neighbours(3,4) )
109  local_grid%halo_size(y_index) = 2
110  local_grid%halo_size(x_index) = 2
111 
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
116 
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
121 
122  retrieve_same_neigh_info = retrieve_same_neighbour_information(local_grid)
123 
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")
127 
128 
129  local_grid%neighbours(y_index,3) = 2
130  local_grid%neighbours(y_index,4) = 2
131 
132  retrieve_same_neigh_info = retrieve_same_neighbour_information(local_grid)
133 
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")
137 
139 
140  ! Test a local data copy for a specific dimension
142  type(local_grid_type) :: local_grid
143  real(kind=DEFAULT_PRECISION), dimension(10,20,30) :: field_data
144  integer :: i,j,k
145 
146  allocate( local_grid%neighbours(3,4) )
147  ! start
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
151  ! end
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
155 
156  do j=1,3
157  do i =1,4
158  local_grid%neighbours(j,i) = i*j
159  enddo
160  enddo
161  do i=1,10
162  do j=1,20
163  do k=1,30
164  field_data(i,j,k) = i*j*k
165  enddo
166  enddo
167  enddo
168 
169  call assert_not_equals(field_data(1,2,21), field_data(1,10,21)&
170  , "Test fields are not equal before calling")
171  call perform_local_data_copy_for_dimension(x_index, 1, 1, local_grid, &
172  field_data)
173 
175 
176  end module test_halo_communication_mod
177 
179  use fruit, only : init_fruit, run_test_case, fruit_summary
184  implicit none
185 
186  call init_fruit
187  call run_test_case(test_get_number_communication_requests, "Test requests")
188  call run_test_case(test_determine_halo_corner, "Test size of halo corners")
189  call run_test_case(test_determine_halo_corner_elements, "Test number elements in halo corners")
190  call run_test_case(test_get_pid_neighbour_location, "Test pid location")
191  call run_test_case(test_pid_been_seen, "Test pid been seen")
192  call run_test_case(test_retrieve_same_neighbour_information, "Test retrieve neighbour info")
193  call run_test_case(test_perform_local_data_copy_for_dimension, "Test local data copies")
194  call fruit_summary
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.
Contains the types used for communication, holding the state of communications and supporting activit...
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
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.
Definition: datadefn.F90:2
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.
Definition: conversions.F90:5
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.
Definition: conversions.F90:36
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, 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...
Converts data types to logical.
Definition: conversions.F90:69
Maintains the state of a halo swap and contains buffers, neighbours etc.
Map data structure that holds string (length 20 maximum) key value pairs.
Definition: collections.F90:86
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...
Defined the local grid, i.e. the grid held on this process after decomposition.
Definition: grids.F90:111
Collection data structures.
Definition: collections.F90:7
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...
Describes the neighbours of a process in a specific dimension and contains the communication buffers ...
Converts data types to real.
Definition: conversions.F90:58
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.
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
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.
program test_halo_communication_driver
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.
Converts data types to integers.
Definition: conversions.F90:47
integer, parameter, public y_index
Definition: grids.F90:14
integer, parameter, public x_index
Definition: grids.F90:14