MONC
instantaneous_manipulation.F90
Go to the documentation of this file.
1 
13  implicit none
14 
15 #ifndef TEST_MODE
16  private
17 #endif
18 
21 
25 contains
26 
30  end subroutine init_instantaneous_manipulation
31 
36 
37  logical function is_instantaneous_time_manipulation_ready_to_write(latest_time, output_frequency, write_time, &
38  latest_timestep, write_timestep)
39  real, intent(in) :: latest_time, output_frequency, write_time
40  integer, intent(in) :: latest_timestep, write_timestep
41 
42  is_instantaneous_time_manipulation_ready_to_write=latest_time + output_frequency .gt. write_time
44 
54  type(data_values_type) function perform_instantaneous_time_manipulation(instant_values, output_frequency, &
55  field_name, timestep, time)
56  real(kind=DEFAULT_PRECISION), dimension(:), intent(in) :: instant_values
57  real, intent(in) :: output_frequency
58  real(kind=DEFAULT_PRECISION), intent(in) :: time
59  character(len=*), intent(in) :: field_name
60  integer, intent(in) :: timestep
61 
62  if (deduce_whether_to_issue_values(field_name, output_frequency, time)) then
63  allocate(perform_instantaneous_time_manipulation%values(size(instant_values)))
64  perform_instantaneous_time_manipulation%values=instant_values
65  end if
67 
73  logical function deduce_whether_to_issue_values(field_name, output_frequency, time)
74  character(len=*), intent(in) :: field_name
75  real, intent(in) :: output_frequency
76  real(kind=DEFAULT_PRECISION), intent(in) :: time
77 
78  real(kind=DEFAULT_PRECISION) :: previous_time_write, time_difference
79 
80  call check_thread_status(forthread_mutex_lock(existing_instantaneous_writes_mutex))
81  if (c_contains(existing_instantaneous_writes, field_name)) then
82  previous_time_write=c_get_real(existing_instantaneous_writes, field_name)
83  time_difference=(aint(time*10000000.0)-aint(previous_time_write*10000000.0))/10000000.0
84  else
85  ! Rethink this as only works if started at time=0
86  time_difference=time
87  end if
88  if (time_difference .ge. output_frequency) then
90  call c_put_real(existing_instantaneous_writes, field_name, time)
91  else
93  end if
94  call check_thread_status(forthread_mutex_unlock(existing_instantaneous_writes_mutex))
96 
99  integer(kind=8) function prepare_to_serialise_instantaneous_state()
100  real(kind=DEFAULT_PRECISION) :: a
101 
102  call check_thread_status(forthread_mutex_lock(existing_instantaneous_writes_mutex))
104  ((string_length + kind(a)) * c_size(existing_instantaneous_writes))
106 
110  subroutine serialise_instantaneous_state(byte_data)
111  character, dimension(:), allocatable, intent(inout) :: byte_data
112 
113  integer :: current_data_point
114  type(mapentry_type) :: map_entry
115  type(iterator_type) :: iterator
116 
117  current_data_point=1
118  current_data_point=pack_scalar_field(byte_data, current_data_point, c_size(existing_instantaneous_writes))
120  do while (c_has_next(iterator))
121  map_entry=c_next_mapentry(iterator)
122  current_data_point=pack_scalar_field(byte_data, current_data_point, string_value=map_entry%key)
123  current_data_point=pack_scalar_field(byte_data, current_data_point, double_real_value=c_get_real(map_entry))
124  end do
125  call check_thread_status(forthread_mutex_unlock(existing_instantaneous_writes_mutex))
126  end subroutine serialise_instantaneous_state
127 
130  subroutine unserialise_instantaneous_state(byte_data)
131  character, dimension(:), allocatable, intent(in) :: byte_data
132 
133  integer :: current_data_point, number_entries, i
134 
135  current_data_point=1
136  number_entries=unpack_scalar_integer_from_bytedata(byte_data, current_data_point)
137  if (number_entries .gt. 0) then
138  do i=1, number_entries
139  call c_put_real(existing_instantaneous_writes, unpack_scalar_string_from_bytedata(byte_data, current_data_point), &
140  unpack_scalar_dp_real_from_bytedata(byte_data, current_data_point))
141  end do
142  end if
143  end subroutine unserialise_instantaneous_state
subroutine, public finalise_instantaneous_manipulation()
Finalises the instantaneous time manipulation.
integer function forthread_mutex_unlock(mutex_id)
Definition: forthread.F90:302
type(data_values_type) function, public perform_instantaneous_time_manipulation(instant_values, output_frequency, field_name, timestep, time)
Performs the instantaneous time manipulation and returns data only if this is to be written to the st...
logical function deduce_whether_to_issue_values(field_name, output_frequency, time)
Determines whether to issue values for write or not. This depends on the time and output frequency...
subroutine, public unserialise_instantaneous_state(byte_data)
Unpacks some serialised byte data to initialise this manipulator to some previous state...
integer function forthread_mutex_destroy(mutex_id)
Definition: forthread.F90:265
Contains functionality for managing and extracting data from the raw data dumps that the IO server re...
Definition: datautils.F90:3
Performs instantaneous time manipulation and only returns a value if the output frequency determines ...
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
A hashmap structure, the same as a map but uses hashing for greatly improved performance when storing...
Definition: collections.F90:94
integer function forthread_mutex_init(mutex_id, attr_id)
Definition: forthread.F90:274
Conversion between common inbuilt FORTRAN data types.
Definition: conversions.F90:5
This defines some constants and procedures that are useful to the IO server and clients that call it...
Definition: ioclient.F90:3
type(hashmap_type), volatile existing_instantaneous_writes
This is a thread pool and the single management "main" thread will spawn out free threads in the pool...
Definition: threadpool.F90:5
integer(kind=8) function, public prepare_to_serialise_instantaneous_state()
Prepares to serialise the instantaneous state, both determines the byte storage size and issues any l...
Returns the number of elements in the collection.
real(kind=double_precision) function, public unpack_scalar_dp_real_from_bytedata(data, start_point)
Unpacks a double precision scalar real from some byte data, this is a very simple unpack routine wrap...
Definition: datautils.F90:89
subroutine, public check_thread_status(ierr)
Checks the error status of any thread operation and reports an error if it failed.
Definition: threadpool.F90:229
subroutine, public init_instantaneous_manipulation()
Initialises the instantaneous time manipulation.
integer function forthread_mutex_lock(mutex_id)
Definition: forthread.F90:284
Collection data structures.
Definition: collections.F90:7
character(len=string_length) function, public unpack_scalar_string_from_bytedata(data, start_point)
Unpacks a string from some byte data with default length, this is a very simple unpack routine wrappi...
Definition: datautils.F90:62
integer, parameter, public string_length
Default length of strings.
Definition: datadefn.F90:10
subroutine, public serialise_instantaneous_state(byte_data)
Will serialise the state of this manipulator so that it can be later restarted. Any locks issued duri...
integer function, public unpack_scalar_integer_from_bytedata(data, start_point)
Unpacks a scalar integer from some byte data, this is a very simple unpack routine wrapping the trans...
Definition: datautils.F90:34
logical function, public is_instantaneous_time_manipulation_ready_to_write(latest_time, output_frequency, write_time, latest_timestep, write_timestep)
Determines whether or not a map contains a specific key.
Gets a specific double precision real element out of the list, stack, queue or map with the correspon...
Parses the XML configuration file to produce the io configuration description which contains the data...
real(kind=double_precision) function, public conv_single_real_to_double(input_real)
Converts from a single to double precision real. This applies some rounding to a certain number of de...
Puts a double precision real key-value pair into the map.
integer function, public pack_scalar_field(buffer, start_offset, int_value, real_value, single_real_value, double_real_value, string_value, logical_value)
Packs the data of a scalar field into a buffer.
Definition: ioclient.F90:312