MOM_ice_shelf_diag_mediator.F90

1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> The subroutines here provide convenient wrappers to the FMS diag_manager
6!! interfaces with additional diagnostic capabilities.
8
9! This file is part of MOM6. See LICENSE.md for the license.
10
11use mom_checksums, only : chksum0, hchksum, uchksum, vchksum, bchksum
12use mom_coms, only : pe_here
13use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
14use mom_cpu_clock, only : clock_module, clock_routine
15use mom_diag_manager_infra, only : mom_diag_manager_init
16use mom_diag_manager_infra, only : mom_diag_axis_init, get_mom_diag_axis_name
17use mom_diag_manager_infra, only : send_data_infra, mom_diag_field_add_attribute, east, north
18use mom_diag_manager_infra, only : register_diag_field_infra, register_static_field_infra
19use mom_diag_manager_infra, only : get_mom_diag_field_id, diag_field_not_found
20use mom_diag_manager_infra, only : diag_send_complete_infra
21use mom_error_handler, only : mom_error, fatal, is_root_pe, assert, calltree_showquery
23use mom_file_parser, only : get_param, log_version, param_file_type
24use mom_grid, only : ocean_grid_type
26use mom_safe_alloc, only : safe_alloc_ptr, safe_alloc_alloc
28use mom_time_manager, only : time_type, get_time
30
31implicit none ; private
32
35public set_is_axes_info, mom_diag_axis_init
36public register_mom_is_diag_field, register_mom_is_static_field, register_mom_is_scalar_field
38public safe_alloc_ptr, safe_alloc_alloc, time_type
43
44!> Make a diagnostic available for averaging or output.
45interface post_is_data
46 module procedure post_is_data_2d, post_is_data_0d
47end interface post_is_data
48
49!> Registers a non-array scalar diagnostic, returning an integer handle
50interface register_mom_is_scalar_field
52end interface register_mom_is_scalar_field
53
54!> A group of 1D axes that comprise a 1D/2D/3D mesh
55type, public :: axes_grp
56 character(len=15) :: id !< The id string for this particular combination of handles.
57 integer :: rank !< Number of dimensions in the list of axes.
58 integer, dimension(:), allocatable :: handles !< Handles to 1D axes.
59 type(diag_ctrl), pointer :: diag_cs => null() !< Circular link back to the main diagnostics control structure
60 !! (Used to avoid passing said structure into every possible call).
61 ! ID's for cell_methods
62 character(len=9) :: x_cell_method = '' !< Default nature of data representation, if axes group
63 !! includes x-direction.
64 character(len=9) :: y_cell_method = '' !< Default nature of data representation, if axes group
65 !! includes y-direction.
66 ! For detecting position on the grid
67 logical :: is_h_point = .false. !< If true, indicates that this axes group is for an h-point located field.
68 logical :: is_q_point = .false. !< If true, indicates that this axes group is for a q-point located field.
69 logical :: is_u_point = .false. !< If true, indicates that this axes group is for a u-point located field.
70 logical :: is_v_point = .false. !< If true, indicates that this axes group is for a v-point located field.
71
72 ! ID's for cell_measures
73 integer :: id_area = -1 !< The diag_manager id for area to be used for cell_measure of variables with this axes_grp.
74 ! For masking
75 real, pointer, dimension(:,:) :: mask2d => null() !< Mask for 2d (x-y) axes [nondim]
76 real, pointer, dimension(:,:) :: mask2d_comp => null() !< Mask for 2-d axes on the computational
77 !! domain for this diagnostic [nondim]
78end type axes_grp
79
80!> This type is used to represent a diagnostic at the diag_mediator level.
81!!
82!! There can be both 'primary' and 'secondary' diagnostics. The primaries
83!! reside in the diag_cs%diags array. They have an id which is an index
84!! into this array. The secondaries are 'variations' on the primary diagnostic.
85!! For example the CMOR diagnostics are secondary. The secondary diagnostics
86!! are kept in a list with the primary diagnostic as the head.
87type, private :: diag_type
88 logical :: in_use !< True if this entry is being used.
89 integer :: fms_diag_id !< Underlying FMS diag_manager id.
90 character(len=64) :: debug_str = '' !< The diagnostic name and module for FATAL errors and debugging.
91 type(axes_grp), pointer :: axes => null() !< The axis group for this diagnostic
92 type(diag_type), pointer :: next => null() !< Pointer to the next diagnostic
93 real :: conversion_factor = 0. !< If non-zero, a factor to multiply data by before posting to FMS,
94 !! often including factors to undo internal scaling in units of [a A-1 ~> 1]
95end type diag_type
96
97!> The diag_ctrl data type contains times to regulate diagnostics along with masks and
98!! axes to use with diagnostics, and a list of structures with data about each diagnostic.
99type, public :: diag_ctrl
100 integer :: available_diag_doc_unit = -1 !< The unit number of a diagnostic documentation file.
101 !! This file is open if available_diag_doc_unit is > 0.
102 integer :: chksum_iounit = -1 !< The unit number of a diagnostic documentation file.
103 !! This file is open if available_diag_doc_unit is > 0.
104 logical :: diag_as_chksum !< If true, log chksums in a text file instead of posting diagnostics
105 logical :: show_call_tree !< Display the call tree while running. Set by VERBOSITY level.
106 logical :: index_space_axes !< If true, diagnostic horizontal coordinates axes are in index space.
107
108 ! The following fields are used for the output of the data.
109 ! These give the computational-domain sizes, and are relative to a start value
110 ! of 1 in memory for the tracer-point arrays.
111 integer :: is !< The start i-index of cell centers within the computational domain
112 integer :: ie !< The end i-index of cell centers within the computational domain
113 integer :: js !< The start j-index of cell centers within the computational domain
114 integer :: je !< The end j-index of cell centers within the computational domain
115 ! These give the memory-domain sizes, and can be start at any value on each PE.
116 integer :: isd !< The start i-index of cell centers within the data domain
117 integer :: ied !< The end i-index of cell centers within the data domain
118 integer :: jsd !< The start j-index of cell centers within the data domain
119 integer :: jed !< The end j-index of cell centers within the data domain
120 real :: time_int !< The time interval for any fields
121 !! that are offered for averaging [s].
122 type(time_type) :: time_end !< The end time of the valid interval for any offered field.
123 logical :: ave_enabled = .false. !< True if averaging is enabled.
124
125 !>@{ The following are 3D and 2D axis groups defined for output. The names indicate
126 !! the horizontal locations (B, T, Cu, or Cv) and vertical locations (here just 1).
127 type(axes_grp) :: axesb1, axest1, axescu1, axescv1
128 !>@}
129 type(axes_grp) :: axesnull !< An axis group for scalars
130
131 ! Mask arrays for 2D diagnostics
132 real, dimension(:,:), pointer :: mask2dt => null() !< 2D mask array for cell-center points [nondim]
133 real, dimension(:,:), pointer :: mask2dbu => null() !< 2D mask array for cell-corner points [nondim]
134 real, dimension(:,:), pointer :: mask2dcu => null() !< 2D mask array for east-face points [nondim]
135 real, dimension(:,:), pointer :: mask2dcv => null() !< 2D mask array for north-face points [nondim]
136 real, dimension(:,:), pointer :: mask2dt_comp => null() !< 2D cell-center mask on the computational domain [nondim]
137
138! Space for diagnostics is dynamically allocated as it is needed.
139! The chunk size is how much the array should grow on each new allocation.
140#define DIAG_ALLOC_CHUNK_SIZE 15
141 type(diag_type), dimension(:), allocatable :: diags !< The list of diagnostics
142 integer :: next_free_diag_id !< The next unused diagnostic ID
143
144 !> default missing value to be sent to ALL diagnostics registrations [various]
145 real :: missing_value = -1.0e34
146
147 type(ocean_grid_type), pointer :: g => null() !< The ocean grid type
148 type(unit_scale_type), pointer :: us => null() !< A dimensional unit scaling type
149
150 !> Number of checksum-only diagnostics
151 integer :: num_chksum_diags
152
153end type diag_ctrl
154
155!>@{ CPU clocks
156integer :: id_clock_diag_mediator
157!>@}
158
159contains
160
161!> Set up the grid and axis information for use by the ice shelf model.
162subroutine set_is_axes_info(G, diag_cs, axes_set_name)
163 type(ocean_grid_type), intent(in) :: g !< The horizontal grid type
164 type(diag_ctrl), intent(inout) :: diag_cs !< A structure that is used to regulate diagnostic output
165 character(len=*), optional, intent(in) :: axes_set_name !< A name to use for this set of axes.
166 !! The default is "ice".
167! This subroutine sets up the grid and axis information for use by the ice shelf model.
168
169 ! Local variables
170 integer :: id_xq, id_yq, id_xh, id_yh, id_null
171 integer :: i, j
172 character(len=80) :: set_name
173 real, allocatable, dimension(:) :: iaxb, iax ! Index-based integer and half-integer i-axis labels [nondim]
174 real, allocatable, dimension(:) :: jaxb, jax ! Index-based integer and half-integer j-axis labels [nondim]
175
176 set_name = "ice_shelf" ; if (present(axes_set_name)) set_name = trim(axes_set_name)
177
178 if (diag_cs%index_space_axes) then
179 allocate(iaxb(g%IsgB:g%IegB))
180 do i=g%IsgB,g%IegB
181 iaxb(i) = real(i)
182 enddo
183 allocate(iax(g%isg:g%ieg))
184 do i=g%isg,g%ieg
185 iax(i) = real(i)-0.5
186 enddo
187 allocate(jaxb(g%JsgB:g%JegB))
188 do j=g%JsgB,g%JegB
189 jaxb(j) = real(j)
190 enddo
191 allocate(jax(g%jsg:g%jeg))
192 do j=g%jsg,g%jeg
193 jax(j) = real(j)-0.5
194 enddo
195 endif
196
197 ! Horizontal axes for the native grids.
198 if (diag_cs%index_space_axes) then
199 if (g%symmetric) then
200 id_xq = mom_diag_axis_init('Iq', iaxb(g%IsgB:g%IegB), 'none', 'x', &
201 'Boundary (q) point grid-space longitude', g%Domain, position=east, set_name=set_name)
202 id_yq = mom_diag_axis_init('Jq', jaxb(g%JsgB:g%JegB), 'none', 'y', &
203 'Boundary (q) point grid-space latitude', g%Domain, position=north, set_name=set_name)
204 else
205 id_xq = mom_diag_axis_init('Iq', iaxb(g%isg:g%ieg), 'none', 'x', &
206 'Boundary (q) point grid-space longitude', g%Domain, position=east, set_name=set_name)
207 id_yq = mom_diag_axis_init('Jq', jaxb(g%jsg:g%jeg), 'none', 'y', &
208 'Boundary (q) point grid-space latitude', g%Domain, position=north, set_name=set_name)
209 endif
210
211 id_xh = mom_diag_axis_init('ih', iax, 'none', 'x', &
212 'Tracer (h) point grid-space longitude', g%Domain, set_name=set_name)
213 id_yh = mom_diag_axis_init('jh', jax, 'none', 'y', &
214 'Tracer (h) point grid-space latitude', g%Domain, set_name=set_name)
215 else
216 if (g%symmetric) then
217 id_xq = mom_diag_axis_init('xB', g%gridLonB(g%isgB:g%iegB), g%x_axis_units, 'x', &
218 'Boundary point nominal longitude', g%Domain, position=east, set_name=set_name)
219 id_yq = mom_diag_axis_init('yB', g%gridLatB(g%jsgB:g%jegB), g%y_axis_units, 'y', &
220 'Boundary point nominal latitude', g%Domain, position=north, set_name=set_name)
221
222 else
223 id_xq = mom_diag_axis_init('xB', g%gridLonB(g%isg:g%ieg), g%x_axis_units, 'x', &
224 'Boundary point nominal longitude', g%Domain, position=east, set_name=set_name)
225 id_yq = mom_diag_axis_init('yB', g%gridLatB(g%jsg:g%jeg), g%y_axis_units, 'y', &
226 'Boundary point nominal latitude', g%Domain, position=north, set_name=set_name)
227
228 endif
229 id_xh = mom_diag_axis_init('xT', g%gridLonT(g%isg:g%ieg), g%x_axis_units, 'x', &
230 'Tracer point nominal longitude', g%Domain, set_name=set_name)
231 id_yh = mom_diag_axis_init('yT', g%gridLatT(g%jsg:g%jeg), g%y_axis_units, 'y', &
232 'Tracer point nominal latitude', g%Domain, set_name=set_name)
233 endif
234
235 ! Axis groupings for 2-D arrays
236 call define_axes_group(diag_cs, (/id_xh, id_yh/), diag_cs%axesT1, &
237 x_cell_method='mean', y_cell_method='mean', is_h_point=.true.)
238 call define_axes_group(diag_cs, (/id_xq, id_yq/), diag_cs%axesB1, &
239 x_cell_method='point', y_cell_method='point', is_q_point=.true.)
240 call define_axes_group(diag_cs, (/id_xq, id_yh/), diag_cs%axesCu1, &
241 x_cell_method='point', y_cell_method='mean', is_u_point=.true.)
242 call define_axes_group(diag_cs, (/id_xh, id_yq/), diag_cs%axesCv1, &
243 x_cell_method='mean', y_cell_method='point', is_v_point=.true.)
244
245 ! Axis group for special null axis for scalars from diag manager.
246 id_null = mom_diag_axis_init('scalar_axis', (/0./), 'none', 'N', 'none', null_axis=.true.)
247 call define_axes_group(diag_cs, (/ id_null /), diag_cs%axesNull)
248
249 if (diag_cs%index_space_axes) then
250 deallocate(iaxb, iax, jaxb, jax)
251 endif
252
253end subroutine set_is_axes_info
254
255!> Attaches the id of cell areas to axes groups for use with cell_measures
256subroutine diag_register_area_ids(diag_cs, id_area_t, id_area_q)
257 type(diag_ctrl), intent(inout) :: diag_cs !< Diagnostics control structure
258 integer, optional, intent(in) :: id_area_t !< Diag_mediator id for area of h-cells
259 integer, optional, intent(in) :: id_area_q !< Diag_mediator id for area of q-cells
260 ! Local variables
261 integer :: fms_id, i
262 if (present(id_area_t)) then
263 fms_id = diag_cs%diags(id_area_t)%fms_diag_id
264 diag_cs%axesT1%id_area = fms_id
265 endif
266 if (present(id_area_q)) then
267 fms_id = diag_cs%diags(id_area_q)%fms_diag_id
268 diag_cs%axesB1%id_area = fms_id
269 endif
270end subroutine diag_register_area_ids
271
272!> Define a group of "axes" from a list of handles and associate a mask with it
273subroutine define_axes_group(diag_cs, handles, axes, &
274 x_cell_method, y_cell_method, &
275 is_h_point, is_q_point, is_u_point, is_v_point)
276 type(diag_ctrl), target, intent(in) :: diag_cs !< Structure used to regulate diagnostic output
277 integer, dimension(:), intent(in) :: handles !< A list of 1D axis handles that define the axis group
278 type(axes_grp), intent(out) :: axes !< The group of axes that is set up here
279 character(len=*), optional, intent(in) :: x_cell_method !< A x-direction cell method used to construct the
280 !! "cell_methods" attribute in CF convention
281 character(len=*), optional, intent(in) :: y_cell_method !< A y-direction cell method used to construct the
282 !! "cell_methods" attribute in CF convention
283 logical, optional, intent(in) :: is_h_point !< If true, indicates this axes group for h-point
284 !! located fields
285 logical, optional, intent(in) :: is_q_point !< If true, indicates this axes group for q-point
286 !! located fields
287 logical, optional, intent(in) :: is_u_point !< If true, indicates this axes group for
288 !! u-point located fields
289 logical, optional, intent(in) :: is_v_point !< If true, indicates this axes group for
290 !! v-point located fields
291
292 ! Local variables
293 integer :: n
294
295 n = size(handles)
296 if (n<1 .or. n>2) call mom_error(fatal, "define_axes_group: wrong size for list of handles!")
297 allocate( axes%handles(n) )
298 axes%id = ints_to_string(handles, max(n,2)) ! Identifying string
299 axes%rank = n
300 axes%handles(:) = handles(:)
301 axes%diag_cs => diag_cs ! A (circular) link back to the diag_ctrl structure
302
303 if ((axes%rank<2) .and. (present(x_cell_method) .or. present(x_cell_method))) &
304 call mom_error(fatal, 'define_axes_group: Can not set x_cell_method or y_cell_method for rank<2.')
305 axes%x_cell_method = '' ; if (present(x_cell_method)) axes%x_cell_method = trim(x_cell_method)
306 axes%y_cell_method = '' ; if (present(y_cell_method)) axes%y_cell_method = trim(y_cell_method)
307
308 if (present(is_h_point)) axes%is_h_point = is_h_point
309 if (present(is_q_point)) axes%is_q_point = is_q_point
310 if (present(is_u_point)) axes%is_u_point = is_u_point
311 if (present(is_v_point)) axes%is_v_point = is_v_point
312
313 ! Setup masks for this axes group
314 axes%mask2d => null()
315 if (axes%rank==2) then
316 if (axes%is_h_point) axes%mask2d => diag_cs%mask2dT
317 if (axes%is_h_point) axes%mask2d_comp => diag_cs%mask2dT_comp
318 if (axes%is_u_point) axes%mask2d => diag_cs%mask2dCu
319 if (axes%is_v_point) axes%mask2d => diag_cs%mask2dCv
320 if (axes%is_q_point) axes%mask2d => diag_cs%mask2dBu
321 endif
322
323end subroutine define_axes_group
324
325!> Set up the array extents for doing diagnostics
326subroutine set_is_diag_mediator_grid(G, diag_cs)
327 type(ocean_grid_type), intent(inout) :: g !< The horizontal grid type
328 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
329
330 diag_cs%is = g%isc - (g%isd-1) ; diag_cs%ie = g%iec - (g%isd-1)
331 diag_cs%js = g%jsc - (g%jsd-1) ; diag_cs%je = g%jec - (g%jsd-1)
332 diag_cs%isd = g%isd ; diag_cs%ied = g%ied
333 diag_cs%jsd = g%jsd ; diag_cs%jed = g%jed
334
335end subroutine set_is_diag_mediator_grid
336
337!> Make a real ice shelf scalar diagnostic available for averaging or output
338subroutine post_is_data_0d(diag_field_id, field, diag_cs, is_static)
339 integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
340 !! previous call to register_MOM_IS_diag_field.
341 real, intent(in) :: field !< real value being offered for output or averaging
342 !! in internally scaled arbitrary units [A ~> a]
343 type(diag_ctrl), target, intent(in) :: diag_CS !< Structure used to regulate diagnostic output
344 logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
345
346 ! Local variables
347 real :: locfield ! The field being offered in arbitrary unscaled units [a]
348 logical :: used, is_stat
349 type(diag_type), pointer :: diag => null()
350
351 integer :: time_days
352 integer :: time_seconds
353 character(len=300) :: debug_mesg
354
355 if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
356 is_stat = .false. ; if (present(is_static)) is_stat = is_static
357
358 ! Iterate over list of diag 'variants', e.g. CMOR aliases, call send_data
359 ! for each one.
360 call assert(diag_field_id < diag_cs%next_free_diag_id, &
361 'post_IS_data_0d: Unregistered diagnostic id')
362 diag => diag_cs%diags(diag_field_id)
363
364 do while (associated(diag))
365 locfield = field
366 if (diag%conversion_factor /= 0.) &
367 locfield = locfield * diag%conversion_factor
368
369 if (diag_cs%diag_as_chksum) then
370 ! Append timestep to mesg
371 call get_time(diag_cs%time_end, time_seconds, days=time_days)
372 write(debug_mesg, '(a, 1x, i0, 1x, i0)') &
373 trim(diag%debug_str), time_days, time_seconds
374
375 call chksum0(locfield, debug_mesg, logunit=diag_cs%chksum_iounit)
376 elseif (is_stat) then
377 used = send_data_infra(diag%fms_diag_id, locfield)
378 elseif (diag_cs%ave_enabled) then
379 used = send_data_infra(diag%fms_diag_id, locfield, diag_cs%time_end)
380 endif
381
382 diag => diag%next
383 enddo
384
385 if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)
386end subroutine post_is_data_0d
387
388
389!> Make a real 2-d array diagnostic available for averaging or output
390subroutine post_is_data_2d(diag_field_id, field, diag_cs, is_static, mask)
391 integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
392 !! previous call to register_MOM_IS_diag_field.
393 real, target, intent(in) :: field(:,:) !< 2-d array being offered for output or averaging
394 !! in internally scaled arbitrary units [A ~> a]
395 type(diag_ctrl), target, intent(in) :: diag_CS !< Structure used to regulate diagnostic output
396 logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
397 real, optional, intent(in) :: mask(:,:) !< If present, use this real array as the data mask [nondim]
398
399 ! Local variables
400 type(diag_type), pointer :: diag => null()
401
402 if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
403
404 ! Iterate over list of diag 'variants' (e.g. CMOR aliases) and post each.
405 call assert(diag_field_id < diag_cs%next_free_diag_id, &
406 'post_IS_data_2d: Unregistered diagnostic id')
407 diag => diag_cs%diags(diag_field_id)
408 do while (associated(diag))
409 call post_data_2d_low(diag, field, diag_cs, is_static, mask)
410 diag => diag%next
411 enddo
412
413 if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)
414end subroutine post_is_data_2d
415
416!> Make a real 2-d array diagnostic available for averaging or output
417!! using a diag_type instead of an integer id.
418subroutine post_data_2d_low(diag, field, diag_cs, is_static, mask)
419 type(diag_type), intent(in) :: diag !< A structure describing the diagnostic to post
420 real, target, intent(in) :: field(:,:) !< 2-d array being offered for output or averaging
421 !! in internally scaled arbitrary units [A ~> a]
422 type(diag_ctrl), intent(in) :: diag_CS !< Structure used to regulate diagnostic output
423 logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
424 real, optional, target, intent(in) :: mask(:,:) !< If present, use this real array as the data mask [nondim]
425
426 ! Local variables
427 real, dimension(:,:), pointer :: locfield ! The field being offered in arbitrary unscaled units [a]
428 real, dimension(:,:), pointer :: locmask ! A pointer to the data mask to use [nondim]
429 logical :: used ! The return value of send_data is not used for anything.
430 logical :: is_stat
431 logical :: i_data, j_data ! True if the field is on the data domain in the i or j directions.
432 integer :: cszi, cszj, dszi, dszj
433 integer :: isv, iev, jsv, jev, i, j
434 integer :: time_days, time_seconds
435 character(len=300) :: mesg
436 character(len=300) :: debug_mesg
437
438 locfield => null()
439 locmask => null()
440 is_stat = .false. ; if (present(is_static)) is_stat = is_static
441
442 ! Determine the proper array indices, noting that because of the (:,:)
443 ! declaration of field, symmetric arrays are using a SW-grid indexing,
444 ! but non-symmetric arrays are using a NE-grid indexing. Send_data
445 ! actually only uses the difference between ie and is to determine
446 ! the output data size and assumes that halos are symmetric.
447 isv = diag_cs%is ; iev = diag_cs%ie ; jsv = diag_cs%js ; jev = diag_cs%je
448
449 cszi = (diag_cs%ie-diag_cs%is) +1 ; dszi = (diag_cs%ied-diag_cs%isd) +1
450 cszj = (diag_cs%je-diag_cs%js) +1 ; dszj = (diag_cs%jed-diag_cs%jsd) +1
451 if ( size(field,1) == dszi ) then
452 isv = diag_cs%is ; iev = diag_cs%ie ; i_data = .true. ! Data domain
453 elseif ( size(field,1) == dszi + 1 ) then
454 isv = diag_cs%is ; iev = diag_cs%ie+1 ; i_data = .true. ! Symmetric data domain
455 elseif ( size(field,1) == cszi ) then
456 isv = 1 ; iev = cszi ; i_data = .false. ! Computational domain
457 elseif ( size(field,1) == cszi + 1 ) then
458 isv = 1 ; iev = cszi+1 ; i_data = .false. ! Symmetric computational domain
459 else
460 write (mesg,*) " peculiar size ",size(field,1)," in i-direction\n"//&
461 "does not match one of ", cszi, cszi+1, dszi, dszi+1
462 call mom_error(fatal,"post_IS_data_2d_low: "//trim(diag%debug_str)//trim(mesg))
463 endif
464
465 if ( size(field,2) == dszj ) then
466 jsv = diag_cs%js ; jev = diag_cs%je ; j_data = .true. ! Data domain
467 elseif ( size(field,2) == dszj + 1 ) then
468 jsv = diag_cs%js ; jev = diag_cs%je+1 ; j_data = .true. ! Symmetric data domain
469 elseif ( size(field,2) == cszj ) then
470 jsv = 1 ; jev = cszj ; j_data = .false. ! Computational domain
471 ! This was: elseif ( size(field,1) == cszj + 1 ) then
472 elseif ( size(field,2) == cszj + 1 ) then
473 jsv = 1 ; jev = cszj+1 ; j_data = .false. ! Symmetric computational domain
474 else
475 write (mesg,*) " peculiar size ",size(field,2)," in j-direction\n"//&
476 "does not match one of ", cszj, cszj+1, dszj, dszj+1
477 call mom_error(fatal,"post_IS_data_2d_low: "//trim(diag%debug_str)//trim(mesg))
478 endif
479
480 if ((diag%conversion_factor /= 0.) .and. (diag%conversion_factor /= 1.)) then
481 allocate( locfield( lbound(field,1):ubound(field,1), lbound(field,2):ubound(field,2) ) )
482 do j=jsv,jev ; do i=isv,iev
483 if (field(i,j) == diag_cs%missing_value) then
484 locfield(i,j) = diag_cs%missing_value
485 else
486 locfield(i,j) = field(i,j) * diag%conversion_factor
487 endif
488 enddo ; enddo
489 locfield(isv:iev,jsv:jev) = field(isv:iev,jsv:jev) * diag%conversion_factor
490 else
491 locfield => field
492 endif
493
494 ! Handle cases where the data and computational domain are the same size.
495 if (diag_cs%ied-diag_cs%isd == diag_cs%ie-diag_cs%is) i_data = j_data
496 if (diag_cs%jed-diag_cs%jsd == diag_cs%je-diag_cs%js) j_data = i_data
497 if ( i_data .NEQV. j_data ) then
498 call mom_error(fatal, "post_IS_data_2d: post_IS_data called for "//&
499 trim(diag%debug_str)//" with mixed computational and data domain array sizes.")
500 endif
501
502 if (present(mask)) then
503 locmask => mask
504 elseif (.not.is_stat) then ! Static fields do not have assigned axes.
505 if (i_data .and. associated(diag%axes%mask2d)) then
506 locmask => diag%axes%mask2d
507 elseif ((.not.i_data) .and. associated(diag%axes%mask2d_comp)) then
508 locmask => diag%axes%mask2d_comp
509 endif
510 endif
511 if (associated(locmask)) call assert(size(locfield) == size(locmask), &
512 'post_data_2d_low: mask size mismatch: '//trim(diag%debug_str))
513
514 if (diag_cs%diag_as_chksum) then
515 ! Append timestep to mesg
516 call get_time(diag_cs%time_end, time_seconds, days=time_days)
517 write(debug_mesg, '(a, 1x, i0, 1x, i0)') &
518 trim(diag%debug_str), time_days, time_seconds
519
520 if (diag%axes%is_h_point) then
521 call hchksum(locfield, debug_mesg, diag_cs%G%HI, &
522 logunit=diag_cs%chksum_iounit)
523 elseif (diag%axes%is_u_point) then
524 call uchksum(locfield, debug_mesg, diag_cs%G%HI, &
525 logunit=diag_cs%chksum_iounit)
526 elseif (diag%axes%is_v_point) then
527 call vchksum(locfield, debug_mesg, diag_cs%G%HI, &
528 logunit=diag_cs%chksum_iounit)
529 elseif (diag%axes%is_q_point) then
530 call bchksum(locfield, debug_mesg, diag_cs%G%HI, &
531 logunit=diag_cs%chksum_iounit)
532 else
533 call mom_error(fatal, "post_data_2d_low: unknown axis type.")
534 endif
535 else
536 if (is_stat) then
537 if (associated(locmask)) then
538 used = send_data_infra(diag%fms_diag_id, locfield, &
539 is_in=isv, ie_in=iev, js_in=jsv, je_in=jev, rmask=locmask)
540 else
541 used = send_data_infra(diag%fms_diag_id, locfield, &
542 is_in=isv, ie_in=iev, js_in=jsv, je_in=jev)
543 endif
544 elseif (diag_cs%ave_enabled) then
545 if (associated(locmask)) then
546 used = send_data_infra(diag%fms_diag_id, locfield, &
547 is_in=isv, ie_in=iev, js_in=jsv, je_in=jev, &
548 time=diag_cs%time_end, weight=diag_cs%time_int, rmask=locmask)
549 else
550 used = send_data_infra(diag%fms_diag_id, locfield, &
551 is_in=isv, ie_in=iev, js_in=jsv, je_in=jev, &
552 time=diag_cs%time_end, weight=diag_cs%time_int)
553 endif
554 endif
555 endif
556
557 if ((diag%conversion_factor /= 0.) .and. (diag%conversion_factor /= 1.)) deallocate( locfield )
558
559end subroutine post_data_2d_low
560
561!> Enable the accumulation of time averages over the specified time interval.
562subroutine enable_averaging(time_int_in, time_end_in, diag_cs)
563 real, intent(in) :: time_int_in !< The time interval [s] over which any
564 !! values that are offered are valid.
565 type(time_type), intent(in) :: time_end_in !< The end time of the valid interval
566 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
567 ! This subroutine enables the accumulation of time averages over the specified time interval.
568
569! if (num_file==0) return
570 diag_cs%time_int = time_int_in
571 diag_cs%time_end = time_end_in
572 diag_cs%ave_enabled = .true.
573end subroutine enable_averaging
574
575!> Enable the accumulation of time averages over the specified time interval in time units.
576subroutine enable_averages(time_int, time_end, diag_CS, T_to_s)
577 real, intent(in) :: time_int !< The time interval over which any values
578 !! that are offered are valid [T ~> s].
579 type(time_type), intent(in) :: time_end !< The end time of the valid interval.
580 type(diag_ctrl), intent(inout) :: diag_cs !< A structure that is used to regulate diagnostic output
581 real, optional, intent(in) :: t_to_s !< A conversion factor for time_int to seconds [s T-1 ~> 1].
582 ! This subroutine enables the accumulation of time averages over the specified time interval.
583
584 if (present(t_to_s)) then
585 diag_cs%time_int = time_int*t_to_s
586 elseif (associated(diag_cs%US)) then
587 diag_cs%time_int = time_int*diag_cs%US%T_to_s
588 else
589 diag_cs%time_int = time_int
590 endif
591 diag_cs%time_end = time_end
592 diag_cs%ave_enabled = .true.
593end subroutine enable_averages
594
595!> Call this subroutine to avoid averaging any offered fields.
596subroutine disable_averaging(diag_cs)
597 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
598
599 diag_cs%time_int = 0.0
600 diag_cs%ave_enabled = .false.
601end subroutine disable_averaging
602
603!> Indicate whether averaging diagnostics is currently enabled
604logical function query_averaging_enabled(diag_cs, time_int, time_end)
605 type(diag_ctrl), intent(in) :: diag_cs !< Structure used to regulate diagnostic output
606 real, optional, intent(out) :: time_int !< Current setting of diag_cs%time_int [s]
607 type(time_type), optional, intent(out) :: time_end !< Current setting of diag_cs%time_end
608
609 if (present(time_int)) time_int = diag_cs%time_int
610 if (present(time_end)) time_end = diag_cs%time_end
611 query_averaging_enabled = diag_cs%ave_enabled
612end function query_averaging_enabled
613
614!> This subroutine initializes the diag_manager via the MOM6 infrastructure
616 character(len=*), optional, intent(out) :: err_msg !< An error message
617
618 call mom_diag_manager_init(err_msg=err_msg)
620
621!> This function returns the valid end time for use with diagnostics that are
622!! handled outside of the MOM6 diagnostics infrastructure.
623function get_diag_time_end(diag_cs)
624 type(diag_ctrl), intent(in) :: diag_cs !< Structure used to regulate diagnostic output
625 type(time_type) :: get_diag_time_end
626 ! This function returns the valid end time for diagnostics that are handled
627 ! outside of the MOM6 infrastructure, such as via the generic tracer code.
628
629 get_diag_time_end = diag_cs%time_end
630end function get_diag_time_end
631
632!> Returns the "diag_mediator" handle for a group (native, CMOR, ...) of diagnostics
633!! derived from one field.
634function register_mom_is_diag_field(module_name, field_name, axes_in, init_time, &
635 long_name, units, missing_value, range, mask_variant, standard_name, &
636 verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, &
637 cmor_long_name, cmor_units, cmor_standard_name, cell_methods, &
638 x_cell_method, y_cell_method, conversion) result (register_diag_field)
639 integer :: register_diag_field !< The returned diagnostic handle
640 character(len=*), intent(in) :: module_name !< Name of this module, usually "ice_model"
641 character(len=*), intent(in) :: field_name !< Name of the diagnostic field
642 type(axes_grp), target, intent(in) :: axes_in !< Container with up to 3 integer handles that
643 !! indicates axes for this field
644 type(time_type), intent(in) :: init_time !< Time at which a field is first available?
645 character(len=*), optional, intent(in) :: long_name !< Long name of a field.
646 character(len=*), optional, intent(in) :: units !< Units of a field.
647 character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
648 real, optional, intent(in) :: missing_value !< A value that indicates missing values in
649 !! output files, in unscaled arbitrary units [a]
650 real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
651 !! in arbitrary units [a]
652 logical, optional, intent(in) :: mask_variant !< If true a logical mask must be provided with
653 !! post_IS_data calls (not used in MOM?)
654 logical, optional, intent(in) :: verbose !< If true, FMS is verbose (not used in MOM?)
655 logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
656 character(len=*), optional, intent(out):: err_msg !< String into which an error message might be
657 !! placed (not used in MOM?)
658 character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should not
659 !! be interpolated as a scalar
660 integer, optional, intent(in) :: tile_count !< no clue (not used in MOM?)
661 character(len=*), optional, intent(in) :: cmor_field_name !< CMOR name of a field
662 character(len=*), optional, intent(in) :: cmor_long_name !< CMOR long name of a field
663 character(len=*), optional, intent(in) :: cmor_units !< CMOR units of a field
664 character(len=*), optional, intent(in) :: cmor_standard_name !< CMOR standardized name associated with a field
665 character(len=*), optional, intent(in) :: cell_methods !< String to append as cell_methods attribute. Use '' to
666 !! have no attribute. If present, this overrides the
667 !! default constructed from the default for
668 !! each individual axis direction.
669 character(len=*), optional, intent(in) :: x_cell_method !< Specifies the cell method for the x-direction.
670 !! Use '' have no method.
671 character(len=*), optional, intent(in) :: y_cell_method !< Specifies the cell method for the y-direction.
672 !! Use '' have no method.
673 real, optional, intent(in) :: conversion !< A value to multiply data by before writing to files,
674 !! often including factors to undo internal scaling and
675 !! in units of [a A-1 ~> 1]
676
677 ! Local variables
678 real :: mom_missing_value ! A value used to indicate missing values in output files, in arbitrary units [a]
679 type(diag_ctrl), pointer :: diag_cs => null() ! A structure that is used to regulate diagnostic output
680 type(axes_grp), pointer :: axes
681 integer :: dm_id
682 character(len=256) :: msg
683 character(len=256) :: cm_string ! A string describing the cell methods returned from attach_cell_methods.
684 character(len=256) :: new_module_name
685 character(len=480) :: module_list, var_list
686 character(len=24) :: dimensions
687 integer :: num_modnm, num_varnm
688 logical :: active
689
690 diag_cs => axes_in%diag_cs
691
692 ! Check if the axes match a standard grid axis.
693 ! If not, allocate the new axis and copy the contents.
694 if (axes_in%id == diag_cs%axesT1%id) then
695 axes => diag_cs%axesT1
696 elseif (axes_in%id == diag_cs%axesB1%id) then
697 axes => diag_cs%axesB1
698 elseif (axes_in%id == diag_cs%axesCu1%id) then
699 axes => diag_cs%axesCu1
700 elseif (axes_in%id == diag_cs%axesCv1%id) then
701 axes => diag_cs%axesCv1
702 else
703 allocate(axes)
704 axes = axes_in
705 endif
706
707 mom_missing_value = axes%diag_cs%missing_value
708 if (present(missing_value)) mom_missing_value = missing_value
709
710 diag_cs => axes%diag_cs
711 dm_id = -1
712
713 module_list = "{"//trim(module_name)
714 num_modnm = 1
715
716 ! Register the native diagnostic
717 active = register_diag_field_expand_cmor(dm_id, module_name, field_name, axes, &
718 init_time, long_name=long_name, units=units, missing_value=mom_missing_value, &
719 range=range, mask_variant=mask_variant, standard_name=standard_name, &
720 verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
721 interp_method=interp_method, tile_count=tile_count, &
722 cmor_field_name=cmor_field_name, cmor_long_name=cmor_long_name, &
723 cmor_units=cmor_units, cmor_standard_name=cmor_standard_name, &
724 cell_methods=cell_methods, x_cell_method=x_cell_method, y_cell_method=y_cell_method, &
725 conversion=conversion)
726 num_varnm = 1 ; var_list = "{"//trim(field_name)
727 if (present(cmor_field_name)) then
728 num_varnm = num_varnm + 1
729 var_list = trim(var_list)//","//trim(cmor_field_name)
730 endif
731 var_list = trim(var_list)//"}"
732
733 dimensions = ""
734 if (axes_in%is_h_point) dimensions = trim(dimensions)//" xh, yh,"
735 if (axes_in%is_q_point) dimensions = trim(dimensions)//" xq, yq,"
736 if (axes_in%is_u_point) dimensions = trim(dimensions)//" xq, yh,"
737 if (axes_in%is_v_point) dimensions = trim(dimensions)//" xh, yq,"
738 if (len_trim(dimensions) > 0) dimensions = trim_trailing_commas(dimensions)
739
740 if (is_root_pe() .and. (diag_cs%available_diag_doc_unit > 0)) then
741 msg = ''
742 if (present(cmor_field_name)) msg = 'CMOR equivalent is "'//trim(cmor_field_name)//'"'
743 call attach_cell_methods(-1, axes, cm_string, cell_methods, x_cell_method, y_cell_method)
744 module_list = trim(module_list)//"}"
745 if (num_modnm <= 1) module_list = module_name
746 if (num_varnm <= 1) var_list = ''
747
748 call log_available_diag(dm_id>0, module_list, field_name, cm_string, msg, diag_cs, &
749 long_name, units, standard_name, variants=var_list, dimensions=dimensions)
750 endif
751
752 register_diag_field = dm_id
753
755
756!> Returns True if either the native or CMOR version of the diagnostic were registered. Updates 'dm_id'
757!! after calling register_diag_field_expand_axes() for both native and CMOR variants of the field.
758logical function register_diag_field_expand_cmor(dm_id, module_name, field_name, axes, init_time, &
759 long_name, units, missing_value, range, mask_variant, standard_name, &
760 verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, &
761 cmor_long_name, cmor_units, cmor_standard_name, cell_methods, &
762 x_cell_method, y_cell_method, conversion)
763 integer, intent(inout) :: dm_id !< The diag_mediator ID for this diagnostic group
764 character(len=*), intent(in) :: module_name !< Name of this module, usually "ice_model" or "ice_model_fast"
765 character(len=*), intent(in) :: field_name !< Name of the diagnostic field
766 type(axes_grp), intent(in) :: axes !< Container with up to 3 integer handles that indicates axes
767 !! for this field
768 type(time_type), intent(in) :: init_time !< Time at which a field is first available?
769 character(len=*), optional, intent(in) :: long_name !< Long name of a field.
770 character(len=*), optional, intent(in) :: units !< Units of a field.
771 character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
772 real, optional, intent(in) :: missing_value !< A value that indicates missing values in
773 !! output files, in unscaled arbitrary units [a]
774 real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
775 !! in arbitrary units [a]
776 logical, optional, intent(in) :: mask_variant !< If true a logical mask must be provided
777 !! with post_data calls (not used in MOM?)
778 logical, optional, intent(in) :: verbose !< If true, FMS is verbose (not used in MOM?)
779 logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
780 character(len=*), optional, intent(out):: err_msg !< String into which an error message might be
781 !! placed (not used in MOM?)
782 character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should
783 !! not be interpolated as a scalar
784 integer, optional, intent(in) :: tile_count !< no clue (not used in MOM?)
785 character(len=*), optional, intent(in) :: cmor_field_name !< CMOR name of a field
786 character(len=*), optional, intent(in) :: cmor_long_name !< CMOR long name of a field
787 character(len=*), optional, intent(in) :: cmor_units !< CMOR units of a field
788 character(len=*), optional, intent(in) :: cmor_standard_name !< CMOR standardized name associated with a field
789 character(len=*), optional, intent(in) :: cell_methods !< String to append as cell_methods attribute.
790 !! Use '' to have no attribute. If present, this
791 !! overrides the default constructed from the default
792 !! for each individual axis direction.
793 character(len=*), optional, intent(in) :: x_cell_method !< Specifies the cell method for the x-direction.
794 !! Use '' have no method.
795 character(len=*), optional, intent(in) :: y_cell_method !< Specifies the cell method for the y-direction.
796 !! Use '' have no method.
797 real, optional, intent(in) :: conversion !< A value to multiply data by before writing to files,
798 !! often including factors to undo internal scaling and
799 !! in units of [a A-1 ~> 1]
800 ! Local variables
801 real :: mom_missing_value ! A value used to indicate missing values in output files, in arbitrary units [a]
802 type(diag_ctrl), pointer :: diag_cs => null()
803 type(diag_type), pointer :: this_diag => null()
804 integer :: fms_id
805 character(len=256) :: posted_cmor_units, posted_cmor_standard_name, posted_cmor_long_name
806 character(len=256) :: cm_string ! A string describing the cell methods returned from attach_cell_methods.
807
808 mom_missing_value = axes%diag_cs%missing_value
809 if (present(missing_value)) mom_missing_value = missing_value
810
812 diag_cs => axes%diag_cs
813
814 ! Set up the 'primary' diagnostic, first get an underlying FMS id
815 fms_id = register_diag_field_expand_axes(module_name, field_name, axes, init_time, &
816 long_name=long_name, units=units, missing_value=mom_missing_value, &
817 range=range, mask_variant=mask_variant, standard_name=standard_name, &
818 verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
819 interp_method=interp_method, tile_count=tile_count)
820 if (.not. diag_cs%diag_as_chksum) &
821 call attach_cell_methods(fms_id, axes, cm_string, cell_methods, x_cell_method, y_cell_method)
822
823 this_diag => null()
824 if (fms_id /= diag_field_not_found) then
825 call add_diag_to_list(diag_cs, dm_id, fms_id, this_diag, axes, module_name, field_name)
826 if (present(conversion)) this_diag%conversion_factor = conversion
828 endif
829
830 ! For the CMOR variation of the above diagnostic
831 if (present(cmor_field_name) .and. .not. diag_cs%diag_as_chksum) then
832 ! Fallback values for strings set to "NULL"
833 posted_cmor_units = "not provided" !
834 posted_cmor_standard_name = "not provided" ! Values might be able to be replaced with a CS%missing field?
835 posted_cmor_long_name = "not provided" !
836
837 ! If attributes are present for MOM variable names, use them first for the register_MOM_IS_diag_field
838 ! call for CMOR verison of the variable
839 if (present(units)) posted_cmor_units = units
840 if (present(standard_name)) posted_cmor_standard_name = standard_name
841 if (present(long_name)) posted_cmor_long_name = long_name
842
843 ! If specified in the call to register_MOM_IS_diag_field, override attributes with the CMOR versions
844 if (present(cmor_units)) posted_cmor_units = cmor_units
845 if (present(cmor_standard_name)) posted_cmor_standard_name = cmor_standard_name
846 if (present(cmor_long_name)) posted_cmor_long_name = cmor_long_name
847
848 fms_id = register_diag_field_expand_axes(module_name, cmor_field_name, axes, init_time, &
849 long_name=trim(posted_cmor_long_name), units=trim(posted_cmor_units), &
850 missing_value=mom_missing_value, range=range, mask_variant=mask_variant, &
851 standard_name=trim(posted_cmor_standard_name), verbose=verbose, do_not_log=do_not_log, &
852 err_msg=err_msg, interp_method=interp_method, tile_count=tile_count)
853 call attach_cell_methods(fms_id, axes, cm_string, cell_methods, x_cell_method, y_cell_method)
854
855 this_diag => null()
856 if (fms_id /= diag_field_not_found) then
857 call add_diag_to_list(diag_cs, dm_id, fms_id, this_diag, axes, module_name, field_name)
858 if (present(conversion)) this_diag%conversion_factor = conversion
860 endif
861 endif
862
864
865!> Returns an FMS id from register_diag_field_fms (the diag_manager routine) after expanding axes
866!! (axes-group) into handles and conditionally adding an FMS area_id for cell_measures.
867integer function register_diag_field_expand_axes(module_name, field_name, axes, init_time, &
868 long_name, units, missing_value, range, mask_variant, standard_name, &
869 verbose, do_not_log, err_msg, interp_method, tile_count)
870 character(len=*), intent(in) :: module_name !< Name of this module, usually "ocean_model"
871 !! or "ice_shelf_model"
872 character(len=*), intent(in) :: field_name !< Name of the diagnostic field
873 type(axes_grp), target, intent(in) :: axes !< Container with up to 3 integer handles that indicates
874 !! axes for this field
875 type(time_type), intent(in) :: init_time !< Time at which a field is first available?
876 character(len=*), optional, intent(in) :: long_name !< Long name of a field.
877 character(len=*), optional, intent(in) :: units !< Units of a field.
878 character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
879 real, optional, intent(in) :: missing_value !< A value that indicates missing values in
880 !! output files, in unscaled arbitrary units [a]
881 real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
882 !! in arbitrary units [a]
883 logical, optional, intent(in) :: mask_variant !< If true a logical mask must be provided
884 !! with post_data calls (not used in MOM?)
885 logical, optional, intent(in) :: verbose !< If true, FMS is verbose (not used in MOM?)
886 logical, optional, intent(in) :: do_not_log !< If true, do not log something
887 !! (not used in MOM?)
888 character(len=*), optional, intent(out):: err_msg !< String into which an error message might be
889 !! placed (not used in MOM?)
890 character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should
891 !! not be interpolated as a scalar
892 integer, optional, intent(in) :: tile_count !< no clue (not used in MOM?)
893 ! Local variables
894 integer :: fms_id, area_id
895
896 ! This gets the cell area associated with the grid location of this variable
897 area_id = axes%id_area
898
899 ! Get the FMS diagnostic id
900 if (axes%diag_cs%diag_as_chksum) then
901 fms_id = axes%diag_cs%num_chksum_diags + 1
902 axes%diag_cs%num_chksum_diags = fms_id
903 elseif (present(interp_method) .or. axes%is_h_point) then
904 ! If interp_method is provided we must use it
905 if (area_id>0) then
906 fms_id = register_diag_field_infra(module_name, field_name, axes%handles, &
907 init_time, long_name=long_name, units=units, missing_value=missing_value, &
908 range=range, mask_variant=mask_variant, standard_name=standard_name, &
909 verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
910 interp_method=interp_method, tile_count=tile_count, area=area_id)
911 else
912 fms_id = register_diag_field_infra(module_name, field_name, axes%handles, &
913 init_time, long_name=long_name, units=units, missing_value=missing_value, &
914 range=range, mask_variant=mask_variant, standard_name=standard_name, &
915 verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
916 interp_method=interp_method, tile_count=tile_count)
917 endif
918 else
919 ! If interp_method is not provided and the field is not at an h-point then interp_method='none'
920 if (area_id>0) then
921 fms_id = register_diag_field_infra(module_name, field_name, axes%handles, &
922 init_time, long_name=long_name, units=units, missing_value=missing_value, &
923 range=range, mask_variant=mask_variant, standard_name=standard_name, &
924 verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
925 interp_method='none', tile_count=tile_count, area=area_id)
926 else
927 fms_id = register_diag_field_infra(module_name, field_name, axes%handles, &
928 init_time, long_name=long_name, units=units, missing_value=missing_value, &
929 range=range, mask_variant=mask_variant, standard_name=standard_name, &
930 verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
931 interp_method='none', tile_count=tile_count)
932 endif
933 endif
934
936
938
939!> Create a diagnostic type and attached to list
940subroutine add_diag_to_list(diag_cs, dm_id, fms_id, this_diag, axes, module_name, field_name)
941 type(diag_ctrl), pointer :: diag_cs !< Diagnostics mediator control structure
942 integer, intent(inout) :: dm_id !< The diag_mediator ID for this diagnostic group
943 integer, intent(in) :: fms_id !< The FMS diag_manager ID for this diagnostic
944 type(diag_type), pointer :: this_diag !< This diagnostic
945 type(axes_grp), target, intent(in) :: axes !< Container with up to 3 integer handles that
946 !! indicates axes for this field
947 character(len=*), intent(in) :: module_name !< Name of this module, usually
948 !! "ocean_model" or "ice_shelf_model"
949 character(len=*), intent(in) :: field_name !< Name of diagnostic
950
951 ! If the diagnostic is needed obtain a diag_mediator ID (if needed)
952 if (dm_id == -1) dm_id = get_new_diag_id(diag_cs)
953 ! Create a new diag_type to store links in
954 call alloc_diag_with_id(dm_id, diag_cs, this_diag)
955 call assert(associated(this_diag), 'add_diag_to_list: allocation failed for '//trim(field_name))
956 ! Record FMS id, masks and conversion factor, in diag_type
957 this_diag%fms_diag_id = fms_id
958 this_diag%debug_str = trim(module_name)//"-"//trim(field_name)
959 this_diag%axes => axes
960
961end subroutine add_diag_to_list
962
963
964!> Attaches "cell_methods" attribute to a variable based on defaults for axes_grp or optional arguments.
965subroutine attach_cell_methods(id, axes, ostring, cell_methods, x_cell_method, y_cell_method)
966 integer, intent(in) :: id !< Handle to diagnostic
967 type(axes_grp), intent(in) :: axes !< Container with up to 3 integer handles that indicates
968 !! axes for this field
969 character(len=*), intent(out) :: ostring !< The cell_methods strings that would appear in the file
970 character(len=*), optional, intent(in) :: cell_methods !< String to append as cell_methods attribute.
971 !! Use '' to have no attribute. If present, this
972 !! overrides the default constructed from the default
973 !! for each individual axis direction.
974 character(len=*), optional, intent(in) :: x_cell_method !< Specifies the cell method for the x-direction.
975 !! Use '' have no method.
976 character(len=*), optional, intent(in) :: y_cell_method !< Specifies the cell method for the y-direction.
977 !! Use '' have no method.
978 ! Local variables
979 character(len=9) :: axis_name
980 logical :: x_mean, y_mean, x_sum, y_sum
981
982 x_mean = .false.
983 y_mean = .false.
984 x_sum = .false.
985 y_sum = .false.
986
987 ostring = ''
988 if (present(cell_methods)) then
989 if (present(x_cell_method) .or. present(y_cell_method)) then
990 call mom_error(fatal, "attach_cell_methods: " // &
991 'Individual direction cell method was specified along with a "cell_methods" string.')
992 endif
993 if (len(trim(cell_methods))>0) then
994 call mom_diag_field_add_attribute(id, 'cell_methods', trim(cell_methods))
995 ostring = trim(cell_methods)
996 endif
997 else
998 if (present(x_cell_method)) then
999 if (len(trim(x_cell_method))>0) then
1000 call get_mom_diag_axis_name(axes%handles(1), axis_name)
1001 call mom_diag_field_add_attribute(id, 'cell_methods', trim(axis_name)//':'//trim(x_cell_method))
1002 ostring = trim(adjustl(ostring))//' '//trim(axis_name)//':'//trim(x_cell_method)
1003 if (trim(x_cell_method)=='mean') x_mean=.true.
1004 if (trim(x_cell_method)=='sum') x_sum=.true.
1005 endif
1006 else
1007 if (len(trim(axes%x_cell_method))>0) then
1008 call get_mom_diag_axis_name(axes%handles(1), axis_name)
1009 call mom_diag_field_add_attribute(id, 'cell_methods', trim(axis_name)//':'//trim(axes%x_cell_method))
1010 ostring = trim(adjustl(ostring))//' '//trim(axis_name)//':'//trim(axes%x_cell_method)
1011 if (trim(axes%x_cell_method)=='mean') x_mean=.true.
1012 if (trim(axes%x_cell_method)=='sum') x_sum=.true.
1013 endif
1014 endif
1015 if (present(y_cell_method)) then
1016 if (len(trim(y_cell_method))>0) then
1017 call get_mom_diag_axis_name(axes%handles(2), axis_name)
1018 call mom_diag_field_add_attribute(id, 'cell_methods', trim(axis_name)//':'//trim(y_cell_method))
1019 ostring = trim(adjustl(ostring))//' '//trim(axis_name)//':'//trim(y_cell_method)
1020 if (trim(y_cell_method)=='mean') y_mean=.true.
1021 if (trim(y_cell_method)=='sum') y_sum=.true.
1022 endif
1023 else
1024 if (len(trim(axes%y_cell_method))>0) then
1025 call get_mom_diag_axis_name(axes%handles(2), axis_name)
1026 call mom_diag_field_add_attribute(id, 'cell_methods', trim(axis_name)//':'//trim(axes%y_cell_method))
1027 ostring = trim(adjustl(ostring))//' '//trim(axis_name)//':'//trim(axes%y_cell_method)
1028 if (trim(axes%y_cell_method)=='mean') y_mean=.true.
1029 if (trim(axes%y_cell_method)=='sum') y_sum=.true.
1030 endif
1031 endif
1032 if (x_mean .and. y_mean) then
1033 call mom_diag_field_add_attribute(id, 'cell_methods', 'area:mean')
1034 ostring = trim(adjustl(ostring))//' area:mean'
1035 elseif (x_sum .and. y_sum) then
1036 call mom_diag_field_add_attribute(id, 'cell_methods', 'area:sum')
1037 ostring = trim(adjustl(ostring))//' area:sum'
1038 endif
1039 endif
1040 ostring = adjustl(ostring)
1041end subroutine attach_cell_methods
1042
1043!> Registers a non-array scalar diagnostic, returning an integer handle
1044function register_scalar_field_axes(module_name, field_name, axes, init_time, &
1045 long_name, units, missing_value, range, standard_name, &
1046 do_not_log, err_msg, interp_method, cmor_field_name, &
1047 cmor_long_name, cmor_units, cmor_standard_name, conversion) result (register_scalar_field)
1048 integer :: register_scalar_field !< An integer handle for a diagnostic array.
1049 character(len=*), intent(in) :: module_name !< Name of this module, usually "ocean_model"
1050 !! or "ice_shelf_model"
1051 character(len=*), intent(in) :: field_name !< Name of the diagnostic field
1052 type(axes_grp), target, intent(in) :: axes !< Container with up to 3 integer handles that
1053 !! indicates axes for this field
1054 type(time_type), intent(in) :: init_time !< Time at which a field is first available?
1055 character(len=*), optional, intent(in) :: long_name !< Long name of a field.
1056 character(len=*), optional, intent(in) :: units !< Units of a field.
1057 character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
1058 real, optional, intent(in) :: missing_value !< A value that indicates missing values in
1059 !! output files, in unscaled arbitrary units [a]
1060 real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
1061 !! in arbitrary units [a]
1062 logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
1063 character(len=*), optional, intent(out):: err_msg !< String into which an error message might be
1064 !! placed (not used in MOM?)
1065 character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should not
1066 !! be interpolated as a scalar
1067 character(len=*), optional, intent(in) :: cmor_field_name !< CMOR name of a field
1068 character(len=*), optional, intent(in) :: cmor_long_name !< CMOR long name of a field
1069 character(len=*), optional, intent(in) :: cmor_units !< CMOR units of a field
1070 character(len=*), optional, intent(in) :: cmor_standard_name !< CMOR standardized name associated with a field
1071 real, optional, intent(in) :: conversion !< A value to multiply data by before writing to files,
1072 !! often including factors to undo internal scaling and
1073 !! in units of [a A-1 ~> 1]
1074
1075 register_scalar_field = register_scalar_field_cs(module_name, field_name, init_time, axes%diag_cs, &
1076 long_name, units, missing_value, range, standard_name, &
1077 do_not_log, err_msg, interp_method, cmor_field_name, &
1078 cmor_long_name, cmor_units, cmor_standard_name, conversion)
1079
1080end function register_scalar_field_axes
1081
1082!> Registers a non-array scalar diagnostic, returning an integer handle
1083function register_scalar_field_cs(module_name, field_name, init_time, diag_cs, &
1084 long_name, units, missing_value, range, standard_name, &
1085 do_not_log, err_msg, interp_method, cmor_field_name, &
1086 cmor_long_name, cmor_units, cmor_standard_name, conversion) result (register_scalar_field)
1087 integer :: register_scalar_field !< An integer handle for a diagnostic array.
1088 character(len=*), intent(in) :: module_name !< Name of this module, usually "ocean_model"
1089 !! or "ice_shelf_model"
1090 character(len=*), intent(in) :: field_name !< Name of the diagnostic field
1091 type(time_type), intent(in) :: init_time !< Time at which a field is first available?
1092 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
1093 character(len=*), optional, intent(in) :: long_name !< Long name of a field.
1094 character(len=*), optional, intent(in) :: units !< Units of a field.
1095 character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
1096 real, optional, intent(in) :: missing_value !< A value that indicates missing values in
1097 !! output files, in unscaled arbitrary units [a]
1098 real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
1099 !! in arbitrary units [a]
1100 logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
1101 character(len=*), optional, intent(out):: err_msg !< String into which an error message might be
1102 !! placed (not used in MOM?)
1103 character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should not
1104 !! be interpolated as a scalar
1105 character(len=*), optional, intent(in) :: cmor_field_name !< CMOR name of a field
1106 character(len=*), optional, intent(in) :: cmor_long_name !< CMOR long name of a field
1107 character(len=*), optional, intent(in) :: cmor_units !< CMOR units of a field
1108 character(len=*), optional, intent(in) :: cmor_standard_name !< CMOR standardized name associated with a field
1109 real, optional, intent(in) :: conversion !< A value to multiply data by before writing to files,
1110 !! often including factors to undo internal scaling and
1111 !! in units of [a A-1 ~> 1]
1112
1113 ! Local variables
1114 real :: mom_missing_value ! A value used to indicate missing values in output files, in arbitrary units [a]
1115 integer :: dm_id, fms_id
1116 type(diag_type), pointer :: diag => null(), cmor_diag => null()
1117 character(len=256) :: posted_cmor_units, posted_cmor_standard_name, posted_cmor_long_name
1118 character(len=16) :: dimensions
1119
1120 mom_missing_value = diag_cs%missing_value
1121 if (present(missing_value)) mom_missing_value = missing_value
1122
1123 dm_id = -1
1124 diag => null()
1125 cmor_diag => null()
1126
1127 if (diag_cs%diag_as_chksum) then
1128 fms_id = diag_cs%num_chksum_diags + 1
1129 diag_cs%num_chksum_diags = fms_id
1130 else
1131 fms_id = register_diag_field_infra(module_name, field_name, init_time, &
1132 long_name=long_name, units=units, missing_value=mom_missing_value, &
1133 range=range, standard_name=standard_name, do_not_log=do_not_log, &
1134 err_msg=err_msg)
1135 endif
1136
1137 if (fms_id /= diag_field_not_found) then
1138 dm_id = get_new_diag_id(diag_cs)
1139 call alloc_diag_with_id(dm_id, diag_cs, diag)
1140 call assert(associated(diag), 'register_scalar_field: diag allocation failed')
1141 diag%fms_diag_id = fms_id
1142 diag%debug_str = trim(module_name)//"-"//trim(field_name)
1143 if (present(conversion)) diag%conversion_factor = conversion
1144 endif
1145
1146 if (present(cmor_field_name)) then
1147 ! Fallback values for strings set to "not provided"
1148 posted_cmor_units = "not provided"
1149 posted_cmor_standard_name = "not provided"
1150 posted_cmor_long_name = "not provided"
1151
1152 ! If attributes are present for MOM variable names, use them as defaults for the
1153 ! register_diag_field_infra call for CMOR verison of the variable
1154 if (present(units)) posted_cmor_units = units
1155 if (present(standard_name)) posted_cmor_standard_name = standard_name
1156 if (present(long_name)) posted_cmor_long_name = long_name
1157
1158 ! If specified in the call to register_MOM_IS_scalar_field, override attributes with the CMOR versions
1159 if (present(cmor_units)) posted_cmor_units = cmor_units
1160 if (present(cmor_standard_name)) posted_cmor_standard_name = cmor_standard_name
1161 if (present(cmor_long_name)) posted_cmor_long_name = cmor_long_name
1162
1163 fms_id = register_diag_field_infra(module_name, cmor_field_name, init_time, &
1164 long_name=trim(posted_cmor_long_name), units=trim(posted_cmor_units), &
1165 missing_value=mom_missing_value, range=range, &
1166 standard_name=trim(posted_cmor_standard_name), do_not_log=do_not_log, err_msg=err_msg)
1167 if (fms_id /= diag_field_not_found) then
1168 if (dm_id == -1) then
1169 dm_id = get_new_diag_id(diag_cs)
1170 endif
1171 call alloc_diag_with_id(dm_id, diag_cs, cmor_diag)
1172 cmor_diag%fms_diag_id = fms_id
1173 cmor_diag%debug_str = trim(module_name)//"-"//trim(cmor_field_name)
1174 if (present(conversion)) cmor_diag%conversion_factor = conversion
1175 endif
1176 endif
1177
1178 dimensions = "scalar"
1179
1180 ! Document diagnostics in list of available diagnostics
1181 if (is_root_pe() .and. diag_cs%available_diag_doc_unit > 0) then
1182 if (present(cmor_field_name)) then
1183 call log_available_diag(associated(diag), module_name, field_name, '', '', diag_cs, &
1184 long_name, units, standard_name, &
1185 variants="{"//trim(field_name)//","//trim(cmor_field_name)//"}", &
1186 dimensions=dimensions)
1187 else
1188 call log_available_diag(associated(diag), module_name, field_name, '', '', diag_cs, &
1189 long_name, units, standard_name, dimensions=dimensions)
1190 endif
1191 endif
1192
1193 register_scalar_field = dm_id
1194
1195end function register_scalar_field_cs
1196
1197!> Registers a static diagnostic, returning an integer handle
1198function register_mom_is_static_field(module_name, field_name, axes, &
1199 long_name, units, missing_value, range, mask_variant, standard_name, &
1200 do_not_log, interp_method, tile_count, &
1201 cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name, area, &
1202 x_cell_method, y_cell_method, area_cell_method, conversion) result(register_static_field)
1203 integer :: register_static_field !< An integer handle for a diagnostic array.
1204 character(len=*), intent(in) :: module_name !< Name of this module, usually "ocean_model"
1205 !! or "ice_shelf_model"
1206 character(len=*), intent(in) :: field_name !< Name of the diagnostic field
1207 type(axes_grp), target, intent(in) :: axes !< Container with up to 3 integer handles that
1208 !! indicates axes for this field
1209 character(len=*), optional, intent(in) :: long_name !< Long name of a field.
1210 character(len=*), optional, intent(in) :: units !< Units of a field.
1211 character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
1212 real, optional, intent(in) :: missing_value !< A value that indicates missing values in
1213 !! output files, in unscaled arbitrary units [a]
1214 real, optional, intent(in) :: range(2) !< Valid range of a variable in arbitrary units [a]
1215 logical, optional, intent(in) :: mask_variant !< If true a logical mask must be provided with
1216 !! post_IS_data calls (not used in MOM?)
1217 logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
1218 character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should not
1219 !! be interpolated as a scalar
1220 integer, optional, intent(in) :: tile_count !< no clue (not used in MOM?)
1221 character(len=*), optional, intent(in) :: cmor_field_name !< CMOR name of a field
1222 character(len=*), optional, intent(in) :: cmor_long_name !< CMOR long name of a field
1223 character(len=*), optional, intent(in) :: cmor_units !< CMOR units of a field
1224 character(len=*), optional, intent(in) :: cmor_standard_name !< CMOR standardized name associated with a field
1225 integer, optional, intent(in) :: area !< fms_id for area_t
1226 character(len=*), optional, intent(in) :: x_cell_method !< Specifies the cell method for the x-direction.
1227 character(len=*), optional, intent(in) :: y_cell_method !< Specifies the cell method for the y-direction.
1228 character(len=*), optional, intent(in) :: area_cell_method !< Specifies the cell method for area
1229 real, optional, intent(in) :: conversion !< A value to multiply data by before writing to files,
1230 !! often including factors to undo internal scaling and
1231 !! in units of [a A-1 ~> 1]
1232
1233 ! Local variables
1234 real :: mom_missing_value ! A value used to indicate missing values in output files, in arbitrary units [a]
1235 type(diag_ctrl), pointer :: diag_cs => null() !< A structure that is used to regulate diagnostic output
1236 type(diag_type), pointer :: diag => null(), cmor_diag => null()
1237 integer :: dm_id, fms_id
1238 character(len=256) :: posted_cmor_units, posted_cmor_standard_name, posted_cmor_long_name
1239 character(len=9) :: axis_name
1240 character(len=24) :: dimensions
1241
1242 mom_missing_value = axes%diag_cs%missing_value
1243 if (present(missing_value)) mom_missing_value = missing_value
1244
1245 diag_cs => axes%diag_cs
1246 dm_id = -1
1247 diag => null()
1248 cmor_diag => null()
1249
1250 if (diag_cs%diag_as_chksum) then
1251 fms_id = diag_cs%num_chksum_diags + 1
1252 diag_cs%num_chksum_diags = fms_id
1253 else
1254 fms_id = register_static_field_infra(module_name, field_name, axes%handles, &
1255 long_name=long_name, units=units, missing_value=mom_missing_value, &
1256 range=range, mask_variant=mask_variant, standard_name=standard_name, &
1257 do_not_log=do_not_log, &
1258 interp_method=interp_method, tile_count=tile_count, area=area)
1259 endif
1260
1261 if (fms_id /= diag_field_not_found) then
1262 dm_id = get_new_diag_id(diag_cs)
1263 call alloc_diag_with_id(dm_id, diag_cs, diag)
1264 call assert(associated(diag), 'register_static_field: diag allocation failed')
1265 diag%fms_diag_id = fms_id
1266 diag%debug_str = trim(module_name)//"-"//trim(field_name)
1267 if (present(conversion)) diag%conversion_factor = conversion
1268
1269 if (diag_cs%diag_as_chksum) then
1270 diag%axes => axes
1271 else
1272 if (present(x_cell_method)) then
1273 call get_mom_diag_axis_name(axes%handles(1), axis_name)
1274 call mom_diag_field_add_attribute(fms_id, 'cell_methods', &
1275 trim(axis_name)//':'//trim(x_cell_method))
1276 endif
1277 if (present(y_cell_method)) then
1278 call get_mom_diag_axis_name(axes%handles(2), axis_name)
1279 call mom_diag_field_add_attribute(fms_id, 'cell_methods', &
1280 trim(axis_name)//':'//trim(y_cell_method))
1281 endif
1282 if (present(area_cell_method)) then
1283 call mom_diag_field_add_attribute(fms_id, 'cell_methods', &
1284 'area:'//trim(area_cell_method))
1285 endif
1286 endif
1287 endif
1288
1289 if (present(cmor_field_name) .and. .not. diag_cs%diag_as_chksum) then
1290 ! Fallback values for strings set to "not provided"
1291 posted_cmor_units = "not provided"
1292 posted_cmor_standard_name = "not provided"
1293 posted_cmor_long_name = "not provided"
1294
1295 ! If attributes are present for MOM variable names, use them first for the register_static_field
1296 ! call for CMOR verison of the variable
1297 if (present(units)) posted_cmor_units = units
1298 if (present(standard_name)) posted_cmor_standard_name = standard_name
1299 if (present(long_name)) posted_cmor_long_name = long_name
1300
1301 ! If specified in the call to register_static_field, override attributes with the CMOR versions
1302 if (present(cmor_units)) posted_cmor_units = cmor_units
1303 if (present(cmor_standard_name)) posted_cmor_standard_name = cmor_standard_name
1304 if (present(cmor_long_name)) posted_cmor_long_name = cmor_long_name
1305
1306 fms_id = register_static_field_infra(module_name, cmor_field_name, axes%handles, &
1307 long_name=trim(posted_cmor_long_name), units=trim(posted_cmor_units), &
1308 missing_value=mom_missing_value, range=range, mask_variant=mask_variant, &
1309 standard_name=trim(posted_cmor_standard_name), do_not_log=do_not_log, &
1310 interp_method=interp_method, tile_count=tile_count, area=area)
1311 if (fms_id /= diag_field_not_found) then
1312 if (dm_id == -1) then
1313 dm_id = get_new_diag_id(diag_cs)
1314 endif
1315 call alloc_diag_with_id(dm_id, diag_cs, cmor_diag)
1316 cmor_diag%fms_diag_id = fms_id
1317 cmor_diag%debug_str = trim(module_name)//"-"//trim(cmor_field_name)
1318 if (present(conversion)) cmor_diag%conversion_factor = conversion
1319 if (present(x_cell_method)) then
1320 call get_mom_diag_axis_name(axes%handles(1), axis_name)
1321 call mom_diag_field_add_attribute(fms_id, 'cell_methods', trim(axis_name)//':'//trim(x_cell_method))
1322 endif
1323 if (present(y_cell_method)) then
1324 call get_mom_diag_axis_name(axes%handles(2), axis_name)
1325 call mom_diag_field_add_attribute(fms_id, 'cell_methods', trim(axis_name)//':'//trim(y_cell_method))
1326 endif
1327 if (present(area_cell_method)) then
1328 call mom_diag_field_add_attribute(fms_id, 'cell_methods', 'area:'//trim(area_cell_method))
1329 endif
1330 endif
1331 endif
1332
1333 dimensions = ""
1334 if (axes%is_h_point) dimensions = trim(dimensions)//" xh, yh,"
1335 if (axes%is_q_point) dimensions = trim(dimensions)//" xq, yq,"
1336 if (axes%is_u_point) dimensions = trim(dimensions)//" xq, yh,"
1337 if (axes%is_v_point) dimensions = trim(dimensions)//" xh, yq,"
1338 if (len_trim(dimensions) > 0) dimensions = trim_trailing_commas(dimensions)
1339
1340 ! Document diagnostics in list of available diagnostics
1341 if (is_root_pe() .and. diag_cs%available_diag_doc_unit > 0) then
1342 if (present(cmor_field_name)) then
1343 call log_available_diag(associated(diag), module_name, field_name, '', '', diag_cs, &
1344 long_name, units, standard_name, &
1345 variants="{"//trim(field_name)//","//trim(cmor_field_name)//"}", &
1346 dimensions=dimensions)
1347 else
1348 call log_available_diag(associated(diag), module_name, field_name, '', '', diag_cs, &
1349 long_name, units, standard_name, dimensions=dimensions)
1350 endif
1351 endif
1352
1353 register_static_field = dm_id
1354
1356
1357!> Add a description of an option to the documentation file
1358subroutine describe_option(opt_name, value, diag_CS)
1359 character(len=*), intent(in) :: opt_name !< The name of the option
1360 character(len=*), intent(in) :: value !< The value of the option
1361 type(diag_ctrl), intent(in) :: diag_CS !< Structure used to regulate diagnostic output
1362
1363 ! Local variables
1364 character(len=480) :: mesg
1365 integer :: len_ind
1366
1367 len_ind = len_trim(value)
1368
1369 mesg = " ! "//trim(opt_name)//": "//trim(value)
1370 write(diag_cs%available_diag_doc_unit, '(a)') trim(mesg)
1371end subroutine describe_option
1372
1373!> Initialize the MOM_IS diag_mediator and opens the available diagnostics file, if appropriate.
1374subroutine mom_is_diag_mediator_init(G, US, param_file, diag_cs, component, err_msg, &
1375 doc_file_dir)
1376 type(ocean_grid_type), target, intent(inout) :: g !< The horizontal grid type
1377 type(unit_scale_type), target, intent(in) :: us !< A dimensional unit scaling type
1378 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1379 type(diag_ctrl), intent(inout) :: diag_cs !< A structure that is used to regulate diagnostic output
1380 character(len=*), optional, intent(in) :: component !< An optional component name
1381 character(len=*), optional, intent(out) :: err_msg !< A string for a returned error message
1382 character(len=*), optional, intent(in) :: doc_file_dir !< A directory in which to create the file
1383
1384 ! This subroutine initializes the diag_mediator and the diag_manager.
1385 ! The grid type should have its dimensions set by this point, but it
1386 ! is not necessary that the metrics and axis labels be set up yet.
1387
1388 ! Local variables
1389 integer :: ios, i, new_unit
1390 logical :: opened, new_file
1391 character(len=8) :: this_pe
1392 character(len=240) :: doc_file, doc_file_dflt, doc_path
1393 character(len=40) :: doc_file_param
1394 ! This include declares and sets the variable "version".
1395# include "version_variable.h"
1396 character(len=40) :: mdl = "MOM_IS_diag_mediator" ! This module's name.
1397 character(len=32) :: filename_appendix = '' !fms appendix to filename for ensemble runs
1398
1399 call mom_diag_manager_init(err_msg=err_msg)
1400
1401 id_clock_diag_mediator = cpu_clock_id('(Ice shelf diagnostics framework)', grain=clock_module)
1402
1403 ! Allocate and initialize list of all diagnostics (and variants)
1404 allocate(diag_cs%diags(diag_alloc_chunk_size))
1405 diag_cs%next_free_diag_id = 1
1406 do i=1, diag_alloc_chunk_size
1407 call initialize_diag_type(diag_cs%diags(i))
1408 enddo
1409
1410 diag_cs%show_call_tree = calltree_showquery()
1411
1412 ! Read all relevant parameters and write them to the model log.
1413 call log_version(param_file, mdl, version, "")
1414
1415 call get_param(param_file, mdl, 'USE_INDEX_DIAGNOSTIC_AXES', diag_cs%index_space_axes, &
1416 'If true, use a grid index coordinate convention for diagnostic axes. ',&
1417 default=.false.)
1418
1419 call get_param(param_file, mdl, 'DIAG_MISVAL', diag_cs%missing_value, &
1420 'Set the default missing value to use for diagnostics.', &
1421 units="various", default=-1.e34)
1422 call get_param(param_file, mdl, 'DIAG_AS_CHKSUM', diag_cs%diag_as_chksum, &
1423 'Instead of writing diagnostics to the diag manager, write '//&
1424 'a text file containing the checksum (bitcount) of the array.', &
1425 default=.false.)
1426
1427 if (diag_cs%diag_as_chksum) &
1428 diag_cs%num_chksum_diags = 0
1429
1430 ! Keep pointers to the grid for diagnostic checksums
1431 diag_cs%G => g
1432 diag_cs%US => us
1433
1434 diag_cs%is = g%isc - (g%isd-1) ; diag_cs%ie = g%iec - (g%isd-1)
1435 diag_cs%js = g%jsc - (g%jsd-1) ; diag_cs%je = g%jec - (g%jsd-1)
1436 diag_cs%isd = g%isd ; diag_cs%ied = g%ied
1437 diag_cs%jsd = g%jsd ; diag_cs%jed = g%jed
1438
1439 ! Initialize available diagnostic log file
1440 if (is_root_pe() .and. (diag_cs%available_diag_doc_unit < 0)) then
1441 if (present(component)) then
1442 doc_file_dflt = trim(component)//".available_diags"
1443 doc_file_param = trim(uppercase(component))//"_AVAILABLE_DIAGS_FILE"
1444 else
1445 write(this_pe,'(i6.6)') pe_here()
1446 doc_file_dflt = "MOM_IS.available_diags."//this_pe
1447 doc_file_param = "AVAILABLE_MOM_IS_DIAGS_FILE"
1448 endif
1449 call get_param(param_file, mdl, trim(doc_file_param), doc_file, &
1450 "A file into which to write a list of all available "//&
1451 "ice shelf diagnostics that can be included in a diag_table.", &
1452 default=doc_file_dflt, do_not_log=(diag_cs%available_diag_doc_unit/=-1))
1453 if (len_trim(doc_file) > 0) then
1454 new_file = .true. ; if (diag_cs%available_diag_doc_unit /= -1) new_file = .false.
1455 ! Find an unused unit number.
1456 do new_unit=512,42,-1
1457 inquire( new_unit, opened=opened)
1458 if (.not.opened) exit
1459 enddo
1460 if (opened) call mom_error(fatal, &
1461 "diag_mediator_init failed to find an unused unit number.")
1462
1463 doc_path = doc_file
1464 if (present(doc_file_dir)) then ; if (len_trim(doc_file_dir) > 0) then
1465 doc_path = trim(slasher(doc_file_dir))//trim(doc_file)
1466 endif ; endif
1467
1468 diag_cs%available_diag_doc_unit = new_unit
1469
1470 if (new_file) then
1471 open(diag_cs%available_diag_doc_unit, file=trim(doc_path), access='SEQUENTIAL', form='FORMATTED', &
1472 action='WRITE', status='REPLACE', iostat=ios)
1473 else ! This file is being reopened, and should be appended.
1474 open(diag_cs%available_diag_doc_unit, file=trim(doc_path), access='SEQUENTIAL', form='FORMATTED', &
1475 action='WRITE', status='OLD', position='APPEND', iostat=ios)
1476 endif
1477 inquire(diag_cs%available_diag_doc_unit, opened=opened)
1478 if ((.not.opened) .or. (ios /= 0)) then
1479 call mom_error(fatal, "Failed to open available diags file "//trim(doc_path)//".")
1480 endif
1481 endif
1482 endif
1483
1484 if (is_root_pe() .and. (diag_cs%chksum_iounit < 0) .and. diag_cs%diag_as_chksum) then
1485 !write(this_pe,'(i6.6)') PE_here()
1486 !doc_file_dflt = "chksum_diag."//this_pe
1487 doc_file_dflt = "chksum_diag"
1488 call get_param(param_file, mdl, "CHKSUM_DIAG_FILE", doc_file, &
1489 "A file into which to write all checksums of the "//&
1490 "diagnostics listed in the diag_table.", &
1491 default=doc_file_dflt, do_not_log=(diag_cs%chksum_iounit/=-1))
1492
1493 call get_filename_appendix(filename_appendix)
1494 if (len_trim(filename_appendix) > 0) then
1495 doc_file = trim(doc_file) //'.'//trim(filename_appendix)
1496 endif
1497#ifdef STATSLABEL
1498 doc_file = trim(doc_file)//"."//trim(adjustl(statslabel))
1499#endif
1500
1501 if (len_trim(doc_file) > 0) then
1502 new_file = .true. ; if (diag_cs%chksum_iounit /= -1) new_file = .false.
1503 ! Find an unused unit number.
1504 do new_unit=512,42,-1
1505 inquire( new_unit, opened=opened)
1506 if (.not.opened) exit
1507 enddo
1508 if (opened) call mom_error(fatal, &
1509 "diag_mediator_init failed to find an unused unit number.")
1510
1511 doc_path = doc_file
1512 if (present(doc_file_dir)) then ; if (len_trim(doc_file_dir) > 0) then
1513 doc_path = trim(slasher(doc_file_dir))//trim(doc_file)
1514 endif ; endif
1515
1516 diag_cs%chksum_iounit = new_unit
1517
1518 if (new_file) then
1519 open(diag_cs%chksum_iounit, file=trim(doc_path), access='SEQUENTIAL', form='FORMATTED', &
1520 action='WRITE', status='REPLACE', iostat=ios)
1521 else ! This file is being reopened, and should be appended.
1522 open(diag_cs%chksum_iounit, file=trim(doc_path), access='SEQUENTIAL', form='FORMATTED', &
1523 action='WRITE', status='OLD', position='APPEND', iostat=ios)
1524 endif
1525 inquire(diag_cs%chksum_iounit, opened=opened)
1526 if ((.not.opened) .or. (ios /= 0)) then
1527 call mom_error(fatal, "Failed to open checksum diags file "//trim(doc_path)//".")
1528 endif
1529 endif
1530 endif
1531
1532 call diag_masks_set(g, diag_cs%missing_value, diag_cs)
1533
1534end subroutine mom_is_diag_mediator_init
1535
1536!> Sets up the 2d masks for native diagnostics
1537subroutine diag_masks_set(G, missing_value, diag_cs)
1538 type(ocean_grid_type), target, intent(in) :: g !< The horizontal grid type
1539 real, intent(in) :: missing_value !< A fill value for missing points
1540 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
1541
1542 ! Local variables
1543 integer :: i, j
1544
1545 ! 2d masks point to the model masks since they are identical
1546 diag_cs%mask2dT => g%mask2dT
1547 diag_cs%mask2dBu => g%mask2dBu
1548 diag_cs%mask2dCu => g%mask2dCu
1549 diag_cs%mask2dCv => g%mask2dCv
1550
1551 allocate(diag_cs%mask2dT_comp(g%isc:g%iec,g%jsc:g%jec))
1552 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1553 diag_cs%mask2dT_comp(i,j) = diag_cs%mask2dT(i,j)
1554 enddo ; enddo
1555
1556 diag_cs%missing_value = missing_value
1557
1558end subroutine diag_masks_set
1559
1560!> Prevent the registration of additional diagnostics, so that the creation of files can occur
1562 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
1563
1564 if (diag_cs%available_diag_doc_unit > -1) then
1565 close(diag_cs%available_diag_doc_unit) ; diag_cs%available_diag_doc_unit = -2
1566 endif
1567
1569
1570!> Deallocate memory associated with the MOM_IS diag mediator
1571subroutine mom_is_diag_mediator_end(diag_CS)
1572 type(diag_ctrl), intent(inout) :: diag_cs !< Structure used to regulate diagnostic output
1573
1574 ! Local variables
1575 type(diag_type), pointer :: diag, next_diag
1576 integer :: i
1577
1578 if (diag_cs%available_diag_doc_unit > -1) then
1579 close(diag_cs%available_diag_doc_unit) ; diag_cs%available_diag_doc_unit = -3
1580 endif
1581 if (diag_cs%chksum_iounit > -1) then
1582 close(diag_cs%chksum_iounit) ; diag_cs%chksum_iounit = -3
1583 endif
1584
1585 do i=1, diag_cs%next_free_diag_id - 1
1586 if (associated(diag_cs%diags(i)%next)) then
1587 next_diag => diag_cs%diags(i)%next
1588 do while (associated(next_diag))
1589 diag => next_diag
1590 next_diag => diag%next
1591 deallocate(diag)
1592 enddo
1593 endif
1594 enddo
1595
1596 deallocate(diag_cs%diags)
1597
1598 ! These points to arrays in the grid type, so they can not be deallocated here.
1599 if (associated(diag_cs%mask2dT)) diag_cs%mask2dT => null()
1600 if (associated(diag_cs%mask2dBu)) diag_cs%mask2dBu => null()
1601 if (associated(diag_cs%mask2dCu)) diag_cs%mask2dCu => null()
1602 if (associated(diag_cs%mask2dCv)) diag_cs%mask2dCv => null()
1603 if (associated(diag_cs%mask2dT_comp)) deallocate(diag_cs%mask2dT_comp)
1604
1605end subroutine mom_is_diag_mediator_end
1606
1607!> Returns a new diagnostic id, it may be necessary to expand the diagnostics array.
1608integer function get_new_diag_id(diag_cs)
1609 type(diag_ctrl), intent(inout) :: diag_cs !< Diagnostics control structure
1610 ! Local variables
1611 type(diag_type), dimension(:), allocatable :: tmp
1612 integer :: i
1613
1614 if (diag_cs%next_free_diag_id > size(diag_cs%diags)) then
1615 call assert(diag_cs%next_free_diag_id - size(diag_cs%diags) == 1, &
1616 'get_new_diag_id: inconsistent diag id')
1617
1618 ! Increase the size of diag_cs%diags and copy data over.
1619 ! Do not use move_alloc() because it is not supported by Fortran 90
1620 allocate(tmp(size(diag_cs%diags)))
1621 tmp(:) = diag_cs%diags(:)
1622 deallocate(diag_cs%diags)
1623 allocate(diag_cs%diags(size(tmp) + diag_alloc_chunk_size))
1624 diag_cs%diags(1:size(tmp)) = tmp(:)
1625 deallocate(tmp)
1626
1627 ! Initialize new part of the diag array.
1628 do i=diag_cs%next_free_diag_id, size(diag_cs%diags)
1629 call initialize_diag_type(diag_cs%diags(i))
1630 enddo
1631 endif
1632
1633 get_new_diag_id = diag_cs%next_free_diag_id
1634 diag_cs%next_free_diag_id = diag_cs%next_free_diag_id + 1
1635
1636end function get_new_diag_id
1637
1638!> Initializes a diag_type (used after allocating new memory)
1639subroutine initialize_diag_type(diag)
1640 type(diag_type), intent(inout) :: diag !< diag_type to be initialized
1641
1642 diag%in_use = .false.
1643 diag%fms_diag_id = -1
1644 diag%axes => null()
1645 diag%next => null()
1646 diag%conversion_factor = 0.
1647
1648end subroutine initialize_diag_type
1649
1650!> Make a new diagnostic. Either use memory which is in the array of 'primary'
1651!! diagnostics, or if that is in use, insert it to the list of secondary diags.
1652subroutine alloc_diag_with_id(diag_id, diag_cs, diag)
1653 integer, intent(in ) :: diag_id !< id for the diagnostic
1654 type(diag_ctrl), target, intent(inout) :: diag_cs !< structure used to regulate diagnostic output
1655 type(diag_type), pointer :: diag !< structure representing a diagnostic (inout)
1656
1657 type(diag_type), pointer :: tmp => null()
1658
1659 if (.not. diag_cs%diags(diag_id)%in_use) then
1660 diag => diag_cs%diags(diag_id)
1661 else
1662 allocate(diag)
1663 tmp => diag_cs%diags(diag_id)%next
1664 diag_cs%diags(diag_id)%next => diag
1665 diag%next => tmp
1666 endif
1667 diag%in_use = .true.
1668
1669end subroutine alloc_diag_with_id
1670
1671!> Log a diagnostic to the available diagnostics file.
1672subroutine log_available_diag(used, module_name, field_name, cell_methods_string, comment, &
1673 diag_CS, long_name, units, standard_name, variants, dimensions)
1674 logical, intent(in) :: used !< Whether this diagnostic was in the diag_table or not
1675 character(len=*), intent(in) :: module_name !< Name of the diagnostic module
1676 character(len=*), intent(in) :: field_name !< Name of this diagnostic field
1677 character(len=*), intent(in) :: cell_methods_string !< The spatial component of the CF cell_methods attribute
1678 character(len=*), intent(in) :: comment !< A comment to append after [Used|Unused]
1679 type(diag_ctrl), intent(in) :: diag_CS !< The diagnotics control structure
1680 character(len=*), optional, intent(in) :: dimensions !< Descriptor of the horizontal and vertical dimensions
1681 character(len=*), optional, intent(in) :: long_name !< CF long name of diagnostic
1682 character(len=*), optional, intent(in) :: units !< Units for diagnostic
1683 character(len=*), optional, intent(in) :: standard_name !< CF standardized name of diagnostic
1684 character(len=*), optional, intent(in) :: variants !< Alternate modules and variable names for
1685 !! this diagnostic and derived diagnostics
1686 ! Local variables
1687 character(len=240) :: mesg
1688
1689 if (used) then
1690 mesg = '"'//trim(field_name)//'" [Used]'
1691 else
1692 mesg = '"'//trim(field_name)//'" [Unused]'
1693 endif
1694 if (len(trim((comment)))>0) then
1695 write(diag_cs%available_diag_doc_unit, '(a,1x,"(",a,")")') trim(mesg),trim(comment)
1696 else
1697 write(diag_cs%available_diag_doc_unit, '(a)') trim(mesg)
1698 endif
1699 call describe_option("modules", module_name, diag_cs)
1700 if (present(dimensions)) then ; if (len(trim(dimensions)) > 0) then
1701 call describe_option("dimensions", dimensions, diag_cs)
1702 endif ; endif
1703 if (present(long_name)) call describe_option("long_name", long_name, diag_cs)
1704 if (present(units)) call describe_option("units", units, diag_cs)
1705 if (present(standard_name)) &
1706 call describe_option("standard_name", standard_name, diag_cs)
1707 if (len(trim((cell_methods_string)))>0) &
1708 call describe_option("cell_methods", trim(cell_methods_string), diag_cs)
1709 if (present(variants)) then ; if (len(trim(variants)) > 0) then
1710 call describe_option("variants", variants, diag_cs)
1711 endif ; endif
1712end subroutine log_available_diag
1713
1714!> Log the diagnostic chksum to the chksum diag file
1715subroutine log_chksum_diag(docunit, description, chksum)
1716 integer, intent(in) :: docunit !< Handle of the log file
1717 character(len=*), intent(in) :: description !< Name of the diagnostic module
1718 integer, intent(in) :: chksum !< chksum of the diagnostic
1719
1720 write(docunit, '(a,1x,i9.8)') description, chksum
1721 flush(docunit)
1722
1723end subroutine log_chksum_diag
1724
1725!> Fakes a register of a diagnostic to find out if an obsolete
1726!! parameter appears in the diag_table.
1727logical function found_in_diagtable(diag, varName)
1728 type(diag_ctrl), intent(in) :: diag !< A structure used to control diagnostics.
1729 character(len=*), intent(in) :: varname !< The obsolete diagnostic name
1730 ! Local
1731 integer :: handle ! Integer handle returned from diag_manager
1732
1733 ! We use register_static_field_fms() instead of register_static_field() so
1734 ! that the diagnostic does not appear in the available diagnostics list.
1735 handle = register_static_field_infra('ice_shelf_model', varname, diag%axesT1%handles)
1736
1737 found_in_diagtable = (handle>0)
1738
1739end function found_in_diagtable
1740
1741!> Finishes the diag manager reduction methods as needed for the time_step
1742subroutine mom_is_diag_send_complete()
1743 call diag_send_complete_infra()
1744end subroutine mom_is_diag_send_complete
1745
1746end module mom_is_diag_mediator