MONC
Functions/Subroutines | Variables
haloswapper_mod Module Reference

Performs halo swapping. In the parallel case this is between neighbouring processes and in the serial case it still needs to wrap the halos around for the boundary conditions. This module determines the policy of halo swapping (i.e. the fields to communicate) and the halo communication module is used to provide the actual mechanism. More...

Functions/Subroutines

type(component_descriptor_type) function, public haloswapper_get_descriptor ()
 Provides a description of this component for the core to register. More...
 
subroutine initialisation_callback (current_state)
 Initialisation callback hook which will set up the halo swapping state and cache some precalculated data for fast(er) halo swapping at each timestep. More...
 
subroutine timestep_callback (current_state)
 Timestep callback hook which performs the halo swapping for each prognostic field. More...
 
subroutine finalisation_callback (current_state)
 The finalisation callback hook which will clean up and free the memory associated with the halo swapping. More...
 
subroutine copy_halo_buffer_to_corners (current_state, neighbour_description, corner_loc, x_target_index, y_target_index, neighbour_location, current_page, source_data)
 Copies the halo buffer to halo location in a corner for a halo cell/column and corner location. The copies are performed for each prognostic field. More...
 
subroutine copy_halo_buffer_to_field (current_state, neighbour_description, dim, target_index, neighbour_location, current_page, source_data)
 Copies the halo buffer to halo location in a field for a specific dimension and halo cell/column. The copies are performed for each prognostic field. More...
 
subroutine copy_fields_to_halo_buffer (current_state, neighbour_description, dim, source_index, pid_location, current_page, source_data)
 Copies the prognostic field data to halo buffers for a specific process in a dimension and halo cell. More...
 
subroutine copy_corners_to_halo_buffer (current_state, neighbour_description, corner_loc, x_source_index, y_source_index, pid_location, current_page, source_data)
 Copies the prognostic corner field data to halo buffers for a specific process in a dimension and halo cell. More...
 
integer function get_fields_per_halo_cell (current_state)
 Deduces the number of fields per halo cell. This depends upon what fields are active in the model. More...
 
subroutine display_debugging_info_if_needed (neighbour_counts, included_fields, rank)
 Displays some debugging information about who is sending what if that logging_mod level is selected. More...
 
subroutine perform_local_data_copy_for_all_prognostics (current_state, halo_depth, involve_corners, source_data)
 Will do any local copying of data required for the boundary conditions. I.e. if all columns in a slice are on a process then it will copy in the y dimension. More...
 

Variables

logical, private first_call = .true.
 
type(halo_communication_type), save halo_swap_state
 

Detailed Description

Performs halo swapping. In the parallel case this is between neighbouring processes and in the serial case it still needs to wrap the halos around for the boundary conditions. This module determines the policy of halo swapping (i.e. the fields to communicate) and the halo communication module is used to provide the actual mechanism.

Function/Subroutine Documentation

◆ copy_corners_to_halo_buffer()

subroutine haloswapper_mod::copy_corners_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  corner_loc,
integer, intent(in)  x_source_index,
integer, intent(in)  y_source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the prognostic corner field data to halo buffers for a specific process in a dimension and halo cell.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
corner_locLocation of the corner
x_source_indexThe X source index of the dimension we are reading from in the prognostic field
y_source_indexThe Y source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 307 of file haloswapper.F90.

