MONC
Functions/Subroutines
fieldslicer_operator_mod Module Reference

Slices a field based upon the selected dimension and index. More...

Functions/Subroutines

subroutine, public perform_fieldslicer_operator (io_configuration, field_values, action_attributes, source_monc_location, source_monc, operator_result_values)
 Performs the actual field slicing. More...
 
subroutine determine_dimension_bounds (corresponding_field_definition, registered_monc_info, dimension_to_slice, index_to_slice, dim_starts, dim_ends, dim_weights, number_dims, sliced_size)
 Determines the bounds for each dimension which is specified via the configuration file. More...
 
type(list_type) function, public fieldslicer_operator_get_required_fields (action_attributes)
 Retrieves a list of the required fields for running this operator. More...
 
integer function get_dimension_to_slice (action_attributes)
 Retrieves the integer index of the dimension to slice. More...
 
integer function convert_dimension_str_to_id (dim_str)
 Converts a dimension string to the corresponding numeric ID. More...
 

Detailed Description

Slices a field based upon the selected dimension and index.

Function/Subroutine Documentation

◆ convert_dimension_str_to_id()

integer function fieldslicer_operator_mod::convert_dimension_str_to_id ( character(len=string_length), intent(in)  dim_str)
private

Converts a dimension string to the corresponding numeric ID.

Parameters
dim_strThe dimension string
Returns
The corresponding numeric dimension ID

Definition at line 150 of file fieldslicer-operator.F90.

150  character(len=STRING_LENGTH), intent(in) :: dim_str
151 
152  if (dim_str .eq. "x") then
153  convert_dimension_str_to_id=x_index
154  else if (dim_str .eq. "y") then
155  convert_dimension_str_to_id=y_index
156  else if (dim_str .eq. "z") then
157  convert_dimension_str_to_id=z_index
158  end if
Here is the caller graph for this function:

◆ determine_dimension_bounds()

subroutine fieldslicer_operator_mod::determine_dimension_bounds ( type(io_configuration_field_type), intent(in)  corresponding_field_definition,
type(io_configuration_registered_monc_type), intent(in)  registered_monc_info,
integer, intent(in)  dimension_to_slice,
integer, intent(in)  index_to_slice,
integer, dimension(:), intent(out), allocatable  dim_starts,
integer, dimension(:), intent(out), allocatable  dim_ends,
integer, dimension(:), intent(out), allocatable  dim_weights,
integer, intent(out)  number_dims,
integer, intent(out)  sliced_size 
)
private

Determines the bounds for each dimension which is specified via the configuration file.

Parameters
corresponding_field_definitionThe definition of the field that is being sliced
registered_monc_infoThis specific MONC process
dimension_to_sliceNumeric code for the dimension that should be sliced
index_to_sliceThe index that we are slicing upon
dim_startsOutput start index for each dimension
dim_endsOutput end index for each dimension
dim_weightsOutput weight for each dimension when dealing with a 1D data
number_dimsOutput number of dimensions that the field comprises of
sliced_sizeOutput overall 1D data size of the sliced array

Definition at line 85 of file fieldslicer-operator.F90.

85  type(io_configuration_field_type), intent(in) :: corresponding_field_definition
86  type(io_configuration_registered_monc_type), intent(in) :: registered_monc_info
87  integer, intent(in) :: dimension_to_slice, index_to_slice
88  integer, dimension(:), allocatable, intent(out) :: dim_starts, dim_ends, dim_weights
89  integer, intent(out) :: number_dims, sliced_size
90 
91  integer :: i, j, dimension_id, amount_to_add
92  integer, dimension(:), allocatable :: dim_sizes
93  logical :: found_slice_field
94 
95  number_dims=corresponding_field_definition%dimensions
96  allocate(dim_sizes(number_dims), dim_starts(number_dims), dim_ends(number_dims), dim_weights(number_dims))
97  found_slice_field=.false.
98 
99  do i=1, number_dims
100  dimension_id=convert_dimension_str_to_id(corresponding_field_definition%dim_size_defns(i))
101  if (dimension_id==dimension_to_slice) then
102  dim_sizes(i)=1
103  dim_starts(i)=index_to_slice
104  dim_ends(i)=index_to_slice
105  found_slice_field=.true.
106  else
107  dim_sizes(i)=registered_monc_info%local_dim_sizes(dimension_id)
108  dim_starts(i)=1
109  dim_ends(i)=registered_monc_info%local_dim_sizes(dimension_id)
110  end if
111  end do
112  if (.not. found_slice_field) call log_log(log_error, "Can not find dimension to slice in provided field")
113  sliced_size=0
114  do i=1, number_dims
115  amount_to_add=1
116  dim_weights(i)=1
117  do j=1, i
118  if (j .lt. i) dim_weights(i)=dim_weights(i)*registered_monc_info%local_dim_sizes(j)
119  amount_to_add=amount_to_add*dim_sizes(j)
120  end do
121  sliced_size=sliced_size+amount_to_add
122  end do
123  deallocate(dim_sizes)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fieldslicer_operator_get_required_fields()

