MONC
Functions/Subroutines | Variables
writer_federator_mod Module Reference

This federates over the writing of diagnostic and prognostic data to the file system. It also manages the time manipulation of fields and groups. More...

Functions/Subroutines

subroutine, public initialise_writer_federator (io_configuration, diagnostic_generation_frequency, continuation_run)
 Initialises the write federator and configures it based on the user configuration. Also initialises the time manipulations. More...
 
subroutine, public finalise_writer_federator ()
 Finalises the write federator and the manipulations. More...
 
subroutine, public inform_writer_federator_time_point (io_configuration, source, data_id, data_dump)
 
subroutine, public inform_writer_federator_fields_present (io_configuration, field_names, diag_field_names_and_roots)
 Informs the writer federator that specific fields are present and should be reflected in the diagnostics output. More...
 
logical function, public is_field_used_by_writer_federator (field_name, field_namespace)
 Determines whether a field is used by the writer federator or not. More...
 
logical function, public is_field_split_on_q (field_name)
 Determines whether a field is split on Q or not. More...
 
subroutine enable_specific_field_by_name (field_name, diagnostics_mode, expected_here)
 Enables a specific field by its name, this will locate all the fields with this name and enable them. More...
 
subroutine, public provide_q_field_names_to_writer_federator (q_provided_field_names)
 Provides the Q field names to the write federator, this is required as on initialisation we don't know what these are and only when MONC register do they inform the IO server of the specifics. More...
 
subroutine, public provide_ordered_field_to_writer_federator (io_configuration, field_name, field_namespace, field_values, timestep, time, source)
 
subroutine provide_ordered_field_to_writer_federator_real_values (io_configuration, field_name, field_namespace, field_values, timestep, time, source)
 Provides fields (either diagnostics or prognostics) to the write federator which will action these as appropriate. This will split Q fields up if appropriate. More...
 
integer function get_size_of_collective_q (io_configuration, field_name, source)
 Retrieves the data size for each Q entry of a collective Q field for the specific source MONC that has sent data. More...
 
subroutine provide_ordered_single_field_to_writer_federator (io_configuration, field_name, field_namespace, field_values, timestep, time, source)
 Provides a single ordered field, i.e. Q fields have been split by this point. More...
 
subroutine write_collective_write_value (result_values, writer_index, contents_index, source, lookup_key)
 Writes the collective values, this is held differently to independent values which are written directly - instead here we need to store the values for each MONC hence a specific type is used instead. More...
 
subroutine determine_if_outstanding_field_can_be_written (io_configuration, writer_entry, specific_field)
 For a specific field wil determine and handle any outstanding fields writes until an outstanding write can not be performed or the outstanding list is empty. More...
 
subroutine determine_if_field_can_be_written (io_configuration, writer_entry, specific_field, timestep, previous_write_timestep, write_time, previous_write_time, field_written)
 Determines if a file can be written to its overarching write representation. If so then a write is issued, otherwise an outstanding write point is registered which will be checked frequency to do a write later on. More...
 
subroutine, public check_writer_for_trigger (io_configuration, source, data_id, data_dump)
 Checks all writer entries for any trigger fires and issues the underlying file storage. More...
 
subroutine check_writer_trigger (io_configuration, writer_entry_index, timestep, time, terminated)
 Checks a writer trigger and issues a file creation along with field write if the conditions (time or timestep) are met. This will either create and write to the file or store a pending state if one is already open (required due to NetCDF/HDF5 limitations with thread safety and parallel access.) More...
 
subroutine, public issue_actual_write (io_configuration, writer_entry, timestep, time, terminated_write)
 Issues the actual file creation, write of available fields and closure if all completed. More...
 
subroutine clean_time_points ()
 Cleans out old timepoints which are no longer going to be of any relavence to the file writing. This ensures that we don't have lots of stale points that need to be processed and searched beyond but are themselves pointless. More...
 
type(map_type) function extract_applicable_time_points (start_time, end_time)
 Extracts the applicable time points from the overall map that lie within a specific range. More...
 
type(map_type) function sort_applicable_time_points (unsorted_timepoints)
 Sorts the time points based upon their timestep, smallest to largest. Note that this is a bubble sort and as such inefficient, so would be good to change to something else but works OK for now. More...
 
subroutine close_diagnostics_file (io_configuration, writer_entry, timestep, time)
 Closes the diagnostics file, this is done via a global callback to issue the closes synchronously (collective operation) More...
 
subroutine handle_close_diagnostics_globalcallback (io_configuration, values, field_name, timestep)
 Call back for the inter IO reduction which actually does the NetCDF file closing which is a collective (synchronous) operation. Calls out to the NetCDF code to do the call and then checks the list of pending file writes to process any others that are waiting in the queue. More...
 
logical function check_for_and_issue_chain_write (io_configuration, writer_entry)
 Will check whether there are any pending writes and if so will issue a chain write for this. More...
 
subroutine register_pending_file_write (writer_entry_index, timestep, time, terminated_write)
 Registers a pending file write which will be actioned later on. More...
 
logical function get_next_applicable_writer_entry (field_name, field_namespace, writer_index_point, contents_index_point)
 Retrieves the index of the next writer which uses a specific field. If none is found then returns false, otherwise true. More...
 
integer function get_total_number_writer_fields (io_configuration, writer_entry_index)
 Determines the total number of fields that make up a writer entry, this is all the fields of the groups that make up this writer and individual fields specified too. More...
 
integer function get_group_number_of_fields (io_configuration, group_members, num_q_fields, namespace)
 Retrieves the number of fields within a group of fields. More...
 
integer function get_field_number_of_fields (io_configuration, field_name, field_namespace, num_q_fields)
 Retrieves the number of fields that make up this field, if it is a Q field then it will be split into many subfields hence it is not a simple 1-1 mapping. More...
 
integer function add_group_of_fields_to_writer_entry (io_configuration, writer_entry_index, facet_index, current_field_index, writer_field_names, duplicate_field_names, diagnostic_generation_frequency)
 Adds a group of fields to a writer entry, groups are expanded out into individual fields, each inherit the properties of the group. More...
 
integer function add_field_to_writer_entry (io_configuration, writer_entry_index, io_config_facet_index, my_facet_index, field_name, field_namespace, writer_field_names, duplicate_field_names, diagnostic_generation_frequency)
 Adds a field to the writer entry, this will split the Q fields. However at initialisation we don't know what the Q fields are called, hence place a marker which will be replaced later on. More...
 
subroutine add_specific_field_to_writer_entry (io_configuration, writer_entry_index, io_config_facet_index, my_facet_index, field_name, field_namespace, writer_field_names, duplicate_field_names, timestep_frequency, diagnostic_field_configuration, prognostic_field_configuration)
 Adds a specific field and its information to a writer entry. More...
 
subroutine handle_duplicate_field_names (writer_entry, duplicate_field_names)
 Marks duplicate field names in a writer entry as duplicates so that the NetCDF layer can then deal with this by issuing unique names. More...
 
integer function get_index_of_group (io_configuration, group_name)
 Searches the IO server configuration for a group with a specific name and returns the index to that group or 0 if no corresponding group is found. More...
 
subroutine determine_collective_type_and_optimise_if_possible (io_configuration, field_to_write_information)
 Determines whether it can optimise a specific collective field. If the field fits into certain limited parameters then it will optimise it. These parameters are very common, hence most fields can be optimised. Basically, it is looking to contiguous blocks of data from different MONCs so that the number of writes to the NetCDF file is limited. More...
 
subroutine initialise_contiguous_data_regions (io_configuration, field_to_write_information)
 Will initialise the collective data regions that form contiguous blocks within the data. This is quite an expensive operation so only done once for each field, but has the potential for very significant performance advantages for the fields that match it. More...
 
subroutine get_common_starts (dim, val, vals, common_starters, num_common)
 Retrieves the number of common starting points that match a specific input value. More...
 
integer function get_dimension_identifier (dim_name, is_auto_dimension)
 Translates a dimension name to its numeric corresponding identifier. More...
 

Variables

type(writer_type), dimension(:), allocatable, volatile writer_entries
 
type(hashset_type), volatile used_field_names
 
type(hashset_type), volatile q_field_names
 
type(hashmap_type), volatile time_points
 
type(hashmap_type), volatile q_field_splits
 
type(hashmap_type), volatile collective_q_field_dims
 
integer, volatile time_points_rwlock
 
integer, volatile collective_contiguous_initialisation_mutex
 
integer, volatile currently_writing_mutex
 
logical, volatile currently_writing
 

Detailed Description

This federates over the writing of diagnostic and prognostic data to the file system. It also manages the time manipulation of fields and groups.

Function/Subroutine Documentation

◆ add_field_to_writer_entry()

integer function writer_federator_mod::add_field_to_writer_entry ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  writer_entry_index,
integer, intent(in)  io_config_facet_index,
integer, intent(in)  my_facet_index,
character(len=*), intent(in)  field_name,
character(len=*), intent(in)  field_namespace,
type(hashset_type), intent(inout)  writer_field_names,
type(hashset_type), intent(inout)  duplicate_field_names,
type(hashmap_type), intent(inout)  diagnostic_generation_frequency 
)
private

Adds a field to the writer entry, this will split the Q fields. However at initialisation we don't know what the Q fields are called, hence place a marker which will be replaced later on.

Parameters
io_configurationThe IO server configuration
writer_entry_indexIndex of the writer entry that we are dealing with
io_config_facet_indexIndex of the facet (group) in the IO server configuration
my_facet_indexThe current field index in this internal module representation of the structure
field_nameThe name of the field that we are constructing
writer_field_namesThe field names in the writer (for duplication checking)
duplicate_field_namesDuplicate field names in the wrier, for duplication checking
diagnostic_generation_frequencyGeneration frequency of the diagnostics
Returns
Location for next field to be written to

Definition at line 1181 of file writer_federator.F90.

1183  type(io_configuration_type), intent(inout) :: io_configuration
1184  integer, intent(in) :: writer_entry_index, io_config_facet_index, my_facet_index
1185  character(len=*), intent(in) :: field_name, field_namespace
1186  type(hashset_type), intent(inout) :: writer_field_names, duplicate_field_names
1187  type(hashmap_type), intent(inout) :: diagnostic_generation_frequency
1188 
1189  integer :: i, number_q_fields, tot_size
1190  type(io_configuration_field_type) :: prognostic_field_configuration
1191  type(io_configuration_data_definition_type) :: prognostic_containing_data_defn
1192  type(io_configuration_diagnostic_field_type) :: diagnostic_field_configuration
1193  type(collective_q_field_representation_type), pointer :: collective_q_field
1194  class(*), pointer :: generic
1195 
1196  if (get_diagnostic_field_configuration(io_configuration, field_name, field_namespace, diagnostic_field_configuration)) then
1197  if (diagnostic_field_configuration%field_type == array_field_type) then
1198  if (diagnostic_field_configuration%dim_size_defns(diagnostic_field_configuration%dimensions) .eq. "qfields") then
1199  number_q_fields=c_get_integer(io_configuration%dimension_sizing, "qfields")
1200  do i=1, number_q_fields
1201  call add_specific_field_to_writer_entry(io_configuration, writer_entry_index, io_config_facet_index, &
1202  my_facet_index+i, trim(field_name)//"_udef"//trim(conv_to_string(i)), field_namespace, writer_field_names, &
1203  duplicate_field_names, c_get_integer(diagnostic_generation_frequency, field_name), diagnostic_field_configuration)
1204  end do
1205  tot_size=1
1206  do i=1, writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dimensions
1207  tot_size=tot_size*writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%actual_dim_size(i)
1208  end do
1209  call c_put_integer(q_field_splits, field_name, tot_size)
1210  add_field_to_writer_entry=number_q_fields
1211  return
1212  end if
1213  end if
1214  call add_specific_field_to_writer_entry(io_configuration, writer_entry_index, io_config_facet_index, &
1215  my_facet_index+1, field_name, field_namespace, writer_field_names, duplicate_field_names, &
1216  c_get_integer(diagnostic_generation_frequency, field_name), diagnostic_field_configuration)
1217  add_field_to_writer_entry=1
1218  else if (get_prognostic_field_configuration(io_configuration, field_name, field_namespace, &
1219  prognostic_field_configuration, prognostic_containing_data_defn)) then
1220  if (prognostic_field_configuration%field_type == array_field_type) then
1221  if (prognostic_field_configuration%dim_size_defns(prognostic_field_configuration%dimensions) .eq. "qfields") then
1222  number_q_fields=c_get_integer(io_configuration%dimension_sizing, "qfields")
1223  do i=1, number_q_fields
1224  call add_specific_field_to_writer_entry(io_configuration, writer_entry_index, io_config_facet_index, &
1225  my_facet_index+i, trim(field_name)//"_udef"//trim(conv_to_string(i)), field_namespace, writer_field_names, &
1226  duplicate_field_names, prognostic_containing_data_defn%frequency, &
1227  prognostic_field_configuration=prognostic_field_configuration)
1228  end do
1229  if (prognostic_field_configuration%collective) then
1230  allocate(collective_q_field)
1231  allocate(collective_q_field%dimensions(&
1232  writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dimensions))
1233  do i=1, writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dimensions
1234  if (trim(writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dim_size_defns(i)) == "z") then
1235  collective_q_field%dimensions(i)=1
1236  else if (trim(writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dim_size_defns(i)) &
1237  == "y") then
1238  collective_q_field%dimensions(i)=2
1239  else if (trim(writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dim_size_defns(i)) &
1240  == "x") then
1241  collective_q_field%dimensions(i)=3
1242  end if
1243  end do
1244  generic=>collective_q_field
1245  call c_put_generic(collective_q_field_dims, field_name, generic, .false.)
1246  else
1247  tot_size=1
1248  do i=1, writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%dimensions
1249  tot_size=tot_size*writer_entries(writer_entry_index)%contents(my_facet_index+number_q_fields)%actual_dim_size(i)
1250  end do
1251  call c_put_integer(q_field_splits, field_name, tot_size)
1252  end if
1253  call c_add_string(q_field_names, field_name)
1254  add_field_to_writer_entry=number_q_fields
1255  return
1256  end if
1257  end if
1258  call add_specific_field_to_writer_entry(io_configuration, writer_entry_index, io_config_facet_index, &
1259  my_facet_index+1, field_name, field_namespace, writer_field_names, duplicate_field_names, &
1260  prognostic_containing_data_defn%frequency, prognostic_field_configuration=prognostic_field_configuration)
1261  add_field_to_writer_entry=1
1262  else
1263  call log_log(log_error, "Field '"//trim(field_name)//&
1264  "' configured for file write but can not find this as a prognostic or diagnostic definition")
1265  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ add_group_of_fields_to_writer_entry()

integer function writer_federator_mod::add_group_of_fields_to_writer_entry ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  writer_entry_index,
integer, intent(in)  facet_index,
integer, intent(in)  current_field_index,
type(hashset_type), intent(inout)  writer_field_names,
type(hashset_type), intent(inout)  duplicate_field_names,
type(hashmap_type), intent(inout)  diagnostic_generation_frequency 
)
private

Adds a group of fields to a writer entry, groups are expanded out into individual fields, each inherit the properties of the group.

Parameters
io_configurationThe IO server configuration
writer_entry_indexIndex of the writer entry that we are dealing with
facet_indexIndex of the facet (group) in the IO server configuration
current_field_indexThe current field index in this internal module representation of the structure
Returns
The next field index to write to

Definition at line 1142 of file writer_federator.F90.

1144  type(io_configuration_type), intent(inout) :: io_configuration
1145  integer, intent(in) :: writer_entry_index, facet_index, current_field_index
1146  type(hashset_type), intent(inout) :: writer_field_names, duplicate_field_names
1147  type(hashmap_type), intent(inout) :: diagnostic_generation_frequency
1148 
1149  integer :: group_index
1150  character(len=STRING_LENGTH) :: field_name
1151  type(iterator_type) :: iterator
1152 
1153  add_group_of_fields_to_writer_entry=current_field_index
1154  group_index=get_index_of_group(io_configuration, &
1155  io_configuration%file_writers(writer_entry_index)%contents(facet_index)%facet_name)
1156  if (group_index == 0) then
1157  call log_log(log_error, "Can not find group '"//&
1158  trim(io_configuration%file_writers(writer_entry_index)%contents(facet_index)%facet_name)//"' in the configuration")
1159  end if
1160  iterator=c_get_iterator(io_configuration%groups(group_index)%members)
1161  do while (c_has_next(iterator))
1162  field_name=c_next_string(iterator)
1163  add_group_of_fields_to_writer_entry=add_group_of_fields_to_writer_entry+add_field_to_writer_entry(io_configuration, &
1164  writer_entry_index, facet_index, add_group_of_fields_to_writer_entry, field_name, &
1165  io_configuration%groups(group_index)%namespace, writer_field_names, duplicate_field_names, &
1166  diagnostic_generation_frequency)
1167  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ add_specific_field_to_writer_entry()

subroutine writer_federator_mod::add_specific_field_to_writer_entry ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  writer_entry_index,
integer, intent(in)  io_config_facet_index,
integer, intent(in)  my_facet_index,
character(len=*), intent(in)  field_name,
character(len=*), intent(in)  field_namespace,
type(hashset_type), intent(inout)  writer_field_names,
type(hashset_type), intent(inout)  duplicate_field_names,
integer, intent(in)  timestep_frequency,
type(io_configuration_diagnostic_field_type), intent(inout), optional  diagnostic_field_configuration,
type(io_configuration_field_type), intent(inout), optional  prognostic_field_configuration 
)
private

Adds a specific field and its information to a writer entry.

Parameters
io_configurationThe IO server configuration
writer_entry_indexIndex of the writer entry that we are dealing with
io_config_facet_indexIndex of the facet (group) in the IO server configuration
my_facet_indexThe current field index in this internal module representation of the structure
field_nameThe name of the field that we are constructing
writer_field_namesThe field names in the writer (for duplication checking)
duplicate_field_namesDuplicate field names in the wrier, for duplication checking
timestep_frequencyTimestepping frequency
diagnostic_field_configurationThe diagnostic field configuration (optional)
prognostic_field_configurationThe prognostic field configuration (optional)

Definition at line 1279 of file writer_federator.F90.

1282  type(io_configuration_type), intent(inout) :: io_configuration
1283  integer, intent(in) :: writer_entry_index, io_config_facet_index, my_facet_index, timestep_frequency
1284  character(len=*), intent(in) :: field_name, field_namespace
1285  type(hashset_type), intent(inout) :: writer_field_names, duplicate_field_names
1286  type(io_configuration_diagnostic_field_type), intent(inout), optional :: diagnostic_field_configuration
1287  type(io_configuration_field_type), intent(inout), optional :: prognostic_field_configuration
1288 
1289  integer :: i
1290 
1291  writer_entries(writer_entry_index)%contents(my_facet_index)%field_name=field_name
1292  writer_entries(writer_entry_index)%contents(my_facet_index)%field_namespace=field_namespace
1293 
1294  call c_add_string(used_field_names, field_name)
1295 
1296  if (.not. c_contains(writer_field_names, field_name)) then
1297  call c_add_string(writer_field_names, writer_entries(writer_entry_index)%contents(my_facet_index)%field_name)
1298  else
1299  call c_add_string(duplicate_field_names, writer_entries(writer_entry_index)%contents(my_facet_index)%field_name)
1300  end if
1301 
1302  if (io_configuration%file_writers(writer_entry_index)%contents(io_config_facet_index)%time_manipulation_type == &
1303  instantaneous_type) then
1304  writer_entries(writer_entry_index)%contents(my_facet_index)%time_manipulation=>perform_instantaneous_time_manipulation
1305  writer_entries(writer_entry_index)%contents(my_facet_index)%ready_to_write=>is_instantaneous_time_manipulation_ready_to_write
1306  else if (io_configuration%file_writers(writer_entry_index)%contents(io_config_facet_index)%time_manipulation_type == &
1307  time_averaged_type) then
1308  writer_entries(writer_entry_index)%contents(my_facet_index)%time_manipulation=>perform_timeaveraged_time_manipulation
1309  writer_entries(writer_entry_index)%contents(my_facet_index)%ready_to_write=>is_time_averaged_time_manipulation_ready_to_write
1310  else if (io_configuration%file_writers(writer_entry_index)%contents(io_config_facet_index)%time_manipulation_type == &
1311  none_type) then
1312  writer_entries(writer_entry_index)%contents(my_facet_index)%time_manipulation=>perform_none_time_manipulation
1313  writer_entries(writer_entry_index)%contents(my_facet_index)%ready_to_write=>is_none_time_manipulation_ready_to_write
1314  end if
1315  writer_entries(writer_entry_index)%contents(my_facet_index)%time_manipulation_type=&
1316  io_configuration%file_writers(writer_entry_index)%contents(io_config_facet_index)%time_manipulation_type
1317  writer_entries(writer_entry_index)%contents(my_facet_index)%output_frequency=&
1318  io_configuration%file_writers(writer_entry_index)%contents(io_config_facet_index)%output_time_frequency
1319  writer_entries(writer_entry_index)%contents(my_facet_index)%previous_write_time=0.0
1320  writer_entries(writer_entry_index)%contents(my_facet_index)%previous_tracked_write_point=0.0
1321  writer_entries(writer_entry_index)%contents(my_facet_index)%duplicate_field_name=.false.
1322  writer_entries(writer_entry_index)%contents(my_facet_index)%pending_to_write=.false.
1323  writer_entries(writer_entry_index)%contents(my_facet_index)%enabled=.false.
1324  writer_entries(writer_entry_index)%contents(my_facet_index)%expected_here=.true.
1325  writer_entries(writer_entry_index)%contents(my_facet_index)%prognostic_field=.false.
1326  writer_entries(writer_entry_index)%contents(my_facet_index)%diagnostic_field=.false.
1327 
1328  if (present(diagnostic_field_configuration)) then
1329  writer_entries(writer_entry_index)%contents(my_facet_index)%timestep_frequency=timestep_frequency
1330  writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions=diagnostic_field_configuration%dimensions
1331  writer_entries(writer_entry_index)%contents(my_facet_index)%data_type=diagnostic_field_configuration%data_type
1332  writer_entries(writer_entry_index)%contents(my_facet_index)%field_type=diagnostic_field_configuration%field_type
1333  writer_entries(writer_entry_index)%contents(my_facet_index)%dim_size_defns=diagnostic_field_configuration%dim_size_defns
1334  writer_entries(writer_entry_index)%contents(my_facet_index)%units=diagnostic_field_configuration%units
1335  writer_entries(writer_entry_index)%contents(my_facet_index)%collective_write=diagnostic_field_configuration%collective
1336  writer_entries(writer_entry_index)%contents(my_facet_index)%collective_initialised=.false.
1337  writer_entries(writer_entry_index)%contents(my_facet_index)%issue_write=.true.
1338  writer_entries(writer_entry_index)%contents(my_facet_index)%diagnostic_field=.true.
1339  else if (present(prognostic_field_configuration)) then
1340  writer_entries(writer_entry_index)%contents(my_facet_index)%timestep_frequency=timestep_frequency
1341  writer_entries(writer_entry_index)%contents(my_facet_index)%data_type=prognostic_field_configuration%data_type
1342  writer_entries(writer_entry_index)%contents(my_facet_index)%field_type=prognostic_field_configuration%field_type
1343  writer_entries(writer_entry_index)%contents(my_facet_index)%units=prognostic_field_configuration%units
1344  writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions=prognostic_field_configuration%dimensions
1345  writer_entries(writer_entry_index)%contents(my_facet_index)%collective_write=prognostic_field_configuration%collective
1346  writer_entries(writer_entry_index)%contents(my_facet_index)%collective_initialised=.false.
1347  writer_entries(writer_entry_index)%contents(my_facet_index)%prognostic_field=.true.
1348  if (.not. prognostic_field_configuration%collective) then
1349  writer_entries(writer_entry_index)%contents(my_facet_index)%issue_write=io_configuration%my_io_rank==0
1350  else
1351  writer_entries(writer_entry_index)%contents(my_facet_index)%issue_write=.true.
1352  end if
1353  if (prognostic_field_configuration%field_type == array_field_type .or. &
1354  prognostic_field_configuration%field_type == map_field_type) then
1355  if (prognostic_field_configuration%dimensions .gt. 0) then
1356  writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions=prognostic_field_configuration%dimensions
1357  writer_entries(writer_entry_index)%contents(my_facet_index)%dim_size_defns=&
1358  prognostic_field_configuration%dim_size_defns
1359  else
1360  call log_log(log_error, "The writing prognostic field '"//trim(field_name)//"' configuration must have dimensions")
1361  end if
1362  end if
1363  else
1364  call log_log(log_error, "A diagnostic or prognostic configuration for the field '"//trim(field_name)//"' was not found")
1365  end if
1366  if (writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions .gt. 0) then
1367  if (writer_entries(writer_entry_index)%contents(my_facet_index)%dim_size_defns(&
1368  writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions) .eq. "qfields") then
1369  writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions=&
1370  writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions-1
1371  end if
1372  do i=1, writer_entries(writer_entry_index)%contents(my_facet_index)%dimensions
1373  writer_entries(writer_entry_index)%contents(my_facet_index)%actual_dim_size(i)=c_get_integer(&
1374  io_configuration%dimension_sizing, writer_entries(writer_entry_index)%contents(my_facet_index)%dim_size_defns(i))
1375  end do
1376  end if
1377  call check_thread_status(forthread_mutex_init(writer_entries(writer_entry_index)%contents(my_facet_index)%values_mutex, -1))
Here is the caller graph for this function:

◆ check_for_and_issue_chain_write()

logical function writer_federator_mod::check_for_and_issue_chain_write ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_type), intent(inout)  writer_entry 
)
private