309 
310  type(model_state_type), intent(inout) :: current_state
311  integer, intent(in) :: corner_loc, pid_location, x_source_index, y_source_index
312  integer, intent(inout) :: current_page(:)
313  type(neighbour_description_type), intent(inout) :: neighbour_description
314  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
315 
316  integer :: page_bookmark, i
317 
318  page_bookmark = current_page(pid_location)
319 
320 #ifdef U_ACTIVE
321  call copy_corner_to_buffer(current_state%local_grid, &
322  neighbour_description%send_corner_buffer,&
323  current_state%u%data, corner_loc, x_source_index, y_source_index, page_bookmark)
324  page_bookmark=page_bookmark+1
325  call copy_corner_to_buffer(current_state%local_grid, &
326  neighbour_description%send_corner_buffer, &
327  current_state%zu%data, corner_loc, x_source_index, y_source_index, page_bookmark)
328  page_bookmark=page_bookmark+1
329 #endif
330 #ifdef V_ACTIVE
331  call copy_corner_to_buffer(current_state%local_grid, &
332  neighbour_description%send_corner_buffer, &
333  current_state%v%data, corner_loc, x_source_index, y_source_index, page_bookmark)
334  page_bookmark=page_bookmark+1
335  call copy_corner_to_buffer(current_state%local_grid, &
336  neighbour_description%send_corner_buffer, &
337  current_state%zv%data, corner_loc, x_source_index, y_source_index, page_bookmark)
338  page_bookmark=page_bookmark+1
339 #endif
340 #ifdef W_ACTIVE
341  call copy_corner_to_buffer(current_state%local_grid, &
342  neighbour_description%send_corner_buffer, &
343  current_state%w%data, corner_loc, x_source_index, y_source_index, page_bookmark)
344  page_bookmark=page_bookmark+1
345  call copy_corner_to_buffer(current_state%local_grid, &
346  neighbour_description%send_corner_buffer, &
347  current_state%zw%data, corner_loc, x_source_index, y_source_index, page_bookmark)
348  page_bookmark=page_bookmark+1
349 #endif
350  if (current_state%th%active) then
351  call copy_corner_to_buffer(current_state%local_grid, &
352  neighbour_description%send_corner_buffer,&
353  current_state%th%data, corner_loc, x_source_index, y_source_index, page_bookmark)
354  page_bookmark=page_bookmark+1
355  call copy_corner_to_buffer(current_state%local_grid, &
356  neighbour_description%send_corner_buffer, &
357  current_state%zth%data, corner_loc, x_source_index, y_source_index, page_bookmark)
358  page_bookmark=page_bookmark+1
359  end if
360  do i=1,current_state%number_q_fields
361  if (current_state%q(i)%active) then
362  call copy_corner_to_buffer(current_state%local_grid, &
363  neighbour_description%send_corner_buffer, &
364  current_state%q(i)%data, corner_loc, x_source_index, y_source_index, page_bookmark)
365  page_bookmark=page_bookmark+1
366  call copy_corner_to_buffer(current_state%local_grid, &
367  neighbour_description%send_corner_buffer, &
368  current_state%zq(i)%data, corner_loc, x_source_index, y_source_index, page_bookmark)
369  page_bookmark=page_bookmark+1
370  end if
371  end do
372  current_page(pid_location)=page_bookmark
Here is the caller graph for this function:

◆ copy_fields_to_halo_buffer()

subroutine haloswapper_mod::copy_fields_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the prognostic field data to halo buffers for a specific process in a dimension and halo cell.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
dimDimension to copy from
source_indexThe source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 231 of file haloswapper.F90.