type(list_type) function, public fieldslicer_operator_mod::fieldslicer_operator_get_required_fields ( type(map_type), intent(inout)  action_attributes)

Retrieves a list of the required fields for running this operator.

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

Definition at line 130 of file fieldslicer-operator.F90.

130  type(map_type), intent(inout) :: action_attributes
131 
132  character(len=STRING_LENGTH) :: field_to_slice
133 
134  field_to_slice=get_action_attribute_string(action_attributes, "field")
135  call c_add_string(fieldslicer_operator_get_required_fields, field_to_slice)
Here is the call graph for this function:

◆ get_dimension_to_slice()

integer function fieldslicer_operator_mod::get_dimension_to_slice ( type(map_type), intent(inout)  action_attributes)
private

Retrieves the integer index of the dimension to slice.

Parameters
action_attributesThe attributes which configure the operator

Definition at line 141 of file fieldslicer-operator.F90.

141  type(map_type), intent(inout) :: action_attributes
142 
143  get_dimension_to_slice=convert_dimension_str_to_id(get_action_attribute_string(action_attributes, "dimension"))
Here is the call graph for this function:
Here is the caller graph for this function:

◆ perform_fieldslicer_operator()

subroutine, public fieldslicer_operator_mod::perform_fieldslicer_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 
)

Performs the actual field slicing.

Parameters
io_configurationConfiguration of the IO server
field_valuesThe field values
action_attributesAttributes associated with the running of this operator
source_monc_locationThe location in configuration data of MONC settings
source_moncProcess ID of the MONC that sent this data
operator_result_valuesThe resulting value or left unallocated if none are appropriate

Definition at line 27 of file fieldslicer-operator.F90.

27  type(io_configuration_type), intent(inout) :: io_configuration
28  type(hashmap_type), intent(inout) :: field_values
29  type(map_type), intent(inout) :: action_attributes
30  integer, intent(in) :: source_monc_location, source_monc
31  real(kind=DEFAULT_PRECISION), dimension(:), allocatable, intent(inout) :: operator_result_values
32 
33  character(len=STRING_LENGTH) :: field_to_slice
34  integer :: i, j, dimension_to_slice, index_to_slice, number_dims, sliced_size, source_dim
35  integer, dimension(:), allocatable :: dim_starts, dim_ends, dim_weights, indexes
36  type(data_values_type), pointer :: field_local_values
37  type(io_configuration_field_type) :: corresponding_field_definition
38 
39  field_to_slice=get_action_attribute_string(action_attributes, "field")
40  ! NSE
41  if (get_prognostic_field_configuration(io_configuration, field_to_slice, "", corresponding_field_definition)) then
42  dimension_to_slice=get_dimension_to_slice(action_attributes)
43  index_to_slice=get_action_attribute_integer(action_attributes, "index")
44 
45  if (io_configuration%registered_moncs(source_monc_location)%local_dim_starts(dimension_to_slice) .le. index_to_slice .and.&
46  io_configuration%registered_moncs(source_monc_location)%local_dim_ends(dimension_to_slice) .ge. index_to_slice) then
47  field_local_values=>get_data_value_by_field_name(field_values, field_to_slice)
48 
49  call determine_dimension_bounds(corresponding_field_definition, io_configuration%registered_moncs(source_monc_location), &
50  dimension_to_slice, index_to_slice, dim_starts, dim_ends, dim_weights, number_dims, sliced_size)
51 
52  allocate(operator_result_values(sliced_size), indexes(number_dims))
53  indexes=dim_starts
54  do i=1, sliced_size
55  source_dim=1
56  do j=1, number_dims
57  source_dim=source_dim+(indexes(j)-1)*dim_weights(j)
58  end do
59  operator_result_values(i)=field_local_values%values(source_dim)
60  indexes(1)=indexes(1)+1
61  do j=1, number_dims
62  if (indexes(j) .gt. dim_ends(j)) then
63  indexes(j)=dim_starts(j)
64  if (j .lt. number_dims) indexes(j+1)=indexes(j+1)+1
65  end if
66  end do
67  end do
68  deallocate(dim_starts, dim_ends, dim_weights, indexes)
69  end if
70  end if
Here is the call graph for this function: