MONC
Data Types | Functions/Subroutines | Variables
iobridge_mod Module Reference

Bridge between MONC and the IO server, this registers the current MONC process, will issue data dumps and deregister MONCs when the model run is completed. More...

Data Types

type  io_configuration_data_definition_type
 
type  io_configuration_field_type
 
type  io_server_sendable_field_sizing
 

Functions/Subroutines

type(component_descriptor_type) function, public iobridge_get_descriptor ()
 
subroutine init_callback (current_state)
 Initialisation call back, called at the start of the model run. More...
 
subroutine timestep_callback (current_state)
 Model dump call back, called for each model dump phase. More...
 
subroutine send_data_to_io_server (current_state, data_index)
 Sends data to the IO server. More...
 
subroutine finalisation_callback (current_state)
 Finalisation call back, called at the end of the model run. More...
 
subroutine build_mpi_data_types ()
 Builds the MPI data types that correspond to the field descriptions and sizings. More...
 
integer function build_mpi_data_type_for_definition (specific_data_definition)
 Builds the MPI data type for a specific definition with sizing information. More...
 
subroutine populate_sendable_fields (current_state)
 Populates the sendable field definitions with the field sizing information. More...
 
subroutine populate_component_public_field (current_state, field_name)
 Populates the field information for a specific publically available field offered by one of the components. More...
 
subroutine populate_globally_visible_sendable_fields (current_state)
 Populates the globally visible sendable fields which is a key value pair mapping between name and description of that field. More...
 
class(*) function, pointer generate_sendable_description (dim1, dim2, dim3, dim4)
 Generates a sendable description based upon the dimension information supplied, missing arguments means that dimension does not exist. More...
 
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 well as meta data such as ZN field and Q field names. More...
 
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. More...
 
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. More...
 
subroutine package_local_monc_decomposition_into_descriptions (current_state, data_description)
 Packages the local MONC decomposition information into descriptions for communication. More...
 
type(io_server_sendable_field_sizing) function get_sendable_field_sizing (field_name, field_found)
 Retrieves the sizing information associated with a specific field. More...
 
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. More...
 
subroutine assemble_individual_description (data_description, index, field_name, field_sizing_description)
 Will assemble an individual description of an array data field. More...
 
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, which is sending the registration command, receiving the data and field definitions from the IO server and then sending back the sizing for the fields that this MONC will contribute. More...
 
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. More...
 
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 representation of these facets, which is a structured tree, data defintions holding their own fields rather than the unstructured data we get from the IO server. More...
 
integer function get_definition_index (name)
 Looks up a specific definition based upon its name and returns the index. More...
 
subroutine pack_send_buffer (current_state, data_definition)
 Packs the current state into the send buffer. This iterates through each field in the data description and adds it to the send buffer. More...
 
integer function pack_scalar_into_send_buffer (current_state, data_definition, field, current_buffer_point)
 Packs scalar fields into the send bufer. More...
 
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 components rather than explicitly available. More...
 
integer function pack_map_into_send_buffer (current_state, data_definition, field, current_buffer_point)
 Packs map fields into the send buffer. More...
 
integer function pack_array_into_send_buffer (current_state, data_definition, field, current_buffer_point)
 Packs array fields into the send bufer. More...
 
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 components rather than explicitly available. More...
 
integer function pack_prognostic_flow_field (buffer, prognostic, start_offset, local_grid)
 Packs the data of a specific prognostic field into a buffer. More...
 
integer function pack_q_fields (buffer, q_fields, number_q_fields, start_offset, local_grid)
 Packs the Q fields into the send buffer. More...
 

Variables

type(io_configuration_data_definition_type), dimension(:), allocatable data_definitions
 
type(map_typeunique_field_names
 
type(map_typesendable_fields
 
type(map_typecomponent_field_descriptions
 
logical io_server_enabled
 
logical in_finalisation_callback
 

Detailed Description

Bridge between MONC and the IO server, this registers the current MONC process, will issue data dumps and deregister MONCs when the model run is completed.

Function/Subroutine Documentation

◆ assemble_individual_description()

subroutine iobridge_mod::assemble_individual_description ( type(data_sizing_description_type), dimension(:), intent(inout)  data_description,
integer, intent(in)  index,
character(len=*), intent(in)  field_name,
type(io_server_sendable_field_sizing), intent(in)  field_sizing_description 
)
private

Will assemble an individual description of an array data field.

Parameters
data_descriptionThe data structure used to describe the size of arrays
indexThe index of this current field
field_nameThe corresponding field name that we are describing
dimensionsThe number of dimensions (zero means the field is inactive)
dim1Optional size of dimension 1
dim2Optional size of dimension 2
dim3Optional size of dimension 3
dim4Optional size of dimension 4

Definition at line 629 of file iobridge.F90.

630  integer, intent(in) :: index
631  character(len=*), intent(in) :: field_name
632  type(io_server_sendable_field_sizing), intent(in) :: field_sizing_description
633  type(data_sizing_description_type), dimension(:), intent(inout) :: data_description
634 
635  data_description(index)%field_name=field_name
636  data_description(index)%dimensions=field_sizing_description%number_dimensions
637  data_description(index)%dim_sizes=field_sizing_description%dimensions
Here is the caller graph for this function:

◆ build_mpi_data_type_for_definition()

integer function iobridge_mod::build_mpi_data_type_for_definition ( type(io_configuration_data_definition_type), intent(inout)  specific_data_definition)
private

Builds the MPI data type for a specific definition with sizing information.

Parameters
specific_data_definitionThe data definition to build the type for
Returns
The size (in bytes) that the send buffer needs to be to store the data for the MPI operation

Definition at line 187 of file iobridge.F90.

188  type(io_configuration_data_definition_type), intent(inout) :: specific_data_definition
189 
190  integer :: type_extents(5), type_counts, i, j, tempsize, field_start, data_type, field_array_sizes, &
191  temp_size, prev_data_type, old_types(20), offsets(20), block_counts(20), ierr, field_ignores
192  logical :: field_found
193  type(io_server_sendable_field_sizing) :: field_size_info
194 
195  type_extents=populate_mpi_type_extents()
196 
197  field_start=1
198  data_type=0
199  type_counts=0
200  field_array_sizes=0
201  field_ignores=0
202  do i=1, specific_data_definition%number_of_data_fields
203  if (data_type == 0) then
204  prev_data_type=data_type
205  data_type=specific_data_definition%fields(i)%data_type
206  else
207  if (data_type .ne. specific_data_definition%fields(i)%data_type) then
208  ! For efficient type packing, combine multiple fields with the same type - therefore when the type changes work the previous one pack
209  call append_mpi_datatype(field_start, i-1-field_ignores, field_array_sizes, data_type, &
210  type_extents, prev_data_type, type_counts+1, old_types, offsets, block_counts)
211  field_start=i
212  field_array_sizes=0
213  field_ignores=0
214  prev_data_type=data_type
215  data_type=specific_data_definition%fields(i)%data_type
216  type_counts=type_counts+1
217  end if
218  end if
219 
220  if (specific_data_definition%fields(i)%field_type .eq. array_field_type .or. &
221  specific_data_definition%fields(i)%field_type .eq. map_field_type) then
222  ! Grab the field info based upon the field name
223  field_size_info=get_sendable_field_sizing(specific_data_definition%fields(i)%name, field_found)
224  specific_data_definition%fields(i)%enabled=field_found
225  if (.not. field_found .or. field_size_info%number_dimensions == 0) then
226  ! If no field info, or the dimension is 0 then this MONC process is not sending that field - check it is optional
227  if (.not. specific_data_definition%fields(i)%optional) then
228  call log_log(log_error, "Non optional field `"//trim(specific_data_definition%fields(i)%name)//&
229  "' omitted from MONC IO server registration")
230  end if
231  field_ignores=field_ignores+1
232  else
233  ! If the field is specified then use the size data to assemble the field size and append to current size
234  temp_size=1
235  do j=1, field_size_info%number_dimensions
236  temp_size=temp_size*field_size_info%dimensions(j)
237  end do
238  if (specific_data_definition%fields(i)%field_type .eq. map_field_type) then
239  field_array_sizes=(field_array_sizes+temp_size*string_length*2)-1
240  else
241  field_array_sizes=(field_array_sizes+temp_size)-1
242  end if
243  end if
244  else
245  if (specific_data_definition%fields(i)%optional) then
246  field_size_info=get_sendable_field_sizing(specific_data_definition%fields(i)%name, field_found)
247  specific_data_definition%fields(i)%enabled=field_found
248  if (.not. field_found) field_ignores=field_ignores+1
249  end if
250  end if
251  end do
252  if (field_start .le. i-1) then
253  ! If there are outstanding fields to process then we do this here
254  call append_mpi_datatype(field_start, i-1, field_array_sizes, data_type, &
255  type_extents, prev_data_type, type_counts+1, old_types, offsets, block_counts)
256  type_counts=type_counts+1
257  end if
258 
259  call mpi_type_struct(type_counts, block_counts, offsets, old_types, specific_data_definition%mpi_datatype, ierr)
260  call mpi_type_commit(specific_data_definition%mpi_datatype, ierr)
261  call mpi_type_size(specific_data_definition%mpi_datatype, tempsize, ierr)
262  build_mpi_data_type_for_definition=tempsize
Here is the call graph for this function:
Here is the caller graph for this function:

◆ build_mpi_data_types()

subroutine iobridge_mod::build_mpi_data_types
private

Builds the MPI data types that correspond to the field descriptions and sizings.

Definition at line 174 of file iobridge.F90.

175  integer :: i, dump_send_buffer_size
176 
177  do i=1, size(data_definitions)
178  dump_send_buffer_size=build_mpi_data_type_for_definition(data_definitions(i))
179  data_definitions(i)%dump_requests=(/mpi_request_null, mpi_request_null/)
180  allocate(data_definitions(i)%send_buffer(dump_send_buffer_size))
181  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

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

Finalisation call back, called at the end of the model run.

Parameters
current_stateThe current model state

Definition at line 150 of file iobridge.F90.

151  type(model_state_type), target, intent(inout) :: current_state
152 
153  integer :: ierr, i
154 
155  if (.not. io_server_enabled) return
156  in_finalisation_callback=.true.
157 
158  do i=1, size(data_definitions)
159  if (data_definitions(i)%send_on_terminate) then
160  call send_data_to_io_server(current_state, i)
161  end if
162  if (data_definitions(i)%dump_requests(1) .ne. mpi_request_null .or. &
163  data_definitions(i)%dump_requests(2) .ne. mpi_request_null) then
164  call mpi_waitall(2, data_definitions(i)%dump_requests, mpi_statuses_ignore, ierr)
165  end if
166  if (allocated(data_definitions(i)%send_buffer)) deallocate(data_definitions(i)%send_buffer)
167  call mpi_type_free(data_definitions(i)%mpi_datatype, ierr)
168  end do
169  call mpi_send(deregister_command, 1, mpi_int, current_state%parallel%corresponding_io_server_process, &
170  command_tag, mpi_comm_world, ierr)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ generate_sendable_description()

class(*) function, pointer iobridge_mod::generate_sendable_description ( integer, intent(in), optional  dim1,
integer, intent(in), optional  dim2,
integer, intent(in), optional  dim3,
integer, intent(in), optional  dim4 
)
private

Generates a sendable description based upon the dimension information supplied, missing arguments means that dimension does not exist.

Parameters
dim1Optional size of dimension one
dim2Optional size of dimension two
dim3Optional size of dimension three
dim4Optional size of dimension four
Returns
The corresponding sendable description for the field, with the number of dimensions and sizing of each of these

Definition at line 420 of file iobridge.F90.

421  integer, intent(in), optional :: dim1, dim2, dim3, dim4
422  class(*), pointer :: generate_sendable_description
423 
424  type(io_server_sendable_field_sizing), pointer :: field
425  integer :: number_dimensions
426 
427  allocate(field)
428  number_dimensions=0
429  field%dimensions=0
430  if (present(dim1)) then
431  field%dimensions(1)=dim1
432  number_dimensions=number_dimensions+1
433  end if
434  if (present(dim2)) then
435  field%dimensions(2)=dim2
436  number_dimensions=number_dimensions+1
437  end if
438  if (present(dim3)) then
439  field%dimensions(3)=dim3
440  number_dimensions=number_dimensions+1
441  end if
442  if (present(dim4)) then
443  field%dimensions(4)=dim4
444  number_dimensions=number_dimensions+1
445  end if
446  field%number_dimensions=number_dimensions
447  generate_sendable_description=>field
Here is the caller graph for this function:

◆ get_component_field_descriptor()

type(component_field_information_type) function iobridge_mod::get_component_field_descriptor ( character(len=*), intent(in)  field_name)
private

Retrieves the descriptor associated with some component's field based upon the field name.

Parameters
field_nameThe field name

Definition at line 607 of file iobridge.F90.

608  character(len=*), intent(in) :: field_name
609 
610  class(*), pointer :: generic
611  if (c_contains(component_field_descriptions, field_name)) then
612  generic=>c_get_generic(component_field_descriptions, field_name)
613  select type(generic)
614  type is (component_field_information_type)
615  get_component_field_descriptor=generic
616  end select
617  end if
Here is the caller graph for this function:

◆ get_definition_index()

integer function iobridge_mod::get_definition_index ( character(len=*), intent(in)  name)
private

Looks up a specific definition based upon its name and returns the index.

Parameters
nameThe defintion name to search for
Returns
The index where the corresponding definition can be found or -1 if no definition is located with that name

Definition at line 729 of file iobridge.F90.

730  character(len=*), intent(in) :: name
731 
732  integer :: i
733  do i=1, size(data_definitions)
734  if (data_definitions(i)%name .eq. name) then
735  get_definition_index=i
736  return
737  end if
738  end do
739  get_definition_index=-1
Here is the caller graph for this function:

◆ get_sendable_field_sizing()

type(io_server_sendable_field_sizing) function iobridge_mod::get_sendable_field_sizing ( character(len=*), intent(in)  field_name,
logical, intent(out), optional  field_found 
)
private

Retrieves the sizing information associated with a specific field.

Parameters
field_nameThe field name to look up
field_foundOptional flag depicting whether the field was found or not

Definition at line 588 of file iobridge.F90.

589  character(len=*), intent(in) :: field_name
590  logical, intent(out), optional :: field_found
591 
592  class(*), pointer :: generic
593 
594  if (present(field_found)) field_found=.false.
595  if (c_contains(sendable_fields, field_name)) then
596  generic=>c_get_generic(sendable_fields, field_name)
597  select type(generic)
598  type is (io_server_sendable_field_sizing)
599  get_sendable_field_sizing=generic
600  if (present(field_found)) field_found=.true.
601  end select
602  end if
Here is the caller graph for this function:

◆ get_total_number_of_fields()

integer function iobridge_mod::get_total_number_of_fields ( type(definition_description_type), dimension(:), intent(inout)  definition_descriptions,
integer, intent(in)  number_defns 
)
private

Retrieve the total number of fields, which is all the fields in all the data definitions.

Parameters
definition_descriptionsData definition descriptions
number_defnsThe number of data definitions
Returns
The total number of fields

Definition at line 676 of file iobridge.F90.

677  type(definition_description_type), dimension(:), intent(inout) :: definition_descriptions
678  integer, intent(in) :: number_defns
679 
680  integer :: i
681 
682  get_total_number_of_fields=0
683  do i=1, number_defns
684  get_total_number_of_fields=get_total_number_of_fields+definition_descriptions(i)%number_fields
685  end do
Here is the caller graph for this function:

◆ handle_component_field_array_packing_into_send_buffer()

integer function iobridge_mod::handle_component_field_array_packing_into_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(io_configuration_data_definition_type), intent(inout)  data_definition,
type(io_configuration_field_type), intent(in)  field,
integer, intent(in)  current_buffer_point 
)
private

Packs a components field array into the send buffer, these are fields that are served up by components rather than explicitly available.

Parameters
current_stateThe current model state
data_definitionThe data definition description
fieldThe specific field we are looking up
current_buffer_pointThe current point in the buffer where this data will be entered
Returns
The new current buffer point which is after the data addition has taken place

Definition at line 1036 of file iobridge.F90.

1038  type(model_state_type), target, intent(inout) :: current_state
1039  type(io_configuration_data_definition_type), intent(inout) :: data_definition
1040  type(io_configuration_field_type), intent(in) :: field
1041  integer, intent(in) :: current_buffer_point
1042 
1043  type(component_field_information_type) :: field_descriptor
1044  type(component_field_value_type) :: published_value
1045 
1046  field_descriptor=get_component_field_descriptor(field%name)
1047  published_value=get_component_field_value(current_state, field%name)
1048  if (field_descriptor%data_type == component_double_data_type) then
1049  if (field_descriptor%number_dimensions == 1) then
1050  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
1051  current_buffer_point, real_array_1d=published_value%real_1d_array)
1052  deallocate(published_value%real_1d_array)
1053  else if (field_descriptor%number_dimensions == 2) then
1054  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
1055  current_buffer_point, real_array_2d=published_value%real_2d_array)
1056  deallocate(published_value%real_2d_array)
1057  else if (field_descriptor%number_dimensions == 3) then
1058  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
1059  current_buffer_point, real_array_3d=published_value%real_3d_array)
1060  deallocate(published_value%real_3d_array)
1061  else if (field_descriptor%number_dimensions == 4) then
1062  handle_component_field_array_packing_into_send_buffer=pack_array_field(data_definition%send_buffer, &
1063  current_buffer_point, real_array_4d=published_value%real_4d_array)
1064  deallocate(published_value%real_4d_array)
1065  end if
1066  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ handle_component_field_scalar_packing_into_send_buffer()

integer function iobridge_mod::handle_component_field_scalar_packing_into_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(io_configuration_data_definition_type), intent(inout)  data_definition,
type(io_configuration_field_type), intent(in)  field,
integer, intent(in)  current_buffer_point 
)
private

Packs a components field scalar into the send buffer, these are fields that are served up by components rather than explicitly available.

Parameters
current_stateThe current model state
data_definitionThe data definition description
fieldThe specific field we are looking up
current_buffer_pointThe current point in the buffer where this data will be entered
Returns
The new current buffer point which is after the data addition has taken place

Definition at line 853 of file iobridge.F90.

855  type(model_state_type), target, intent(inout) :: current_state
856  type(io_configuration_data_definition_type), intent(inout) :: data_definition
857  type(io_configuration_field_type), intent(in) :: field
858  integer, intent(in) :: current_buffer_point
859 
860  type(component_field_information_type) :: field_descriptor
861  type(component_field_value_type) :: published_value
862 
863  field_descriptor=get_component_field_descriptor(field%name)
864  published_value=get_component_field_value(current_state, field%name)
865  if (field_descriptor%data_type == component_double_data_type) then
866  handle_component_field_scalar_packing_into_send_buffer=pack_scalar_field(data_definition%send_buffer, &
867  current_buffer_point, real_value=published_value%scalar_real)
868  else if (field_descriptor%data_type == component_integer_data_type) then
869  handle_component_field_scalar_packing_into_send_buffer=pack_scalar_field(data_definition%send_buffer, &
870  current_buffer_point, int_value=published_value%scalar_int)
871  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_callback()

subroutine iobridge_mod::init_callback ( type(model_state_type), intent(inout), target  current_state)
private

Initialisation call back, called at the start of the model run.

Parameters
current_stateThe current model state

Definition at line 74 of file iobridge.F90.

75  type(model_state_type), target, intent(inout) :: current_state
76 
77  integer :: mpi_type_data_sizing_description, mpi_type_definition_description, mpi_type_field_description, ierr
78 
79  if (.not. options_get_logical(current_state%options_database, "enable_io_server")) then
80  io_server_enabled=.false.
81  call log_master_log(log_warn, "Enabled IO bridge but missing IO server compilation, therefore ignoring IO bridge component")
82  return
83  end if
84 
85  io_server_enabled=.true.
86  in_finalisation_callback=.false.
87 
88  call populate_sendable_fields(current_state)
89 
90  mpi_type_data_sizing_description=build_mpi_type_data_sizing_description()
91  mpi_type_definition_description=build_mpi_type_definition_description()
92  mpi_type_field_description=build_mpi_type_field_description()
93 
94  call register_with_io_server(current_state, mpi_type_definition_description, mpi_type_field_description)
95  call send_monc_specific_data_to_server(current_state, mpi_type_data_sizing_description)
96 
97  call mpi_type_free(mpi_type_data_sizing_description, ierr)
98  call mpi_type_free(mpi_type_definition_description, ierr)
99  call mpi_type_free(mpi_type_field_description, ierr)
100 
101  call build_mpi_data_types()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ iobridge_get_descriptor()

type(component_descriptor_type) function, public iobridge_mod::iobridge_get_descriptor

Definition at line 64 of file iobridge.F90.

65  iobridge_get_descriptor%name="iobridge"
66  iobridge_get_descriptor%version=0.1
67  iobridge_get_descriptor%initialisation=>init_callback
68  iobridge_get_descriptor%timestep=>timestep_callback
69  iobridge_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ pack_array_into_send_buffer()

integer function iobridge_mod::pack_array_into_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(io_configuration_data_definition_type), intent(inout)  data_definition,
type(io_configuration_field_type), intent(in)  field,
integer, intent(in)  current_buffer_point 
)
private

Packs array fields into the send bufer.

Parameters
current_stateThe current model state
data_definitionThe data definition description
fieldThe specific field we are looking up
current_buffer_pointThe current point in the buffer where this data will be entered
Returns
The new current buffer point which is after the data addition has taken place

Definition at line 910 of file iobridge.F90.