Will check whether there are any pending writes and if so will issue a chain write for this.

Parameters
io_configurationThe IO server configuration
writer_entryThe specific writer entry
Returns
Whether a chain write was issued or not

Definition at line 962 of file writer_federator.F90.

963  type(io_configuration_type), intent(inout) :: io_configuration
964  type(writer_type), intent(inout) :: writer_entry
965 
966  class(*), pointer :: generic
967 
968  call check_thread_status(forthread_mutex_lock(writer_entry%pending_writes_mutex))
969  if (.not. c_is_empty(writer_entry%pending_writes)) then
970  check_for_and_issue_chain_write=.true.
971  generic=>c_pop_generic(writer_entry%pending_writes)
972  call check_thread_status(forthread_mutex_unlock(writer_entry%pending_writes_mutex))
973  select type(generic)
974  type is (pending_write_type)
975  if (log_get_logging_level() .ge. log_debug) then
976  call log_log(log_debug, "Chain to next pending entry ts= "//trim(conv_to_string(generic%timestep)))
977  end if
978  call issue_actual_write(io_configuration, writer_entry, generic%timestep, &
979  generic%write_time, generic%terminated_write)
980  deallocate(generic)
981  end select
982  else
983  check_for_and_issue_chain_write=.false.
984  call check_thread_status(forthread_mutex_unlock(writer_entry%pending_writes_mutex))
985  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_writer_for_trigger()

subroutine, public writer_federator_mod::check_writer_for_trigger ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  source,
integer, intent(in)  data_id,
character, dimension(:), intent(in), allocatable  data_dump 
)

Checks all writer entries for any trigger fires and issues the underlying file storage.

Parameters
io_configurationConfiguration of the IO server
sourceThe source PID of the MONC process
data_idThe ID of the data definition that is represented by the dump
data_dumpThe data we have received from MONC

Definition at line 642 of file writer_federator.F90.

643  type(io_configuration_type), intent(inout) :: io_configuration
644  integer, intent(in) :: source, data_id
645  character, dimension(:), allocatable, intent(in) :: data_dump
646 
647  integer :: i, timestep
648  real(kind=default_precision) :: time
649  logical :: terminated
650 
651  if (is_field_present(io_configuration, source, data_id, "timestep") .and. &
652  is_field_present(io_configuration, source, data_id, "time")) then
653  timestep=get_scalar_integer_from_monc(io_configuration, source, data_id, data_dump, "timestep")
654  time=get_scalar_real_from_monc(io_configuration, source, data_id, data_dump, "time")
655 
656  if (is_field_present(io_configuration, source, data_id, "terminated")) then
657  terminated=get_scalar_logical_from_monc(io_configuration, source, data_id, data_dump, "terminated")
658  else
659  terminated=.false.
660  end if
661  do i=1, size(writer_entries)
662  call check_writer_trigger(io_configuration, i, timestep, real(time, kind=4), terminated)
663  end do
664  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_writer_trigger()

subroutine writer_federator_mod::check_writer_trigger ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  writer_entry_index,
integer, intent(in)  timestep,
real, intent(in)  time,
logical, intent(in)  terminated 
)
private

Checks a writer trigger and issues a file creation along with field write if the conditions (time or timestep) are met. This will either create and write to the file or store a pending state if one is already open (required due to NetCDF/HDF5 limitations with thread safety and parallel access.)

Parameters
io_configurationThe IO server configuration
writer_entry_indexIndex of the writer we are concerned with
timestepThe corresponding timestep
timeThe corresponding model time

Definition at line 674 of file writer_federator.F90.

675  type(io_configuration_type), intent(inout) :: io_configuration
676  integer, intent(in) :: writer_entry_index, timestep
677  real, intent(in) :: time
678  logical, intent(in) :: terminated
679 
680  real :: time_difference
681  integer :: i
682  logical :: issue_write, issue_terminated_write
683 
684  call check_thread_status(forthread_mutex_lock(writer_entries(writer_entry_index)%trigger_and_write_mutex))
685  issue_terminated_write=writer_entries(writer_entry_index)%write_on_terminate .and. terminated
686  if (writer_entries(writer_entry_index)%write_on_model_time) then
687  time_difference=time-writer_entries(writer_entry_index)%latest_pending_write_time
688  issue_write=time_difference .ge. writer_entries(writer_entry_index)%write_time_frequency
689  else
690  if (writer_entries(writer_entry_index)%write_timestep_frequency .gt. 0) then
691  issue_write=writer_entries(writer_entry_index)%latest_pending_write_timestep .ne. timestep .and. &
692  mod(timestep, writer_entries(writer_entry_index)%write_timestep_frequency) == 0
693  else
694  issue_write=.false.
695  end if
696  issue_terminated_write=issue_terminated_write .and. &
697  writer_entries(writer_entry_index)%latest_pending_write_timestep .ne. timestep
698  end if
699 
700  if (issue_write .or. issue_terminated_write) then
701  writer_entries(writer_entry_index)%latest_pending_write_time=time
702  writer_entries(writer_entry_index)%latest_pending_write_timestep=timestep
703 
704  call check_thread_status(forthread_mutex_lock(currently_writing_mutex))
705 
706  if (currently_writing) then
707  call check_thread_status(forthread_mutex_unlock(currently_writing_mutex))
708  call check_thread_status(forthread_mutex_unlock(writer_entries(writer_entry_index)%trigger_and_write_mutex))
709  call register_pending_file_write(writer_entry_index, timestep, time, &
710  writer_entries(writer_entry_index)%write_on_terminate .and. terminated)
711  else
712  currently_writing=.true.
713  call check_thread_status(forthread_mutex_unlock(currently_writing_mutex))
714  call check_thread_status(forthread_mutex_unlock(writer_entries(writer_entry_index)%trigger_and_write_mutex))
715  call issue_actual_write(io_configuration, writer_entries(writer_entry_index), timestep, time, &
716  writer_entries(writer_entry_index)%write_on_terminate .and. terminated)
717  end if
718  else
719  call check_thread_status(forthread_mutex_unlock(writer_entries(writer_entry_index)%trigger_and_write_mutex))
720  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ clean_time_points()

subroutine writer_federator_mod::clean_time_points
private

Cleans out old timepoints which are no longer going to be of any relavence to the file writing. This ensures that we don't have lots of stale points that need to be processed and searched beyond but are themselves pointless.

Definition at line 784 of file writer_federator.F90.

785  real :: time_entry
786  type(iterator_type) :: iterator
787  type(mapentry_type) :: map_entry
788  type(list_type) :: removed_entries
789  integer :: i
790  logical :: remove_timepoint
791 
792  call check_thread_status(forthread_rwlock_rdlock(time_points_rwlock))
793  iterator=c_get_iterator(time_points)
794  do while (c_has_next(iterator))
795  map_entry=c_next_mapentry(iterator)
796  time_entry=real(c_get_real(map_entry))
797  remove_timepoint=.true.
798  do i=1, size(writer_entries)
799  if (writer_entries(i)%previous_write_time .le. time_entry) then
800  remove_timepoint=.false.
801  if (writer_entries(i)%write_on_terminate) then
802  !print *, "Ignore CTP ", time_entry, writer_entries(i)%previous_write_time
803  end if
804  exit
805  end if
806  end do
807  if (remove_timepoint) call c_add_string(removed_entries, map_entry%key)
808  end do
809  iterator=c_get_iterator(removed_entries)
810  do while (c_has_next(iterator))
811  call c_remove(time_points, c_next_string(iterator))
812  end do
813  call check_thread_status(forthread_rwlock_unlock(time_points_rwlock))
814  call c_free(removed_entries)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ close_diagnostics_file()

subroutine writer_federator_mod::close_diagnostics_file ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_type), intent(inout)  writer_entry,
integer, intent(in)  timestep,
real, intent(in)  time 
)
private

Closes the diagnostics file, this is done via a global callback to issue the closes synchronously (collective operation)

Parameters
io_configurationThe IO server configuration
writer_entryThe writer entry
timestepThe file write timestep
timeThe file write time

Definition at line 880 of file writer_federator.F90.

881  type(io_configuration_type), intent(inout) :: io_configuration
882  type(writer_type), intent(inout) :: writer_entry
883  integer, intent(in) :: timestep
884  real, intent(in) :: time
885 
886  if (log_get_logging_level() .ge. log_debug) then
887  call log_log(log_debug, "Issue close for NetCDF file at timestep "//trim(conv_to_string(timestep)))
888  end if
889  call perform_global_callback(io_configuration, writer_entry%filename, timestep, handle_close_diagnostics_globalcallback)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ determine_collective_type_and_optimise_if_possible()

subroutine writer_federator_mod::determine_collective_type_and_optimise_if_possible ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_field_type), intent(inout)  field_to_write_information 
)
private

Determines whether it can optimise a specific collective field. If the field fits into certain limited parameters then it will optimise it. These parameters are very common, hence most fields can be optimised. Basically, it is looking to contiguous blocks of data from different MONCs so that the number of writes to the NetCDF file is limited.

Parameters
io_configurationThe IO server configuration
field_to_write_informationDescription of the the specific field that we are looking at here

Definition at line 1423 of file writer_federator.F90.

1424  type(io_configuration_type), intent(inout) :: io_configuration
1425  type(writer_field_type), intent(inout) :: field_to_write_information
1426 
1427  if (field_to_write_information%dimensions .eq. 3 .and. &
1428  get_dimension_identifier(field_to_write_information%dim_size_defns(1)) == z_index .and. &
1429  get_dimension_identifier(field_to_write_information%dim_size_defns(2)) == y_index .and. &
1430  get_dimension_identifier(field_to_write_information%dim_size_defns(3)) == x_index) then
1431  field_to_write_information%collective_contiguous_optimisation=.true.
1432  call initialise_contiguous_data_regions(io_configuration, field_to_write_information)
1433  else
1434  field_to_write_information%collective_contiguous_optimisation=.false.
1435  end if
1436  field_to_write_information%collective_initialised=.true.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ determine_if_field_can_be_written()

subroutine writer_federator_mod::determine_if_field_can_be_written ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_type), intent(inout)  writer_entry,
type(writer_field_type), intent(inout)  specific_field,
integer, intent(in)  timestep,
integer, intent(in)  previous_write_timestep,
real, intent(in)  write_time,
real, intent(in)  previous_write_time,
logical, intent(out), optional  field_written 
)
private

Determines if a file can be written to its overarching write representation. If so then a write is issued, otherwise an outstanding write point is registered which will be checked frequency to do a write later on.

Parameters
specific_fieldThe specific field we are checking and going to write if possible
timestepThe current timestep that we are at for this write
write_timeThe current time that we are at for this write
field_writtenAn optional output logical representing whether a write was performed or not

Definition at line 573 of file writer_federator.F90.

