MOM6
MOM_transform_FMS.F90
1 !> Support functions and interfaces to permit transformed model domains to
2 !! interact with FMS operations registered on the non-transformed domains.
3 
5 
6 use horiz_interp_mod, only : horiz_interp_type
7 use mom_error_handler, only : mom_error, fatal
8 use mom_io, only : fieldtype, write_field
9 use mpp_domains_mod, only : domain2d
10 use fms_mod, only : mpp_chksum
11 use time_manager_mod, only : time_type
12 use time_interp_external_mod, only : time_interp_external
13 
15 
16 implicit none
17 
18 private
19 public rotated_mpp_chksum
22 
23 !> Rotate and compute the FMS (mpp) checksum of a field
25  module procedure rotated_mpp_chksum_real_0d
26  module procedure rotated_mpp_chksum_real_1d
27  module procedure rotated_mpp_chksum_real_2d
28  module procedure rotated_mpp_chksum_real_3d
29  module procedure rotated_mpp_chksum_real_4d
30 end interface rotated_mpp_chksum
31 
32 !> Rotate and write a registered field to an FMS output file
34  module procedure rotated_write_field_real_0d
35  module procedure rotated_write_field_real_1d
36  module procedure rotated_write_field_real_2d
37  module procedure rotated_write_field_real_3d
38  module procedure rotated_write_field_real_4d
39 end interface rotated_write_field
40 
41 !> Read a field based on model time, and rotate to the model domain
43  module procedure rotated_time_interp_external_0d
44  module procedure rotated_time_interp_external_2d
45  module procedure rotated_time_interp_external_3d
46 end interface rotated_time_interp_external
47 
48 contains
49 
50 ! NOTE: No transformations are applied to the 0d and 1d field implementations,
51 ! but are provided to maintain compatibility with the FMS interfaces.
52 
53 
54 !> Compute the FMS (mpp) checksum of a scalar.
55 !! This function is provided to support the full FMS mpp_chksum interface.
56 function rotated_mpp_chksum_real_0d(field, pelist, mask_val, turns) &
57  result(chksum)
58  real, intent(in) :: field !> Input scalar
59  integer, optional, intent(in) :: pelist(:) !> PE list of ranks to checksum
60  real, optional, intent(in) :: mask_val !> FMS mask value
61  integer, optional, intent(in) :: turns !> Number of quarter turns
62  integer :: chksum !> FMS checksum of scalar
63 
64  if (present(turns)) &
65  call mom_error(fatal, "Rotation not supported for 0d fields.")
66 
67  chksum = mpp_chksum(field, pelist=pelist, mask_val=mask_val)
68 end function rotated_mpp_chksum_real_0d
69 
70 
71 !> Compute the FMS (mpp) checksum of a 1d field.
72 !! This function is provided to support the full FMS mpp_chksum interface.
73 function rotated_mpp_chksum_real_1d(field, pelist, mask_val, turns) &
74  result(chksum)
75  real, intent(in) :: field(:) !> Input field
76  integer, optional, intent(in) :: pelist(:) !> PE list of ranks to checksum
77  real, optional, intent(in) :: mask_val !> FMS mask value
78  integer, optional, intent(in) :: turns !> Number of quarter-turns
79  integer :: chksum !> FMS checksum of field
80 
81  if (present(turns)) &
82  call mom_error(fatal, "Rotation not supported for 1d fields.")
83 
84  chksum = mpp_chksum(field, pelist=pelist, mask_val=mask_val)
85 end function rotated_mpp_chksum_real_1d
86 
87 
88 !> Compute the FMS (mpp) checksum of a rotated 2d field.
89 function rotated_mpp_chksum_real_2d(field, pelist, mask_val, turns) &
90  result(chksum)
91  real, intent(in) :: field(:,:) !> Unrotated input field
92  integer, optional, intent(in) :: pelist(:) !> PE list of ranks to checksum
93  real, optional, intent(in) :: mask_val !> FMS mask value
94  integer, optional, intent(in) :: turns !> Number of quarter-turns
95  integer :: chksum !> FMS checksum of field
96 
97  real, allocatable :: field_rot(:,:)
98  integer :: qturns
99 
100  qturns = 0
101  if (present(turns)) &
102  qturns = modulo(turns, 4)
103 
104  if (qturns == 0) then
105  chksum = mpp_chksum(field, pelist=pelist, mask_val=mask_val)
106  else
107  call allocate_rotated_array(field, [1,1], qturns, field_rot)
108  call rotate_array(field, qturns, field_rot)
109  chksum = mpp_chksum(field_rot, pelist=pelist, mask_val=mask_val)
110  deallocate(field_rot)
111  endif
112 end function rotated_mpp_chksum_real_2d
113 
114 
115 !> Compute the FMS (mpp) checksum of a rotated 3d field.
116 function rotated_mpp_chksum_real_3d(field, pelist, mask_val, turns) &
117  result(chksum)
118  real, intent(in) :: field(:,:,:) !> Unrotated input field
119  integer, optional, intent(in) :: pelist(:) !> PE list of ranks to checksum
120  real, optional, intent(in) :: mask_val !> FMS mask value
121  integer, optional, intent(in) :: turns !> Number of quarter-turns
122  integer :: chksum !> FMS checksum of field
123 
124  real, allocatable :: field_rot(:,:,:)
125  integer :: qturns
126 
127  qturns = 0
128  if (present(turns)) &
129  qturns = modulo(turns, 4)
130 
131  if (qturns == 0) then
132  chksum = mpp_chksum(field, pelist=pelist, mask_val=mask_val)
133  else
134  call allocate_rotated_array(field, [1,1,1], qturns, field_rot)
135  call rotate_array(field, qturns, field_rot)
136  chksum = mpp_chksum(field_rot, pelist=pelist, mask_val=mask_val)
137  deallocate(field_rot)
138  endif
139 end function rotated_mpp_chksum_real_3d
140 
141 
142 !> Compute the FMS (mpp) checksum of a rotated 4d field.
143 function rotated_mpp_chksum_real_4d(field, pelist, mask_val, turns) &
144  result(chksum)
145  real, intent(in) :: field(:,:,:,:) !> Unrotated input field
146  integer, optional, intent(in) :: pelist(:) !> PE list of ranks to checksum
147  real, optional, intent(in) :: mask_val !> FMS mask value
148  integer, optional, intent(in) :: turns !> Number of quarter-turns
149  integer :: chksum !> FMS checksum of field
150 
151  real, allocatable :: field_rot(:,:,:,:)
152  integer :: qturns
153 
154  qturns = 0
155  if (present(turns)) &
156  qturns = modulo(turns, 4)
157 
158  if (qturns == 0) then
159  chksum = mpp_chksum(field, pelist=pelist, mask_val=mask_val)
160  else
161  call allocate_rotated_array(field, [1,1,1,1], qturns, field_rot)
162  call rotate_array(field, qturns, field_rot)
163  chksum = mpp_chksum(field_rot, pelist=pelist, mask_val=mask_val)
164  deallocate(field_rot)
165  endif
166 end function rotated_mpp_chksum_real_4d
167 
168 
169 ! NOTE: In MOM_io, write_field points to mpp_write, which supports a very broad
170 ! range of interfaces. Here, we only support the much more narrow family of
171 ! mpp_write_2ddecomp functions used to write tiled data.
172 
173 
174 !> Write the rotation of a 1d field to an FMS output file
175 !! This function is provided to support the full FMS write_field interface.
176 subroutine rotated_write_field_real_0d(io_unit, field_md, field, tstamp, turns)
177  integer, intent(in) :: io_unit !> File I/O unit handle
178  type(fieldtype), intent(in) :: field_md !> FMS field metadata
179  real, intent(inout) :: field !> Unrotated field array
180  real, optional, intent(in) :: tstamp !> Model timestamp
181  integer, optional, intent(in) :: turns !> Number of quarter-turns
182 
183  if (present(turns)) &
184  call mom_error(fatal, "Rotation not supported for 0d fields.")
185 
186  call write_field(io_unit, field_md, field, tstamp=tstamp)
187 end subroutine rotated_write_field_real_0d
188 
189 
190 !> Write the rotation of a 1d field to an FMS output file
191 !! This function is provided to support the full FMS write_field interface.
192 subroutine rotated_write_field_real_1d(io_unit, field_md, field, tstamp, turns)
193  integer, intent(in) :: io_unit !> File I/O unit handle
194  type(fieldtype), intent(in) :: field_md !> FMS field metadata
195  real, intent(inout) :: field(:) !> Unrotated field array
196  real, optional, intent(in) :: tstamp !> Model timestamp
197  integer, optional, intent(in) :: turns !> Number of quarter-turns
198 
199  if (present(turns)) &
200  call mom_error(fatal, "Rotation not supported for 0d fields.")
201 
202  call write_field(io_unit, field_md, field, tstamp=tstamp)
203 end subroutine rotated_write_field_real_1d
204 
205 
206 !> Write the rotation of a 2d field to an FMS output file
207 subroutine rotated_write_field_real_2d(io_unit, field_md, domain, field, &
208  tstamp, tile_count, default_data, turns)
209  integer, intent(in) :: io_unit !> File I/O unit handle
210  type(fieldtype), intent(in) :: field_md !> FMS field metadata
211  type(domain2d), intent(inout) :: domain !> FMS MPP domain
212  real, intent(inout) :: field(:,:) !> Unrotated field array
213  real, optional, intent(in) :: tstamp !> Model timestamp
214  integer, optional, intent(in) :: tile_count !> PEs per tile (default: 1)
215  real, optional, intent(in) :: default_data !> Default fill value
216  integer, optional, intent(in) :: turns !> Number of quarter-turns
217 
218  real, allocatable :: field_rot(:,:)
219  integer :: qturns
220 
221  qturns = 0
222  if (present(turns)) &
223  qturns = modulo(turns, 4)
224 
225  if (qturns == 0) then
226  call write_field(io_unit, field_md, domain, field, tstamp=tstamp, &
227  tile_count=tile_count, default_data=default_data)
228  else
229  call allocate_rotated_array(field, [1,1], qturns, field_rot)
230  call rotate_array(field, qturns, field_rot)
231  call write_field(io_unit, field_md, domain, field_rot, tstamp=tstamp, &
232  tile_count=tile_count, default_data=default_data)
233  deallocate(field_rot)
234  endif
235 end subroutine rotated_write_field_real_2d
236 
237 
238 !> Write the rotation of a 3d field to an FMS output file
239 subroutine rotated_write_field_real_3d(io_unit, field_md, domain, field, &
240  tstamp, tile_count, default_data, turns)
241  integer, intent(in) :: io_unit !> File I/O unit handle
242  type(fieldtype), intent(in) :: field_md !> FMS field metadata
243  type(domain2d), intent(inout) :: domain !> FMS MPP domain
244  real, intent(inout) :: field(:,:,:) !> Unrotated field array
245  real, optional, intent(in) :: tstamp !> Model timestamp
246  integer, optional, intent(in) :: tile_count !> PEs per tile (default: 1)
247  real, optional, intent(in) :: default_data !> Default fill value
248  integer, optional, intent(in) :: turns !> Number of quarter-turns
249 
250  real, allocatable :: field_rot(:,:,:)
251  integer :: qturns
252 
253  qturns = 0
254  if (present(turns)) &
255  qturns = modulo(turns, 4)
256 
257  if (qturns == 0) then
258  call write_field(io_unit, field_md, domain, field, tstamp=tstamp, &
259  tile_count=tile_count, default_data=default_data)
260  else
261  call allocate_rotated_array(field, [1,1,1], qturns, field_rot)
262  call rotate_array(field, qturns, field_rot)
263  call write_field(io_unit, field_md, domain, field_rot, tstamp=tstamp, &
264  tile_count=tile_count, default_data=default_data)
265  deallocate(field_rot)
266  endif
267 end subroutine rotated_write_field_real_3d
268 
269 
270 !> Write the rotation of a 4d field to an FMS output file
271 subroutine rotated_write_field_real_4d(io_unit, field_md, domain, field, &
272  tstamp, tile_count, default_data, turns)
273  integer, intent(in) :: io_unit !> File I/O unit handle
274  type(fieldtype), intent(in) :: field_md !> FMS field metadata
275  type(domain2d), intent(inout) :: domain !> FMS MPP domain
276  real, intent(inout) :: field(:,:,:,:) !> Unrotated field array
277  real, optional, intent(in) :: tstamp !> Model timestamp
278  integer, optional, intent(in) :: tile_count !> PEs per tile (default: 1)
279  real, optional, intent(in) :: default_data !> Default fill value
280  integer, optional, intent(in) :: turns !> Number of quarter-turns
281 
282  real, allocatable :: field_rot(:,:,:,:)
283  integer :: qturns
284 
285  qturns = 0
286  if (present(turns)) &
287  qturns = modulo(turns, 4)
288 
289  if (qturns == 0) then
290  call write_field(io_unit, field_md, domain, field, tstamp=tstamp, &
291  tile_count=tile_count, default_data=default_data)
292  else
293  call allocate_rotated_array(field, [1,1,1,1], qturns, field_rot)
294  call rotate_array(field, qturns, field_rot)
295  call write_field(io_unit, field_md, domain, field_rot, tstamp=tstamp, &
296  tile_count=tile_count, default_data=default_data)
297  deallocate(field_rot)
298  endif
299 end subroutine rotated_write_field_real_4d
300 
301 
302 !> Read a scalar field based on model time
303 !! This function is provided to support the full FMS time_interp_external
304 !! interface.
305 subroutine rotated_time_interp_external_0d(fms_id, time, data_in, verbose, &
306  turns)
307  integer, intent(in) :: fms_id !< FMS field ID
308  type(time_type), intent(in) :: time !< Model time
309  real, intent(inout) :: data_in !< field to write data
310  logical, intent(in), optional :: verbose !< Verbose output
311  integer, intent(in), optional :: turns !< Number of quarter turns
312 
313  if (present(turns)) &
314  call mom_error(fatal, "Rotation not supported for 0d fields.")
315 
316  call time_interp_external(fms_id, time, data_in, verbose=verbose)
317 end subroutine rotated_time_interp_external_0d
318 
319 !> Read a 2d field based on model time, and rotate to the model grid
320 subroutine rotated_time_interp_external_2d(fms_id, time, data_in, interp, &
321  verbose, horz_interp, mask_out, is_in, ie_in, js_in, je_in, window_id, &
322  turns)
323  integer, intent(in) :: fms_id
324  type(time_type), intent(in) :: time
325  real, dimension(:,:), intent(inout) :: data_in
326  integer, intent(in), optional :: interp
327  logical, intent(in), optional :: verbose
328  type(horiz_interp_type),intent(in), optional :: horz_interp
329  logical, dimension(:,:), intent(out), optional :: mask_out
330  integer, intent(in), optional :: is_in, ie_in, js_in, je_in
331  integer, intent(in), optional :: window_id
332  integer, intent(in), optional :: turns
333 
334  real, allocatable :: data_pre(:,:)
335  integer :: qturns
336 
337  ! TODO: Mask rotation requires logical array rotation support
338  if (present(mask_out)) &
339  call mom_error(fatal, "Rotation of masked output not yet support")
340 
341  qturns = 0
342  if (present(turns)) &
343  qturns = modulo(turns, 4)
344 
345 
346  if (qturns == 0) then
347  call time_interp_external(fms_id, time, data_in, interp=interp, &
348  verbose=verbose, horz_interp=horz_interp, mask_out=mask_out, &
349  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in, &
350  window_id=window_id)
351  else
352  call allocate_rotated_array(data_in, [1,1], -qturns, data_pre)
353  call time_interp_external(fms_id, time, data_pre, interp=interp, &
354  verbose=verbose, horz_interp=horz_interp, mask_out=mask_out, &
355  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in, &
356  window_id=window_id)
357  call rotate_array(data_pre, turns, data_in)
358  deallocate(data_pre)
359  endif
360 end subroutine rotated_time_interp_external_2d
361 
362 
363 !> Read a 3d field based on model time, and rotate to the model grid
364 subroutine rotated_time_interp_external_3d(fms_id, time, data_in, interp, &
365  verbose, horz_interp, mask_out, is_in, ie_in, js_in, je_in, window_id, &
366  turns)
367  integer, intent(in) :: fms_id
368  type(time_type), intent(in) :: time
369  real, dimension(:,:,:), intent(inout) :: data_in
370  integer, intent(in), optional :: interp
371  logical, intent(in), optional :: verbose
372  type(horiz_interp_type),intent(in), optional :: horz_interp
373  logical, dimension(:,:,:), intent(out), optional :: mask_out
374  integer, intent(in), optional :: is_in, ie_in, js_in, je_in
375  integer, intent(in), optional :: window_id
376  integer, intent(in), optional :: turns
377 
378  real, allocatable :: data_pre(:,:,:)
379  integer :: qturns
380 
381  ! TODO: Mask rotation requires logical array rotation support
382  if (present(mask_out)) &
383  call mom_error(fatal, "Rotation of masked output not yet support")
384 
385  qturns = 0
386  if (present(turns)) &
387  qturns = modulo(turns, 4)
388 
389  if (qturns == 0) then
390  call time_interp_external(fms_id, time, data_in, interp=interp, &
391  verbose=verbose, horz_interp=horz_interp, mask_out=mask_out, &
392  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in, &
393  window_id=window_id)
394  else
395  call allocate_rotated_array(data_in, [1,1,1], -qturns, data_pre)
396  call time_interp_external(fms_id, time, data_pre, interp=interp, &
397  verbose=verbose, horz_interp=horz_interp, mask_out=mask_out, &
398  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in, &
399  window_id=window_id)
400  call rotate_array(data_pre, turns, data_in)
401  deallocate(data_pre)
402  endif
403 end subroutine rotated_time_interp_external_3d
404 
405 end module mom_transform_fms
mom_array_transform
Module for supporting the rotation of a field's index map. The implementation of each angle is descri...
Definition: MOM_array_transform.F90:14
mom_array_transform::rotate_array
Rotate the elements of an array to the rotated set of indices. Rotation is applied across the first a...
Definition: MOM_array_transform.F90:27
mom_transform_fms
Support functions and interfaces to permit transformed model domains to interact with FMS operations ...
Definition: MOM_transform_FMS.F90:4
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_transform_fms::rotated_write_field
Rotate and write a registered field to an FMS output file.
Definition: MOM_transform_FMS.F90:33
mom_array_transform::allocate_rotated_array
Allocate an array based on the rotated index map of an unrotated reference array.
Definition: MOM_array_transform.F90:64
mom_transform_fms::rotated_time_interp_external
Read a field based on model time, and rotate to the model domain.
Definition: MOM_transform_FMS.F90:42
mom_transform_fms::rotated_mpp_chksum
Rotate and compute the FMS (mpp) checksum of a field.
Definition: MOM_transform_FMS.F90:24
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2