911  type(model_state_type), target, intent(inout) :: current_state
912  type(io_configuration_data_definition_type), intent(inout) :: data_definition
913  type(io_configuration_field_type), intent(in) :: field
914  integer, intent(in) :: current_buffer_point
915 
916  if (field%name .eq. "local_grid_size") then
917  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
918  int_array=current_state%local_grid%size)
919  else if (field%name .eq. "local_grid_start") then
920  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
921  int_array=current_state%local_grid%start)
922  else if (field%name .eq. "z") then
923  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
924  real_array_1d=current_state%global_grid%configuration%vertical%z)
925  else if (field%name .eq. "olubar") then
926  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
927  real_array_1d=current_state%global_grid%configuration%vertical%olubar)
928  else if (field%name .eq. "olzubar") then
929  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
930  real_array_1d=current_state%global_grid%configuration%vertical%olzubar)
931  else if (field%name .eq. "olvbar") then
932  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
933  real_array_1d=current_state%global_grid%configuration%vertical%olvbar)
934  else if (field%name .eq. "olzvbar") then
935  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
936  real_array_1d=current_state%global_grid%configuration%vertical%olzvbar)
937  else if (field%name .eq. "olthbar") then
938  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
939  real_array_1d=current_state%global_grid%configuration%vertical%olthbar)
940  else if (field%name .eq. "olzthbar") then
941  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
942  real_array_1d=current_state%global_grid%configuration%vertical%olzthbar)
943  else if (field%name .eq. "olqbar") then
944  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
945  real_array_2d=current_state%global_grid%configuration%vertical%olqbar)
946  else if (field%name .eq. "olzqbar") then
947  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
948  real_array_2d=current_state%global_grid%configuration%vertical%olzqbar)
949  else if (field%name .eq. "thref") then
950  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
951  real_array_1d=current_state%global_grid%configuration%vertical%thref)
952  else if (field%name .eq. "prefn") then
953  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
954  real_array_1d=current_state%global_grid%configuration%vertical%prefn)
955  else if (field%name .eq. "rhon") then
956  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
957  real_array_1d=current_state%global_grid%configuration%vertical%rhon)
958  else if (field%name .eq. "rho") then
959  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
960  real_array_1d=current_state%global_grid%configuration%vertical%rho)
961  else if (field%name .eq. "u") then
962  current_state%u%data=current_state%u%data+current_state%ugal
963  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%u, &
964  current_buffer_point, current_state%local_grid)
965  current_state%u%data=current_state%u%data-current_state%ugal
966  else if (field%name .eq. "u_nogal") then
967  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%u, current_buffer_point, &
968  current_state%local_grid)
969  else if (field%name .eq. "zu") then
970  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zu, current_buffer_point,&
971  current_state%local_grid)
972  else if (field%name .eq. "v") then
973  current_state%v%data=current_state%v%data+current_state%vgal
974  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%v, current_buffer_point, &
975  current_state%local_grid)
976  current_state%v%data=current_state%v%data-current_state%vgal
977  else if (field%name .eq. "v_nogal") then
978  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%v, current_buffer_point, &
979  current_state%local_grid)
980  else if (field%name .eq. "zv") then
981  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zv, current_buffer_point,&
982  current_state%local_grid)
983  else if (field%name .eq. "w") then
984  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%w, current_buffer_point, &
985  current_state%local_grid)
986  else if (field%name .eq. "zw") then
987  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zw, current_buffer_point,&
988  current_state%local_grid)
989  else if (field%name .eq. "q") then
990  pack_array_into_send_buffer=pack_q_fields(data_definition%send_buffer, current_state%q, current_state%number_q_fields, &
991  current_buffer_point, current_state%local_grid)
992  else if (field%name .eq. "zq") then
993  pack_array_into_send_buffer=pack_q_fields(data_definition%send_buffer, current_state%zq, current_state%number_q_fields, &
994  current_buffer_point, current_state%local_grid)
995  else if (field%name .eq. "th") then
996  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%th, current_buffer_point,&
997  current_state%local_grid)
998  else if (field%name .eq. "zth") then
999  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%zth, current_buffer_point,&
1000  current_state%local_grid)
1001  else if (field%name .eq. "p") then
1002  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, current_state%p, current_buffer_point, &
1003  current_state%local_grid)
1004  else if (field%name .eq. "sth_lw") then
1005  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, &
1006  current_state%sth_lw, current_buffer_point, current_state%local_grid)
1007  else if (field%name .eq. "sth_sw") then
1008  pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, &
1009  current_state%sth_sw, current_buffer_point, current_state%local_grid)
1010  else if (field%name .eq. "w_up") then
1011  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
1012  real_array_1d=current_state%global_grid%configuration%vertical%w_up)
1013  else if (field%name .eq. "w_dwn") then
1014  pack_array_into_send_buffer=pack_array_field(data_definition%send_buffer, current_buffer_point, &
1015  real_array_1d=current_state%global_grid%configuration%vertical%w_dwn)
1016 !!$ else if (field%name .eq. "sw_down_surf") then
1017 !!$ pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, &
1018 !!$ current_state%sth_sw, current_buffer_point, current_state%local_grid)
1019 !!$ else if (field%name .eq. "lww_down_surf") then
1020 !!$ pack_array_into_send_buffer=pack_prognostic_flow_field(data_definition%send_buffer, &
1021 !!$ current_state%sth_sw, current_buffer_point, current_state%local_grid)
1022  else
1023  ! Handle component field here
1024  pack_array_into_send_buffer=handle_component_field_array_packing_into_send_buffer(current_state, &
1025  data_definition, field, current_buffer_point)
1026  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ pack_map_into_send_buffer()

integer function iobridge_mod::pack_map_into_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(io_configuration_data_definition_type), intent(inout)  data_definition,
type(io_configuration_field_type), intent(in)  field,
integer, intent(in)  current_buffer_point 
)
private

Packs map fields into the send buffer.

Parameters
current_stateThe current model state
data_definitionThe data definition description
fieldThe specific field we are looking up
current_buffer_pointThe current point in the buffer where this data will be entered
Returns
The new current buffer point which is after the data addition has taken place

Definition at line 880 of file iobridge.F90.

881  type(model_state_type), target, intent(inout) :: current_state
882  type(io_configuration_data_definition_type), intent(inout) :: data_definition
883  type(io_configuration_field_type), intent(in) :: field
884  integer, intent(in) :: current_buffer_point
885 
886  integer :: i
887  type(q_metadata_type) :: specific_q_data
888  type(hashmap_type) :: q_indicies_map
889 
890  if (field%name .eq. "options_database") then
891  pack_map_into_send_buffer=pack_map_field(data_definition%send_buffer, current_buffer_point, current_state%options_database)
892  else if (field%name .eq. "q_indicies") then
893  do i=1, get_max_number_q_indices()
894  specific_q_data=get_indices_descriptor(i)
895  if (specific_q_data%l_used) then
896  call c_put_integer(q_indicies_map, specific_q_data%name, i)
897  end if
898  end do
899  pack_map_into_send_buffer=pack_map_field(data_definition%send_buffer, current_buffer_point, q_indicies_map)
900  call c_free(q_indicies_map)
901  end if
Here is the caller graph for this function:

◆ pack_prognostic_flow_field()

integer function iobridge_mod::pack_prognostic_flow_field ( character, dimension(:), intent(inout), allocatable  buffer,
type(prognostic_field_type), intent(inout)  prognostic,
integer, intent(in)  start_offset,
type(local_grid_type), intent(inout)  local_grid 
)
private

Packs the data of a specific prognostic field into a buffer.

Parameters
bufferThe buffer to pack the field into
prognosticThe prognostic field
start_offsetThe starting offset to write into the buffer
local_gridDescription of the local grid
Returns
The next location in the buffer to write to (next start offset)

Definition at line 1075 of file iobridge.F90.

1076  character, dimension(:), allocatable, intent(inout) :: buffer
1077  type(prognostic_field_type), intent(inout) :: prognostic
1078  integer, intent(in) :: start_offset
1079  type(local_grid_type), intent(inout) :: local_grid
1080 
1081  integer :: target_end
1082 
1083  target_end=start_offset + (local_grid%size(z_index)*local_grid%size(y_index)*local_grid%size(x_index)*kind(prognostic%data)-1)
1084 
1085  buffer(start_offset : target_end) = transfer(prognostic%data(&
1086  local_grid%local_domain_start_index(z_index): local_grid%local_domain_end_index(z_index),&
1087  local_grid%local_domain_start_index(y_index): local_grid%local_domain_end_index(y_index), &
1088  local_grid%local_domain_start_index(x_index): local_grid%local_domain_end_index(x_index)), &
1089  buffer(start_offset : target_end))
1090  pack_prognostic_flow_field=target_end+1
Here is the caller graph for this function:

◆ pack_q_fields()

integer function iobridge_mod::pack_q_fields ( character, dimension(:), intent(inout), allocatable  buffer,
type(prognostic_field_type), dimension(:), intent(inout)  q_fields,
integer, intent(in)  number_q_fields,
integer, intent(in)  start_offset,
type(local_grid_type), intent(inout)  local_grid 
)
private

Packs the Q fields into the send buffer.

Parameters
bufferThe send buffer to pack into
q_fieldsQ prognostic fields
number_q_fieldsThe number of Q fields
start_offsetStarting offset in the buffer to pack into
local_gridLocal grid description
Returns
Updated write location, which is the next location in the buffer to write to

Definition at line 1100 of file iobridge.F90.

1101  character, dimension(:), allocatable, intent(inout) :: buffer
1102  type(prognostic_field_type), dimension(:), intent(inout) :: q_fields
1103  integer, intent(in) :: start_offset, number_q_fields
1104  type(local_grid_type), intent(inout) :: local_grid
1105 
1106  integer :: target_end, i, current_starting_index
1107 
1108  current_starting_index=start_offset
1109 
1110  do i=1,number_q_fields
1111  target_end=current_starting_index + (local_grid%size(z_index)*local_grid%size(y_index)*&
1112  local_grid%size(x_index)*kind(q_fields(i)%data)-1)
1113  buffer(current_starting_index : target_end) = transfer(q_fields(i)%data(&
1114  local_grid%local_domain_start_index(z_index): local_grid%local_domain_end_index(z_index),&
1115  local_grid%local_domain_start_index(y_index): local_grid%local_domain_end_index(y_index), &
1116  local_grid%local_domain_start_index(x_index): local_grid%local_domain_end_index(x_index)), &
1117  buffer(current_starting_index : target_end))
1118  current_starting_index=target_end+1
1119  end do
1120  pack_q_fields=target_end+1
Here is the caller graph for this function:

◆ pack_scalar_into_send_buffer()

integer function iobridge_mod::pack_scalar_into_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(io_configuration_data_definition_type), intent(inout)  data_definition,
type(io_configuration_field_type), intent(in)  field,
integer, intent(in)  current_buffer_point 
)
private

Packs scalar fields into the send bufer.

Parameters
current_stateThe current model state
data_definitionThe data definition description
fieldThe specific field we are looking up
current_buffer_pointThe current point in the buffer where this data will be entered
Returns
The new current buffer point which is after the data addition has taken place

Definition at line 775 of file iobridge.F90.

776  type(model_state_type), target, intent(inout) :: current_state
777  type(io_configuration_data_definition_type), intent(inout) :: data_definition
778  type(io_configuration_field_type), intent(in) :: field
779  integer, intent(in) :: current_buffer_point
780 
781  if (field%name .eq. "timestep") then
782  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
783  int_value=current_state%timestep)
784  else if (field%name .eq. "terminated") then
785  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
786  logical_value=.not. current_state%continue_timestep .and. in_finalisation_callback)
787  else if (field%name .eq. "z_size") then
788  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
789  int_value=current_state%global_grid%size(z_index))
790  else if (field%name .eq. "y_size") then
791  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
792  int_value=current_state%global_grid%size(y_index))
793  else if (field%name .eq. "y_bottom") then
794  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
795  real_value=current_state%global_grid%bottom(y_index))
796  else if (field%name .eq. "y_top") then
797  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
798  real_value=current_state%global_grid%top(y_index))
799  else if (field%name .eq. "y_resolution") then
800  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
801  real_value=current_state%global_grid%resolution(y_index))
802  else if (field%name .eq. "x_size") then
803  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
804  int_value=current_state%global_grid%size(x_index))
805  else if (field%name .eq. "x_bottom") then
806  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
807  real_value=current_state%global_grid%bottom(x_index))
808  else if (field%name .eq. "x_top") then
809  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
810  real_value=current_state%global_grid%top(x_index))
811  else if (field%name .eq. "x_resolution") then
812  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
813  real_value=current_state%global_grid%resolution(x_index))
814  else if (field%name .eq. "time") then
815  ! The time is incremented with dtm as the model was about to increment for the next step and this is needed for diagnostics
816  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
817  real_value=current_state%time+current_state%dtm)
818  else if (field%name .eq. "ugal") then
819  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
820  real_value=current_state%ugal)
821  else if (field%name .eq. "vgal") then
822  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
823  real_value=current_state%vgal)
824  else if (field%name .eq. "nqfields") then
825  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
826  int_value=current_state%number_q_fields)
827  else if (field%name .eq. "dtm") then
828  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
829  real_value=current_state%dtm)
830  else if (field%name .eq. "dtm_new") then
831  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
832  real_value=current_state%dtm_new)
833  else if (field%name .eq. "absolute_new_dtm") then
834  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
835  real_value=current_state%absolute_new_dtm)
836  else if (field%name .eq. "rad_last_time") then
837  pack_scalar_into_send_buffer=pack_scalar_field(data_definition%send_buffer, current_buffer_point, &
838  real_value=current_state%rad_last_time)
839  else
840  ! Handle component field here
841  pack_scalar_into_send_buffer=handle_component_field_scalar_packing_into_send_buffer(current_state, &
842  data_definition, field, current_buffer_point)
843  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ pack_send_buffer()

subroutine iobridge_mod::pack_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(io_configuration_data_definition_type), intent(inout)  data_definition 
)
private

Packs the current state into the send buffer. This iterates through each field in the data description and adds it to the send buffer.

Parameters
current_stateThe current model state
data_definitionThe definition of the data which hold the send buffer and the fields

Definition at line 746 of file iobridge.F90.

747  type(model_state_type), target, intent(inout) :: current_state
748  type(io_configuration_data_definition_type), intent(inout) :: data_definition
749 
750  integer :: current_buffer_point, i
751 
752  current_buffer_point=1
753  do i=1, data_definition%number_of_data_fields
754  if (data_definition%fields(i)%enabled) then
755  if (data_definition%fields(i)%field_type == array_field_type) then
756  current_buffer_point=pack_array_into_send_buffer(current_state, data_definition, data_definition%fields(i), &
757  current_buffer_point)
758  else if (data_definition%fields(i)%field_type == map_field_type) then
759  current_buffer_point=pack_map_into_send_buffer(current_state, data_definition, data_definition%fields(i), &
760  current_buffer_point)
761  else if (data_definition%fields(i)%field_type == scalar_field_type) then
762  current_buffer_point=pack_scalar_into_send_buffer(current_state, data_definition, data_definition%fields(i), &
763  current_buffer_point)
764  end if
765  end if
766  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ package_local_monc_decomposition_into_descriptions()

subroutine iobridge_mod::package_local_monc_decomposition_into_descriptions ( type(model_state_type), intent(inout), target  current_state,
type(data_sizing_description_type), dimension(:), intent(inout)  data_description 
)
private

Packages the local MONC decomposition information into descriptions for communication.

Parameters
current_stateThe current model state
data_descriptionTHe data description to pack into

Definition at line 561 of file iobridge.F90.