575  type(io_configuration_type), intent(inout) :: io_configuration
576  type(writer_type), intent(inout) :: writer_entry
577  type(writer_field_type), intent(inout) :: specific_field
578  integer, intent(in) :: timestep, previous_write_timestep
579  real, intent(in) :: write_time, previous_write_time
580  logical, intent(out), optional :: field_written
581 
582  real :: value_to_test, largest_value_found
583  integer :: num_matching
584  logical :: entry_beyond_this_write
585  type(iterator_type) :: iterator
586  type(mapentry_type) :: map_entry
587  type(write_field_collective_values_type), pointer :: multi_monc_entries
588  class(*), pointer :: generic
589 
590  num_matching=0
591  largest_value_found=0.0
592  entry_beyond_this_write=.false.
593  call check_thread_status(forthread_mutex_lock(specific_field%values_mutex))
594  if (.not. c_is_empty(specific_field%values_to_write)) then
595  iterator=c_get_iterator(specific_field%values_to_write)
596  do while (c_has_next(iterator))
597  map_entry=c_next_mapentry(iterator)
598  value_to_test=conv_to_real(map_entry%key)
599  if (specific_field%collective_write) then
600  generic=>c_get_generic(map_entry)
601  select type(generic)
602  type is(write_field_collective_values_type)
603  multi_monc_entries=>generic
604  end select
605  if (c_size(multi_monc_entries%monc_values) .ne. io_configuration%number_of_moncs) cycle
606  end if
607  if (value_to_test .gt. write_time) entry_beyond_this_write=.true.
608  if (value_to_test .le. write_time .and. value_to_test .gt. previous_write_time) then
609  num_matching=num_matching+1
610  if (largest_value_found .lt. value_to_test) largest_value_found=value_to_test
611  end if
612  end do
613  end if
614 
615  if (num_matching .gt. 0 .and. (specific_field%ready_to_write(largest_value_found, specific_field%output_frequency, write_time, &
616  specific_field%latest_timestep_values, timestep) .or. entry_beyond_this_write)) then
617  if (.not. specific_field%collective_write .or. .not. specific_field%collective_contiguous_optimisation) then
618  if (specific_field%issue_write) then
619  call write_variable(io_configuration, specific_field, writer_entry%filename, timestep, write_time)
620  end if
621  specific_field%previous_write_time=writer_entry%write_time
622  end if
623  specific_field%pending_to_write=.false.
624  if (present(field_written)) field_written=.true.
625  else
626  if (log_get_logging_level() .ge. log_debug) then
627  call log_log(log_debug, "Setting outstanding field ts="//conv_to_string(writer_entry%write_timestep)//&
628  " write time="//conv_to_string(writer_entry%write_time)//" prev="//conv_to_string(previous_write_time)//&
629  " largest entry="//conv_to_string(largest_value_found)//" num matching="//conv_to_string(num_matching))
630  end if
631  specific_field%pending_to_write=.true.
632  if (present(field_written)) field_written=.false.
633  end if
634  call check_thread_status(forthread_mutex_unlock(specific_field%values_mutex))
Here is the call graph for this function:
Here is the caller graph for this function:

◆ determine_if_outstanding_field_can_be_written()

subroutine writer_federator_mod::determine_if_outstanding_field_can_be_written ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_type), intent(inout)  writer_entry,
type(writer_field_type), intent(inout)  specific_field 
)
private

For a specific field wil determine and handle any outstanding fields writes until an outstanding write can not be performed or the outstanding list is empty.

Parameters
specific_fieldThe specific field that we are concerned with

Definition at line 541 of file writer_federator.F90.