233 
234  type(model_state_type), intent(inout) :: current_state
235  integer, intent(in) :: dim, pid_location, source_index
236  integer, intent(inout) :: current_page(:)
237  type(neighbour_description_type), intent(inout) :: neighbour_description
238  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
239 
240  integer :: page_bookmark, i
241 
242  page_bookmark = current_page(pid_location)
243 
244 #ifdef U_ACTIVE
245  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer,&
246  current_state%u%data, dim, source_index, page_bookmark)
247  page_bookmark = page_bookmark + 1
248  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
249  current_state%zu%data, dim, source_index, page_bookmark)
250  page_bookmark = page_bookmark + 1
251 #endif
252 #ifdef V_ACTIVE
253  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
254  current_state%v%data, dim, source_index, page_bookmark)
255  page_bookmark = page_bookmark + 1
256  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
257  current_state%zv%data, dim, source_index, page_bookmark)
258  page_bookmark = page_bookmark + 1
259 #endif
260 #ifdef W_ACTIVE
261  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
262  current_state%w%data, dim, source_index, page_bookmark)
263  page_bookmark = page_bookmark+1
264  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
265  current_state%zw%data, dim, source_index, page_bookmark)
266  page_bookmark = page_bookmark+1
267 #endif
268  if (current_state%th%active) then
269  call copy_field_to_buffer(current_state%local_grid, &
270  neighbour_description%send_halo_buffer, current_state%th%data, dim, source_index,&
271  page_bookmark)
272  page_bookmark = page_bookmark+1
273  call copy_field_to_buffer(current_state%local_grid, &
274  neighbour_description%send_halo_buffer, current_state%zth%data, dim, source_index, &
275  page_bookmark)
276  page_bookmark = page_bookmark+1
277  end if
278 
279  do i = 1,current_state%number_q_fields
280  if (current_state%q(i)%active) then
281  call copy_field_to_buffer(current_state%local_grid, &
282  neighbour_description%send_halo_buffer, current_state%q(i)%data, dim, &
283  source_index, page_bookmark)
284  page_bookmark = page_bookmark + 1
285  call copy_field_to_buffer(current_state%local_grid, &
286  neighbour_description%send_halo_buffer, &
287  current_state%zq(i)%data, dim, source_index, page_bookmark)
288  page_bookmark = page_bookmark + 1
289  end if
290  end do
291  current_page(pid_location) = page_bookmark
Here is the caller graph for this function:

◆ copy_halo_buffer_to_corners()

subroutine haloswapper_mod::copy_halo_buffer_to_corners ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  corner_loc,
integer, intent(in)  x_target_index,
integer, intent(in)  y_target_index,
integer, intent(in)  neighbour_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the halo buffer to halo location in a corner for a halo cell/column and corner location. The copies are performed for each prognostic field.

Parameters
current_stateThe current model state
neighbour_descriptionThe halo swapping description of the neighbour we are accessing the buffer of
corner_locThe location of the corner
x_target_indexThe target index for the x dimension we are receiving for
y_target_indexThe target index for the y dimension we are receiving for
neighbour_locationThe location in the local neighbour data stores of this neighbour
current_pageThe current, next, halo swap page to read from (all previous have been read and copied already)
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 96 of file haloswapper.F90.

98 
99  type(model_state_type), intent(inout) :: current_state
100  integer, intent(in) :: corner_loc, x_target_index, y_target_index, neighbour_location
101  integer, intent(inout) :: current_page(:)
102  type(neighbour_description_type), intent(inout) :: neighbour_description
103  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
104 
105  integer :: page_bookmark, i
106 
107  page_bookmark = current_page(neighbour_location)
108 #ifdef U_ACTIVE
109  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer,&
110  current_state%u%data, corner_loc, x_target_index, y_target_index, page_bookmark)
111  page_bookmark=page_bookmark+1
112  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
113  current_state%zu%data, corner_loc, x_target_index, y_target_index, page_bookmark)
114  page_bookmark=page_bookmark+1
115 #endif
116 #ifdef V_ACTIVE
117  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
118  current_state%v%data, corner_loc, x_target_index, y_target_index, page_bookmark)
119  page_bookmark=page_bookmark+1
120  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
121  current_state%zv%data, corner_loc, x_target_index, y_target_index, page_bookmark)
122  page_bookmark=page_bookmark+1
123 #endif
124 #ifdef W_ACTIVE
125  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
126  current_state%w%data, corner_loc, x_target_index, y_target_index, page_bookmark)
127  page_bookmark=page_bookmark+1
128  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
129  current_state%zw%data, corner_loc, x_target_index, y_target_index, page_bookmark)
130  page_bookmark=page_bookmark+1
131 #endif
132  if (current_state%th%active) then
133  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
134  current_state%th%data, corner_loc, x_target_index, y_target_index, page_bookmark)
135  page_bookmark=page_bookmark+1
136  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
137  current_state%zth%data, corner_loc, x_target_index, y_target_index, page_bookmark)
138  page_bookmark=page_bookmark+1
139  end if
140  do i=1,current_state%number_q_fields
141  if (current_state%q(i)%active) then
142  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
143  current_state%q(i)%data, corner_loc, x_target_index, y_target_index, page_bookmark)
144  page_bookmark=page_bookmark + 1
145  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
146  current_state%zq(i)%data, corner_loc, x_target_index, y_target_index, page_bookmark)
147  page_bookmark=page_bookmark + 1
148  end if
149  end do
150  current_page(neighbour_location)=page_bookmark
Here is the caller graph for this function:

◆ copy_halo_buffer_to_field()

subroutine haloswapper_mod::copy_halo_buffer_to_field ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  target_index,
integer, intent(in)  neighbour_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the halo buffer to halo location in a field for a specific dimension and halo cell/column. The copies are performed for each prognostic field.

Parameters
current_stateThe current model state
neighbour_descriptionThe halo swapping description of the neighbour we are accessing the buffer of
dimThe dimension we receive for
target_indexThe target index for the dimension we are receiving for
neighbour_locationThe location in the local neighbour data stores of this neighbour
current_pageThe current, next, halo swap page to read from (all previous have been read and copied already)
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 164 of file haloswapper.F90.

166  type(model_state_type), intent(inout) :: current_state
167  integer, intent(in) :: dim, target_index, neighbour_location
168  integer, intent(inout) :: current_page(:)
169  type(neighbour_description_type), intent(inout) :: neighbour_description
170  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
171 
172  integer :: page_bookmark, i
173 
174  page_bookmark = current_page(neighbour_location)
175 #ifdef U_ACTIVE
176  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
177  current_state%u%data, dim, target_index, page_bookmark)
178  page_bookmark=page_bookmark+1
179  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
180  current_state%zu%data, dim, target_index, page_bookmark)
181  page_bookmark=page_bookmark+1
182 #endif
183 #ifdef V_ACTIVE
184  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
185  current_state%v%data, dim, target_index, page_bookmark)
186  page_bookmark=page_bookmark+1
187  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
188  current_state%zv%data, dim, target_index, page_bookmark)
189  page_bookmark=page_bookmark+1
190 #endif
191 #ifdef W_ACTIVE
192  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
193  current_state%w%data, dim, target_index, page_bookmark)
194  page_bookmark=page_bookmark+1
195  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
196  current_state%zw%data, dim, target_index, page_bookmark)
197  page_bookmark=page_bookmark+1
198 #endif
199  if (current_state%th%active) then
200  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
201  current_state%th%data, dim, target_index, page_bookmark)
202  page_bookmark=page_bookmark+1
203  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
204  current_state%zth%data, dim, target_index, page_bookmark)
205  page_bookmark=page_bookmark+1
206  end if
207  do i=1,current_state%number_q_fields
208  if (current_state%q(i)%active) then
209  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
210  current_state%q(i)%data, dim, target_index, page_bookmark)
211  page_bookmark=page_bookmark+1
212  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
213  current_state%zq(i)%data, dim, target_index, page_bookmark)
214  page_bookmark=page_bookmark+1
215  end if
216  end do
217  current_page(neighbour_location)=page_bookmark
Here is the caller graph for this function:

◆ display_debugging_info_if_needed()

subroutine haloswapper_mod::display_debugging_info_if_needed ( integer, intent(in)  neighbour_counts,
integer, intent(in)  included_fields,
integer, intent(in)  rank 
)
private

Displays some debugging information about who is sending what if that logging_mod level is selected.