562  type(model_state_type), target, intent(inout) :: current_state
563  type(data_sizing_description_type), dimension(:), intent(inout) :: data_description
564 
565  type(io_server_sendable_field_sizing) :: sizing_info
566 
567  sizing_info%number_dimensions=3
568  sizing_info%dimensions(z_index)=current_state%local_grid%size(z_index)
569  sizing_info%dimensions(y_index)=current_state%local_grid%size(y_index)
570  sizing_info%dimensions(x_index)=current_state%local_grid%size(x_index)
571  call assemble_individual_description(data_description, 1, local_sizes_key, sizing_info)
572  sizing_info%dimensions(z_index)=current_state%local_grid%start(z_index)
573  sizing_info%dimensions(y_index)=current_state%local_grid%start(y_index)
574  sizing_info%dimensions(x_index)=current_state%local_grid%start(x_index)
575  call assemble_individual_description(data_description, 2, local_start_points_key, sizing_info)
576  sizing_info%dimensions(z_index)=current_state%local_grid%end(z_index)
577  sizing_info%dimensions(y_index)=current_state%local_grid%end(y_index)
578  sizing_info%dimensions(x_index)=current_state%local_grid%end(x_index)
579  call assemble_individual_description(data_description, 3, local_end_points_key, sizing_info)
580  sizing_info%number_dimensions=1
581  sizing_info%dimensions(1)=get_number_active_q_indices()
582  call assemble_individual_description(data_description, 4, number_q_indicies_key, sizing_info)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ populate_component_public_field()

subroutine iobridge_mod::populate_component_public_field ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  field_name 
)
private

Populates the field information for a specific publically available field offered by one of the components.

Parameters
field_visibility_descriptorThe field descriptor which contains sizing information

Definition at line 282 of file iobridge.F90.

283  type(model_state_type), target, intent(inout) :: current_state
284  character(len=*), intent(in) :: field_name
285 
286  class(*), pointer :: generic_data
287  type(io_server_sendable_field_sizing), pointer :: field_sizing
288  type(component_field_information_type) :: field_information
289  type(component_field_information_type), pointer :: field_information_info_alloc
290 
291 
292  field_information=get_component_field_information(current_state, field_name)
293  if (field_information%enabled) then
294  allocate(field_information_info_alloc, source=field_information)
295  generic_data=>field_information_info_alloc
296  call c_put_generic(component_field_descriptions, field_name, generic_data, .false.)
297 
298  allocate(field_sizing)
299  field_sizing%number_dimensions=field_information%number_dimensions
300  field_sizing%dimensions=field_information%dimension_sizes
301  generic_data=>field_sizing
302  call c_put_generic(sendable_fields, field_name, generic_data, .false.)
303  end if
Here is the caller graph for this function:

◆ populate_data_definition_configuration()

subroutine iobridge_mod::populate_data_definition_configuration ( type(definition_description_type), dimension(:), intent(inout)  definition_descriptions,
integer, intent(in)  number_defns,
type(field_description_type), dimension(:), intent(inout)  field_descriptions,
integer, intent(in)  number_fields 
)
private

Based upon the received data and field definitions this will configure the IO bridge internal representation of these facets, which is a structured tree, data defintions holding their own fields rather than the unstructured data we get from the IO server.

Parameters
definition_descriptionsdata definitions received from the IO server
number_defnsThe number of data definitions
field_descriptionsThe field descriptions that we received from the IO server
number_fieldsThe total number of field descriptions received

Definition at line 695 of file iobridge.F90.

696  type(definition_description_type), dimension(:), intent(inout) :: definition_descriptions
697  type(field_description_type), dimension(:), intent(inout) :: field_descriptions
698  integer, intent(in) :: number_defns, number_fields
699 
700  integer :: i, definition_index, field_index
701 
702  allocate(data_definitions(number_defns))
703  do i=1, number_defns
704  data_definitions(i)%name=definition_descriptions(i)%definition_name
705  data_definitions(i)%send_on_terminate=definition_descriptions(i)%send_on_terminate
706  data_definitions(i)%number_of_data_fields=0 ! Will increment this for each field
707  data_definitions(i)%frequency=definition_descriptions(i)%frequency
708  allocate(data_definitions(i)%fields(definition_descriptions(i)%number_fields))
709  end do
710  do i=1, number_fields
711  definition_index=get_definition_index(field_descriptions(i)%definition_name)
712  field_index=data_definitions(definition_index)%number_of_data_fields+1
713  data_definitions(definition_index)%number_of_data_fields=field_index
714  data_definitions(definition_index)%fields(field_index)%name=field_descriptions(i)%field_name
715  data_definitions(definition_index)%fields(field_index)%field_type=field_descriptions(i)%field_type
716  data_definitions(definition_index)%fields(field_index)%data_type=field_descriptions(i)%data_type
717  data_definitions(definition_index)%fields(field_index)%optional=field_descriptions(i)%optional
718  if (field_descriptions(i)%optional .or. field_descriptions(i)%field_type == array_field_type .or. &
719  field_descriptions(i)%field_type == map_field_type) then
720  call c_put_integer(unique_field_names, field_descriptions(i)%field_name, 1)
721  end if
722  if (.not. field_descriptions(i)%optional) data_definitions(definition_index)%fields(field_index)%enabled=.true.
723  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ populate_globally_visible_sendable_fields()

subroutine iobridge_mod::populate_globally_visible_sendable_fields ( type(model_state_type), intent(inout), target  current_state)
private

Populates the globally visible sendable fields which is a key value pair mapping between name and description of that field.

Parameters
current_stateThe current model state

Definition at line 308 of file iobridge.F90.