542  type(io_configuration_type), intent(inout) :: io_configuration
543  type(writer_type), intent(inout) :: writer_entry
544  type(writer_field_type), intent(inout) :: specific_field
545 
546  logical :: field_write_success, do_close_num_fields
547 
548  if (specific_field%pending_to_write) then
549  call determine_if_field_can_be_written(io_configuration, writer_entry, specific_field, writer_entry%write_timestep, &
550  writer_entry%previous_write_timestep, writer_entry%write_time, writer_entry%previous_write_time, field_write_success)
551  if (field_write_success) then
552  if (log_get_logging_level() .ge. log_debug) then
553  call log_log(log_debug, "Flushed outstanding field ts="//conv_to_string(writer_entry%write_timestep)//&
554  " write time="//conv_to_string(writer_entry%write_time))
555  end if
556  call check_thread_status(forthread_mutex_lock(writer_entry%num_fields_to_write_mutex))
557  writer_entry%num_fields_to_write=writer_entry%num_fields_to_write-1
558  do_close_num_fields=writer_entry%num_fields_to_write == 0
559  call check_thread_status(forthread_mutex_unlock(writer_entry%num_fields_to_write_mutex))
560  if (do_close_num_fields) then
561  call close_diagnostics_file(io_configuration, writer_entry, writer_entry%write_timestep, writer_entry%write_time)
562  end if
563  end if
564  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ enable_specific_field_by_name()

subroutine writer_federator_mod::enable_specific_field_by_name ( character(len=*), intent(in)  field_name,
logical, intent(in)  diagnostics_mode,
logical, intent(in), optional  expected_here 
)
private

Enables a specific field by its name, this will locate all the fields with this name and enable them.

Parameters
field_nameThe name of the field to enable

Definition at line 247 of file writer_federator.F90.

248  character(len=*), intent(in) :: field_name
249  logical, intent(in) :: diagnostics_mode
250  logical, intent(in), optional :: expected_here
251 
252  logical :: continue_search
253  integer :: writer_index, contents_index
254 
255  continue_search=.true.
256  writer_index=1
257  contents_index=0
258  do while (continue_search)
259  contents_index=contents_index+1
260  continue_search=get_next_applicable_writer_entry(field_name, writer_index_point=writer_index, &
261  contents_index_point=contents_index)
262  if (continue_search) then
263  if ((writer_entries(writer_index)%contents(contents_index)%diagnostic_field .and. diagnostics_mode) .or. &
264  (writer_entries(writer_index)%contents(contents_index)%prognostic_field .and. .not. diagnostics_mode)) then
265  writer_entries(writer_index)%contents(contents_index)%enabled=.true.
266  if (present(expected_here)) then
267  writer_entries(writer_index)%contents(contents_index)%expected_here=expected_here
268  end if
269  end if
270  end if
271  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ extract_applicable_time_points()

type(map_type) function writer_federator_mod::extract_applicable_time_points ( real, intent(in)  start_time,
real, intent(in)  end_time 
)
private

Extracts the applicable time points from the overall map that lie within a specific range.

Parameters
start_timeThe start time where values must be greater than
end_timeThe end time where values must be equal or less than
Returns
The sorted (based on timestep) list of MONC communication time points that correspond to this

Definition at line 822 of file writer_federator.F90.

823  real, intent(in) :: start_time, end_time
824 
825  real :: time_entry
826  type(iterator_type) :: iterator
827  type(mapentry_type) :: map_entry
828 
829  call check_thread_status(forthread_rwlock_rdlock(time_points_rwlock))
830  iterator=c_get_iterator(time_points)
831  do while (c_has_next(iterator))
832  map_entry=c_next_mapentry(iterator)
833  time_entry=real(c_get_real(map_entry))
834  if (time_entry .gt. start_time .and. time_entry .le. end_time) then
835  call c_put_real(extract_applicable_time_points, map_entry%key, conv_single_real_to_double(time_entry))
836  end if
837  end do
838  call check_thread_status(forthread_rwlock_unlock(time_points_rwlock))
839  extract_applicable_time_points=sort_applicable_time_points(extract_applicable_time_points)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalise_writer_federator()

subroutine, public writer_federator_mod::finalise_writer_federator

Finalises the write federator and the manipulations.

Definition at line 122 of file writer_federator.F90.

123  call check_thread_status(forthread_rwlock_destroy(time_points_rwlock))
124  call check_thread_status(forthread_mutex_destroy(collective_contiguous_initialisation_mutex))
125  call check_thread_status(forthread_mutex_destroy(currently_writing_mutex))
126  call finalise_time_averaged_manipulation()
127  call finalise_instantaneous_manipulation()
128  call finalise_netcdf_filetype()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_common_starts()

subroutine writer_federator_mod::get_common_starts ( integer, intent(in)  dim,
integer, intent(in)  val,
integer, dimension(:,:), intent(in)  vals,
integer, dimension(:), intent(out)  common_starters,
integer, intent(out)  num_common 
)
private

Retrieves the number of common starting points that match a specific input value.

Parameters
dimThe dimension we are searching in
valThe value that we are looking to match against
valsAll the other dimensions information to search against
common_startersThe common starting point locations
num_commonThe number of common starting points identified

Definition at line 1549 of file writer_federator.F90.

1550  integer, intent(in) :: dim, val
1551  integer, dimension(:,:), intent(in) :: vals
1552  integer, dimension(:), intent(out) :: common_starters
1553  integer, intent(out) :: num_common
1554 
1555  integer :: i
1556 
1557  num_common=0
1558  do i=1, size(vals, 2)
1559  if (vals(dim, i) == val) then
1560  num_common=num_common+1
1561  common_starters(num_common)=i
1562  end if
1563  end do
Here is the caller graph for this function:

◆ get_dimension_identifier()

integer function writer_federator_mod::get_dimension_identifier ( character(len=*), intent(in)  dim_name,
logical, intent(out), optional  is_auto_dimension 
)
private

Translates a dimension name to its numeric corresponding identifier.

Parameters
dim_nameThe name of the dimension to look up
is_auto_dimensionOptional parameter determining whether dimension is auto or not
Returns
Corresponding identifier

Definition at line 1570 of file writer_federator.F90.

1571  character(len=*), intent(in) :: dim_name
1572  logical, intent(out), optional :: is_auto_dimension
1573 
1574  integer :: dash_idx
1575  logical :: is_modified_size
1576 
1577  dash_idx=index(dim_name, "_")
1578  dash_idx=dash_idx-1
1579  is_modified_size=dash_idx .ne. -1
1580  if (.not. is_modified_size) dash_idx=len_trim(dim_name)
1581 
1582  if (dim_name(:dash_idx) .eq. "z" .or. dim_name(:dash_idx) .eq. "zn") then
1583  get_dimension_identifier=z_index
1584  else if (dim_name(:dash_idx) .eq. "y") then
1585  get_dimension_identifier=y_index
1586  else if (dim_name(:dash_idx) .eq. "x") then
1587  get_dimension_identifier=x_index
1588  else
1589  get_dimension_identifier=-1
1590  end if
1591 
1592  if (present(is_auto_dimension)) is_auto_dimension=is_modified_size
Here is the caller graph for this function:

◆ get_field_number_of_fields()

integer function writer_federator_mod::get_field_number_of_fields ( type(io_configuration_type), intent(inout)  io_configuration,
character(len=string_length), intent(in)  field_name,
character(len=string_length), intent(in)  field_namespace,
integer, intent(in)  num_q_fields 
)
private

Retrieves the number of fields that make up this field, if it is a Q field then it will be split into many subfields hence it is not a simple 1-1 mapping.

Parameters
io_configurationThe IO server configuration
field_nameThe name of the field
num_q_fieldsThe number of Q fields
Returns
The number of fields that make up this field

Definition at line 1106 of file writer_federator.F90.

1107  type(io_configuration_type), intent(inout) :: io_configuration
1108  character(len=STRING_LENGTH), intent(in) :: field_name, field_namespace
1109  integer, intent(in) :: num_q_fields
1110 
1111  type(io_configuration_field_type) :: prognostic_field_configuration
1112  type(io_configuration_data_definition_type) :: prognostic_containing_data_defn
1113  type(io_configuration_diagnostic_field_type) :: diagnostic_field_configuration
1114 
1115  if (get_diagnostic_field_configuration(io_configuration, field_name, field_namespace, diagnostic_field_configuration)) then
1116  if (diagnostic_field_configuration%field_type == array_field_type) then
1117  if (diagnostic_field_configuration%dim_size_defns(diagnostic_field_configuration%dimensions) .eq. "qfields") then
1118  get_field_number_of_fields=num_q_fields
1119  return
1120  end if
1121  end if
1122  get_field_number_of_fields=1
1123  else if (get_prognostic_field_configuration(io_configuration, field_name, field_namespace, &
1124  prognostic_field_configuration, prognostic_containing_data_defn)) then
1125  if (prognostic_field_configuration%field_type == array_field_type) then
1126  if (prognostic_field_configuration%dim_size_defns(prognostic_field_configuration%dimensions) .eq. "qfields") then
1127  get_field_number_of_fields=num_q_fields
1128  return
1129  end if
1130  end if
1131  get_field_number_of_fields=1
1132  end if
Here is the caller graph for this function:

◆ get_group_number_of_fields()

integer function writer_federator_mod::get_group_number_of_fields ( type(io_configuration_type), intent(inout)  io_configuration,
type(list_type group_members,
integer, intent(in)  num_q_fields,
character(len=string_length), intent(in)  namespace 
)
private

Retrieves the number of fields within a group of fields.

Parameters
io_configurationThe IO server configuration
group_membersThe members of the group
num_q_fieldsThe number of Q fields
Returns
The number of fields that make up this group

Definition at line 1082 of file writer_federator.F90.

1083  type(io_configuration_type), intent(inout) :: io_configuration
1084  type(list_type) :: group_members
1085  integer, intent(in) :: num_q_fields
1086  character(len=STRING_LENGTH), intent(in) :: namespace
1087 
1088  type(iterator_type) :: iterator
1089  character(len=STRING_LENGTH) :: field_name
1090 
1091  get_group_number_of_fields=0
1092  iterator=c_get_iterator(group_members)
1093  do while (c_has_next(iterator))
1094  field_name=c_next_string(iterator)
1095  get_group_number_of_fields=get_group_number_of_fields+get_field_number_of_fields(io_configuration, field_name, namespace, &
1096  num_q_fields)
1097  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_index_of_group()

integer function writer_federator_mod::get_index_of_group ( type(io_configuration_type), intent(inout)  io_configuration,
character(len=*), intent(in)  group_name 
)
private

Searches the IO server configuration for a group with a specific name and returns the index to that group or 0 if no corresponding group is found.

Parameters
io_configurationThe IO server configuration
group_nameThe group name to look up
Returns
The index to this group or 0 if no group is found

Definition at line 1402 of file writer_federator.F90.

1403  type(io_configuration_type), intent(inout) :: io_configuration
1404  character(len=*), intent(in) :: group_name
1405 
1406  integer :: i, entries
1407 
1408  entries=io_configuration%number_of_groups
1409  do i=1, entries
1410  if (io_configuration%groups(i)%name == group_name) then
1411  get_index_of_group=i
1412  return
1413  end if
1414  end do
1415  get_index_of_group=0
Here is the caller graph for this function:

◆ get_next_applicable_writer_entry()

logical function writer_federator_mod::get_next_applicable_writer_entry ( character(len=*), intent(in)  field_name,
character(len=*), intent(in), optional  field_namespace,
integer, intent(inout)  writer_index_point,
integer, intent(inout)  contents_index_point 
)
private

Retrieves the index of the next writer which uses a specific field. If none is found then returns false, otherwise true.

Parameters
field_nameThe field name to search for
writer_start_pointThe index that we start searching from
result_writerThe index of the (next) writer that requires this field is written in here
result_contentsThe index of the contents in the next writer that corresponds to this field
Returns
Whether or not a next entry has been found

Definition at line 1017 of file writer_federator.F90.

1018  character(len=*), intent(in) :: field_name
1019  character(len=*), intent(in), optional :: field_namespace
1020  integer, intent(inout) :: writer_index_point, contents_index_point
1021 
1022  integer :: i, j
1023 
1024  if (writer_index_point .le. size(writer_entries)) then
1025  do i=writer_index_point, size(writer_entries)
1026  if (contents_index_point .le. size(writer_entries(i)%contents)) then
1027  do j=contents_index_point, size(writer_entries(i)%contents)
1028  if (writer_entries(i)%contents(j)%field_name==field_name) then
1029  if (present(field_namespace)) then
1030  if (writer_entries(i)%contents(j)%field_namespace .ne. field_namespace) cycle
1031  end if
1032  writer_index_point=i
1033  contents_index_point=j
1034  get_next_applicable_writer_entry=.true.
1035  return
1036  end if
1037  end do
1038  end if
1039  contents_index_point=1
1040  end do
1041  end if
1042  get_next_applicable_writer_entry=.false.
Here is the caller graph for this function:

◆ get_size_of_collective_q()

integer function writer_federator_mod::get_size_of_collective_q ( type(io_configuration_type), intent(inout)  io_configuration,
character(len=*), intent(in)  field_name,
integer, intent(in)  source 
)
private

Retrieves the data size for each Q entry of a collective Q field for the specific source MONC that has sent data.

Parameters
io_configurationThe IO server configuration
field_nameThe field name to write (if appropriate)
sourceMONC source for the communicated fields
Returns
The size (elements) per Q split field

Definition at line 410 of file writer_federator.F90.

411  type(io_configuration_type), intent(inout) :: io_configuration
412  character(len=*), intent(in) :: field_name
413  integer, intent(in) :: source
414 
415  class(*), pointer :: generic
416  integer :: i, monc_index
417 
418  get_size_of_collective_q=1
419  monc_index=get_monc_location(io_configuration, source)
420  generic=>c_get_generic(collective_q_field_dims, field_name)
421  select type(generic)
422  type is(collective_q_field_representation_type)
423  do i=1, size(generic%dimensions)
424  get_size_of_collective_q=&
425  get_size_of_collective_q*io_configuration%registered_moncs(monc_index)%local_dim_sizes(generic%dimensions(i))
426  end do
427  end select
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_total_number_writer_fields()

integer function writer_federator_mod::get_total_number_writer_fields ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  writer_entry_index 
)
private

Determines the total number of fields that make up a writer entry, this is all the fields of the groups that make up this writer and individual fields specified too.

Parameters
io_configurationThe IO server configuration
writer_entry_indexThe index of the writer that we are inquiring about
Returns
The total number of fields that make up this writer

Definition at line 1050 of file writer_federator.F90.

1051  type(io_configuration_type), intent(inout) :: io_configuration
1052  integer, intent(in) :: writer_entry_index
1053 
1054  integer :: i, number_contents, group_index, number_q_fields
1055 
1056  get_total_number_writer_fields=0
1057  number_q_fields=c_get_integer(io_configuration%dimension_sizing, "qfields")
1058 
1059  number_contents=io_configuration%file_writers(writer_entry_index)%number_of_contents
1060  do i=1, number_contents
1061  if (io_configuration%file_writers(writer_entry_index)%contents(i)%facet_type == group_type) then
1062  group_index=get_index_of_group(io_configuration, &
1063  io_configuration%file_writers(writer_entry_index)%contents(i)%facet_name)
1064  if (group_index == 0) call log_log(log_error, "Can not find group '"//trim(&
1065  io_configuration%file_writers(writer_entry_index)%contents(i)%facet_name)//"'")
1066  get_total_number_writer_fields=get_total_number_writer_fields+&
1067  get_group_number_of_fields(io_configuration, io_configuration%groups(group_index)%members, number_q_fields, &
1068  io_configuration%groups(group_index)%namespace)
1069  else if (io_configuration%file_writers(writer_entry_index)%contents(i)%facet_type == field_type) then
1070  ! NSE
1071  get_total_number_writer_fields=get_total_number_writer_fields+get_field_number_of_fields(io_configuration, &
1072  io_configuration%file_writers(writer_entry_index)%contents(i)%facet_name, "", number_q_fields)
1073  end if
1074  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ handle_close_diagnostics_globalcallback()

subroutine writer_federator_mod::handle_close_diagnostics_globalcallback ( type(io_configuration_type), intent(inout)  io_configuration,
real(default_precision), dimension(:)  values,
character(len=string_length)  field_name,
integer  timestep 
)
private

Call back for the inter IO reduction which actually does the NetCDF file closing which is a collective (synchronous) operation. Calls out to the NetCDF code to do the call and then checks the list of pending file writes to process any others that are waiting in the queue.

Parameters
io_configurationThe IO server configuration
valuesThe inter IO resulting values, we don't care about these
field_nameThe field name that is being communicated
timestepThe write timestep

Definition at line 899 of file writer_federator.F90.

900  type(io_configuration_type), intent(inout) :: io_configuration
901  real(DEFAULT_PRECISION), dimension(:) :: values
902  character(len=STRING_LENGTH) :: field_name
903  integer :: timestep
904 
905  type(writer_type), pointer :: writer_entry
906  integer :: i
907  logical :: terminated, done_chain_run
908 
909  writer_entry=>get_writer_entry_from_netcdf(field_name, timestep, terminated)
910 
911  do i=1, size(writer_entry%contents)
912  if (writer_entry%contents(i)%enabled .and. writer_entry%contents(i)%collective_write .and. &
913  writer_entry%contents(i)%collective_contiguous_optimisation) then
914  call check_thread_status(forthread_mutex_lock(writer_entry%contents(i)%values_mutex))
915  call write_variable(io_configuration, writer_entry%contents(i), writer_entry%filename, timestep, writer_entry%write_time)
916  writer_entry%contents(i)%previous_write_time=writer_entry%write_time
917  call check_thread_status(forthread_mutex_unlock(writer_entry%contents(i)%values_mutex))
918  end if
919  end do
920 
921  writer_entry%previous_write_time=writer_entry%write_time
922  writer_entry%previous_write_timestep=writer_entry%write_timestep
923  writer_entry%defined_write_time=writer_entry%defined_write_time+writer_entry%write_time_frequency
924 
925  call clean_time_points()
926 
927  if (writer_entry%contains_io_status_dump) then
928  if (.not. terminated) then
929  do while (.not. is_io_server_state_writer_ready(timestep))
930  end do
931  end if
932  call check_thread_status(forthread_rwlock_rdlock(time_points_rwlock))
933  call store_io_server_state(io_configuration, writer_entries, time_points, writer_entry, timestep)
934  call check_thread_status(forthread_rwlock_unlock(time_points_rwlock))
935  end if
936 
937  call close_netcdf_file(io_configuration, field_name, timestep)
938 
939  done_chain_run=.false.
940  do i=1, size(writer_entries)
941  if (writer_entries(i)%filename .ne. writer_entry%filename) then
942  done_chain_run=check_for_and_issue_chain_write(io_configuration, writer_entries(i))
943  if (done_chain_run) exit
944  end if
945  end do
946  if (.not. done_chain_run) done_chain_run=check_for_and_issue_chain_write(io_configuration, writer_entry)
947 
948  if (.not. done_chain_run) then
949  call check_thread_status(forthread_mutex_lock(currently_writing_mutex))
950  currently_writing=.false.
951  call check_thread_status(forthread_mutex_unlock(currently_writing_mutex))
952  if (log_get_logging_level() .ge. log_debug) then
953  call log_log(log_debug, "No more pending entries to chain to at ts= "//trim(conv_to_string(timestep)))
954  end if
955  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ handle_duplicate_field_names()

subroutine writer_federator_mod::handle_duplicate_field_names ( type(writer_type), intent(inout)  writer_entry,
type(hashset_type), intent(inout)  duplicate_field_names 
)
private

Marks duplicate field names in a writer entry as duplicates so that the NetCDF layer can then deal with this by issuing unique names.

Parameters
writer_entryThe writer state
duplicate_field_namesA hashset of duplicate field names which are marked

Definition at line 1384 of file writer_federator.F90.

1385  type(writer_type), intent(inout) :: writer_entry
1386  type(hashset_type), intent(inout) :: duplicate_field_names
1387 
1388  integer :: i
1389 
1390  do i=1, size(writer_entry%contents)
1391  if (c_contains(duplicate_field_names, writer_entry%contents(i)%field_name)) then
1392  writer_entry%contents(i)%duplicate_field_name=.true.
1393  end if
1394  end do
Here is the caller graph for this function:

◆ inform_writer_federator_fields_present()

subroutine, public writer_federator_mod::inform_writer_federator_fields_present ( type(io_configuration_type), intent(inout)  io_configuration,
type(hashset_type), intent(inout), optional  field_names,
type(hashmap_type), intent(inout), optional  diag_field_names_and_roots 
)

Informs the writer federator that specific fields are present and should be reflected in the diagnostics output.

Parameters
field_namesThe set of field names that are present

Definition at line 161 of file writer_federator.F90.

162  type(io_configuration_type), intent(inout) :: io_configuration
163  type(hashset_type), intent(inout), optional :: field_names
164  type(hashmap_type), intent(inout), optional :: diag_field_names_and_roots
165 
166  type(iterator_type) :: iterator
167  character(len=STRING_LENGTH) :: specific_name
168  integer :: i, number_q_fields, expected_io
169  logical :: field_found, expected_here, diagnostics_mode
170 
171  iterator=c_get_iterator(used_field_names)
172  expected_io=-1
173  do while (c_has_next(iterator))
174  specific_name=c_next_string(iterator)
175  if (present(field_names)) then
176  field_found=c_contains(field_names, specific_name)
177  diagnostics_mode=.false.
178  else if (present(diag_field_names_and_roots)) then
179  field_found=c_contains(diag_field_names_and_roots, specific_name)
180  if (field_found) expected_io=c_get_integer(diag_field_names_and_roots, specific_name)
181  diagnostics_mode=.true.
182  else
183  field_found=.false.
184  end if
185  if (field_found) then
186  expected_here=expected_io == -1 .or. expected_io == io_configuration%my_io_rank
187  call enable_specific_field_by_name(specific_name, diagnostics_mode, expected_here)
188  end if
189  end do
190  iterator=c_get_iterator(q_field_names)
191  do while (c_has_next(iterator))
192  specific_name=c_next_string(iterator)
193  if (present(field_names)) then
194  field_found=c_contains(field_names, specific_name)
195  diagnostics_mode=.false.
196  else if (present(diag_field_names_and_roots)) then
197  field_found=c_contains(diag_field_names_and_roots, specific_name)
198  if (field_found) expected_io=c_get_integer(diag_field_names_and_roots, specific_name)
199  diagnostics_mode=.true.
200  else
201  field_found=.false.
202  end if
203  if (field_found) then
204  expected_here=expected_io == -1 .or. expected_io == io_configuration%my_io_rank
205  number_q_fields=c_get_integer(io_configuration%dimension_sizing, "qfields")
206  do i=1, number_q_fields
207  if (c_size(io_configuration%q_field_names) .ge. i) then
208  call enable_specific_field_by_name(trim(specific_name)//"_"//trim(c_get_string(io_configuration%q_field_names, i)), &
209  diagnostics_mode, expected_here)
210  else
211  call enable_specific_field_by_name(trim(specific_name)//"_udef"//trim(conv_to_string(i)), &
212  diagnostics_mode, expected_here)
213  end if
214  end do
215  end if
216  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ inform_writer_federator_time_point()

subroutine, public writer_federator_mod::inform_writer_federator_time_point ( type(io_configuration_type), intent(inout)  io_configuration,
integer, intent(in)  source,
integer, intent(in)  data_id,
character, dimension(:), intent(in), allocatable  data_dump 
)

Definition at line 131 of file writer_federator.F90.

132  type(io_configuration_type), intent(inout) :: io_configuration
133  integer, intent(in) :: source, data_id
134  character, dimension(:), allocatable, intent(in) :: data_dump
135 
136  real(kind=default_precision) :: time
137  integer :: timestep
138  character(len=STRING_LENGTH) :: timestep_key
139 
140  if (is_field_present(io_configuration, source, data_id, "time") .and. &
141  is_field_present(io_configuration, source, data_id, "timestep")) then
142  time=get_scalar_real_from_monc(io_configuration, source, data_id, data_dump, "time")
143  timestep=get_scalar_integer_from_monc(io_configuration, source, data_id, data_dump, "timestep")
144 
145  timestep_key=conv_to_string(timestep)
146 
147  call check_thread_status(forthread_rwlock_rdlock(time_points_rwlock))
148  if (.not. c_contains(time_points, timestep_key)) then
149  call check_thread_status(forthread_rwlock_unlock(time_points_rwlock))
150  call check_thread_status(forthread_rwlock_wrlock(time_points_rwlock))
151  if (.not. c_contains(time_points, timestep_key)) then
152  call c_put_real(time_points, timestep_key, time)
153  end if
154  end if
155  call check_thread_status(forthread_rwlock_unlock(time_points_rwlock))
156  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initialise_contiguous_data_regions()

subroutine writer_federator_mod::initialise_contiguous_data_regions ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_field_type), intent(inout)  field_to_write_information 
)
private