Parameters
neighbour_countsNumber of neighbours that I have
included_fieldsNumber of prognostic fields that we include in this halo swap
rankMy current rank

Definition at line 407 of file haloswapper.F90.

408  integer, intent(in) :: neighbour_counts, rank, included_fields
409 
410  if (.not. first_call .or. log_get_logging_level() .lt. log_debug) return
411  call log_log(log_debug, "Rank "//trim(conv_to_string(rank))//": "//trim(conv_to_string(neighbour_counts))//&
412  " neighbours per timestep over "//trim(conv_to_string(included_fields))//" fields")
413  first_call = .false.
Here is the caller graph for this function:

◆ finalisation_callback()

subroutine haloswapper_mod::finalisation_callback ( type(model_state_type), intent(inout), target  current_state)
private

The finalisation callback hook which will clean up and free the memory associated with the halo swapping.

Parameters
current_stateThe current model state

Definition at line 76 of file haloswapper.F90.

77  type(model_state_type), target, intent(inout) :: current_state
78 
79  call finalise_halo_communication(halo_swap_state)
Here is the caller graph for this function:

◆ get_fields_per_halo_cell()

integer function haloswapper_mod::get_fields_per_halo_cell ( type(model_state_type), intent(inout)  current_state)
private

Deduces the number of fields per halo cell. This depends upon what fields are active in the model.

Parameters
current_stateThe current model state

Definition at line 377 of file haloswapper.F90.

378  type(model_state_type), intent(inout) :: current_state
379 
380  integer :: i
381 
382  get_fields_per_halo_cell=0
383 
384 #ifdef U_ACTIVE
385  get_fields_per_halo_cell=get_fields_per_halo_cell+2
386 #endif
387 #ifdef V_ACTIVE
388  get_fields_per_halo_cell=get_fields_per_halo_cell+2
389 #endif
390 #ifdef W_ACTIVE
391  get_fields_per_halo_cell=get_fields_per_halo_cell+2
392 #endif
393  if (current_state%th%active) then
394  get_fields_per_halo_cell=get_fields_per_halo_cell+2
395  end if
396  do i=1,current_state%number_q_fields
397  if (current_state%q(i)%active) then
398  get_fields_per_halo_cell=get_fields_per_halo_cell+2
399  end if
400  end do
Here is the caller graph for this function:

◆ haloswapper_get_descriptor()

type(component_descriptor_type) function, public haloswapper_mod::haloswapper_get_descriptor

Provides a description of this component for the core to register.

Returns
The descriptor containing registration information for this component

Definition at line 36 of file haloswapper.F90.

37  haloswapper_get_descriptor%name="halo_swapper"
38  haloswapper_get_descriptor%version=0.1
39  haloswapper_get_descriptor%initialisation=>initialisation_callback
40  haloswapper_get_descriptor%timestep=>timestep_callback
41  haloswapper_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ initialisation_callback()

subroutine haloswapper_mod::initialisation_callback ( type(model_state_type), intent(inout), target  current_state)
private

Initialisation callback hook which will set up the halo swapping state and cache some precalculated data for fast(er) halo swapping at each timestep.

Parameters
current_stateThe current model state

Definition at line 47 of file haloswapper.F90.

48  type(model_state_type), target, intent(inout) :: current_state
49  integer :: halo_depth
50 
51  ! get halo_depth and pass it to the halo_swapping routines
52  halo_depth = options_get_integer(current_state%options_database, "halo_depth")
53  call init_halo_communication(current_state, get_fields_per_halo_cell, halo_swap_state, &
54  halo_depth, .true.)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ perform_local_data_copy_for_all_prognostics()

subroutine haloswapper_mod::perform_local_data_copy_for_all_prognostics ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  halo_depth,
logical, intent(in)  involve_corners,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Will do any local copying of data required for the boundary conditions. I.e. if all columns in a slice are on a process then it will copy in the y dimension.

Parameters
current_stateThe current model state_mod
copy_countsNumber of local copies performed
source_dataOptional source data which is read from into send buffers and written into by receieve buffers

Definition at line 422 of file haloswapper.F90.

423  type(model_state_type), intent(inout) :: current_state
424  integer, intent(in) :: halo_depth
425  logical, intent(in) :: involve_corners
426  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
427 
428  integer :: i
429 
430 #ifdef U_ACTIVE
431  call perform_local_data_copy_for_field(current_state%u%data, current_state%local_grid, &
432  current_state%parallel%my_rank, halo_depth, involve_corners)
433  call perform_local_data_copy_for_field(current_state%zu%data, current_state%local_grid, &
434  current_state%parallel%my_rank, halo_depth, involve_corners)
435 #endif
436 #ifdef V_ACTIVE
437  call perform_local_data_copy_for_field(current_state%v%data, current_state%local_grid, &
438  current_state%parallel%my_rank, halo_depth, involve_corners)
439  call perform_local_data_copy_for_field(current_state%zv%data, current_state%local_grid, &
440  current_state%parallel%my_rank, halo_depth, involve_corners)
441 #endif
442 #ifdef W_ACTIVE
443  call perform_local_data_copy_for_field(current_state%w%data, current_state%local_grid, &
444  current_state%parallel%my_rank, halo_depth, involve_corners)
445  call perform_local_data_copy_for_field(current_state%zw%data, current_state%local_grid, &
446  current_state%parallel%my_rank, halo_depth, involve_corners)
447 #endif
448  if (current_state%th%active) then
449  call perform_local_data_copy_for_field(current_state%th%data, current_state%local_grid, &
450  current_state%parallel%my_rank, halo_depth, involve_corners)
451  call perform_local_data_copy_for_field(current_state%zth%data, current_state%local_grid, &
452  current_state%parallel%my_rank, halo_depth, involve_corners)
453  end if
454  do i=1,current_state%number_q_fields
455  if (current_state%q(i)%active) then
456  call perform_local_data_copy_for_field(current_state%q(i)%data, current_state%local_grid, &
457  current_state%parallel%my_rank, halo_depth, involve_corners)
458  call perform_local_data_copy_for_field(current_state%zq(i)%data, current_state%local_grid, &
459  current_state%parallel%my_rank, halo_depth, involve_corners)
460  end if
461  end do
Here is the caller graph for this function:

◆ timestep_callback()

subroutine haloswapper_mod::timestep_callback ( type(model_state_type), intent(inout), target  current_state)
private

Timestep callback hook which performs the halo swapping for each prognostic field.

In parallel this is performed with MPI communication calls and wrapping around. In serial still need to wrap data around

Parameters
current_stateThe current model state_mod

Definition at line 62 of file haloswapper.F90.

63  type(model_state_type), target, intent(inout) :: current_state
64 
65  call blocking_halo_swap(current_state, halo_swap_state, copy_fields_to_halo_buffer, &
66  perform_local_data_copy_for_all_prognostics, copy_halo_buffer_to_field, &
67  copy_corners_to_halo_buffer, copy_halo_buffer_to_corners)
68  call display_debugging_info_if_needed(halo_swap_state%number_distinct_neighbours, &
69  halo_swap_state%fields_per_cell, current_state%parallel%my_rank)
70 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ first_call

logical, private haloswapper_mod::first_call = .true.
private

Definition at line 27 of file haloswapper.F90.

27  logical, private :: first_call = .true.

◆ halo_swap_state

type(halo_communication_type), save haloswapper_mod::halo_swap_state
private

Definition at line 28 of file haloswapper.F90.

28  type(halo_communication_type), save :: halo_swap_state
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
logging_mod::log_get_logging_level
integer function, public log_get_logging_level()
Retrieves the current logging level.
Definition: logging.F90:122
logging_mod::log_debug
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
Definition: logging.F90:14