309  type(model_state_type), target, intent(inout) :: current_state
310 
311  integer :: x_size, y_size, z_size
312  class(*), pointer :: raw_generic
313 
314  z_size=current_state%local_grid%size(z_index)
315  y_size=current_state%local_grid%size(y_index)
316  x_size=current_state%local_grid%size(x_index)
317 
318  raw_generic=>generate_sendable_description(options_size(current_state%options_database))
319  call c_put_generic(sendable_fields, "options_database", raw_generic, .false.)
320  if (get_number_active_q_indices() .gt. 0) then
321  raw_generic=>generate_sendable_description(get_number_active_q_indices())
322  call c_put_generic(sendable_fields, "q_indicies", raw_generic, .false.)
323  end if
324  raw_generic=>generate_sendable_description(3)
325  call c_put_generic(sendable_fields, "local_grid_size", raw_generic, .false.)
326  call c_put_generic(sendable_fields, "local_grid_start", raw_generic, .false.)
327 #ifdef U_ACTIVE
328  raw_generic=>generate_sendable_description()
329  call c_put_generic(sendable_fields, "x_size", raw_generic, .false.)
330  call c_put_generic(sendable_fields, "x_bottom", raw_generic, .false.)
331  call c_put_generic(sendable_fields, "x_resolution", raw_generic, .false.)
332  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
333  call c_put_generic(sendable_fields, "u", raw_generic, .false.)
334  call c_put_generic(sendable_fields, "u_nogal", raw_generic, .false.)
335  call c_put_generic(sendable_fields, "zu", raw_generic, .false.)
336  if (allocated(current_state%global_grid%configuration%vertical%olubar)) then
337  raw_generic=>generate_sendable_description(z_size)
338  call c_put_generic(sendable_fields, "olubar", raw_generic, .false.)
339  call c_put_generic(sendable_fields, "olzubar", raw_generic, .false.)
340  end if
341 #endif
342 #ifdef V_ACTIVE
343  raw_generic=>generate_sendable_description()
344  call c_put_generic(sendable_fields, "y_size", raw_generic, .false.)
345  call c_put_generic(sendable_fields, "y_bottom", raw_generic, .false.)
346  call c_put_generic(sendable_fields, "y_resolution", raw_generic, .false.)
347  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
348  call c_put_generic(sendable_fields, "v", raw_generic, .false.)
349  call c_put_generic(sendable_fields, "v_nogal", raw_generic, .false.)
350  call c_put_generic(sendable_fields, "zv", raw_generic, .false.)
351  if (allocated(current_state%global_grid%configuration%vertical%olvbar)) then
352  raw_generic=>generate_sendable_description(z_size)
353  call c_put_generic(sendable_fields, "olvbar", raw_generic, .false.)
354  call c_put_generic(sendable_fields, "olzvbar", raw_generic, .false.)
355  end if
356 #endif
357 #ifdef W_ACTIVE
358  raw_generic=>generate_sendable_description(z_size)
359  call c_put_generic(sendable_fields, "z", raw_generic, .false.)
360  call c_put_generic(sendable_fields, "thref", raw_generic, .false.)
361  call c_put_generic(sendable_fields, "prefn", raw_generic, .false.)
362  call c_put_generic(sendable_fields, "rhon", raw_generic, .false.)
363  call c_put_generic(sendable_fields, "rho", raw_generic, .false.)
364  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
365  call c_put_generic(sendable_fields, "w", raw_generic, .false.)
366  call c_put_generic(sendable_fields, "zw", raw_generic, .false.)
367 #endif
368  if (current_state%number_q_fields .gt. 0) then
369  raw_generic=>generate_sendable_description(z_size, y_size, x_size, current_state%number_q_fields)
370  call c_put_generic(sendable_fields, "q", raw_generic, .false.)
371  call c_put_generic(sendable_fields, "zq", raw_generic, .false.)
372  if (allocated(current_state%global_grid%configuration%vertical%olqbar)) then
373  raw_generic=>generate_sendable_description(z_size, current_state%number_q_fields)
374  call c_put_generic(sendable_fields, "olqbar", raw_generic, .false.)
375  call c_put_generic(sendable_fields, "olzqbar", raw_generic, .false.)
376  end if
377  end if
378  if (current_state%th%active) then
379  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
380  call c_put_generic(sendable_fields, "th", raw_generic, .false.)
381  call c_put_generic(sendable_fields, "zth", raw_generic, .false.)
382  if (allocated(current_state%global_grid%configuration%vertical%olthbar)) then
383  raw_generic=>generate_sendable_description(z_size)
384  call c_put_generic(sendable_fields, "olthbar", raw_generic, .false.)
385  call c_put_generic(sendable_fields, "olzthbar", raw_generic, .false.)
386  end if
387  end if
388  if (current_state%p%active) then
389  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
390  call c_put_generic(sendable_fields, "p", raw_generic, .false.)
391  end if
392  ! need to dump heating rate tendency from socrates radiation
393  if (is_component_enabled(current_state%options_database, "socrates_couple")) then
394  raw_generic=>generate_sendable_description(z_size, y_size, x_size)
395  call c_put_generic(sendable_fields, "sth_lw", raw_generic, .false.)
396  call c_put_generic(sendable_fields, "sth_sw", raw_generic, .false.)
397  endif
398 
399  if (is_component_enabled(current_state%options_database, "pdf_analysis")) then
400  if (allocated(current_state%global_grid%configuration%vertical%w_up)) then
401  raw_generic=>generate_sendable_description(z_size)
402  call c_put_generic(sendable_fields, "w_up", raw_generic, .false.)
403  end if
404 
405  if (allocated(current_state%global_grid%configuration%vertical%w_dwn)) then
406  raw_generic=>generate_sendable_description(z_size)
407  call c_put_generic(sendable_fields, "w_dwn", raw_generic, .false.)
408  end if
409  end if
410 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ populate_sendable_fields()

subroutine iobridge_mod::populate_sendable_fields ( type(model_state_type), intent(inout), target  current_state)
private

Populates the sendable field definitions with the field sizing information.

Parameters
current_stateThe current model state

Definition at line 267 of file iobridge.F90.

268  type(model_state_type), target, intent(inout) :: current_state
269 
270  type(list_type) :: published_field_descriptors
271  integer :: i
272 
273  call populate_globally_visible_sendable_fields(current_state)
274  published_field_descriptors=get_all_component_published_fields()
275  do i=1, c_size(published_field_descriptors)
276  call populate_component_public_field(current_state, c_get_string(published_field_descriptors, i))
277  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ register_with_io_server()

subroutine iobridge_mod::register_with_io_server ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  mpi_type_definition_description,
integer, intent(in)  mpi_type_field_description 
)
private

Registers this MONC with the corresponding IO server. This will encapsulate the entire protocol, which is sending the registration command, receiving the data and field definitions from the IO server and then sending back the sizing for the fields that this MONC will contribute.

Parameters
current_stateThe current model state
mpi_type_definition_descriptionMPI data type for data definition message
mpi_type_field_descriptionMPI data type for field definition message

Definition at line 646 of file iobridge.F90.

647  type(model_state_type), target, intent(inout) :: current_state
648  integer, intent(in) :: mpi_type_definition_description, mpi_type_field_description
649 
650  type(definition_description_type), dimension(:), allocatable :: definition_descriptions
651  type(field_description_type), dimension(:), allocatable :: field_descriptions
652  integer :: number_defns, number_fields, status(MPI_STATUS_SIZE), ierr
653 
654  call mpi_send(register_command, 1, mpi_int, current_state%parallel%corresponding_io_server_process, &
655  command_tag, mpi_comm_world, ierr)
656 
657  call mpi_probe(current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, status, ierr)
658  call mpi_get_count(status, mpi_type_definition_description, number_defns, ierr)
659  allocate(definition_descriptions(number_defns))
660 
661  call mpi_recv(definition_descriptions, number_defns, mpi_type_definition_description, &
662  current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, mpi_status_ignore, ierr)
663  number_fields=get_total_number_of_fields(definition_descriptions, number_defns)
664 
665  allocate(field_descriptions(number_fields))
666  call mpi_recv(field_descriptions, number_fields, mpi_type_field_description, &
667  current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, mpi_status_ignore, ierr)
668  call populate_data_definition_configuration(definition_descriptions, number_defns, field_descriptions, number_fields)
669  deallocate(definition_descriptions)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ send_data_field_sizes_to_server()

integer function iobridge_mod::send_data_field_sizes_to_server ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  mpi_type_data_sizing_description,
type(data_sizing_description_type), dimension(:), intent(inout)  data_description,
integer, intent(in)  number_unique_fields 
)
private

Assembles all the data field sizing information and sends this to the IO server.

Parameters
current_stateThe current model state
mpi_type_data_sizing_descriptionMPI data type representing the sizing message
data_descriptionData descriptions which will be populated and then sent
number_unique_fieldsThe number of unique fields that we are sending over

Definition at line 481 of file iobridge.F90.

483  type(model_state_type), target, intent(inout) :: current_state
484  integer, intent(in) :: mpi_type_data_sizing_description, number_unique_fields
485  type(data_sizing_description_type), dimension(:), intent(inout) :: data_description
486 
487  integer :: ierr, i, next_index, request_handle
488  character(len=STRING_LENGTH) :: field_name
489 
490  call package_local_monc_decomposition_into_descriptions(current_state, data_description)
491  next_index=5
492  do i=1, number_unique_fields
493  field_name=c_key_at(unique_field_names, i)
494  if (c_contains(sendable_fields, field_name)) then
495  call assemble_individual_description(data_description, next_index, field_name, get_sendable_field_sizing(field_name))
496  next_index=next_index+1
497  end if
498  end do
499 
500  call mpi_isend(data_description, next_index-1, mpi_type_data_sizing_description, &
501  current_state%parallel%corresponding_io_server_process, data_tag, mpi_comm_world, request_handle, ierr)
502  send_data_field_sizes_to_server=request_handle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ send_data_to_io_server()