Will initialise the collective data regions that form contiguous blocks within the data. This is quite an expensive operation so only done once for each field, but has the potential for very significant performance advantages for the fields that match it.

Parameters
io_configurationThe IO server configuration
field_to_write_informationThe specific field that is being looked at and identified for contiguous data

Definition at line 1444 of file writer_federator.F90.

1445  type(io_configuration_type), intent(inout) :: io_configuration
1446  type(writer_field_type), intent(inout) :: field_to_write_information
1447 
1448  integer :: start(field_to_write_information%dimensions, io_configuration%number_of_moncs), &
1449  count(field_to_write_information%dimensions, io_configuration%number_of_moncs), &
1450  common_starters(io_configuration%number_of_moncs), num_common, num_current_contents, active_dim, other_dim, &
1451  j, k, i, dim_identifier, number_distinct_writes, start_blocks(io_configuration%number_of_moncs), ierr, &
1452  count_blocks(io_configuration%number_of_moncs), current_contents(io_configuration%number_of_moncs), &
1453  monc_write_start_offset_per_dim(field_to_write_information%dimensions,io_configuration%number_of_moncs)
1454  logical :: processed(io_configuration%number_of_moncs)
1455 
1456  type(write_field_collective_descriptor_type), pointer :: collective_descriptor
1457  type(write_field_collective_monc_info_type), pointer :: specific_monc_collective
1458  class(*), pointer :: generic
1459 
1460  processed=.false.
1461  number_distinct_writes=0
1462  do j=1, io_configuration%number_of_moncs
1463  do k=1, field_to_write_information%dimensions
1464  dim_identifier=get_dimension_identifier(field_to_write_information%dim_size_defns(k))
1465  start(k, j)=io_configuration%registered_moncs(j)%local_dim_starts(dim_identifier)
1466  count(k, j)=io_configuration%registered_moncs(j)%local_dim_sizes(dim_identifier)
1467  end do
1468  end do
1469 
1470  do j=1, io_configuration%number_of_moncs
1471  if (.not. processed(j)) then
1472  call get_common_starts(y_index, start(y_index, j), start, common_starters, num_common)
1473  if (num_common == 0) then
1474  call get_common_starts(x_index, start(x_index, j), start, common_starters, num_common)
1475  if (num_common .gt. 0) then
1476  active_dim=x_index
1477  other_dim=y_index
1478  end if
1479  else
1480  active_dim=y_index
1481  other_dim=x_index
1482  end if
1483  number_distinct_writes=number_distinct_writes+1
1484  allocate(collective_descriptor)
1485  allocate(collective_descriptor%absolute_start(field_to_write_information%dimensions), &
1486  collective_descriptor%count(field_to_write_information%dimensions))
1487  start_blocks(number_distinct_writes)=start(other_dim, j)
1488  count_blocks(number_distinct_writes)=count(other_dim, j)
1489  num_current_contents=1
1490  current_contents(num_current_contents)=j
1491  processed(j)=.true.
1492  if (num_common .gt. 0) then
1493  do k=1, num_common
1494  do i=1, num_common
1495  if (.not. processed(common_starters(i)) .and. count(active_dim, j) == count(active_dim, i)) then
1496  if (start(other_dim, common_starters(i)) .lt. start_blocks(number_distinct_writes) .and. &
1497  start(other_dim, common_starters(i)) + count(other_dim, common_starters(i)) &
1498  == start_blocks(number_distinct_writes)) then
1499  start_blocks(number_distinct_writes)=start(other_dim, common_starters(i))
1500  count_blocks(number_distinct_writes)=count_blocks(number_distinct_writes)+count(other_dim, common_starters(i))
1501  processed(common_starters(i))=.true.
1502  num_current_contents=num_current_contents+1
1503  current_contents(num_current_contents)=common_starters(i)
1504  else if (start(other_dim, common_starters(i)) .gt. start_blocks(number_distinct_writes) .and. &
1505  start_blocks(number_distinct_writes) + count_blocks(number_distinct_writes) &
1506  == start(other_dim, common_starters(i))) then
1507  count_blocks(number_distinct_writes)=count_blocks(number_distinct_writes)+count(other_dim, common_starters(i))
1508  processed(common_starters(i))=.true.
1509  num_current_contents=num_current_contents+1
1510  current_contents(num_current_contents)=common_starters(i)
1511  end if
1512  end if
1513  end do
1514  end do
1515  end if
1516  collective_descriptor%absolute_start=start(:,j)
1517  collective_descriptor%count=count(:,j)
1518  collective_descriptor%absolute_start(other_dim)=start_blocks(number_distinct_writes)
1519  collective_descriptor%count(other_dim)=count_blocks(number_distinct_writes)
1520  collective_descriptor%split_dim=other_dim
1521  if (num_current_contents .gt. 0) then
1522  do k=1, num_current_contents
1523  allocate(specific_monc_collective)
1524  specific_monc_collective%relative_dimension_start=(start(other_dim,current_contents(k))-&
1525  start_blocks(number_distinct_writes)) + 1
1526  specific_monc_collective%counts=count(:, current_contents(k))
1527  specific_monc_collective%monc_location=current_contents(k)
1528  specific_monc_collective%monc_source=io_configuration%registered_moncs(current_contents(k))%source_id
1529  generic=>specific_monc_collective
1530  call c_add_generic(collective_descriptor%specific_monc_info, generic, .false.)
1531  end do
1532  end if
1533  generic=>collective_descriptor
1534  call c_add_generic(field_to_write_information%collective_descriptors, generic, .false.)
1535  end if
1536  end do
1537  call lock_mpi()
1538  call mpi_iallreduce(number_distinct_writes, field_to_write_information%max_num_collective_writes, 1, mpi_int, mpi_max, &
1539  io_configuration%io_communicator, field_to_write_information%max_num_collective_writes_request_handle, ierr)
1540  call unlock_mpi()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initialise_writer_federator()

subroutine, public writer_federator_mod::initialise_writer_federator ( type(io_configuration_type), intent(inout)  io_configuration,
type(hashmap_type), intent(inout)  diagnostic_generation_frequency,
logical, intent(in)  continuation_run 
)

Initialises the write federator and configures it based on the user configuration. Also initialises the time manipulations.

Parameters
io_configurationThe IO server configuration

Definition at line 57 of file writer_federator.F90.

