MONC
iobridge.F90
Go to the documentation of this file.
1 
9  use state_mod, only : model_state_type
24  use mpi, only : mpi_comm_world, mpi_int, mpi_byte, mpi_request_null, mpi_statuses_ignore, mpi_status_ignore, mpi_status_size
26  implicit none
27 
28 #ifndef TEST_MODE
29  private
30 #endif
31 
33  integer :: number_dimensions, dimensions(4)
35 
37  character(len=STRING_LENGTH) :: name
38  integer :: field_type, data_type
39  logical :: optional, enabled
41 
43  character(len=STRING_LENGTH) :: name
44  logical :: send_on_terminate
45  integer :: number_of_data_fields, frequency, mpi_datatype
46  type(io_configuration_field_type), dimension(:), allocatable :: fields
47  integer :: dump_requests(2)
48  character, dimension(:), allocatable :: send_buffer
50 
51  type(io_configuration_data_definition_type), dimension(:), allocatable :: data_definitions
54 
56 
57 contains
58 
60  iobridge_get_descriptor%name="iobridge"
61  iobridge_get_descriptor%version=0.1
65  end function iobridge_get_descriptor
66 
69  subroutine init_callback(current_state)
70  type(model_state_type), target, intent(inout) :: current_state
71 
72  integer :: mpi_type_data_sizing_description, mpi_type_definition_description, mpi_type_field_description, ierr
73 
74  if (.not. options_get_logical(current_state%options_database, "enable_io_server")) then
75  io_server_enabled=.false.
76  call log_master_log(log_warn, "Enabled IO bridge but missing IO server compilation, therefore ignoring IO bridge component")
77  return
78  end if
79 
80  io_server_enabled=.true.
82 
83  call populate_sendable_fields(current_state)
84 
85  mpi_type_data_sizing_description=build_mpi_type_data_sizing_description()
86  mpi_type_definition_description=build_mpi_type_definition_description()
87  mpi_type_field_description=build_mpi_type_field_description()
88 
89  call register_with_io_server(current_state, mpi_type_definition_description, mpi_type_field_description)
90  call send_monc_specific_data_to_server(current_state, mpi_type_data_sizing_description)
91 
92  call mpi_type_free(mpi_type_data_sizing_description, ierr)
93  call mpi_type_free(mpi_type_definition_description, ierr)
94  call mpi_type_free(mpi_type_field_description, ierr)
95 
97  end subroutine init_callback
98 
101  subroutine timestep_callback(current_state)
102  type(model_state_type), target, intent(inout) :: current_state
103 
104  integer :: i
105 
106  if (.not. io_server_enabled) return
107 
108  do i=1, size(data_definitions)
109  if (data_definitions(i)%frequency .gt. 0) then
110  if (mod(current_state%timestep, data_definitions(i)%frequency) == 0) then
111  call send_data_to_io_server(current_state, i)
112  end if
113  end if
114  end do
115  end subroutine timestep_callback
116 
120  subroutine send_data_to_io_server(current_state, data_index)
121  type(model_state_type), target, intent(inout) :: current_state
122  integer, intent(in) :: data_index
123 
124  integer :: command_to_send, ierr
125 
126  if (data_definitions(data_index)%dump_requests(1) .ne. mpi_request_null .or. &
127  data_definitions(data_index)%dump_requests(2) .ne. mpi_request_null) then
128  ! Here wait for previous data dump to complete (consider extending to using buffers for performance)
129  call mpi_waitall(2, data_definitions(data_index)%dump_requests, mpi_statuses_ignore, ierr)
130  end if
131 
132  ! Pack the send buffer and send it to the IO server
133  call pack_send_buffer(current_state, data_definitions(data_index))
134 
135  command_to_send=data_command_start+data_index
136  call mpi_issend(command_to_send, 1, mpi_int, current_state%parallel%corresponding_io_server_process, &
137  command_tag, mpi_comm_world, data_definitions(data_index)%dump_requests(1), ierr)
138  call mpi_issend(data_definitions(data_index)%send_buffer, 1, data_definitions(data_index)%mpi_datatype, &
139  current_state%parallel%corresponding_io_server_process, data_tag+data_index, mpi_comm_world, &
140  data_definitions(data_index)%dump_requests(2), ierr)
141  end subroutine send_data_to_io_server
142 
145  subroutine finalisation_callback(current_state)
146  type(model_state_type), target, intent(inout) :: current_state
147 
148  integer :: ierr, i
149 
150  if (.not. io_server_enabled) return
152 
153  do i=1, size(data_definitions)
154  if (data_definitions(i)%send_on_terminate) then
155  call send_data_to_io_server(current_state, i)
156  end if
157  if (data_definitions(i)%dump_requests(1) .ne. mpi_request_null .or. &
158  data_definitions(i)%dump_requests(2) .ne. mpi_request_null) then
159  call mpi_waitall(2, data_definitions(i)%dump_requests, mpi_statuses_ignore, ierr)
160  end if
161  if (allocated(data_definitions(i)%send_buffer)) deallocate(data_definitions(i)%send_buffer)
162  call mpi_type_free(data_definitions(i)%mpi_datatype, ierr)
163  end do
164  call mpi_send(deregister_command, 1, mpi_int, current_state%parallel%corresponding_io_server_process, &
165  command_tag, mpi_comm_world, ierr)
166  end subroutine finalisation_callback
167 
169  subroutine build_mpi_data_types()
170  integer :: i, dump_send_buffer_size
171 
172  do i=1, size(data_definitions)
173  dump_send_buffer_size=build_mpi_data_type_for_definition(data_definitions(i))
174  data_definitions(i)%dump_requests=(/mpi_request_null, mpi_request_null/)
175  allocate(data_definitions(i)%send_buffer(dump_send_buffer_size))
176  end do
177  end subroutine build_mpi_data_types
178 
182  integer function build_mpi_data_type_for_definition(specific_data_definition)
183  type(io_configuration_data_definition_type), intent(inout) :: specific_data_definition
184 
185  integer :: type_extents(5), type_counts, i, j, tempsize, field_start, data_type, field_array_sizes, &
186  temp_size, prev_data_type, old_types(20), offsets(20), block_counts(20), ierr, field_ignores
187  logical :: field_found
188  type(io_server_sendable_field_sizing) :: field_size_info
189 
190  type_extents=populate_mpi_type_extents()
191 
192  field_start=1
193  data_type=0
194  type_counts=0
195  field_array_sizes=0
196  field_ignores=0
197  do i=1, specific_data_definition%number_of_data_fields
198  if (data_type == 0) then
199  prev_data_type=data_type
200  data_type=specific_data_definition%fields(i)%data_type
201  else
202  if (data_type .ne. specific_data_definition%fields(i)%data_type) then
203  ! For efficient type packing, combine multiple fields with the same type - therefore when the type changes work the previous one pack
204  call append_mpi_datatype(field_start, i-1-field_ignores, field_array_sizes, data_type, &
205  type_extents, prev_data_type, type_counts+1, old_types, offsets, block_counts)
206  field_start=i
207  field_array_sizes=0
208  field_ignores=0
209  prev_data_type=data_type
210  data_type=specific_data_definition%fields(i)%data_type
211  type_counts=type_counts+1
212  end if
213  end if
214 
215  if (specific_data_definition%fields(i)%field_type .eq. array_field_type .or. &
216  specific_data_definition%fields(i)%field_type .eq. map_field_type) then
217  ! Grab the field info based upon the field name
218  field_size_info=get_sendable_field_sizing(specific_data_definition%fields(i)%name, field_found)
219  specific_data_definition%fields(i)%enabled=field_found
220  if (.not. field_found .or. field_size_info%number_dimensions == 0) then
221  ! If no field info, or the dimension is 0 then this MONC process is not sending that field - check it is optional
222  if (.not. specific_data_definition%fields(i)%optional) then
223  call log_log(log_error, "Non optional field `"//trim(specific_data_definition%fields(i)%name)//&
224  "' omitted from MONC IO server registration")
225  end if
226  field_ignores=field_ignores+1
227  else
228  ! If the field is specified then use the size data to assemble the field size and append to current size
229  temp_size=1
230  do j=1, field_size_info%number_dimensions
231  temp_size=temp_size*field_size_info%dimensions(j)
232  end do
233  if (specific_data_definition%fields(i)%field_type .eq. map_field_type) then
234  field_array_sizes=(field_array_sizes+temp_size*string_length*2)-1
235  else
236  field_array_sizes=(field_array_sizes+temp_size)-1
237  end if
238  end if
239  else
240  if (specific_data_definition%fields(i)%optional) then
241  field_size_info=get_sendable_field_sizing(specific_data_definition%fields(i)%name, field_found)
242  specific_data_definition%fields(i)%enabled=field_found
243  if (.not. field_found) field_ignores=field_ignores+1
244  end if
245  end if
246  end do
247  if (field_start .le. i-1) then
248  ! If there are outstanding fields to process then we do this here
249  call append_mpi_datatype(field_start, i-1, field_array_sizes, data_type, &
250  type_extents, prev_data_type, type_counts+1, old_types, offsets, block_counts)
251  type_counts=type_counts+1
252  end if
253 
254  call mpi_type_struct(type_counts, block_counts, offsets, old_types, specific_data_definition%mpi_datatype, ierr)
255  call mpi_type_commit(specific_data_definition%mpi_datatype, ierr)
256  call mpi_type_size(specific_data_definition%mpi_datatype, tempsize, ierr)
259 
262  subroutine populate_sendable_fields(current_state)
263  type(model_state_type), target, intent(inout) :: current_state
264 
265  type(list_type) :: published_field_descriptors
266  integer :: i
267 
269  published_field_descriptors=get_all_component_published_fields()
270  do i=1, c_size(published_field_descriptors)
271  call populate_component_public_field(current_state, c_get_string(published_field_descriptors, i))
272  end do
273  end subroutine populate_sendable_fields
274 
277  subroutine populate_component_public_field(current_state, field_name)
278  type(model_state_type), target, intent(inout) :: current_state
279  character(len=*), intent(in) :: field_name
280 
281  class(*), pointer :: generic_data
282  type(io_server_sendable_field_sizing), pointer :: field_sizing
283  type(component_field_information_type) :: field_information
284  type(component_field_information_type), pointer :: field_information_info_alloc
285 
286 
287  field_information=get_component_field_information(current_state, field_name)
288  if (field_information%enabled) then
289  allocate(field_information_info_alloc, source=field_information)
290  generic_data=>field_information_info_alloc
291  call c_put_generic(component_field_descriptions, field_name, generic_data, .false.)
292 
293  allocate(field_sizing)
294  field_sizing%number_dimensions=field_information%number_dimensions
295  field_sizing%dimensions=field_information%dimension_sizes
296  generic_data=>field_sizing
297  call c_put_generic(sendable_fields, field_name, generic_data, .false.)
298  end if
299  end subroutine populate_component_public_field
300 
303  subroutine populate_globally_visible_sendable_fields(current_state)
304  type(model_state_type), target, intent(inout) :: current_state
305 
306  integer :: x_size, y_size, z_size
307  class(*), pointer :: raw_generic
308 
309  z_size=current_state%local_grid%size(z_index)
310  y_size=current_state%local_grid%size(y_index)
311  x_size=current_state%local_grid%size(x_index)
312 
313  raw_generic=>generate_sendable_description(options_size(current_state%options_database))
314  call c_put_generic(sendable_fields, "options_database", raw_generic, .false.)
315  if (get_number_active_q_indices() .gt. 0) then
316  raw_generic=>generate_sendable_description(get_number_active_q_indices())
317  call c_put_generic(sendable_fields, "q_indicies", raw_generic, .false.)
318  end if
319  raw_generic=>generate_sendable_description(3)
320  call c_put_generic(sendable_fields, "local_grid_size", raw_generic, .false.)
321  call c_put_generic(sendable_fields, "local_grid_start", raw_generic, .false.)
322 #ifdef U_ACTIVE
323  raw_generic=>generate_sendable_description()
324  call c_put_generic(sendable_fields, "x_size", raw_generic, .false.)
325  call c_put_generic(sendable_fields, "x_bottom", raw_generic, .false.)
326  call c_put_generic(sendable_fields, "x_resolution", raw_generic, .false.)
327  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
328  call c_put_generic(sendable_fields, "u", raw_generic, .false.)
329  call c_put_generic(sendable_fields, "u_nogal", raw_generic, .false.)
330  call c_put_generic(sendable_fields, "zu", raw_generic, .false.)
331  if (allocated(current_state%global_grid%configuration%vertical%olubar)) then
332  raw_generic=>generate_sendable_description(z_size)
333  call c_put_generic(sendable_fields, "olubar", raw_generic, .false.)
334  call c_put_generic(sendable_fields, "olzubar", raw_generic, .false.)
335  end if
336 #endif
337 #ifdef V_ACTIVE
338  raw_generic=>generate_sendable_description()
339  call c_put_generic(sendable_fields, "y_size", raw_generic, .false.)
340  call c_put_generic(sendable_fields, "y_bottom", raw_generic, .false.)
341  call c_put_generic(sendable_fields, "y_resolution", raw_generic, .false.)
342  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
343  call c_put_generic(sendable_fields, "v", raw_generic, .false.)
344  call c_put_generic(sendable_fields, "v_nogal", raw_generic, .false.)
345  call c_put_generic(sendable_fields, "zv", raw_generic, .false.)
346  if (allocated(current_state%global_grid%configuration%vertical%olvbar)) then
347  raw_generic=>generate_sendable_description(z_size)
348  call c_put_generic(sendable_fields, "olvbar", raw_generic, .false.)
349  call c_put_generic(sendable_fields, "olzvbar", raw_generic, .false.)
350  end if
351 #endif
352 #ifdef W_ACTIVE
353  raw_generic=>generate_sendable_description(z_size)
354  call c_put_generic(sendable_fields, "z", raw_generic, .false.)
355  call c_put_generic(sendable_fields, "thref", raw_generic, .false.)
356  call c_put_generic(sendable_fields, "prefn", raw_generic, .false.)
357  call c_put_generic(sendable_fields, "rhon", raw_generic, .false.)
358  call c_put_generic(sendable_fields, "rho", raw_generic, .false.)
359  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
360  call c_put_generic(sendable_fields, "w", raw_generic, .false.)
361  call c_put_generic(sendable_fields, "zw", raw_generic, .false.)
362 #endif
363  if (current_state%number_q_fields .gt. 0) then
364  raw_generic=>generate_sendable_description(z_size, y_size, x_size, current_state%number_q_fields)
365  call c_put_generic(sendable_fields, "q", raw_generic, .false.)
366  call c_put_generic(sendable_fields, "zq", raw_generic, .false.)
367  if (allocated(current_state%global_grid%configuration%vertical%olqbar)) then
368  raw_generic=>generate_sendable_description(z_size, current_state%number_q_fields)
369  call c_put_generic(sendable_fields, "olqbar", raw_generic, .false.)
370  call c_put_generic(sendable_fields, "olzqbar", raw_generic, .false.)
371  end if
372  end if
373  if (current_state%th%active) then
374  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
375  call c_put_generic(sendable_fields, "th", raw_generic, .false.)
376  call c_put_generic(sendable_fields, "zth", raw_generic, .false.)
377  if (allocated(current_state%global_grid%configuration%vertical%olthbar)) then
378  raw_generic=>generate_sendable_description(z_size)
379  call c_put_generic(sendable_fields, "olthbar", raw_generic, .false.)
380  call c_put_generic(sendable_fields, "olzthbar", raw_generic, .false.)
381  end if
382  end if
383  if (current_state%p%active) then
384  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
385  call c_put_generic(sendable_fields, "p", raw_generic, .false.)
386  end if
388 
396  function generate_sendable_description(dim1, dim2, dim3, dim4)
397  integer, intent(in), optional :: dim1, dim2, dim3, dim4
398  class(*), pointer :: generate_sendable_description
399 
400  type(io_server_sendable_field_sizing), pointer :: field
401  integer :: number_dimensions
402 
403  allocate(field)
404  number_dimensions=0
405  field%dimensions=0
406  if (present(dim1)) then
407  field%dimensions(1)=dim1
408  number_dimensions=number_dimensions+1
409  end if
410  if (present(dim2)) then
411  field%dimensions(2)=dim2
412  number_dimensions=number_dimensions+1
413  end if
414  if (present(dim3)) then
415  field%dimensions(3)=dim3
416  number_dimensions=number_dimensions+1
417  end if
418  if (present(dim4)) then
419  field%dimensions(4)=dim4
420  number_dimensions=number_dimensions+1
421  end if
422  field%number_dimensions=number_dimensions
423  generate_sendable_description=>field
424  end function generate_sendable_description
425 
430  subroutine send_monc_specific_data_to_server(current_state, mpi_type_data_sizing_description)
431  type(model_state_type), target, intent(inout) :: current_state
432  integer, intent(in) :: mpi_type_data_sizing_description
433 
434  type(data_sizing_description_type), dimension(:), allocatable :: data_description
435  character, dimension(:), allocatable :: buffer
436  integer :: number_unique_fields, buffer_size, request_handles(2), ierr
437  real(kind=DEFAULT_PRECISION) :: dreal
438 
439  number_unique_fields=c_size(unique_field_names)
440  allocate(data_description(number_unique_fields+4))
441  request_handles(1)=send_data_field_sizes_to_server(current_state, mpi_type_data_sizing_description, &
442  data_description, number_unique_fields)
443  buffer_size=(kind(dreal)*current_state%local_grid%size(z_index)) + (string_length * current_state%number_q_fields)
444  allocate(buffer(buffer_size))
445  request_handles(2)=send_general_monc_information_to_server(current_state, buffer)
446  call mpi_waitall(2, request_handles, mpi_statuses_ignore, ierr)
447  deallocate(data_description)
448  deallocate(buffer)
449  end subroutine send_monc_specific_data_to_server
450 
456  integer function send_data_field_sizes_to_server(current_state, mpi_type_data_sizing_description, &
457  data_description, number_unique_fields)
458  type(model_state_type), target, intent(inout) :: current_state
459  integer, intent(in) :: mpi_type_data_sizing_description, number_unique_fields
460  type(data_sizing_description_type), dimension(:), intent(inout) :: data_description
461 
462  integer :: ierr, i, next_index, request_handle
463  character(len=STRING_LENGTH) :: field_name
464 
465  call package_local_monc_decomposition_into_descriptions(current_state, data_description)
466  next_index=5
467  do i=1, number_unique_fields
468  field_name=c_key_at(unique_field_names, i)
469  if (c_contains(sendable_fields, field_name)) then
470  call assemble_individual_description(data_description, next_index, field_name, get_sendable_field_sizing(field_name))
471  next_index=next_index+1
472  end if
473  end do
474 
475  call mpi_isend(data_description, next_index-1, mpi_type_data_sizing_description, &
476  current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, request_handle, ierr)
477  send_data_field_sizes_to_server=request_handle
479 
484  integer function send_general_monc_information_to_server(current_state, buffer)
485  type(model_state_type), target, intent(inout) :: current_state
486  character, dimension(:), intent(inout) :: buffer
487 
488  character(len=STRING_LENGTH) :: q_field_name
489  type(q_metadata_type) :: q_meta_data
490  integer :: current_loc, n, ierr, request_handle
491 
492  current_loc=1
493  current_loc=pack_array_field(buffer, current_loc, real_array_1d=current_state%global_grid%configuration%vertical%zn)
494  if (current_state%number_q_fields .gt. 0) then
495  do n=1, current_state%number_q_fields
496  q_meta_data=get_indices_descriptor(n)
497  if (q_meta_data%l_used) then
498  q_field_name=q_meta_data%name
499  else
500  q_field_name="qfield_"//trim(conv_to_string(n))
501  end if
502  current_loc=pack_scalar_field(buffer, current_loc, string_value=q_field_name)
503  end do
504  end if
505  call mpi_isend(buffer, current_loc-1, mpi_byte, current_state%parallel%corresponding_io_server_process, &
506  data_tag, mpi_comm_world, request_handle, ierr)
509 
513  subroutine package_local_monc_decomposition_into_descriptions(current_state, data_description)
514  type(model_state_type), target, intent(inout) :: current_state
515  type(data_sizing_description_type), dimension(:), intent(inout) :: data_description
516 
517  type(io_server_sendable_field_sizing) :: sizing_info
518 
519  sizing_info%number_dimensions=3
520  sizing_info%dimensions(z_index)=current_state%local_grid%size(z_index)
521  sizing_info%dimensions(y_index)=current_state%local_grid%size(y_index)
522  sizing_info%dimensions(x_index)=current_state%local_grid%size(x_index)
523  call assemble_individual_description(data_description, 1, local_sizes_key, sizing_info)
524  sizing_info%dimensions(z_index)=current_state%local_grid%start(z_index)
525  sizing_info%dimensions(y_index)=current_state%local_grid%start(y_index)
526  sizing_info%dimensions(x_index)=current_state%local_grid%start(x_index)
527  call assemble_individual_description(data_description, 2, local_start_points_key, sizing_info)
528  sizing_info%dimensions(z_index)=current_state%local_grid%end(z_index)
529  sizing_info%dimensions(y_index)=current_state%local_grid%end(y_index)
530  sizing_info%dimensions(x_index)=current_state%local_grid%end(x_index)
531  call assemble_individual_description(data_description, 3, local_end_points_key, sizing_info)
532  sizing_info%number_dimensions=1
533  sizing_info%dimensions(1)=get_number_active_q_indices()
534  call assemble_individual_description(data_description, 4, number_q_indicies_key, sizing_info)
536 
540  type(io_server_sendable_field_sizing) function get_sendable_field_sizing(field_name, field_found)
541  character(len=*), intent(in) :: field_name
542  logical, intent(out), optional :: field_found
543 
544  class(*), pointer :: generic
545 
546  if (present(field_found)) field_found=.false.
547  if (c_contains(sendable_fields, field_name)) then
548  generic=>c_get_generic(sendable_fields, field_name)
549  select type(generic)
552  if (present(field_found)) field_found=.true.
553  end select
554  end if
555  end function get_sendable_field_sizing
556 
559  type(component_field_information_type) function get_component_field_descriptor(field_name)
560  character(len=*), intent(in) :: field_name
561 
562  class(*), pointer :: generic
563  if (c_contains(component_field_descriptions, field_name)) then
564  generic=>c_get_generic(component_field_descriptions, field_name)
565  select type(generic)
566  type is (component_field_information_type)
568  end select
569  end if
570  end function get_component_field_descriptor
571 
581  subroutine assemble_individual_description(data_description, index, field_name, field_sizing_description)
582  integer, intent(in) :: index
583  character(len=*), intent(in) :: field_name
584  type(io_server_sendable_field_sizing), intent(in) :: field_sizing_description
585  type(data_sizing_description_type), dimension(:), intent(inout) :: data_description
586 
587  data_description(index)%field_name=field_name
588  data_description(index)%dimensions=field_sizing_description%number_dimensions
589  data_description(index)%dim_sizes=field_sizing_description%dimensions
590  end subroutine assemble_individual_description
591 
598  subroutine register_with_io_server(current_state, mpi_type_definition_description, mpi_type_field_description)
599  type(model_state_type), target, intent(inout) :: current_state
600  integer, intent(in) :: mpi_type_definition_description, mpi_type_field_description
601 
602  type(definition_description_type), dimension(:), allocatable :: definition_descriptions
603  type(field_description_type), dimension(:), allocatable :: field_descriptions
604  integer :: number_defns, number_fields, status(mpi_status_size), ierr
605 
606  call mpi_send(register_command, 1, mpi_int, current_state%parallel%corresponding_io_server_process, &
607  command_tag, mpi_comm_world, ierr)
608 
609  call mpi_probe(current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, status, ierr)
610  call mpi_get_count(status, mpi_type_definition_description, number_defns, ierr)
611  allocate(definition_descriptions(number_defns))
612 
613  call mpi_recv(definition_descriptions, number_defns, mpi_type_definition_description, &
614  current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, mpi_status_ignore, ierr)
615  number_fields=get_total_number_of_fields(definition_descriptions, number_defns)
616 
617  allocate(field_descriptions(number_fields))
618  call mpi_recv(field_descriptions, number_fields, mpi_type_field_description, &
619  current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, mpi_status_ignore, ierr)
620  call populate_data_definition_configuration(definition_descriptions, number_defns, field_descriptions, number_fields)
621  deallocate(definition_descriptions)
622  end subroutine register_with_io_server
623 
628  integer function get_total_number_of_fields(definition_descriptions, number_defns)
629  type(definition_description_type), dimension(:), intent(inout) :: definition_descriptions
630  integer, intent(in) :: number_defns
631 
632  integer :: i
633 
635  do i=1, number_defns
636  get_total_number_of_fields=get_total_number_of_fields+definition_descriptions(i)%number_fields
637  end do
638  end function get_total_number_of_fields
639 
647  subroutine populate_data_definition_configuration(definition_descriptions, number_defns, field_descriptions, number_fields)
648  type(definition_description_type), dimension(:), intent(inout) :: definition_descriptions
649  type(field_description_type), dimension(:), intent(inout) :: field_descriptions
650  integer, intent(in) :: number_defns, number_fields
651 
652  integer :: i, definition_index, field_index
653 
654  allocate(data_definitions(number_defns))
655  do i=1, number_defns
656  data_definitions(i)%name=definition_descriptions(i)%definition_name
657  data_definitions(i)%send_on_terminate=definition_descriptions(i)%send_on_terminate
658  data_definitions(i)%number_of_data_fields=0 ! Will increment this for each field
659  data_definitions(i)%frequency=definition_descriptions(i)%frequency
660  allocate(data_definitions(i)%fields(definition_descriptions(i)%number_fields))
661  end do
662  do i=1, number_fields
663  definition_index=get_definition_index(field_descriptions(i)%definition_name)
664  field_index=data_definitions(definition_index)%number_of_data_fields+1
665  data_definitions(definition_index)%number_of_data_fields=field_index
666  data_definitions(definition_index)%fields(field_index)%name=field_descriptions(i)%field_name
667  data_definitions(definition_index)%fields(field_index)%field_type=field_descriptions(i)%field_type
668  data_definitions(definition_index)%fields(field_index)%data_type=field_descriptions(i)%data_type
669  data_definitions(definition_index)%fields(field_index)%optional=field_descriptions(i)%optional
670  if (field_descriptions(i)%optional .or. field_descriptions(i)%field_type == array_field_type .or. &
671  field_descriptions(i)%field_type == map_field_type) then
672  call c_put_integer(unique_field_names, field_descriptions(i)%field_name, 1)
673  end if
674  if (.not. field_descriptions(i)%optional) data_definitions(definition_index)%fields(field_index)%enabled=.true.
675  end do
677 
681  integer function get_definition_index(name)
682  character(len=*), intent(in) :: name
683 
684  integer :: i
685  do i=1, size(data_definitions)
686  if (data_definitions(i)%name .eq. name) then
688  return
689  end if
690  end do
692  end function get_definition_index
693 
698  subroutine pack_send_buffer(current_state, data_definition)
699  type(model_state_type), target, intent(inout) :: current_state
700  type(io_configuration_data_definition_type), intent(inout) :: data_definition
701 
702  integer :: current_buffer_point, i
703 
704  current_buffer_point=1
705  do i=1, data_definition%number_of_data_fields
706  if (data_definition%fields(i)%enabled) then
707  if (data_definition%fields(i)%field_type == array_field_type) then
708  current_buffer_point=pack_array_into_send_buffer(current_state, data_definition, data_definition%fields(i), &
709  current_buffer_point)
710  else if (data_definition%fields(i)%field_type == map_field_type) then
711  current_buffer_point=pack_map_into_send_buffer(current_state, data_definition, data_definition%fields(i), &
712  current_buffer_point)
713  else if (data_definition%fields(i)%field_type == scalar_field_type) then
714  current_buffer_point=pack_scalar_into_send_buffer(current_state, data_definition, data_definition%fields(i), &
715  current_buffer_point)
716  end if
717  end if
718  end do
719  end subroutine pack_send_buffer
720 
727  integer function pack_scalar_into_send_buffer(current_state, data_definition, field, current_buffer_point)
728  type(model_state_type), target, intent(inout) :: current_state
729  type(io_configuration_data_definition_type), intent(inout) :: data_definition
730  type(io_configuration_field_type), intent(in) :: field
731  integer, intent(in) :: current_buffer_point
732 
733  if (field%name .eq. "timestep") then
734  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
735  int_value=current_state%timestep)
736  else if (field%name .eq. "terminated") then
737  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
738  logical_value=.not. current_state%continue_timestep .and. in_finalisation_callback)
739  else if (field%name .eq. "z_size") then
740  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
741  int_value=current_state%global_grid%size(z_index))
742  else if (field%name .eq. "y_size") then
743  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
744  int_value=current_state%global_grid%size(y_index))
745  else if (field%name .eq. "y_bottom") then
746  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
747  real_value=current_state%global_grid%bottom(y_index))
748  else if (field%name .eq. "y_top") then
749  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
750  real_value=current_state%global_grid%top(y_index))
751  else if (field%name .eq. "y_resolution") then
752  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
753  real_value=current_state%global_grid%resolution(y_index))
754  else if (field%name .eq. "x_size") then
755  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
756  int_value=current_state%global_grid%size(x_index))
757  else if (field%name .eq. "x_bottom") then
758  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
759  real_value=current_state%global_grid%bottom(x_index))
760  else if (field%name .eq. "x_top") then
761  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
762  real_value=current_state%global_grid%top(x_index))
763  else if (field%name .eq. "x_resolution") then
764  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
765  real_value=current_state%global_grid%resolution(x_index))
766  else if (field%name .eq. "time") then
767  ! The time is incremented with dtm as the model was about to increment for the next step and this is needed for diagnostics
768  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
769  real_value=current_state%time+current_state%dtm)
770  else if (field%name .eq. "ugal") then
771  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
772  real_value=current_state%ugal)
773  else if (field%name .eq. "vgal") then
774  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
775  real_value=current_state%vgal)
776  else if (field%name .eq. "nqfields") then
777  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
778  int_value=current_state%number_q_fields)
779  else if (field%name .eq. "dtm") then
780  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
781  real_value=current_state%dtm)
782  else if (field%name .eq. "dtm_new") then
783  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
784  real_value=current_state%dtm_new)
785  else if (field%name .eq. "absolute_new_dtm") then
786  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
787  real_value=current_state%absolute_new_dtm)
788  else
789  ! Handle component field here
791  data_definition, field, current_buffer_point)
792  end if
793  end function pack_scalar_into_send_buffer
794 
802  integer function handle_component_field_scalar_packing_into_send_buffer(current_state, data_definition, &
803  field, current_buffer_point)
804  type(model_state_type), target, intent(inout) :: current_state
805  type(io_configuration_data_definition_type), intent(inout) :: data_definition
806  type(io_configuration_field_type), intent(in) :: field
807  integer, intent(in) :: current_buffer_point
808 
809  type(component_field_information_type) :: field_descriptor
810  type(component_field_value_type) :: published_value
811 
812  field_descriptor=get_component_field_descriptor(field%name)
813  published_value=get_component_field_value(current_state, field%name)
814  if (field_descriptor%data_type == component_double_data_type) then
815  handle_component_field_scalar_packing_into_send_buffer=pack_scalar_field(data_definition%send_buffer, &
816  current_buffer_point, real_value=published_value%scalar_real)
817  else if (field_descriptor%data_type == component_integer_data_type) then
818  handle_component_field_scalar_packing_into_send_buffer=pack_scalar_field(data_definition%send_buffer, &
819  current_buffer_point, int_value=published_value%scalar_int)
820  end if
822 
829  integer function pack_map_into_send_buffer(current_state, data_definition, field, current_buffer_point)
830  type(model_state_type), target, intent(inout) :: current_state
831  type(io_configuration_data_definition_type), intent(inout) :: data_definition
832  type(io_configuration_field_type), intent(in) :: field
833  integer, intent(in) :: current_buffer_point
834 
835  integer :: i
836  type(q_metadata_type) :: specific_q_data
837  type(hashmap_type) :: q_indicies_map
838 
839  if (field%name .eq. "options_database") then
840  pack_map_into_send_buffer=pack_map_field(data_definition%send_buffer, current_buffer_point, current_state%options_database)
841  else if (field%name .eq. "q_indicies") then
842  do i=1, get_max_number_q_indices()
843  specific_q_data=get_indices_descriptor(i)
844  if (specific_q_data%l_used) then
845  call c_put_integer(q_indicies_map, specific_q_data%name, i)
846  end if
847  end do
848  pack_map_into_send_buffer=pack_map_field(data_definition%send_buffer, current_buffer_point, q_indicies_map)
849  call c_free(q_indicies_map)
850  end if
851  end function pack_map_into_send_buffer
852 
859  integer function pack_array_into_send_buffer(current_state, data_definition, field, current_buffer_point)
860  type(model_state_type), target, intent(inout) :: current_state
861  type(io_configuration_data_definition_type), intent(inout) :: data_definition
862  type(io_configuration_field_type), intent(in) :: field
863  integer, intent(in) :: current_buffer_point
864 
865  if (field%name .eq. "local_grid_size") then
866  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
867  int_array=current_state%local_grid%size)
868  else if (field%name .eq. "local_grid_start") then
869  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
870  int_array=current_state%local_grid%start)
871  else if (field%name .eq. "z") then
872  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
873  real_array_1d=current_state%global_grid%configuration%vertical%z)
874  else if (field%name .eq. "olubar") then
875  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
876  real_array_1d=current_state%global_grid%configuration%vertical%olubar)
877  else if (field%name .eq. "olzubar") then
878  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
879  real_array_1d=current_state%global_grid%configuration%vertical%olzubar)
880  else if (field%name .eq. "olvbar") then
881  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
882  real_array_1d=current_state%global_grid%configuration%vertical%olvbar)
883  else if (field%name .eq. "olzvbar") then
884  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
885  real_array_1d=current_state%global_grid%configuration%vertical%olzvbar)
886  else if (field%name .eq. "olthbar") then
887  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
888  real_array_1d=current_state%global_grid%configuration%vertical%olthbar)
889  else if (field%name .eq. "olzthbar") then
890  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
891  real_array_1d=current_state%global_grid%configuration%vertical%olzthbar)
892  else if (field%name .eq. "olqbar") then
893  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
894  real_array_2d=current_state%global_grid%configuration%vertical%olqbar)
895  else if (field%name .eq. "olzqbar") then
896  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
897  real_array_2d=current_state%global_grid%configuration%vertical%olzqbar)
898  else if (field%name .eq. "thref") then
899  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
900  real_array_1d=current_state%global_grid%configuration%vertical%thref)
901  else if (field%name .eq. "prefn") then
902  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
903  real_array_1d=current_state%global_grid%configuration%vertical%prefn)
904  else if (field%name .eq. "rhon") then
905  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
906  real_array_1d=current_state%global_grid%configuration%vertical%rhon)
907  else if (field%name .eq. "rho") then
908  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
909  real_array_1d=current_state%global_grid%configuration%vertical%rho)
910  else if (field%name .eq. "u") then
911  current_state%u%data=current_state%u%data+current_state%ugal
912  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%u, &
913  current_buffer_point, current_state%local_grid)
914  current_state%u%data=current_state%u%data-current_state%ugal
915  else if (field%name .eq. "u_nogal") then
916  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%u, current_buffer_point, &
917  current_state%local_grid)
918  else if (field%name .eq. "zu") then
919  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zu, current_buffer_point,&
920  current_state%local_grid)
921  else if (field%name .eq. "v") then
922  current_state%v%data=current_state%v%data+current_state%vgal
923  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%v, current_buffer_point, &
924  current_state%local_grid)
925  current_state%v%data=current_state%v%data-current_state%vgal
926  else if (field%name .eq. "v_nogal") then
927  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%v, current_buffer_point, &
928  current_state%local_grid)
929  else if (field%name .eq. "zv") then
930  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zv, current_buffer_point,&
931  current_state%local_grid)
932  else if (field%name .eq. "w") then
933  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%w, current_buffer_point, &
934  current_state%local_grid)
935  else if (field%name .eq. "zw") then
936  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zw, current_buffer_point,&
937  current_state%local_grid)
938  else if (field%name .eq. "q") then
939  pack_array_into_send_buffer=pack_q_fields(data_definition%send_buffer, current_state%q, current_state%number_q_fields, &
940  current_buffer_point, current_state%local_grid)
941  else if (field%name .eq. "zq") then
942  pack_array_into_send_buffer=pack_q_fields(data_definition%send_buffer, current_state%zq, current_state%number_q_fields, &
943  current_buffer_point, current_state%local_grid)
944  else if (field%name .eq. "th") then
945  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%th, current_buffer_point,&
946  current_state%local_grid)
947  else if (field%name .eq. "zth") then
948  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zth, current_buffer_point,&
949  current_state%local_grid)
950  else if (field%name .eq. "p") then
951  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%p, current_buffer_point, &
952  current_state%local_grid)
953  else
954  ! Handle component field here
956  data_definition, field, current_buffer_point)
957  end if
958  end function pack_array_into_send_buffer
959 
967  integer function handle_component_field_array_packing_into_send_buffer(current_state, data_definition, &
968  field, current_buffer_point)
969  type(model_state_type), target, intent(inout) :: current_state
970  type(io_configuration_data_definition_type), intent(inout) :: data_definition
971  type(io_configuration_field_type), intent(in) :: field
972  integer, intent(in) :: current_buffer_point
973 
974  type(component_field_information_type) :: field_descriptor
975  type(component_field_value_type) :: published_value
976 
977  field_descriptor=get_component_field_descriptor(field%name)
978  published_value=get_component_field_value(current_state, field%name)
979  if (field_descriptor%data_type == component_double_data_type) then
980  if (field_descriptor%number_dimensions == 1) then
981  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
982  current_buffer_point, real_array_1d=published_value%real_1d_array)
983  deallocate(published_value%real_1d_array)
984  else if (field_descriptor%number_dimensions == 2) then
985  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
986  current_buffer_point, real_array_2d=published_value%real_2d_array)
987  deallocate(published_value%real_2d_array)
988  else if (field_descriptor%number_dimensions == 3) then
989  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
990  current_buffer_point, real_array_3d=published_value%real_3d_array)
991  deallocate(published_value%real_3d_array)
992  else if (field_descriptor%number_dimensions == 4) then
993  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
994  current_buffer_point, real_array_4d=published_value%real_4d_array)
995  deallocate(published_value%real_4d_array)
996  end if
997  end if
999 
1006  integer function pack_prognostic_flow_field(buffer, prognostic, start_offset, local_grid)
1007  character, dimension(:), allocatable, intent(inout) :: buffer
1008  type(prognostic_field_type), intent(inout) :: prognostic
1009  integer, intent(in) :: start_offset
1010  type(local_grid_type), intent(inout) :: local_grid
1011 
1012  integer :: target_end
1013 
1014  target_end=start_offset + (local_grid%size(z_index)*local_grid%size(y_index)*local_grid%size(x_index)*kind(prognostic%data)-1)
1015 
1016  buffer(start_offset : target_end) = transfer(prognostic%data(&
1017  local_grid%local_domain_start_index(z_index): local_grid%local_domain_end_index(z_index),&
1018  local_grid%local_domain_start_index(y_index): local_grid%local_domain_end_index(y_index), &
1019  local_grid%local_domain_start_index(x_index): local_grid%local_domain_end_index(x_index)), &
1020  buffer(start_offset : target_end))
1021  pack_prognostic_flow_field=target_end+1
1022  end function pack_prognostic_flow_field
1023 
1031  integer function pack_q_fields(buffer, q_fields, number_q_fields, start_offset, local_grid)
1032  character, dimension(:), allocatable, intent(inout) :: buffer
1033  type(prognostic_field_type), dimension(:), intent(inout) :: q_fields
1034  integer, intent(in) :: start_offset, number_q_fields
1035  type(local_grid_type), intent(inout) :: local_grid
1036 
1037  integer :: target_end, i, current_starting_index
1038 
1039  current_starting_index=start_offset
1040 
1041  do i=1,number_q_fields
1042  target_end=current_starting_index + (local_grid%size(z_index)*local_grid%size(y_index)*&
1043  local_grid%size(x_index)*kind(q_fields(i)%data)-1)
1044  buffer(current_starting_index : target_end) = transfer(q_fields(i)%data(&
1045  local_grid%local_domain_start_index(z_index): local_grid%local_domain_end_index(z_index),&
1046  local_grid%local_domain_start_index(y_index): local_grid%local_domain_end_index(y_index), &
1047  local_grid%local_domain_start_index(x_index): local_grid%local_domain_end_index(x_index)), &
1048  buffer(current_starting_index : target_end))
1049  current_starting_index=target_end+1
1050  end do
1051  pack_q_fields=target_end+1
1052  end function pack_q_fields
1053 end module iobridge_mod
integer function get_total_number_of_fields(definition_descriptions, number_defns)
Retrieve the total number of fields, which is all the fields in all the data definitions.
Definition: iobridge.F90:629
logical io_server_enabled
Definition: iobridge.F90:53
Retrieves the key currently being held at a specific index in the map or "" if the index > map elemen...
Gets a specific generic element out of the list, stack, queue or map with the corresponding key...
Wrapper type for the value returned for a published field from a component.
Puts an integer key-value pair into the map.
integer function handle_component_field_array_packing_into_send_buffer(current_state, data_definition, field, current_buffer_point)
Packs a components field array into the send buffer, these are fields that are served up by component...
Definition: iobridge.F90:969
integer, parameter, public float_data_type
Definition: ioclient.F90:40
integer function send_general_monc_information_to_server(current_state, buffer)
Sends the general MONC information (ZN field and Q field names) to the IO server. ...
Definition: iobridge.F90:485
integer function, public get_mpi_datatype_from_internal_representation(type_code)
Gets the MPI datatype from out internal representation of the field data type (as in the configuratio...
Definition: ioclient.F90:203
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
integer, parameter, public array_field_type
Definition: ioclient.F90:38
integer function, public get_max_number_q_indices()
Gets the maximum number of Q indicies.
Definition: q_indices.F90:81
Contains prognostic field definitions and functions.
Definition: prognostics.F90:2
integer, parameter, public boolean_data_type
Definition: ioclient.F90:40
character(len=string_length), parameter, public local_end_points_key
Definition: ioclient.F90:43
integer, parameter, public register_command
Definition: ioclient.F90:34
A prognostic field which is assumed to be 3D.
Definition: prognostics.F90:13
Logging utility.
Definition: logging.F90:2
integer function send_data_field_sizes_to_server(current_state, mpi_type_data_sizing_description, data_description, number_unique_fields)
Assembles all the data field sizing information and sends this to the IO server.
Definition: iobridge.F90:458
Bridge between MONC and the IO server, this registers the current MONC process, will issue data dumps...
Definition: iobridge.F90:3
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
integer function build_mpi_data_type_for_definition(specific_data_definition)
Builds the MPI data type for a specific definition with sizing information.
Definition: iobridge.F90:183
integer, parameter, public command_tag
Definition: ioclient.F90:34
integer function pack_prognostic_flow_field(buffer, prognostic, start_offset, local_grid)
Packs the data of a specific prognostic field into a buffer.
Definition: iobridge.F90:1007
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
type(map_type) component_field_descriptions
Definition: iobridge.F90:52
The ModelState which represents the current state of a run.
Definition: state.F90:39
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
A hashmap structure, the same as a map but uses hashing for greatly improved performance when storing...
Definition: collections.F90:94
Conversion between common inbuilt FORTRAN data types.
Definition: conversions.F90:5
integer, parameter, public double_data_type
Definition: ioclient.F90:40
character(len=string_length), parameter, public local_sizes_key
Definition: ioclient.F90:43
Converts data types to strings.
Definition: conversions.F90:36
integer, parameter, public single_precision
Single precision (32 bit) kind.
Definition: datadefn.F90:13
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
subroutine timestep_callback(current_state)
Model dump call back, called for each model dump phase.
Definition: iobridge.F90:102
type(map_type) sendable_fields
Definition: iobridge.F90:52
type(component_field_value_type) function, public get_component_field_value(current_state, name)
Retrieves the value wrapper of a components published field.
Definition: registry.F90:144
type(q_metadata_type) function, public get_indices_descriptor(i)
Retrieves the indicies descriptor at a specific location.
Definition: q_indices.F90:100
integer, parameter, public double_precision
Double precision (64 bit) kind.
Definition: datadefn.F90:14
subroutine package_local_monc_decomposition_into_descriptions(current_state, data_description)
Packages the local MONC decomposition information into descriptions for communication.
Definition: iobridge.F90:514
integer, parameter, public string_data_type
Definition: ioclient.F90:40
subroutine build_mpi_data_types()
Builds the MPI data types that correspond to the field descriptions and sizings.
Definition: iobridge.F90:170
This defines some constants and procedures that are useful to the IO server and clients that call it...
Definition: ioclient.F90:3
integer, parameter, public integer_data_type
Definition: ioclient.F90:40
Map data structure that holds string (length 20 maximum) key value pairs.
Definition: collections.F90:86
subroutine init_callback(current_state)
Initialisation call back, called at the start of the model run.
Definition: iobridge.F90:70
integer function, public pack_map_field(buffer, start_offset, map_to_pack)
Packs a map into the send buffer.
Definition: ioclient.F90:224
This manages the Q variables and specifically the mapping between names and the index that they are s...
Definition: q_indices.F90:2
Returns the number of elements in the collection.
type(component_field_information_type) function, public get_component_field_information(current_state, name)
Retrieves information about a components published field which includes its type and size...
Definition: registry.F90:165
subroutine pack_send_buffer(current_state, data_definition)
Packs the current state into the send buffer. This iterates through each field in the data descriptio...
Definition: iobridge.F90:699
integer function, public get_number_active_q_indices()
Gets the number of active Q indicies (i.e. those allocated to specific uses)
Definition: q_indices.F90:87
Interfaces and types that MONC components must specify.
integer function pack_scalar_into_send_buffer(current_state, data_definition, field, current_buffer_point)
Packs scalar fields into the send bufer.
Definition: iobridge.F90:728
subroutine register_with_io_server(current_state, mpi_type_definition_description, mpi_type_field_description)
Registers this MONC with the corresponding IO server. This will encapsulate the entire protocol...
Definition: iobridge.F90:599
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
type(component_descriptor_type) function, public iobridge_get_descriptor()
Definition: iobridge.F90:60
integer, parameter, public data_command_start
Definition: ioclient.F90:34
integer function pack_q_fields(buffer, q_fields, number_q_fields, start_offset, local_grid)
Packs the Q fields into the send buffer.
Definition: iobridge.F90:1032
integer, parameter, public component_integer_data_type
integer, parameter, public log_warn
Log WARNING and ERROR messages.
Definition: logging.F90:12
type(map_type) unique_field_names
Definition: iobridge.F90:52
integer function handle_component_field_scalar_packing_into_send_buffer(current_state, data_definition, field, current_buffer_point)
Packs a components field scalar into the send buffer, these are fields that are served up by componen...
Definition: iobridge.F90:804
integer, parameter, public string_length
Default length of strings.
Definition: datadefn.F90:10
integer function, dimension(5), public populate_mpi_type_extents()
Provides the type extents of the types that we are using in construction of the MPI data type...
Definition: ioclient.F90:60
subroutine populate_globally_visible_sendable_fields(current_state)
Populates the globally visible sendable fields which is a key value pair mapping between name and des...
Definition: iobridge.F90:304
subroutine populate_component_public_field(current_state, field_name)
Populates the field information for a specific publically available field offered by one of the compo...
Definition: iobridge.F90:278
subroutine populate_sendable_fields(current_state)
Populates the sendable field definitions with the field sizing information.
Definition: iobridge.F90:263
List data structure which implements a doubly linked list. This list will preserve its order...
Definition: collections.F90:60
class(*) function, pointer generate_sendable_description(dim1, dim2, dim3, dim4)
Generates a sendable description based upon the dimension information supplied, missing arguments mea...
Definition: iobridge.F90:397
subroutine finalisation_callback(current_state)
Finalisation call back, called at the end of the model run.
Definition: iobridge.F90:146
subroutine send_monc_specific_data_to_server(current_state, mpi_type_data_sizing_description)
Sends this MONC specific information to the IO server, which is field info (sizing & availability) as...
Definition: iobridge.F90:431
integer function, public options_size(options_database)
Returns the number of entries in the options database.
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
type(list_type) function, public get_all_component_published_fields()
Retrieves all of the published field information.
Definition: registry.F90:185
integer function, public build_mpi_type_definition_description()
Builds the MPI data type for sending data descriptions to registree MONCs.
Definition: ioclient.F90:76
subroutine, public append_mpi_datatype(field_start, field_end, field_array_sizes, data_type, type_extents, prev_data_type, type_index, old_types, offsets, block_counts)
Appends the MPI datatype details to the block counts, old types and offsets arrays. This will lump together multiple concurrent fields with the same type.
Definition: ioclient.F90:186
integer function pack_map_into_send_buffer(current_state, data_definition, field, current_buffer_point)
Packs map fields into the send buffer.
Definition: iobridge.F90:830
integer function, public options_get_integer(options_database, key, index)
Retrieves an integer value from the database that matches the provided key.
Manages the options database. Contains administration functions and deduce runtime options from the c...
integer function, public pack_array_field(buffer, start_offset, int_array, real_array_1d, real_array_2d, real_array_3d, real_array_4d)
Packs an array field into the sending buffer.
Definition: ioclient.F90:273
logical function, public options_get_logical(options_database, key, index)
Retrieves a logical value from the database that matches the provided key.
Frees up all the allocatable, heap, memory associated with a list, stack, queue or map...
Puts a generic key-value pair into the map.
integer function pack_array_into_send_buffer(current_state, data_definition, field, current_buffer_point)
Packs array fields into the send bufer.
Definition: iobridge.F90:860
type(component_field_information_type) function get_component_field_descriptor(field_name)
Retrieves the descriptor associated with some component's field based upon the field name...
Definition: iobridge.F90:560
character(len=string_length), parameter, public number_q_indicies_key
Definition: ioclient.F90:43
subroutine populate_data_definition_configuration(definition_descriptions, number_defns, field_descriptions, number_fields)
Based upon the received data and field definitions this will configure the IO bridge internal represe...
Definition: iobridge.F90:648
Determines whether or not a map contains a specific key.
integer, parameter, public scalar_field_type
Definition: ioclient.F90:38
type(io_configuration_data_definition_type), dimension(:), allocatable data_definitions
Definition: iobridge.F90:51
type(io_server_sendable_field_sizing) function get_sendable_field_sizing(field_name, field_found)
Retrieves the sizing information associated with a specific field.
Definition: iobridge.F90:541
integer function, public build_mpi_type_data_sizing_description()
Builds the MPI type used for sending to the IO server a description of the data, namely the size of t...
Definition: ioclient.F90:147
integer function, public build_mpi_type_field_description()
Builds the MPI data type for sending field descriptions to registree MONCs.
Definition: ioclient.F90:108
The model state which represents the current state of a run.
Definition: state.F90:2
integer, parameter, public map_field_type
Field data type identifiers.
Definition: ioclient.F90:38
logical in_finalisation_callback
Definition: iobridge.F90:53
integer, parameter, public y_index
Definition: grids.F90:14
integer, parameter, public data_tag
Definition: ioclient.F90:34
Gets a specific string element out of the list, stack, queue or map with the corresponding key...
subroutine assemble_individual_description(data_description, index, field_name, field_sizing_description)
Will assemble an individual description of an array data field.
Definition: iobridge.F90:582
integer function get_definition_index(name)
Looks up a specific definition based upon its name and returns the index.
Definition: iobridge.F90:682
integer, parameter, public x_index
Definition: grids.F90:14
subroutine send_data_to_io_server(current_state, data_index)
Sends data to the IO server.
Definition: iobridge.F90:121
integer, parameter, public deregister_command
Definition: ioclient.F90:34
integer function, public pack_scalar_field(buffer, start_offset, int_value, real_value, single_real_value, double_real_value, string_value, logical_value)
Packs the data of a scalar field into a buffer.
Definition: ioclient.F90:312
MONC component registry.
Definition: registry.F90:5
character(len=string_length), parameter, public local_start_points_key
Definition: ioclient.F90:43
integer, parameter, public component_double_data_type