subroutine iobridge_mod::send_data_to_io_server ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  data_index 
)
private

Sends data to the IO server.

Parameters
current_stateThe current model state
data_indexThe specific data index to send over

Definition at line 125 of file iobridge.F90.

126  type(model_state_type), target, intent(inout) :: current_state
127  integer, intent(in) :: data_index
128 
129  integer :: command_to_send, ierr
130 
131  if (data_definitions(data_index)%dump_requests(1) .ne. mpi_request_null .or. &
132  data_definitions(data_index)%dump_requests(2) .ne. mpi_request_null) then
133  ! Here wait for previous data dump to complete (consider extending to using buffers for performance)
134  call mpi_waitall(2, data_definitions(data_index)%dump_requests, mpi_statuses_ignore, ierr)
135  end if
136 
137  ! Pack the send buffer and send it to the IO server
138  call pack_send_buffer(current_state, data_definitions(data_index))
139 
140  command_to_send=data_command_start+data_index
141  call mpi_issend(command_to_send, 1, mpi_int, current_state%parallel%corresponding_io_server_process, &
142  command_tag, mpi_comm_world, data_definitions(data_index)%dump_requests(1), ierr)
143  call mpi_issend(data_definitions(data_index)%send_buffer, 1, data_definitions(data_index)%mpi_datatype, &
144  current_state%parallel%corresponding_io_server_process, data_tag+data_index, mpi_comm_world, &
145  data_definitions(data_index)%dump_requests(2), ierr)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ send_general_monc_information_to_server()

integer function iobridge_mod::send_general_monc_information_to_server ( type(model_state_type), intent(inout), target  current_state,
character, dimension(:), intent(inout)  buffer 
)
private

Sends the general MONC information (ZN field and Q field names) to the IO server.

Parameters
current_stateThe current model state
bufferThe communication buffer to use
Returns
Handle to nonblocking send

Definition at line 509 of file iobridge.F90.

510  type(model_state_type), target, intent(inout) :: current_state
511  character, dimension(:), intent(inout) :: buffer
512 
513  character(len=STRING_LENGTH) :: q_field_name, cd_field_name
514  type(q_metadata_type) :: q_meta_data
515  integer :: current_loc, n, ierr, request_handle
516 
517  current_loc=1
518  current_loc=pack_array_field(buffer, current_loc, real_array_1d=current_state%global_grid%configuration%vertical%zn)
519  if (current_state%number_q_fields .gt. 0) then
520  do n=1, current_state%number_q_fields
521  q_meta_data=get_indices_descriptor(n)
522  if (q_meta_data%l_used) then
523  q_field_name=q_meta_data%name
524  else
525  q_field_name="qfield_"//trim(conv_to_string(n))
526  end if
527  current_loc=pack_scalar_field(buffer, current_loc, string_value=q_field_name)
528  end do
529  end if
530  current_loc=pack_array_field(buffer, current_loc, real_array_1d=current_state%global_grid%configuration%vertical%z)
531 
532  do n=1,ncond*2
533  if (n .le. ncond) then
534  cd_field_name = cond_request(n)
535  current_loc=pack_scalar_field(buffer, current_loc, string_value=cd_field_name)
536  cd_field_name = cond_long(n)
537  current_loc=pack_scalar_field(buffer, current_loc, string_value=cd_field_name)
538  else
539  cd_field_name = .not." "//trim(cond_request(n-ncond))
540  current_loc=pack_scalar_field(buffer, current_loc, string_value=cd_field_name)
541  cd_field_name = .not." "//trim(cond_long(n-ncond))
542  current_loc=pack_scalar_field(buffer, current_loc, string_value=cd_field_name)
543  end if
544  end do
545  do n=1,ndiag
546  cd_field_name = diag_request(n)
547  current_loc=pack_scalar_field(buffer, current_loc, string_value=cd_field_name)
548  cd_field_name = diag_long(n)
549  current_loc=pack_scalar_field(buffer, current_loc, string_value=cd_field_name)
550  end do
551 
552 
553  call mpi_isend(buffer, current_loc-1, mpi_byte, current_state%parallel%corresponding_io_server_process, &
554  data_tag, mpi_comm_world, request_handle, ierr)
555  send_general_monc_information_to_server=request_handle
Here is the caller graph for this function:

◆ send_monc_specific_data_to_server()

subroutine iobridge_mod::send_monc_specific_data_to_server ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  mpi_type_data_sizing_description 
)
private

Sends this MONC specific information to the IO server, which is field info (sizing & availability) as well as meta data such as ZN field and Q field names.

Parameters
current_stateThe current model state
mpi_type_data_sizing_descriptionMPI data type representing the sizing message

Definition at line 454 of file iobridge.F90.

455  type(model_state_type), target, intent(inout) :: current_state
456  integer, intent(in) :: mpi_type_data_sizing_description
457 
458  type(data_sizing_description_type), dimension(:), allocatable :: data_description
459  character, dimension(:), allocatable :: buffer
460  integer :: number_unique_fields, buffer_size, request_handles(2), ierr
461  real(kind=default_precision) :: dreal
462 
463  number_unique_fields=c_size(unique_field_names)
464  allocate(data_description(number_unique_fields+4))
465  request_handles(1)=send_data_field_sizes_to_server(current_state, mpi_type_data_sizing_description, &
466  data_description, number_unique_fields)
467  buffer_size=(kind(dreal)*current_state%local_grid%size(z_index))*2 + (string_length * current_state%number_q_fields &
468  + 4*ncond*string_length + 2*ndiag*string_length )
469  allocate(buffer(buffer_size))
470  request_handles(2)=send_general_monc_information_to_server(current_state, buffer)
471  call mpi_waitall(2, request_handles, mpi_statuses_ignore, ierr)
472  deallocate(data_description)
473  deallocate(buffer)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ timestep_callback()

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

Model dump call back, called for each model dump phase.

Parameters
current_stateThe current model state

Definition at line 106 of file iobridge.F90.

107  type(model_state_type), target, intent(inout) :: current_state
108 
109  integer :: i
110 
111  if (.not. io_server_enabled) return
112 
113  do i=1, size(data_definitions)
114  if (data_definitions(i)%frequency .gt. 0) then
115  if (mod(current_state%timestep, data_definitions(i)%frequency) == 0) then
116  call send_data_to_io_server(current_state, i)
117  end if
118  end if
119  end do
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ component_field_descriptions

type(map_type) iobridge_mod::component_field_descriptions
private

Definition at line 57 of file iobridge.F90.

◆ data_definitions

type(io_configuration_data_definition_type), dimension(:), allocatable iobridge_mod::data_definitions
private

Definition at line 56 of file iobridge.F90.

56  type(io_configuration_data_definition_type), dimension(:), allocatable :: data_definitions

◆ in_finalisation_callback

logical iobridge_mod::in_finalisation_callback
private

Definition at line 58 of file iobridge.F90.

◆ io_server_enabled

logical iobridge_mod::io_server_enabled
private

Definition at line 58 of file iobridge.F90.

58  logical :: io_server_enabled, in_finalisation_callback

◆ sendable_fields

type(map_type) iobridge_mod::sendable_fields
private

Definition at line 57 of file iobridge.F90.

◆ unique_field_names

type(map_type) iobridge_mod::unique_field_names
private

Definition at line 57 of file iobridge.F90.

57  type(map_type) :: unique_field_names, sendable_fields, component_field_descriptions
logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
logging_mod::log_warn
integer, parameter, public log_warn
Log WARNING and ERROR messages.
Definition: logging.F90:12
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_master_log
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17