58  type(io_configuration_type), intent(inout) :: io_configuration
59  type(hashmap_type), intent(inout) :: diagnostic_generation_frequency
60  logical, intent(in) :: continuation_run
61 
62  integer :: i, j, number_contents, current_field_index
63  type(hashset_type) :: writer_field_names, duplicate_field_names
64 
65  call check_thread_status(forthread_rwlock_init(time_points_rwlock, -1))
66  call check_thread_status(forthread_mutex_init(collective_contiguous_initialisation_mutex, -1))
67  call check_thread_status(forthread_mutex_init(currently_writing_mutex, -1))
68 
69  currently_writing=.false.
70 
71  call init_time_averaged_manipulation()
72  call init_instantaneous_manipulation()
73  call initialise_netcdf_filetype()
74 
75  allocate(writer_entries(io_configuration%number_of_writers))
76  do i=1, io_configuration%number_of_writers
77  current_field_index=0
78  number_contents=io_configuration%file_writers(i)%number_of_contents
79  allocate(writer_entries(i)%contents(get_total_number_writer_fields(io_configuration, i)))
80  writer_entries(i)%filename=io_configuration%file_writers(i)%file_name
81  writer_entries(i)%title=io_configuration%file_writers(i)%title
82  writer_entries(i)%write_on_terminate=io_configuration%file_writers(i)%write_on_terminate
83  writer_entries(i)%include_in_io_state_write=io_configuration%file_writers(i)%include_in_io_state_write
84  call check_thread_status(forthread_mutex_init(writer_entries(i)%trigger_and_write_mutex, -1))
85  call check_thread_status(forthread_mutex_init(writer_entries(i)%num_fields_to_write_mutex, -1))
86  call check_thread_status(forthread_mutex_init(writer_entries(i)%pending_writes_mutex, -1))
87  writer_entries(i)%write_on_model_time=io_configuration%file_writers(i)%write_on_model_time
88  if (writer_entries(i)%write_on_model_time) then
89  writer_entries(i)%write_timestep_frequency=0
90  writer_entries(i)%write_time_frequency=io_configuration%file_writers(i)%write_time_frequency
91  else
92  writer_entries(i)%write_time_frequency=0
93  writer_entries(i)%write_timestep_frequency=io_configuration%file_writers(i)%write_timestep_frequency
94  end if
95  writer_entries(i)%previous_write_time=0
96  writer_entries(i)%defined_write_time=io_configuration%file_writers(i)%write_time_frequency
97  writer_entries(i)%latest_pending_write_time=0
98  writer_entries(i)%latest_pending_write_timestep=0
99  writer_entries(i)%contains_io_status_dump=.false.
100  do j=1, number_contents
101  if (io_configuration%file_writers(i)%contents(j)%facet_type == group_type) then
102  current_field_index=add_group_of_fields_to_writer_entry(io_configuration, i, j, current_field_index, &
103  writer_field_names, duplicate_field_names, diagnostic_generation_frequency)
104  else if (io_configuration%file_writers(i)%contents(j)%facet_type == field_type) then
105  current_field_index=current_field_index+add_field_to_writer_entry(io_configuration, &
106  i, j, current_field_index, io_configuration%file_writers(i)%contents(j)%facet_name, "", writer_field_names, &
107  duplicate_field_names, diagnostic_generation_frequency)
108  else if (io_configuration%file_writers(i)%contents(j)%facet_type == io_state_type) then
109  writer_entries(i)%contains_io_status_dump=.true.
110  end if
111  end do
112  if (.not. c_is_empty(duplicate_field_names)) call handle_duplicate_field_names(writer_entries(i), duplicate_field_names)
113  call c_free(writer_field_names)
114  call c_free(duplicate_field_names)
115  end do
116  if (continuation_run) then
117  call reactivate_writer_federator_state(io_configuration, writer_entries, time_points)
118  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ is_field_split_on_q()

logical function, public writer_federator_mod::is_field_split_on_q ( character(len=*), intent(in)  field_name)

Determines whether a field is split on Q or not.

Parameters
field_nameThe field name to check whether it is being used or not
Returns
Whether this field is used or not (and then further split to be constituient parts of Q)

Definition at line 239 of file writer_federator.F90.

240  character(len=*), intent(in) :: field_name
241 
242  is_field_split_on_q=c_contains(q_field_names, field_name)
Here is the caller graph for this function:

◆ is_field_used_by_writer_federator()

logical function, public writer_federator_mod::is_field_used_by_writer_federator ( character(len=*), intent(in)  field_name,
character(len=*), intent(in)  field_namespace 
)

Determines whether a field is used by the writer federator or not.

Parameters
field_nameThe field name to check whether it is being used or not
Returns
Whether this field is used or not

Definition at line 222 of file writer_federator.F90.

223  character(len=*), intent(in) :: field_name, field_namespace
224 
225  integer :: writer_index, contents_index
226 
227  writer_index=1
228  contents_index=1
229  if (c_contains(used_field_names, field_name)) then
230  is_field_used_by_writer_federator=get_next_applicable_writer_entry(field_name, field_namespace, writer_index, contents_index)
231  else
232  is_field_used_by_writer_federator=.false.
233  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ issue_actual_write()

subroutine, public writer_federator_mod::issue_actual_write ( type(io_configuration_type), intent(inout)  io_configuration,
type(writer_type), intent(inout)  writer_entry,
integer, intent(in)  timestep,
real, intent(in)  time,
logical, intent(in)  terminated_write 
)

Issues the actual file creation, write of available fields and closure if all completed.

Parameters
io_configurationThe IO server configuration
writer_entry_indexIndex of the writer we are concerned with
timestepThe corresponding timestep
timeThe corresponding model time

Definition at line 728 of file writer_federator.F90.

