MONC
Functions/Subroutines
localreduce_operator_mod Module Reference

Performs a local reduction, reducing a local array into a single scalar value. More...

Functions/Subroutines

subroutine, public perform_localreduce_operator (io_configuration, field_values, action_attributes, source_monc_location, source_monc, operator_result_values)
 Executes this local reduction operator. More...
 
real(kind=default_precision) function do_local_reduction (data, reduction_operator)
 Does the actual local reduction, translating the array into a vector based upon the operator. More...
 
type(list_type) function, public localreduce_operator_get_required_fields (action_attributes)
 Retrieves the list of fields needed by this operator for a specific configuration. More...
 

Detailed Description

Performs a local reduction, reducing a local array into a single scalar value.

Function/Subroutine Documentation

◆ do_local_reduction()

real(kind=default_precision) function localreduce_operator_mod::do_local_reduction ( real(kind=default_precision), dimension(:), intent(in)  data,
character(len=string_length), intent(in)  reduction_operator 
)
private

Does the actual local reduction, translating the array into a vector based upon the operator.

Parameters
dataThe array data to reduce
reduction_operatorThe operator to apply for this reduction
Returns
The resulting scalar value

Definition at line 46 of file localreduce-operator.F90.

47  real(kind=default_precision), dimension(:), intent(in) :: data
48  character(len=STRING_LENGTH), intent(in) :: reduction_operator
49 
50  if (reduction_operator .eq. "max") then
51  do_local_reduction=maxval(data)
52  else if (reduction_operator .eq. "min") then
53  do_local_reduction=minval(data)
54  else if (reduction_operator .eq. "sum") then
55  do_local_reduction=sum(data)
56  else if (reduction_operator .eq. "mean") then
57  do_local_reduction=sum(data)/size(data)
58  end if
Here is the caller graph for this function:

◆ localreduce_operator_get_required_fields()

type(list_type) function, public localreduce_operator_mod::localreduce_operator_get_required_fields ( type(map_type), intent(inout)  action_attributes)

Retrieves the list of fields needed by this operator for a specific configuration.

Parameters
action_attributesThe attributes which configure the operator
Returns
A list of required fields before the operator can run

Definition at line 64 of file localreduce-operator.F90.

65  type(map_type), intent(inout) :: action_attributes
66 
67  character(len=STRING_LENGTH) :: field_to_reduce
68 
69  field_to_reduce=get_action_attribute_string(action_attributes, "field")
70  call c_add_string(localreduce_operator_get_required_fields, field_to_reduce)
Here is the call graph for this function:

◆ perform_localreduce_operator()

subroutine, public localreduce_operator_mod::perform_localreduce_operator ( type(io_configuration_type), intent(inout)  io_configuration,
type(hashmap_type), intent(inout)  field_values,
type(map_type), intent(inout)  action_attributes,
integer, intent(in)  source_monc_location,
integer, intent(in)  source_monc,
real(kind=default_precision), dimension(:), intent(inout), allocatable  operator_result_values 
)

Executes this local reduction operator.

Parameters
io_configurationConfiguration of the IO server
field_valuesThe field values
action_attributesAttributes associated with the running of this operator
source_monc_locationLocation of the source MONC
source_moncThe source MONC
operator_result_valuesThe operator resulting (scalar) value

Definition at line 23 of file localreduce-operator.F90.

25  type(io_configuration_type), intent(inout) :: io_configuration
26  type(hashmap_type), intent(inout) :: field_values
27  type(map_type), intent(inout) :: action_attributes
28  integer, intent(in) :: source_monc_location, source_monc
29  real(kind=default_precision), dimension(:), allocatable, intent(inout) :: operator_result_values
30 
31  character(len=STRING_LENGTH) :: field_to_reduce, reduction_operator
32  type(data_values_type), pointer :: field_local_values
33 
34  allocate(operator_result_values(1))
35 
36  field_to_reduce=get_action_attribute_string(action_attributes, "field")
37  reduction_operator=get_action_attribute_string(action_attributes, "operator")
38  field_local_values=>get_data_value_by_field_name(field_values, field_to_reduce)
39  operator_result_values=do_local_reduction(field_local_values%values, reduction_operator)
Here is the call graph for this function:
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