729  type(io_configuration_type), intent(inout) :: io_configuration
730  type(writer_type), intent(inout) :: writer_entry
731  integer, intent(in) :: timestep
732  real, intent(in) :: time
733  logical, intent(in) :: terminated_write
734 
735  integer :: i, j, total_outstanding, num_written, total_flds
736  logical :: field_written
737  type(map_type) :: applicable_time_points
738 
739  call check_thread_status(forthread_mutex_lock(collective_contiguous_initialisation_mutex))
740  do i=1, size(writer_entry%contents)
741  if (writer_entry%contents(i)%enabled .and. writer_entry%contents(i)%collective_write) then
742  if (.not. writer_entry%contents(i)%collective_initialised) then
743  call determine_collective_type_and_optimise_if_possible(io_configuration, writer_entry%contents(i))
744  end if
745  end if
746  end do
747  call check_thread_status(forthread_mutex_unlock(collective_contiguous_initialisation_mutex))
748 
749  writer_entry%write_time=time
750  writer_entry%write_timestep=timestep
751  applicable_time_points=extract_applicable_time_points(writer_entry%previous_write_time, time)
752  call define_netcdf_file(io_configuration, writer_entry, timestep, time, applicable_time_points, terminated_write)
753  call c_free(applicable_time_points)
754  total_outstanding=0
755  total_flds=0
756  num_written=0
757  call check_thread_status(forthread_mutex_lock(writer_entry%num_fields_to_write_mutex))
758  do j=1, size(writer_entry%contents)
759  if (writer_entry%contents(j)%enabled .and. writer_entry%contents(j)%expected_here) then
760  total_flds=total_flds+1
761  call determine_if_field_can_be_written(io_configuration, writer_entry, writer_entry%contents(j), timestep, &
762  writer_entry%previous_write_timestep, time, writer_entry%contents(j)%previous_write_time, field_written)
763  if (.not. field_written) then
764  total_outstanding=total_outstanding+1
765  else
766  num_written=num_written+1
767  end if
768  end if
769  end do
770  writer_entry%num_fields_to_write=total_outstanding
771  call check_thread_status(forthread_mutex_unlock(writer_entry%num_fields_to_write_mutex))
772  if (log_get_logging_level() .ge. log_debug) then
773  call log_log(log_debug, "Started write for NetCDF file, timestep= "//trim(conv_to_string(timestep))&
774  //" total="//trim(conv_to_string(total_flds))//" written="//trim(conv_to_string(num_written))//&
775  " outstanding="//trim(conv_to_string(total_outstanding)))
776  end if
777  if (total_outstanding == 0) then
778  call close_diagnostics_file(io_configuration, writer_entry, timestep, time)
779  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ provide_ordered_field_to_writer_federator()

subroutine, public writer_federator_mod::provide_ordered_field_to_writer_federator ( type(io_configuration_type), intent(inout)  io_configuration,
character(len=*), intent(in)  field_name,
character(len=*), intent(in)  field_namespace,
type(data_values_type), target  field_values,
integer, intent(in)  timestep,
real(kind=default_precision), intent(in)  time,
integer, intent(in)  source 
)

Definition at line 311 of file writer_federator.F90.

313  type(io_configuration_type), intent(inout) :: io_configuration
314  character(len=*), intent(in) :: field_name, field_namespace
315  integer, intent(in) :: timestep, source
316  type(data_values_type), target :: field_values
317  real(kind=default_precision), intent(in) :: time
318 
319  integer :: writer_index, contents_index
320  logical :: continue_search
321  type(data_values_type), pointer :: result_values
322  type(hashmap_type) :: typed_result_values
323  class(*), pointer :: generic
324 
325  if (field_values%data_type == double_data_type) then
326  call provide_ordered_field_to_writer_federator_real_values(io_configuration, field_name, field_namespace, &
327  field_values%values, timestep, time, source)
328  else if (field_values%data_type == string_data_type) then
329  continue_search=.true.
330  writer_index=1
331  contents_index=0
332  allocate(result_values, source=field_values)
333  generic=>result_values
334  if (c_contains(used_field_names, field_name)) then
335  do while (continue_search)
336  contents_index=contents_index+1
337  continue_search=get_next_applicable_writer_entry(field_name, field_namespace, writer_index, contents_index)
338  if (continue_search) then
339  if (.not. writer_entries(writer_index)%contents(contents_index)%enabled) then
340  call log_log(log_warn, "Received data for previously un-enabled field '"//&
341  writer_entries(writer_index)%contents(contents_index)%field_name//"'")
342  end if
343  writer_entries(writer_index)%contents(contents_index)%enabled=.true.
344  writer_entries(writer_index)%contents(contents_index)%latest_timestep_values=timestep
345  if (log_get_logging_level() .ge. log_debug) then
346  call log_log(log_debug, "[WRITE FED VALUE STORE] Storing value for field "//trim(field_name)//" ts="//&
347  trim(conv_to_string(timestep))// " t="//trim(conv_to_string(time)))
348  end if
349  call check_thread_status(forthread_mutex_lock(writer_entries(writer_index)%contents(contents_index)%values_mutex))
350  call c_put_generic(writer_entries(writer_index)%contents(contents_index)%values_to_write, conv_to_string(time), &
351  generic, .false.)
352  call check_thread_status(forthread_mutex_unlock(writer_entries(writer_index)%contents(contents_index)%values_mutex))
353  if (writer_entries(writer_index)%contents(contents_index)%pending_to_write) then
354  call determine_if_outstanding_field_can_be_written(io_configuration, writer_entries(writer_index), &
355  writer_entries(writer_index)%contents(contents_index))
356  end if
357  end if
358  end do
359  end if
360  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ provide_ordered_field_to_writer_federator_real_values()

subroutine writer_federator_mod::provide_ordered_field_to_writer_federator_real_values ( type(io_configuration_type), intent(inout)  io_configuration,
character(len=*), intent(in)  field_name,
character(len=*), intent(in)  field_namespace,
real(kind=default_precision), dimension(:), intent(in)  field_values,
integer, intent(in)  timestep,
real(kind=default_precision), intent(in)  time,
integer, intent(in)  source 
)
private

Provides fields (either diagnostics or prognostics) to the write federator which will action these as appropriate. This will split Q fields up if appropriate.

Parameters
io_configurationThe IO server configuration
field_nameThe field name to write (if appropriate)
field_valuesThe field values to write (if appropriate)
timestepCorresponding MONC timestep
timeCorresponding MONC model time
sourceOptional MONC source for the communicated fields

Definition at line 371 of file writer_federator.F90.

373  type(io_configuration_type), intent(inout) :: io_configuration
374  character(len=*), intent(in) :: field_name, field_namespace
375  integer, intent(in) :: timestep, source
376  real(kind=default_precision), dimension(:), intent(in) :: field_values
377  real(kind=default_precision), intent(in) :: time
378 
379  type(iterator_type) :: iterator
380  integer :: individual_size, index
381 
382  if (c_contains(used_field_names, field_name)) then
383  call provide_ordered_single_field_to_writer_federator(io_configuration, field_name, field_namespace, field_values, &
384  timestep, time, source)
385  else if (c_contains(q_field_names, field_name)) then
386  if (c_contains(q_field_splits, field_name)) then
387  individual_size=c_get_integer(q_field_splits, field_name)
388  else if (source .gt. -1) then
389  individual_size=get_size_of_collective_q(io_configuration, field_name, source)
390  else
391  call log_log(log_warn, "Can not find Q split field in Q field names or collective field names with source, ignoring")
392  return
393  end if
394  iterator=c_get_iterator(io_configuration%q_field_names)
395  index=1
396  do while (c_has_next(iterator))
397  call provide_ordered_single_field_to_writer_federator(io_configuration, &
398  trim(field_name)//"_"//trim(c_next_string(iterator)), field_namespace, field_values(index:index+individual_size-1), &
399  timestep, time, source)
400  index=index+individual_size
401  end do
402  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ provide_ordered_single_field_to_writer_federator()

subroutine writer_federator_mod::provide_ordered_single_field_to_writer_federator ( type(io_configuration_type), intent(inout)  io_configuration,
character(len=*), intent(in)  field_name,
character(len=*), intent(in)  field_namespace,
real(kind=default_precision), dimension(:), intent(in)  field_values,
integer, intent(in)  timestep,
real(kind=default_precision), intent(in)  time,
integer, intent(in)  source 
)
private

Provides a single ordered field, i.e. Q fields have been split by this point.

Parameters
io_configurationThe IO server configuration
field_nameThe field name to write (if appropriate)
field_valuesThe field values to write (if appropriate)
timestepCorresponding MONC timestep
timeCorresponding MONC model time
sourceOptional MONC source for the communicated fields

Definition at line 437 of file writer_federator.F90.

439  type(io_configuration_type), intent(inout) :: io_configuration
440  character(len=*), intent(in) :: field_name, field_namespace
441  integer, intent(in) :: timestep, source
442  real(kind=default_precision), dimension(:), intent(in) :: field_values
443  real(kind=default_precision), intent(in) :: time
444 
445  integer :: writer_index, contents_index
446  logical :: continue_search
447  type(data_values_type), pointer :: result_values
448  type(hashmap_type) :: typed_result_values
449  class(*), pointer :: generic
450 
451  continue_search=.true.
452  writer_index=1
453  contents_index=0
454  if (c_contains(used_field_names, field_name)) then
455  do while (continue_search)
456  contents_index=contents_index+1
457  continue_search=get_next_applicable_writer_entry(field_name, field_namespace, writer_index, contents_index)
458  if (continue_search) then
459  if (.not. writer_entries(writer_index)%contents(contents_index)%enabled) then
460  call log_log(log_warn, "Received data for previously un-enabled field '"//&
461  writer_entries(writer_index)%contents(contents_index)%field_name//"'")
462  end if
463  writer_entries(writer_index)%contents(contents_index)%enabled=.true.
464  if (.not. c_contains(typed_result_values, conv_to_string(&
465  writer_entries(writer_index)%contents(contents_index)%time_manipulation_type))) then
466  allocate(result_values)
467  if (writer_entries(writer_index)%contents(contents_index)%collective_write .and. source .gt. -1) then
468  result_values=writer_entries(writer_index)%contents(contents_index)%time_manipulation(field_values, &
469  writer_entries(writer_index)%contents(contents_index)%output_frequency, &
470  trim(field_name)//"#"//conv_to_string(source), timestep, time)
471  else
472  result_values=writer_entries(writer_index)%contents(contents_index)%time_manipulation(field_values, &
473  writer_entries(writer_index)%contents(contents_index)%output_frequency, &
474  field_name, timestep, time)
475  end if
476  generic=>result_values
477  call c_put_generic(typed_result_values, conv_to_string(&
478  writer_entries(writer_index)%contents(contents_index)%time_manipulation_type), generic, .false.)
479  else
480  result_values=>get_data_value_by_field_name(typed_result_values, conv_to_string(&
481  writer_entries(writer_index)%contents(contents_index)%time_manipulation_type))
482  end if
483  if (allocated(result_values%values)) then
484  writer_entries(writer_index)%contents(contents_index)%latest_timestep_values=timestep
485  if (log_get_logging_level() .ge. log_debug) then
486  call log_log(log_debug, "[WRITE FED VALUE STORE] Storing value for field "//trim(field_name)//" ts="//&
487  trim(conv_to_string(timestep))// " t="//trim(conv_to_string(time)))
488  end if
489  call check_thread_status(forthread_mutex_lock(writer_entries(writer_index)%contents(contents_index)%values_mutex))
490  if (writer_entries(writer_index)%contents(contents_index)%collective_write .and. source .gt. -1) then
491  call write_collective_write_value(result_values, writer_index, contents_index, source, conv_to_string(time))
492  else
493  call c_put_generic(writer_entries(writer_index)%contents(contents_index)%values_to_write, conv_to_string(time), &
494  generic, .false.)
495  end if
496  call check_thread_status(forthread_mutex_unlock(writer_entries(writer_index)%contents(contents_index)%values_mutex))
497  if (writer_entries(writer_index)%contents(contents_index)%pending_to_write) then
498  call determine_if_outstanding_field_can_be_written(io_configuration, writer_entries(writer_index), &
499  writer_entries(writer_index)%contents(contents_index))
500  end if
501  end if
502  end if
503  end do
504  end if
505  call c_free(typed_result_values)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ provide_q_field_names_to_writer_federator()

subroutine, public writer_federator_mod::provide_q_field_names_to_writer_federator ( type(list_type), intent(inout)  q_provided_field_names)

Provides the Q field names to the write federator, this is required as on initialisation we don't know what these are and only when MONC register do they inform the IO server of the specifics.

Parameters
q_field_namesAn ordered list of Q field names

Definition at line 277 of file writer_federator.F90.

278  type(list_type), intent(inout) :: q_provided_field_names
279 
280  type(iterator_type) :: iterator, q_field_iterator
281  logical :: continue_search
282  integer :: writer_index, contents_index, i
283  character(len=STRING_LENGTH) :: search_field, field_name, specific_name
284 
285  iterator=c_get_iterator(q_field_names)
286  do while (c_has_next(iterator))
287  specific_name=c_next_string(iterator)
288  q_field_iterator=c_get_iterator(q_provided_field_names)
289  i=1
290  do while (c_has_next(q_field_iterator))
291  search_field=trim(specific_name)//"_udef"//trim(conv_to_string(i))
292  field_name=trim(specific_name)//"_"//trim(c_next_string(q_field_iterator))
293  continue_search=.true.
294  writer_index=1
295  contents_index=0
296  do while (continue_search)
297  contents_index=contents_index+1
298  continue_search=get_next_applicable_writer_entry(search_field, writer_index_point=writer_index, &
299  contents_index_point=contents_index)
300  if (continue_search) then
301  writer_entries(writer_index)%contents(contents_index)%field_name=field_name
302  end if
303  end do
304  i=i+1
305  call c_add_string(used_field_names, field_name)
306  call c_remove(used_field_names, search_field)
307  end do
308  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ register_pending_file_write()

subroutine writer_federator_mod::register_pending_file_write ( integer, intent(in)  writer_entry_index,
integer, intent(in)  timestep,
real, intent(in)  time,
logical, intent(in)  terminated_write 
)
private

Registers a pending file write which will be actioned later on.

Parameters
writer_entry_indexIndex of the writer entry
timestepThe timestep that the pending write represents
timeThe time of the pending write

Definition at line 992 of file writer_federator.F90.

993  integer, intent(in) :: writer_entry_index, timestep
994  real, intent(in) :: time
995  logical, intent(in) :: terminated_write
996 
997  type(pending_write_type), pointer :: pending_write
998  class(*), pointer :: generic
999 
1000  allocate(pending_write)
1001  pending_write%write_time=time
1002  pending_write%timestep=timestep
1003  pending_write%terminated_write=terminated_write
1004 
1005  generic=>pending_write
1006  call check_thread_status(forthread_mutex_lock(writer_entries(writer_entry_index)%pending_writes_mutex))
1007  call c_push_generic(writer_entries(writer_entry_index)%pending_writes, generic, .false.)
1008  call check_thread_status(forthread_mutex_unlock(writer_entries(writer_entry_index)%pending_writes_mutex))
Here is the caller graph for this function:

◆ sort_applicable_time_points()

type(map_type) function writer_federator_mod::sort_applicable_time_points ( type(map_type), intent(inout)  unsorted_timepoints)
private

Sorts the time points based upon their timestep, smallest to largest. Note that this is a bubble sort and as such inefficient, so would be good to change to something else but works OK for now.

Parameters
unsorted_timepointsThe unsorted timepoints, which is changed by destroying the map
Returns
The sorted version of the input map based upon timestep

Definition at line 846 of file writer_federator.F90.

847  type(map_type), intent(inout) :: unsorted_timepoints
848 
849  integer :: i, entries, specific_ts, smallest_ts
850  character(len=STRING_LENGTH) :: smallest_key
851  real(kind=default_precision) :: rvalue
852  type(iterator_type) :: iterator
853  type(mapentry_type) :: map_entry
854 
855  entries=c_size(unsorted_timepoints)
856  do i=1, entries
857  smallest_key=""
858  iterator=c_get_iterator(unsorted_timepoints)
859  do while (c_has_next(iterator))
860  map_entry=c_next_mapentry(iterator)
861  specific_ts=conv_to_integer(map_entry%key)
862  if (len_trim(smallest_key) == 0 .or. smallest_ts .gt. specific_ts) then
863  smallest_ts=specific_ts
864  smallest_key=map_entry%key
865  rvalue=c_get_real(map_entry)
866  end if
867  end do
868  call c_put_real(sort_applicable_time_points, smallest_key, rvalue)
869  call c_remove(unsorted_timepoints, smallest_key)
870  end do
871  call c_free(unsorted_timepoints)
Here is the caller graph for this function:

◆ write_collective_write_value()

subroutine writer_federator_mod::write_collective_write_value ( type(data_values_type), pointer  result_values,
integer, intent(in)  writer_index,
integer, intent(in)  contents_index,
integer, intent(in)  source,
character(len=*), intent(in)  lookup_key 
)
private

Writes the collective values, this is held differently to independent values which are written directly - instead here we need to store the values for each MONC hence a specific type is used instead.

Parameters
result_valuesThe data values to store
writer_indexThe file writer index
contents_indexThe contents index
sourceThe MONC process id
lookup_keyThe values lookup key

Definition at line 515 of file writer_federator.F90.

516  integer, intent(in) :: writer_index, contents_index, source
517  type(data_values_type), pointer :: result_values
518  character(len=*), intent(in) :: lookup_key
519 
520  class(*), pointer :: generic
521  type(write_field_collective_values_type), pointer :: stored_monc_values
522 
523  if (c_contains(writer_entries(writer_index)%contents(contents_index)%values_to_write, lookup_key)) then
524  generic=>c_get_generic(writer_entries(writer_index)%contents(contents_index)%values_to_write, lookup_key)
525  select type(generic)
526  type is(write_field_collective_values_type)
527  stored_monc_values=>generic
528  end select
529  else
530  allocate(stored_monc_values)
531  generic=>stored_monc_values
532  call c_put_generic(writer_entries(writer_index)%contents(contents_index)%values_to_write, lookup_key, generic, .false.)
533  end if
534  generic=>result_values
535  call c_put_generic(stored_monc_values%monc_values, conv_to_string(source), generic, .false.)
Here is the caller graph for this function:

Variable Documentation

◆ collective_contiguous_initialisation_mutex

integer, volatile writer_federator_mod::collective_contiguous_initialisation_mutex
private

Definition at line 47 of file writer_federator.F90.

◆ collective_q_field_dims

type(hashmap_type), volatile writer_federator_mod::collective_q_field_dims
private

Definition at line 45 of file writer_federator.F90.

◆ currently_writing

logical, volatile writer_federator_mod::currently_writing
private

Definition at line 48 of file writer_federator.F90.

48  logical, volatile :: currently_writing

◆ currently_writing_mutex

integer, volatile writer_federator_mod::currently_writing_mutex
private

Definition at line 47 of file writer_federator.F90.

◆ q_field_names

type(hashset_type), volatile writer_federator_mod::q_field_names
private

Definition at line 44 of file writer_federator.F90.

◆ q_field_splits

type(hashmap_type), volatile writer_federator_mod::q_field_splits
private

Definition at line 45 of file writer_federator.F90.

◆ time_points

type(hashmap_type), volatile writer_federator_mod::time_points
private

Definition at line 45 of file writer_federator.F90.

45  type(hashmap_type), volatile :: time_points, q_field_splits, collective_q_field_dims

◆ time_points_rwlock

integer, volatile writer_federator_mod::time_points_rwlock
private

Definition at line 47 of file writer_federator.F90.

47  integer, volatile :: time_points_rwlock, collective_contiguous_initialisation_mutex, currently_writing_mutex

◆ used_field_names

type(hashset_type), volatile writer_federator_mod::used_field_names
private

Definition at line 44 of file writer_federator.F90.

44  type(hashset_type), volatile :: used_field_names, q_field_names

◆ writer_entries

type(writer_type), dimension(:), allocatable, volatile writer_federator_mod::writer_entries
private

Definition at line 43 of file writer_federator.F90.

43  type(writer_type), volatile, dimension(:), allocatable :: writer_entries
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
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
logging_mod::log_get_logging_level
integer function, public log_get_logging_level()
Retrieves the current logging level.
Definition: logging.F90:122
logging_mod::log_debug
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
Definition: logging.F90:14
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
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