FMS  2024.03
Flexible Modeling System
diag_axis.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup diag_axis_mod diag_axis_mod
20 !> @ingroup diag_manager
21 !> @brief An integral part of @ref diag_manager_mod. It helps to create axis IDs
22 !! that are used in @ref register_diag_field.
23 !!
24 !> @author Seth Underwood
25 !!
26 !! Users first create axis ID by calling diag_axis_init, then use this axis ID in
27 !! register_diag_field.
28 
29 !> @addtogroup diag_axis_mod
30 !> @{
31 MODULE diag_axis_mod
32 use platform_mod
33 
34  USE mpp_domains_mod, ONLY: domainug, domain1d, domain2d, mpp_get_compute_domain,&
35  & mpp_get_domain_components, null_domain1d, null_domain2d, null_domainug,&
36  & north, east, center, &
38  USE fms_mod, ONLY: error_mesg, write_version_number, lowercase, uppercase,&
39  & fms_error_handler, fatal, note
40  USE diag_data_mod, ONLY: diag_axis_type, max_subaxes, max_axes,&
41  & max_num_axis_sets, max_axis_attributes, debug_diag_manager,&
42  & first_send_data_call, diag_atttype, use_modern_diag
43  use fms_diag_object_mod, only:fms_diag_object
44  USE netcdf, ONLY: nf90_int, nf90_float, nf90_char
45 
46  IMPLICIT NONE
47 
48  PRIVATE
56  & north, east, center, diag_axis_type
57 
58  ! Include variable "version" to be written to log file
59 #include<file_version.h>
60 
61 !----------
62  integer(I4_KIND),parameter,public :: DIAG_AXIS_NODOMAIN = 0 !< For unstructured grid support
63  integer(I4_KIND),parameter,public :: diag_axis_2ddomain = 1 !< For unstructured grid support
64  integer(I4_KIND),parameter,public :: diag_axis_ugdomain = 2 !< For unstructured grid support
65 !----------
66 
67  INTEGER, DIMENSION(:), ALLOCATABLE :: num_subaxes !< counter of number of axes defined
68  INTEGER :: num_def_axes = 0
69 
70  CHARACTER(len=128), DIMENSION(:), ALLOCATABLE, SAVE :: axis_sets !< storage for axis set names
71  INTEGER :: num_axis_sets = 0
72 
73  TYPE(diag_axis_type), ALLOCATABLE, SAVE :: axes(:) !< global storage for all defined axes
74  LOGICAL :: module_is_initialized = .false.
75 
76  !> @}
77 
78  !> @brief Add an arbitrary attribute and value to the diagnostic axis.
79  !!
80  !> Any number of attributes can be added to a given axis. All attribute addition must
81  !! be done before first <TT>send_data</TT> call.<br>
82  !!
83  !! If a real or integer attribute is already defined, a FATAL error will be called.
84  !! If a character attribute is already defined, then it will be prepended to the
85  !! existing attribute value.
86  !! <br>Example usage:
87  !! @code{.F90} call diag_axis_add_attribute(diag_axis_id, att_name, att_value) @endcode
88  !> @ingroup diag_axis_mod
90  MODULE PROCEDURE diag_axis_add_attribute_scalar_r
91  MODULE PROCEDURE diag_axis_add_attribute_scalar_i
92  MODULE PROCEDURE diag_axis_add_attribute_scalar_c
93  MODULE PROCEDURE diag_axis_add_attribute_r1d
94  MODULE PROCEDURE diag_axis_add_attribute_i1d
95  END INTERFACE diag_axis_add_attribute
96 
97  !> @addtogroup diag_axis_mod
98  !> @{
99 
100 CONTAINS
101 
102  !> @brief Initialize the axis, and return the axis ID.
103  !!
104  !> <TT>diag_axis_init</TT> initializes an axis and returns the axis ID that
105  !! is to be used with <TT>register_diag_field</TT>. This function also
106  !! increments the axis counter and fills in the axes
107  !!
108  !! @return integer axis ID
109  INTEGER FUNCTION diag_axis_init(name, array_data, units, cart_name, long_name, direction,&
110  & set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count, domain_position )
111  CHARACTER(len=*), INTENT(in) :: name !< Short name for axis
112  CLASS(*), DIMENSION(:), INTENT(in) :: array_data !< Array of coordinate values
113  CHARACTER(len=*), INTENT(in) :: units !< Units for the axis
114  CHARACTER(len=*), INTENT(in) :: cart_name !< Cartesian axis ("X", "Y", "Z", "T")
115  CHARACTER(len=*), INTENT(in), OPTIONAL :: long_name !< Long name for the axis.
116  CHARACTER(len=*), INTENT(in), OPTIONAL :: set_name
117  INTEGER, INTENT(in), OPTIONAL :: direction !< Indicates the direction of the axis
118  INTEGER, INTENT(in), OPTIONAL :: edges !< Axis ID for the previously defined "edges axis"
119  TYPE(domain1d), INTENT(in), OPTIONAL :: domain
120  TYPE(domain2d), INTENT(in), OPTIONAL :: domain2
121  TYPE(domainug), INTENT(in), OPTIONAL :: domainu
122  CHARACTER(len=*), INTENT(in), OPTIONAL :: aux !< Auxiliary name, can only be <TT>geolon_t</TT> or <TT>geolat_t</TT>
123  CHARACTER(len=*), INTENT(in), OPTIONAL :: req !< Required field names.
124  INTEGER, INTENT(in), OPTIONAL :: tile_count
125  INTEGER, INTENT(in), OPTIONAL :: domain_position
126 
127 
128  TYPE(domain1d) :: domain_x, domain_y
129  INTEGER :: ierr, axlen
130  INTEGER :: i, set, tile
131  INTEGER :: isc, iec, isg, ieg
132  CHARACTER(len=128) :: emsg
133 
134  IF ( .NOT.module_is_initialized ) THEN
135  CALL write_version_number("DIAG_AXIS_MOD", version)
136  ENDIF
137 
138  if (use_modern_diag) then
139  !TODO Passing in the axis_length because of a gnu issue where inside fms_diag_axis_init, the size of DATA
140  !was 2 which was causing the axis_data to not be written correctly...
141  diag_axis_init = fms_diag_object%fms_diag_axis_init(name, array_data, units, cart_name, size(array_data(:)), &
142  & long_name=long_name, direction=direction, set_name=set_name, edges=edges, domain=domain, domain2=domain2, &
143  & domainu=domainu, aux=aux, req=req, tile_count=tile_count, domain_position=domain_position)
144  return
145  endif
146  IF ( PRESENT(tile_count)) THEN
147  tile = tile_count
148  ELSE
149  tile = 1
150  END IF
151 
152  ! Allocate the axes
153  IF (.NOT. ALLOCATED(axis_sets)) ALLOCATE(axis_sets(max_num_axis_sets))
154  IF (.NOT. ALLOCATED(axes)) ALLOCATE(axes(max_axes))
155  IF (.NOT. ALLOCATED(num_subaxes)) THEN
156  ALLOCATE(num_subaxes(max_axes))
157  num_subaxes = 0
158  END IF
159 
160  !---- is there an axis set? ----
161  IF ( PRESENT(set_name) ) THEN
162  set = get_axis_set_num(set_name)
163  !---- add new set name ----
164  IF (set == 0) THEN
165  num_axis_sets = num_axis_sets + 1
166  IF ( num_axis_sets > max_num_axis_sets ) THEN
167  WRITE (emsg, fmt='("num_axis_sets (",I2,") exceeds max_num_axis_sets (",I2,"). ")')&
168  & num_axis_sets, max_num_axis_sets
169  ! <ERROR STATUS="FATAL">
170  ! num_axis_sets (<num_axis_sets>) exceeds max_num_axis_sets(<num_axis_sets>).
171  ! Increase max_num_axis_sets via diag_manager_nml.
172  ! </ERROR>
173  CALL error_mesg('diag_axis_mod::diag_axis_init',&
174  & trim(emsg)//' Increase max_num_axis_sets via diag_manager_nml.', fatal)
175  END IF
176  set = num_axis_sets
177  axis_sets(set) = set_name
178  END IF
179  ELSE
180  set = 0
181  END IF
182 
183  !---- see if axis already exists --
184  ! if this is time axis, return the ID of a previously defined
185  ! if this is spatial axis, FATAL error
186  DO i = 1, num_def_axes
187  IF ( trim(name) == axes(i)%name ) THEN
188  IF ( trim(name) == 'Stations' .OR. trim(name) == 'Levels') THEN
189  diag_axis_init = i
190  RETURN
191  ELSE IF ( set == axes(i)%set ) THEN
192  IF ( trim(lowercase(name)) == 'time' .OR.&
193  & trim(lowercase(cart_name)) == 't' .OR.&
194  & trim(lowercase(name)) == 'nv' .OR.&
195  & trim(lowercase(cart_name)) == 'n' ) THEN
196  diag_axis_init = i
197  RETURN
198  ELSE IF ( (lowercase(cart_name) /= 'x' .AND. lowercase(cart_name) /= 'y')&
199  & .OR. tile /= axes(i)%tile_count) THEN
200  ! <ERROR STATUS="FATAL">axis_name <NAME> and axis_set already exist.</ERROR>
201  CALL error_mesg('diag_axis_mod::diag_axis_init',&
202  & 'axis_name '//trim(name)//' and axis_set already exist.', fatal)
203  END IF
204  END IF
205  END IF
206  END DO
207 
208  !---- register axis ----
209  num_def_axes = num_def_axes + 1
210  ! <ERROR STATUS="FATAL">max_axes exceeded, increase it via diag_manager_nml</ERROR>
211  IF (num_def_axes > max_axes) CALL error_mesg ('diag_axis_mod::diag_axis_init',&
212  & 'max_axes exceeded, increase via diag_manager_nml', fatal)
213  diag_axis_init = num_def_axes
214 
215  !---- check and then save cart_name name ----
216  IF ( trim(uppercase(cart_name)) == 'X' .OR.&
217  & trim(uppercase(cart_name)) == 'Y' .OR.&
218  & trim(uppercase(cart_name)) == 'Z' .OR.&
219  & trim(uppercase(cart_name)) == 'T' .OR.&
220  & trim(uppercase(cart_name)) == 'U' .OR.&
221  & trim(uppercase(cart_name)) == 'N' ) THEN
222  axes(diag_axis_init)%cart_name = trim(uppercase(cart_name))
223  ELSE
224  ! <ERROR STATUS="FATAL">Invalid cart_name name.</ERROR>
225  CALL error_mesg('diag_axis_mod::diag_axis_init', 'Invalid cart_name name. '//trim(uppercase(cart_name)), fatal)
226  END IF
227 
228  !---- allocate storage for coordinate values of axis ----
229  IF ( axes(diag_axis_init)%cart_name == 'T' ) THEN
230  axlen = 0
231  ELSE
232  axlen = SIZE(array_data(:))
233  END IF
234  ALLOCATE ( axes(diag_axis_init)%diag_type_data(1:axlen) )
235 
236  ! Initialize Axes(diag_axis_init)
237  axes(diag_axis_init)%name = trim(name)
238  SELECT TYPE (array_data)
239  TYPE IS (real(kind=r4_kind))
240  axes(diag_axis_init)%diag_type_data = array_data(1:axlen)
241  TYPE IS (real(kind=r8_kind))
242  axes(diag_axis_init)%diag_type_data = real(array_data(1:axlen))
243  CLASS DEFAULT
244  CALL error_mesg('diag_axis_mod::diag_axis_init',&
245  & 'The axis data is not one of the supported types of real(kind=4) or real(kind=8)', fatal)
246  END SELECT
247  axes(diag_axis_init)%units = units
248  axes(diag_axis_init)%length = axlen
249  axes(diag_axis_init)%set = set
250  ! start and end are used in subaxes information only
251  axes(diag_axis_init)%start = -1
252  axes(diag_axis_init)%end = -1
253  axes(diag_axis_init)%subaxis_name = ""
254  axes(diag_axis_init)%shift = 0
255  axes(diag_axis_init)%num_attributes = 0
256 
257  IF ( PRESENT(long_name) ) THEN
258  axes(diag_axis_init)%long_name = long_name
259  ELSE
260  axes(diag_axis_init)%long_name = name
261  END IF
262 
263  IF ( PRESENT(aux) ) THEN
264  axes(diag_axis_init)%aux = trim(aux)
265  ELSE
266  axes(diag_axis_init)%aux = 'none'
267  END IF
268 
269  IF ( PRESENT(req) ) THEN
270  axes(diag_axis_init)%req = trim(req)
271  ELSE
272  axes(diag_axis_init)%req = 'none'
273  END IF
274  IF ( PRESENT(domain_position) ) THEN
275  if (domain_position == north .or. domain_position == east .or. domain_position == center) then
276  axes(diag_axis_init)%domain_position = domain_position
277  else
278  CALL error_mesg('diag_axis_mod::diag_axis_init', "Position must be NORTH, EAST, or CENTER" ,&
279  fatal)
280  endif
281  ELSE
282  axes(diag_axis_init)%domain_position = center
283  END IF
284 
285  !---- axis direction (-1, 0, or +1) ----
286  IF ( PRESENT(direction) )THEN
287  IF ( abs(direction) /= 1 .AND. direction /= 0 )&
288  ! <ERROR STATUS="FATAL">direction must be 0, +1, or -1</ERROR>
289  & CALL error_mesg('diag_axis_mod::diag_axis_init', 'direction must be 0, +1 or -1', fatal)
290  axes(diag_axis_init)%direction = direction
291  ELSE
292  axes(diag_axis_init)%direction = 0
293  END IF
294 
295  !---- Handle the DomainU check
296  IF (present(domainu) .AND. (PRESENT(domain2) .OR. PRESENT(domain)) ) THEN
297  ! <ERROR STATUS="FATAL">Presence of DomainU and another Domain at the same time is prohibited</ERROR>
298  CALL error_mesg('diag_axis_mod::diag_axis_init',&
299  & 'Presence of DomainU and another Domain at the same time is prohibited', fatal)
300  !---- domain2d type ----
301  ELSE IF ( PRESENT(domain2) .AND. PRESENT(domain)) THEN
302  ! <ERROR STATUS="FATAL">Presence of both Domain and Domain2 at the same time is prohibited</ERROR>
303  CALL error_mesg('diag_axis_mod::diag_axis_init',&
304  & 'Presence of both Domain and Domain2 at the same time is prohibited', fatal)
305  ELSE IF ( PRESENT(domain2) .OR. PRESENT(domain)) THEN
306  IF ( axes(diag_axis_init)%cart_name /= 'X' .AND. axes(diag_axis_init)%cart_name /= 'Y') THEN
307  ! <ERROR STATUS="FATAL">Domain must not be present for an axis which is not in the X or Y direction.</ERROR>
308  CALL error_mesg('diag_axis_mod::diag_axis_init',&
309  & 'A Structured Domain must not be present for an axis which is not in the X or Y direction', fatal)
310  END IF
311  ELSE IF (present(domainu) .AND. axes(diag_axis_init)%cart_name /= 'U') THEN
312  CALL error_mesg('diag_axis_mod::diag_axis_init',&
313  & 'In the unstructured domain, the axis cart_name must be U', fatal)
314  END IF
315 
316  axes(diag_axis_init)%tile_count = tile
317 
318  IF ( PRESENT(domain2) ) THEN
319  axes(diag_axis_init)%Domain2 = domain2
320  CALL mpp_get_domain_components(domain2, domain_x, domain_y, tile_count=tile_count)
321  IF ( axes(diag_axis_init)%cart_name == 'X' ) axes(diag_axis_init)%Domain = domain_x
322  IF ( axes(diag_axis_init)%cart_name == 'Y' ) axes(diag_axis_init)%Domain = domain_y
323  axes(diag_axis_init)%DomainUG = null_domainug
324  ELSE IF ( PRESENT(domain)) THEN
325  !---- domain1d type ----
326  axes(diag_axis_init)%Domain2 = null_domain2d ! needed since not 2-D domain
327  axes(diag_axis_init)%Domain = domain
328  axes(diag_axis_init)%DomainUG = null_domainug
329  ELSE IF (present(domainu)) THEN
330  axes(diag_axis_init)%Domain2 = null_domain2d
331  axes(diag_axis_init)%Domain = null_domain1d
332  axes(diag_axis_init)%DomainUG = domainu
333  ELSE
334  axes(diag_axis_init)%Domain2 = null_domain2d
335  axes(diag_axis_init)%Domain = null_domain1d
336  axes(diag_axis_init)%DomainUG = null_domainug
337  END IF
338 
339  !--- set up the shift value for x-y axis
340  IF ( axes(diag_axis_init)%Domain .NE. null_domain1d ) THEN
341  CALL mpp_get_compute_domain(axes(diag_axis_init)%Domain, isc, iec)
342  CALL mpp_get_global_domain(axes(diag_axis_init)%Domain, isg, ieg)
343  IF ( axes(diag_axis_init)%length == ieg - isg + 2 ) THEN
344  axes(diag_axis_init)%shift = 1
345  END IF
346  END IF
347 
348  !---- have axis edges been defined ? ----
349  axes(diag_axis_init)%edges = 0
350  IF (PRESENT(edges) ) THEN
351  IF ( edges > 0 .AND. edges < num_def_axes ) THEN
352  ierr=0
353  IF ( axes(edges)%cart_name /= axes(diag_axis_init)%cart_name) ierr=1
354  IF ( axes(edges)%length /= axes(diag_axis_init)%length+1 ) ierr=ierr+2
355  IF ( axes(edges)%set /= axes(diag_axis_init)%set ) ierr=ierr+4
356  IF ( ierr > 0 ) THEN
357  ! <ERROR STATUS="FATAL">Edges axis does not match axis (code <CODE>).</ERROR>
358  WRITE (emsg,'("Edges axis does not match axis (code ",I1,").")') ierr
359  CALL error_mesg('diag_axis_mod::diag_axis_init', emsg, fatal)
360  END IF
361  axes(diag_axis_init)%edges = edges
362  ELSE
363  ! <ERROR STATUS="FATAL">Edges axis is not defined.</ERROR>
364  CALL error_mesg('diag_axis_mod::diag_axis_init', 'Edges axis is not defined', fatal)
365  END IF
366  END IF
367 
368  ! Module is now initialized
369  module_is_initialized = .true.
370 
371  END FUNCTION diag_axis_init
372 
373  !> @brief Create a subaxis on a parent axis.
374  !!
375  !> Given the ID of a parent axis, create a subaxis and fill it with data,
376  !! and return the ID of the corresponding subaxis.
377  !!
378  !! The subaxis is defined on the parent axis from <TT>start_indx</TT>
379  !! to <TT>end_indx</TT>.
380  !!
381  !! @return Integer ID of the corresponding subaxis.
382  INTEGER FUNCTION diag_subaxes_init(axis, subdata, start_indx, end_indx, domain_2d)
383  INTEGER, INTENT(in) :: axis !< ID of the parent axis
384  REAL, DIMENSION(:), INTENT(in) :: subdata !< Data of the subaxis
385  INTEGER, INTENT(in) :: start_indx !< Start index of the subaxis
386  INTEGER, INTENT(in) :: end_indx !< End index of the subaxis
387  TYPE(domain2d), INTENT(in), OPTIONAL :: domain_2d
388 
389  INTEGER :: i, nsub_axis, direction
390  INTEGER :: xbegin, xend, ybegin, yend
391  INTEGER :: ad_xbegin, ad_xend, ad_ybegin, ad_yend
392  CHARACTER(len=128) :: name, nsub_name
393  CHARACTER(len=128) :: units
394  CHARACTER(len=128) :: cart_name
395  CHARACTER(len=128) :: long_name
396  CHARACTER(len=128) :: emsg
397  LOGICAL :: subaxis_set, hasdomain
398 
399  ! there may be more than 1 subaxis on a parent axis, check for redundancy
400  nsub_axis = 0
401  subaxis_set = .false.
402 
403  IF ( PRESENT(domain_2d) ) THEN
404  hasdomain = .true.
405  CALL mpp_get_compute_domain(domain_2d, xbegin, xend, ybegin, yend)
406  ELSE
407  hasdomain = .false.
408  END IF
409  sa_search: DO i = 1, num_subaxes(axis)
410  IF ( start_indx == axes(axis)%start(i) .AND. end_indx == axes(axis)%end(i) ) THEN
411  IF ( hasdomain ) THEN
412  CALL mpp_get_compute_domain(axes(axis)%subaxis_domain2(i), ad_xbegin, ad_xend, ad_ybegin, ad_yend)
413  IF ( .NOT.((xbegin == ad_xbegin .AND. xend == ad_xend) .AND.&
414  & (ybegin == ad_ybegin .AND. yend == ad_yend)) ) THEN
415  cycle sa_search
416  END IF
417  END IF
418  nsub_axis = i
419  subaxis_set = .true. !subaxis already exists
420  name = trim(axes(axis)%subaxis_name(nsub_axis))
421  EXIT sa_search
422  END IF
423  END DO sa_search
424 
425  IF ( nsub_axis == 0 ) THEN ! create new subaxis
426  num_subaxes(axis) = num_subaxes(axis) + 1
427  IF (num_subaxes(axis) > max_subaxes) THEN
428  ! <ERROR STATUS="FATAL">max_subaxes (value <VALUE>) is too small. Consider increasing max_subaxes.</ERROR>
429  WRITE (emsg,'("max_subaxes (value ",I4,") is too small. Consider increasing max_subaxes.")') max_subaxes
430  CALL error_mesg('diag_axis_mod::diag_subaxes_init', emsg, fatal)
431  END IF
432  nsub_axis = num_subaxes(axis)
433  axes(axis)%start(nsub_axis) = start_indx
434  axes(axis)%end(nsub_axis) = end_indx
435  if ( hasdomain ) axes(axis)%subaxis_domain2(nsub_axis) = domain_2d
436  END IF
437 
438  ! Create new name for the subaxis from name of parent axis
439  ! If subaxis already exists, get the index and return
440  IF(subaxis_set) THEN
441  IF ( axes(axis)%set > 0 ) THEN
442  diag_subaxes_init = get_axis_num(name, set_name=trim(axis_sets(axes(axis)%set)))
443  ELSE
445  END IF
446  ELSE
447  ! get a new index for subaxis
448  !::sdu:: Need a check to allow larger numbers in the index number.
449  WRITE (nsub_name,'(I2.2)') nsub_axis
450  name = trim(axes(axis)%name)//'_sub'//trim(nsub_name)
451  axes(axis)%subaxis_name(nsub_axis) = name
452  long_name = trim(axes(axis)%long_name)
453  units = trim(axes(axis)%units)
454  cart_name = trim(axes(axis)%cart_name)
455  direction = axes(axis)%direction
456  IF (axes(axis)%set > 0) THEN
457  diag_subaxes_init = diag_axis_init(trim(name), subdata, trim(units), trim(cart_name), trim(long_name),&
458  & set_name=trim(axis_sets(axes(axis)%set)), direction=direction, domain2=domain_2d)
459  ELSE
460  diag_subaxes_init = diag_axis_init(trim(name), subdata, trim(units), trim(cart_name), trim(long_name),&
461  & direction=direction, domain2=domain_2d)
462  END IF
463  END IF
464  END FUNCTION diag_subaxes_init
465  !> @brief Return information about the axis with index ID
466  SUBROUTINE get_diag_axis(id, name, units, long_name, cart_name,&
467  & direction, edges, Domain, DomainU, array_data, num_attributes, attributes, domain_position)
468  CHARACTER(len=*), INTENT(out) :: name, units, long_name, cart_name
469  INTEGER, INTENT(in) :: id !< Axis ID
470  TYPE(domain1d), INTENT(out) :: domain
471  TYPE(domainug), INTENT(out) :: domainu
472  INTEGER, INTENT(out) :: direction !< Direction of data. (See <TT>@ref diag_axis_init</TT> for a description of
473  !! allowed values)
474  INTEGER, INTENT(out) :: edges !< Axis ID for the previously defined "edges axis".
475  CLASS(*), DIMENSION(:), INTENT(out) :: array_data !< Array of coordinate values for this axis.
476  INTEGER, INTENT(out), OPTIONAL :: num_attributes
477  TYPE(diag_atttype), ALLOCATABLE, DIMENSION(:), INTENT(out), OPTIONAL :: attributes
478  INTEGER, INTENT(out), OPTIONAL :: domain_position
479 
480  INTEGER :: i, j, istat
481 
482  CALL valid_id_check(id, 'get_diag_axis')
483  name = axes(id)%name
484  units = axes(id)%units
485  long_name = axes(id)%long_name
486  cart_name = axes(id)%cart_name
487  direction = axes(id)%direction
488  edges = axes(id)%edges
489  domain = axes(id)%Domain
490  domainu = axes(id)%DomainUG
491  if (present(domain_position)) domain_position = axes(id)%domain_position
492  IF ( axes(id)%length > SIZE(array_data(:)) ) THEN
493  ! <ERROR STATUS="FATAL">array data is too small.</ERROR>
494  CALL error_mesg('diag_axis_mod::get_diag_axis', 'array data is too small', fatal)
495  ELSE
496  SELECT TYPE (array_data)
497  TYPE IS (real(kind=r4_kind))
498  array_data(1:axes(id)%length) = real(axes(id)%diag_type_data(1:axes(id)%length), kind=r4_kind)
499  TYPE IS (real(kind=r8_kind))
500  array_data(1:axes(id)%length) = axes(id)%diag_type_data(1:axes(id)%length)
501  CLASS DEFAULT
502  CALL error_mesg('diag_axis_mod::get_diag_axis',&
503  & 'The axis data is not one of the supported types of real(kind=4) or real(kind=8)', fatal)
504  END SELECT
505  END IF
506  IF ( PRESENT(num_attributes) ) THEN
507  num_attributes = axes(id)%num_attributes
508  END IF
509  IF ( PRESENT(attributes) ) THEN
510  IF ( allocated(axes(id)%attributes) ) THEN
511  IF ( ALLOCATED(attributes) ) THEN
512  ! If allocate, make sure attributes is large enough to hold Axis(id)%attributes
513  IF ( axes(id)%num_attributes .GT. SIZE(attributes(:)) ) THEN
514  CALL error_mesg('diag_axis_mod::get_diag_axis', 'array attribute is too small', fatal)
515  END IF
516  ELSE
517  ! Allocate attributes
518  ALLOCATE(attributes(axes(id)%num_attributes), stat=istat)
519  IF ( istat .NE. 0 ) THEN
520  CALL error_mesg('diag_axis_mod::get_diag_axis', 'Unable to allocate memory for attribute', fatal)
521  END IF
522  END IF
523  DO i=1, axes(id)%num_attributes
524  ! Unallocate all att arrays in preparation for new data
525  IF ( allocated(attributes(i)%fatt) ) THEN
526  DEALLOCATE(attributes(i)%fatt)
527  END IF
528  IF ( allocated(attributes(i)%iatt) ) THEN
529  DEALLOCATE(attributes(i)%iatt)
530  END IF
531 
532  ! Copy in attribute data
533  attributes(i)%type = axes(id)%attributes(i)%type
534  attributes(i)%len = axes(id)%attributes(i)%len
535  attributes(i)%name = axes(id)%attributes(i)%name
536  attributes(i)%catt = axes(id)%attributes(i)%catt
537  ! Allocate fatt arrays (if needed), and copy in data
538  IF ( allocated(axes(id)%attributes(i)%fatt) ) THEN
539  ALLOCATE(attributes(i)%fatt(SIZE(axes(id)%attributes(i)%fatt(:))), stat=istat)
540  IF ( istat .NE. 0 ) THEN
541  CALL error_mesg('diag_axis_mod::get_diag_axis', &
542  & 'Unable to allocate memory for attribute%fatt', fatal)
543  END IF
544  DO j=1, SIZE(attributes(i)%fatt(:))
545  attributes(i)%fatt(j) = axes(id)%attributes(i)%fatt(j)
546  END DO
547  END IF
548  ! Allocate iatt arrays (if needed), and copy in data
549  IF ( allocated(axes(id)%attributes(i)%iatt) ) THEN
550  ALLOCATE(attributes(i)%iatt(SIZE(axes(id)%attributes(i)%iatt(:))), stat=istat)
551  IF ( istat .NE. 0 ) THEN
552  CALL error_mesg('diag_axis_mod::get_diag_axis', &
553  & 'Unable to allocate memory for attribute%iatt', fatal)
554  END IF
555  DO j=1, SIZE(attributes(i)%iatt(:))
556  attributes(i)%iatt(j) = axes(id)%attributes(i)%iatt(j)
557  END DO
558  END IF
559  END DO
560  END IF
561  END IF
562  END SUBROUTINE get_diag_axis
563 
564  !> @brief Return the axis cartesian.
565  SUBROUTINE get_diag_axis_cart(id, cart_name)
566  INTEGER, INTENT(in) :: id !< Axis ID
567  CHARACTER(len=*), INTENT(out) :: cart_name !< Cartesian axis
568 
569  CALL valid_id_check(id, 'get_diag_axis_cart')
570  cart_name = axes(id)%cart_name
571  END SUBROUTINE get_diag_axis_cart
572 
573  !> @brief Return the axis data.
574  SUBROUTINE get_diag_axis_data(id, axis_data)
575  INTEGER, INTENT(in) :: id !< Axis ID
576  REAL, DIMENSION(:), INTENT(out) :: axis_data !< Axis data
577 
578  CALL valid_id_check(id, 'get_diag_axis_data')
579  IF (axes(id)%length > SIZE(axis_data(:))) THEN
580  ! <ERROR STATUS="FATAL">array data is too small</ERROR>
581  CALL error_mesg('diag_axis_mod::get_diag_axis_data', 'array data is too small', fatal)
582  ELSE
583  axis_data(1:axes(id)%length) = axes(id)%diag_type_data
584  END IF
585  END SUBROUTINE get_diag_axis_data
586 
587  !> @brief Return the short name of the axis.
588  SUBROUTINE get_diag_axis_name(id, axis_name)
589  INTEGER , INTENT(in) :: id !< Axis ID
590  CHARACTER(len=*), INTENT(out) :: axis_name !< Axis short name
591 
592  if (use_modern_diag) then
593  axis_name = fms_diag_object%fms_get_axis_name_from_id(id)
594  else
595  CALL valid_id_check(id, 'get_diag_axis_name')
596  axis_name = axes(id)%name
597  endif
598  END SUBROUTINE get_diag_axis_name
599 
600  !> @brief Return the name of the axis' domain
601  SUBROUTINE get_diag_axis_domain_name(id, name)
602  INTEGER, INTENT(in) :: id !< Axis ID
603  CHARACTER(len=*), INTENT(out) :: name !< Axis' domain name
604 
605  CALL valid_id_check(id, 'get_diag_axis_domain_name')
606  name = mpp_get_domain_name(axes(id)%domain2)
607  END SUBROUTINE get_diag_axis_domain_name
608 
609  !> @brief Return the length of the axis.
610  !> @return length of axis as an integer
611  INTEGER FUNCTION get_axis_length(id)
612  INTEGER, INTENT(in) :: id !< Axis ID
613  INTEGER :: length
614 
615  if (use_modern_diag) then
616  get_axis_length = fms_diag_object%fms_get_axis_length(id)
617  else
618  CALL valid_id_check(id, 'get_axis_length')
619  IF ( axes(id)%Domain .NE. null_domain1d ) THEN
620  CALL mpp_get_compute_domain(axes(id)%Domain,size=length)
621  !---one extra point is needed for some case. ( like symmetry domain )
622  get_axis_length = length + axes(id)%shift
623  ELSE
624  get_axis_length = axes(id)%length
625  END IF
626  endif
627  END FUNCTION get_axis_length
628 
629  !> @brief Return the auxiliary name for the axis.
630  !! @return auxiliary name for the axis
631  CHARACTER(len=128) FUNCTION get_axis_aux(id)
632  INTEGER, INTENT(in) :: id !< Axis ID
633 
634  CALL valid_id_check(id, 'get_axis_aux')
635  get_axis_aux = axes(id)%aux
636  END FUNCTION get_axis_aux
637 
638  !> @brief Return the required field names for the axis.
639  !! @return required field names for the axis
640  CHARACTER(len=128) FUNCTION get_axis_reqfld(id)
641  INTEGER, INTENT(in) :: id !< Axis ID
642 
643  CALL valid_id_check(id, 'get_axis_reqfld')
644  get_axis_reqfld = axes(id)%req
645  END FUNCTION get_axis_reqfld
646 
647  !> @brief Return the global length of the axis.
648  !! @return global length of the axis
649  INTEGER FUNCTION get_axis_global_length(id)
650  INTEGER, INTENT(in) :: id !< Axis ID
651 
652  CALL valid_id_check(id, 'get_axis_global_length')
653  get_axis_global_length = axes(id)%length
654  END FUNCTION get_axis_global_length
655 
656  !> @brief Return the tile count for the axis.
657  !! @return tile count for the axis
658  INTEGER FUNCTION get_tile_count(ids)
659  INTEGER, DIMENSION(:), INTENT(in) :: ids !< Axis IDs.
660  !! Possible dimensions: 1 <= <TT>size(ids(:))</TT> <= 4.
661 
662  INTEGER :: i, id, flag
663 
664  IF ( SIZE(ids(:)) < 1 ) THEN
665  ! <ERROR STATUS="FATAL">input argument has incorrect size.</ERROR>
666  CALL error_mesg('diag_axis_mod::get_tile_count', 'input argument has incorrect size', fatal)
667  END IF
668  get_tile_count = 1
669  flag = 0
670  DO i = 1, SIZE(ids(:))
671  id = ids(i)
672  CALL valid_id_check(id, 'get_tile_count')
673  IF ( axes(id)%cart_name == 'X' .OR. &
674  axes(id)%cart_name == 'Y' ) flag = flag + 1
675  ! --- both x/y axes found ---
676  IF ( flag == 2 ) THEN
677  get_tile_count = axes(id)%tile_count
678  EXIT
679  END IF
680  END DO
681  END FUNCTION get_tile_count
682 
683  !> @brief Retrun the 1D domain for the axis ID given.
684  !! @return 1D domain for the axis ID given
685  TYPE(domain1d) function get_domain1d(id)
686  INTEGER, INTENT(in) :: id !< Axis ID
687 
688  CALL valid_id_check(id, 'get_domain1d')
689  IF (axes(id)%Domain .NE. null_domain1d) THEN
690  get_domain1d = axes(id)%Domain
691  ELSE
692  get_domain1d = null_domain1d
693  ENDIF
694  END FUNCTION get_domain1d
695 
696  !> @brief Return the 2D domain for the axis IDs given.
697  !! @return 2D domain for the axis IDs given
698  TYPE(domain2d) function get_domain2d(ids)
699  INTEGER, DIMENSION(:), INTENT(in) :: ids !< Axis IDs.
700  !! Possible dimensions: 1 <= <TT>size(ids(:))</TT> <= 4.
701 
702  INTEGER :: i, id, flag
703 
704  IF ( SIZE(ids(:)) < 1 ) THEN
705  ! <ERROR STATUS="FATAL">input argument has incorrect size.</ERROR>
706  CALL error_mesg('diag_axis_mod::get_domain2d', 'input argument has incorrect size', fatal)
707  END IF
708 
709  if (use_modern_diag) then
710  get_domain2d = fms_diag_object%fms_get_domain2d(ids)
711  return
712  endif
713 
714  get_domain2d = null_domain2d
715  flag = 0
716  DO i = 1, SIZE(ids(:))
717  id = ids(i)
718  CALL valid_id_check(id, 'get_domain2d')
719  IF ( axes(id)%cart_name == 'X' .OR. axes(id)%cart_name == 'Y' ) flag = flag + 1
720  ! --- both x/y axes found ---
721  IF ( flag == 2 ) THEN
722  IF (axes(id)%Domain2 .NE. null_domain2d) get_domain2d = axes(id)%Domain2
723  EXIT
724  END IF
725  END DO
726  END FUNCTION get_domain2d
727 
728  !> @brief Retrun the 1D domain for the axis ID given.
729  !! @return 1D domain for the axis ID given
730  TYPE(domainug) function get_domainug(id)
731  INTEGER, INTENT(in) :: id !< Axis ID
732 
733  CALL valid_id_check(id, 'get_domainUG')
734  IF (axes(id)%DomainUG .NE. null_domainug) THEN
735  get_domainug = axes(id)%DomainUG
736  ELSE
737  get_domainug = null_domainug
738  ENDIF
739  END FUNCTION get_domainug
740 
741 !ug support
742  !> @brief Checks if the axes are compatible
743  !! @return integer domain_type
744  function axis_compatible_check(id,varname) result(domain_type)
745 
746  !Inputs/Outputs
747  integer,dimension(:),intent(in) :: id !<The array of axis IDs
748  character(*),intent(in),optional :: varname !<The name of the variable
749  integer(I4_KIND) :: domain_type !<DIAG_AXIS_NODOMAIN = no domain.
750  !<DIAG_AXIS_2DDOMAIN = structured domain.
751  !<DIAG_AXIS_UGDOMAIN = unstructured domain.
752 
753  !Local variables
754  logical :: xory !<XorY set to true if X or Y is found as a cart_name.
755  logical :: ug !<UG set to true if U is found as a cart_name.
756  integer :: n !<Looping index.
757  logical :: uses_domain2d !<True if an axis is associated with a 2D domain.
758  logical :: uses_domainug !<True if an axis is associated with an unstructured domain.
759 
760  !Initialize flags.
761  xory = .false.
762  ug = .false.
763  uses_domain2d = .false.
764  uses_domainug = .false.
765 
766  !Make sure that valid set of axes was passed, and determine the domain-type
767  !associated with the axes.
768  do n = 1,size(id)
769  call valid_id_check(id(n), &
770  "axis_compatible_check")
771  if (axes(id(n))%cart_name .eq. "X" .or. &
772  axes(id(n))%cart_name .eq. "Y") then
773  xory = .true.
774  elseif (axes(id(n))%cart_name .eq. "U") then
775  ug = .true.
776  endif
777  if (axes(id(n))%Domain2 .ne. null_domain2d) then
778  uses_domain2d = .true.
779  elseif (axes(id(n))%DomainUG .ne. null_domainug) then
780  uses_domainug = .true.
781  endif
782  enddo
783  if (ug .and. xory) then
784  if (present(varname)) then
785  call error_mesg("axis_compatible_check", &
786  "Can not use an unstructured grid with a "// &
787  "horizontal cartesian coordinate for the field " &
788  //trim(varname), &
789  fatal)
790  else
791  call error_mesg("axis_compatible_check", &
792  "Can not use an unstructured grid with a horizontal "// &
793  "cartesian coordinate", &
794  fatal)
795  endif
796  endif
797  if (uses_domain2d .and. uses_domainug) then
798  if (present(varname)) then
799  call error_mesg("axis_compatible_check", &
800  "Can not use an unstructured grid with a"// &
801  "structured grid for the field "//trim(varname), &
802  fatal)
803  else
804  call error_mesg("axis_compatible_check", &
805  "Can not use an unstructured grid with a"// &
806  "structured grid.", &
807  fatal)
808  endif
809  endif
810  if (uses_domain2d) then
811  domain_type = diag_axis_2ddomain
812  elseif (uses_domainug) then
813  domain_type = diag_axis_ugdomain
814  else
815  domain_type = diag_axis_nodomain
816  endif
817 
818  return
819  end function axis_compatible_check
820 
821  !> @brief Return the value of the shift for the axis IDs given.
822  SUBROUTINE get_axes_shift(ids, ishift, jshift)
823  INTEGER, DIMENSION(:), INTENT(in) :: ids
824  INTEGER, INTENT(out) :: ishift !< X shift value.
825  INTEGER, INTENT(out) :: jshift !< Y shift value.
826 
827  INTEGER :: i, id
828 
829  !-- get the value of the shift.
830  ishift = 0
831  jshift = 0
832  DO i = 1, SIZE(ids(:))
833  id = ids(i)
834  CALL valid_id_check(id, 'get_axes_shift')
835  SELECT CASE (axes(id)%cart_name)
836  CASE ( 'X' )
837  ishift = axes(id)%shift
838  CASE ( 'Y' )
839  jshift = axes(id)%shift
840  END SELECT
841  END DO
842  END SUBROUTINE get_axes_shift
843 
844  !> @brief Returns index into axis table corresponding to a given axis name.
845  !! @return Returns index into axis table corresponding to a given axis name.
846  INTEGER FUNCTION get_axis_num(axis_name, set_name)
847  CHARACTER(len=*), INTENT(in) :: axis_name !< Axis name
848  CHARACTER(len=*), INTENT(in), OPTIONAL :: set_name !< Set name
849 
850  INTEGER :: set, n
851 
852  IF ( PRESENT(set_name) ) THEN
853  set = get_axis_set_num(trim(set_name))
854  ELSE
855  set = 0
856  END IF
857  get_axis_num = 0
858  DO n = 1, num_def_axes
859  IF ( trim(axis_name) == trim(axes(n)%name) .AND. axes(n)%set == set ) THEN
860  get_axis_num = n
861  RETURN
862  END IF
863  END DO
864  END FUNCTION get_axis_num
865 
866  !> @brief Returns index in axis set table corresponding to a given axis set name
867  !! @return Returns index in axis set table corresponding to a given axis set name
868  INTEGER FUNCTION get_axis_set_num (set_name)
869  CHARACTER(len=*), INTENT(in) :: set_name !< Set name
870 
871  INTEGER :: iset
872 
873  get_axis_set_num = 0
874  DO iset = 1, num_axis_sets
875  IF ( set_name == axis_sets(iset) ) THEN
876  get_axis_set_num = iset
877  RETURN
878  END IF
879  END DO
880  END FUNCTION get_axis_set_num
881 
882  !> @brief Check to see if the given axis id is a valid id. If the axis id is invalid,
883  !! call a FATAL error. If the ID is valid, just return.
884  SUBROUTINE valid_id_check(id, routine_name)
885  INTEGER, INTENT(in) :: id !< Axis is to check for validity
886  CHARACTER(len=*), INTENT(in) :: routine_name !< Name of the subroutine checking for a valid axis id
887 
888  CHARACTER(len=5) :: emsg
889 
890  IF ( id < 1 .OR. id > num_def_axes) THEN
891  ! <ERROR STATUS="FATAL">
892  ! Illegal value for axis used (value <VALUE>).
893  ! </ERROR>
894  WRITE (emsg, '(I2)') id
895  CALL error_mesg('diag_axis_mod::'//trim(routine_name),&
896  & 'Illegal value for axis_id used (value '//trim(emsg)//').', fatal)
897  END IF
898  END SUBROUTINE valid_id_check
899 
900  SUBROUTINE diag_axis_attribute_init(diag_axis_id, name, type, cval, ival, rval)
901  INTEGER, INTENT(in) :: diag_axis_id !< input field ID, obtained from diag_axis_mod::diag_axis_init.
902  CHARACTER(len=*) :: name !< Name of the attribute
903  INTEGER, INTENT(in) :: type !< NetCDF type (NF_FLOAT, NF_INT, NF_CHAR)
904  CHARACTER(len=*), INTENT(in), OPTIONAL :: cval !< Character string attribute value
905  INTEGER, DIMENSION(:), INTENT(in), OPTIONAL :: ival !< Integer attribute value(s)
906  REAL, DIMENSION(:), INTENT(in), OPTIONAL :: rval !< Real attribute value(s)
907 
908  INTEGER :: istat, length, i, this_attribute
909  CHARACTER(len=1024) :: err_msg
910 
911  IF ( .NOT.first_send_data_call ) THEN
912  ! Call error due to unable to add attribute after send_data called
913  ! <ERROR STATUS="FATAL">
914  ! Attempting to add attribute <name> to axis <axis_name>
915  ! after first send_data call. Too late.
916  ! </ERROR>
917  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute', 'Attempting to add attribute "'&
918  &//trim(name)//'" to axis after first send_data call. Too late.', fatal)
919  END IF
920 
921  ! Simply return if diag_axis_id <= 0 --- not registered
922  IF ( diag_axis_id .LE. 0 ) THEN
923  RETURN
924  ELSE IF ( diag_axis_id .GT. num_def_axes ) THEN
925  ! Call error axis_id greater than num_def_axes. Axis is unknown
926  ! <ERROR STATUS="FATAL">
927  ! Attempting to add attribute <name> to axis ID <axis_ID>, however ID unknown.
928  ! </ERROR>
929  WRITE(err_msg, '(I5)') diag_axis_id
930  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute', 'Attempting to add attribute "'&
931  &//trim(name)//'" to axis ID "'//trim(err_msg)//'", however ID unknown.', fatal)
932 
933  ELSE
934  ! Allocate memory for the attributes
935  CALL attribute_init_axis(axes(diag_axis_id))
936 
937  ! Check if attribute already exists
938  this_attribute = 0
939  DO i=1, axes(diag_axis_id)%num_attributes
940  IF ( trim(axes(diag_axis_id)%attributes(i)%name) .EQ. trim(name) ) THEN
941  this_attribute = i
942  EXIT
943  END IF
944  END DO
945 
946  IF ( this_attribute.NE.0 .AND. (type.EQ.nf90_int .OR. type.EQ.nf90_float) ) THEN
947  ! <ERROR STATUS="FATAL">
948  ! Attribute <name> already defined for axis <axis_name>.
949  ! Contact the developers
950  ! </ERROR>
951  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute',&
952  & 'Attribute "'//trim(name)//'" already defined for axis "'&
953  &//trim(axes(diag_axis_id)%name)//'". Contact the developers.', fatal)
954  ELSE IF ( this_attribute.NE.0 .AND. type.EQ.nf90_char .AND. debug_diag_manager ) THEN
955  ! <ERROR STATUS="NOTE">
956  ! Attribute <name> already defined for axis <axis_name>.
957  ! Prepending.
958  ! </ERROR>
959  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute',&
960  & 'Attribute "'//trim(name)//'" already defined for axis"'&
961  &//trim(axes(diag_axis_id)%name)//'". Prepending.', note)
962  ELSE
963  ! Defining a new attribute
964  ! Increase the number of field attributes
965  this_attribute = axes(diag_axis_id)%num_attributes + 1
966  ! Checking to see if num_attributes == max_axis_attributes, and return error message
967  IF ( this_attribute .GT. max_axis_attributes ) THEN
968  ! <ERROR STATUS="FATAL">
969  ! Number of attributes exceeds max_axis_attributes for attribute <name> for axis <axis_name>.
970  ! Increase diag_manager_nml:max_axis_attributes.
971  ! </ERROR>
972  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute',&
973  & 'Number of attributes exceeds max_axis_attributes for attribute "'&
974  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)&
975  & //'". Increase diag_manager_nml:max_axis_attributes.',&
976  & fatal)
977  ELSE
978  axes(diag_axis_id)%num_attributes = this_attribute
979  ! Set name and type
980  axes(diag_axis_id)%attributes(this_attribute)%name = name
981  axes(diag_axis_id)%attributes(this_attribute)%type = type
982  ! Initialize catt to a blank string, as len_trim doesn't always work on an uninitialized string
983  axes(diag_axis_id)%attributes(this_attribute)%catt = ''
984  END IF
985  END IF
986 
987  SELECT CASE (type)
988  CASE (nf90_int)
989  IF ( .NOT.PRESENT(ival) ) THEN
990  ! <ERROR STATUS="FATAL">
991  ! Number type claims INTEGER, but ival not present for attribute <name> for axis <axis_name>.
992  ! Contact the developers.
993  ! </ERROR>
994  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute',&
995  & 'Attribute type claims INTEGER, but ival not present for attribute "'&
996  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)&
997  & //'". Contact the developers.', fatal)
998  END IF
999  length = SIZE(ival)
1000  ! Allocate iatt(:) to size of ival
1001  ALLOCATE(axes(diag_axis_id)%attributes(this_attribute)%iatt(length), stat=istat)
1002  IF ( istat.NE.0 ) THEN
1003  ! <ERROR STATUS="FATAL">
1004  ! Unable to allocate iatt for attribute <name> for axis <axis_name>
1005  ! </ERROR>
1006  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute', 'Unable to allocate iatt for attribute "'&
1007  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)//'"', fatal)
1008  END IF
1009  ! Set remaining fields
1010  axes(diag_axis_id)%attributes(this_attribute)%len = length
1011  axes(diag_axis_id)%attributes(this_attribute)%iatt = ival
1012  CASE (nf90_float)
1013  IF ( .NOT.PRESENT(rval) ) THEN
1014  ! <ERROR STATUS="FATAL">
1015  ! Attribute type claims REAL, but rval not present for attribute <name> for axis <axis_name>.
1016  ! Contact the developers.
1017  ! </ERROR>
1018  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute',&
1019  & 'Attribute type claims REAL, but rval not present for attribute "'&
1020  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)&
1021  & //'". Contact the developers.', fatal)
1022  END IF
1023  length = SIZE(rval)
1024  ! Allocate iatt(:) to size of rval
1025  ALLOCATE(axes(diag_axis_id)%attributes(this_attribute)%fatt(length), stat=istat)
1026  IF ( istat.NE.0 ) THEN
1027  ! <ERROR STATUS="FATAL">
1028  ! Unable to allocate fatt for attribute <name> for axis <axis_name>
1029  ! </ERROR>
1030  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute', 'Unable to allocate fatt for attribute "'&
1031  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)&
1032  & //'"', fatal)
1033  END IF
1034  ! Set remaining fields
1035  axes(diag_axis_id)%attributes(this_attribute)%len = length
1036  axes(diag_axis_id)%attributes(this_attribute)%fatt = rval
1037  CASE (nf90_char)
1038  IF ( .NOT.PRESENT(cval) ) THEN
1039  ! <ERROR STATUS="FATAL">
1040  ! Attribute type claims CHARACTER, but cval not present for attribute <name> for axis <axis_name>.
1041  ! Contact the developers.
1042  ! </ERROR>
1043  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute',&
1044  & 'Attribute type claims CHARACTER, but cval not present for attribute "'&
1045  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)&
1046  & //'". Contact the developers.', fatal)
1047  END IF
1048  CALL prepend_attribute_axis(axes(diag_axis_id), trim(name), trim(cval))
1049  CASE default
1050  ! <ERROR STATUS="FATAL">
1051  ! Unknown attribute type for attribute <name> for axis <axis_name>.
1052  ! Contact the developers.
1053  ! </ERROR>
1054  CALL error_mesg('diag_manager_mod::diag_axis_add_attribute', 'Unknown attribute type for attribute "'&
1055  & //trim(name)//'" for axis "'//trim(axes(diag_axis_id)%name)&
1056  & //'". Contact the developers.', fatal)
1057  END SELECT
1058  END IF
1059  END SUBROUTINE diag_axis_attribute_init
1060 
1061  SUBROUTINE diag_axis_add_attribute_scalar_r(diag_axis_id, att_name, att_value)
1062  INTEGER, INTENT(in) :: diag_axis_id
1063  CHARACTER(len=*), INTENT(in) :: att_name
1064  REAL, INTENT(in) :: att_value
1065 
1066  if (use_modern_diag) then
1067  call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, (/ att_value /))
1068  else
1069  CALL diag_axis_add_attribute_r1d(diag_axis_id, att_name, (/ att_value /))
1070  endif
1071  END SUBROUTINE diag_axis_add_attribute_scalar_r
1072 
1073  SUBROUTINE diag_axis_add_attribute_scalar_i(diag_axis_id, att_name, att_value)
1074  INTEGER, INTENT(in) :: diag_axis_id
1075  CHARACTER(len=*), INTENT(in) :: att_name
1076  INTEGER, INTENT(in) :: att_value
1077 
1078  if (use_modern_diag) then
1079  call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, (/ att_value /))
1080  else
1081  CALL diag_axis_add_attribute_i1d(diag_axis_id, att_name, (/ att_value /))
1082  endif
1083  END SUBROUTINE diag_axis_add_attribute_scalar_i
1084 
1085  SUBROUTINE diag_axis_add_attribute_scalar_c(diag_axis_id, att_name, att_value)
1086  INTEGER, INTENT(in) :: diag_axis_id
1087  CHARACTER(len=*), INTENT(in) :: att_name
1088  CHARACTER(len=*), INTENT(in) :: att_value
1089 
1090  if (use_modern_diag) then
1091  call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, (/ att_value /))
1092  else
1093  CALL diag_axis_attribute_init(diag_axis_id, att_name, nf90_char, cval=att_value)
1094  endif
1095  END SUBROUTINE diag_axis_add_attribute_scalar_c
1096 
1097  SUBROUTINE diag_axis_add_attribute_r1d(diag_axis_id, att_name, att_value)
1098  INTEGER, INTENT(in) :: diag_axis_id
1099  CHARACTER(len=*), INTENT(in) :: att_name
1100  REAL, DIMENSION(:), INTENT(in) :: att_value
1101 
1102  if (use_modern_diag) then
1103  call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, att_value)
1104  else
1105  CALL diag_axis_attribute_init(diag_axis_id, att_name, nf90_float, rval=att_value)
1106  endif
1107  END SUBROUTINE diag_axis_add_attribute_r1d
1108 
1109  SUBROUTINE diag_axis_add_attribute_i1d(diag_axis_id, att_name, att_value)
1110  INTEGER, INTENT(in) :: diag_axis_id
1111  CHARACTER(len=*), INTENT(in) :: att_name
1112  INTEGER, DIMENSION(:), INTENT(in) :: att_value
1113  if (use_modern_diag) then
1114  call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, att_value)
1115  else
1116  CALL diag_axis_attribute_init(diag_axis_id, att_name, nf90_int, ival=att_value)
1117  endif
1118  END SUBROUTINE diag_axis_add_attribute_i1d
1119 
1120  !> @brief Allocates memory in out_file for the attributes. Will <TT>FATAL</TT> if err_msg is not included
1121  !! in the subroutine call.
1122  SUBROUTINE attribute_init_axis(out_axis, err_msg)
1123  TYPE(diag_axis_type), INTENT(inout) :: out_axis !< output file to allocate memory for attribute
1124  CHARACTER(LEN=*), INTENT(out), OPTIONAL :: err_msg !< Error message, passed back to calling function
1125 
1126  INTEGER :: istat
1127 
1128  ! Initialize err_msg
1129  IF ( PRESENT(err_msg) ) err_msg = ''
1130 
1131  ! Allocate memory for the attributes
1132  IF ( .NOT.allocated(out_axis%attributes) ) THEN
1133  ALLOCATE(out_axis%attributes(max_axis_attributes), stat=istat)
1134  IF ( istat.NE.0 ) THEN
1135  ! <ERROR STATUS="FATAL">
1136  ! Unable to allocate memory for diag axis attributes
1137  ! </ERROR>
1138  IF ( fms_error_handler('diag_util_mod::attribute_init_axis', &
1139  & 'Unable to allocate memory for diag axis attributes', err_msg) ) THEN
1140  RETURN
1141  END IF
1142  ELSE
1143  ! Set equal to 0. It will be increased when attributes added
1144  out_axis%num_attributes = 0
1145  END IF
1146  END IF
1147  END SUBROUTINE attribute_init_axis
1148 
1149  !> @brief Prepends the attribute value to an already existing attribute. If the
1150  !! attribute isn't yet defined, then creates a new attribute
1151  SUBROUTINE prepend_attribute_axis(out_axis, att_name, prepend_value, err_msg)
1152  TYPE(diag_axis_type), INTENT(inout) :: out_axis !< diagnostic axis that will get the attribute
1153  CHARACTER(len=*), INTENT(in) :: att_name !< Name of the attribute
1154  CHARACTER(len=*), INTENT(in) :: prepend_value !< Value to prepend
1155  CHARACTER(len=*), INTENT(out) , OPTIONAL :: err_msg !< Error message, passed back to calling routine
1156 
1157  INTEGER :: length, i, this_attribute
1158  CHARACTER(len=512) :: err_msg_local
1159 
1160  ! Initialize string variables
1161  err_msg_local = ''
1162  IF ( PRESENT(err_msg) ) err_msg = ''
1163 
1164  ! Make sure the attributes for this out file have been initialized
1165  CALL attribute_init_axis(out_axis, err_msg_local)
1166  IF ( trim(err_msg_local) .NE. '' ) THEN
1167  IF ( fms_error_handler('diag_util_mod::prepend_attribute_axis', trim(err_msg_local), err_msg) ) THEN
1168  RETURN
1169  END IF
1170  END IF
1171 
1172  ! Find if attribute exists
1173  this_attribute = 0
1174  DO i=1, out_axis%num_attributes
1175  IF ( trim(out_axis%attributes(i)%name) .EQ. trim(att_name) ) THEN
1176  this_attribute = i
1177  EXIT
1178  END IF
1179  END DO
1180 
1181  IF ( this_attribute > 0 ) THEN
1182  IF ( out_axis%attributes(this_attribute)%type .NE. nf90_char ) THEN
1183  ! <ERROR STATUS="FATAL">
1184  ! Attribute <name> is not a character attribute.
1185  ! </ERROR>
1186  IF ( fms_error_handler('diag_util_mod::prepend_attribute_axis',&
1187  & 'Attribute "'//trim(att_name)//'" is not a character attribute.',&
1188  & err_msg) ) THEN
1189  RETURN
1190  END IF
1191  END IF
1192  ELSE
1193  ! Defining a new attribute
1194  ! Increase the number of file attributes
1195  this_attribute = out_axis%num_attributes + 1
1196  IF ( this_attribute .GT. max_axis_attributes ) THEN
1197  ! <ERROR STATUS="FATAL">
1198  ! Number of attributes exceeds max_axis_attributes for attribute <name>.
1199  ! Increase diag_manager_nml:max_axis_attributes.
1200  ! </ERROR>
1201  IF ( fms_error_handler('diag_util_mod::prepend_attribute_axis',&
1202  & 'Number of attributes exceeds max_axis_attributes for attribute "'&
1203  &//trim(att_name)//'". Increase diag_manager_nml:max_axis_attributes.',&
1204  & err_msg) ) THEN
1205  RETURN
1206  END IF
1207  ELSE
1208  out_axis%num_attributes = this_attribute
1209  ! Set name and type
1210  out_axis%attributes(this_attribute)%name = att_name
1211  out_axis%attributes(this_attribute)%type = nf90_char
1212  ! Initialize catt to a blank string, as len_trim doesn't always work on an uninitialized string
1213  out_axis%attributes(this_attribute)%catt = ''
1214  END IF
1215  END IF
1216 
1217  ! Only add string only if not already defined
1218  IF ( index(trim(out_axis%attributes(this_attribute)%catt), trim(prepend_value)).EQ.0 ) THEN
1219  ! Check if new string length goes beyond the length of catt
1220  length = len_trim(trim(prepend_value)//" "//trim(out_axis%attributes(this_attribute)%catt))
1221  IF ( length.GT.len(out_axis%attributes(this_attribute)%catt) ) THEN
1222  ! <ERROR STATUS="FATAL">
1223  ! Prepend length for attribute <name> is longer than allowed.
1224  ! </ERROR>
1225  IF ( fms_error_handler('diag_util_mod::prepend_attribute_file',&
1226  & 'Prepend length for attribute "'//trim(att_name)//'" is longer than allowed.',&
1227  & err_msg) ) THEN
1228  RETURN
1229  END IF
1230  END IF
1231  ! Set files
1232  out_axis%attributes(this_attribute)%catt =&
1233  & trim(prepend_value)//' '//trim(out_axis%attributes(this_attribute)%catt)
1234  out_axis%attributes(this_attribute)%len = length
1235  END IF
1236  END SUBROUTINE prepend_attribute_axis
1237 
1238  !> @brief given an axis, returns TRUE if the axis uses compression-by-gathering: that is, if
1239  !! this is an axis for fields on unstructured grid
1240  !! @return logical whether or not the axis uses compression-by-gathering
1241  logical function axis_is_compressed(id)
1242  integer, intent(in) :: id
1243 
1244  integer :: i
1245 
1246  CALL valid_id_check(id, 'axis_is_compressed')
1247 
1248  axis_is_compressed = .false.
1249  if (.not.allocated(axes(id)%attributes)) return
1250  do i = 1, axes(id)%num_attributes
1251  if (trim(axes(id)%attributes(i)%name)=='compress') then
1252  axis_is_compressed = .true.
1253  return
1254  endif
1255  enddo
1256  end function axis_is_compressed
1257 
1258 
1259  !> @brief given an index of compressed-by-gathering axis, return an array of axes used in
1260  !! compression. It is a fatal error to call it on axis that is not compressed
1261  subroutine get_compressed_axes_ids(id, r)
1262  integer, intent(in) :: id
1263  integer, intent(out), allocatable :: r(:)
1264 
1265  integer iatt, k, k1, k2, n
1266  logical :: space
1267 
1268  character(*), parameter :: tag = 'get_compressed_axes_ids'
1269 
1270  CALL valid_id_check(id, tag)
1271 
1272  associate(axis=>axes(id))
1273  if (.not.allocated(axis%attributes)) call error_mesg(tag, &
1274  'attempt to get compression dimensions from axis "'//trim(axis%name)// &
1275  & '" which is not compressed (does not have any attributes)', fatal)
1276 
1277  iatt = 0
1278  do k = 1,axis%num_attributes
1279  if (trim(axis%attributes(k)%name)=='compress') then
1280  iatt = k; exit ! from loop
1281  endif
1282  enddo
1283 
1284  if (iatt == 0) call error_mesg(tag, &
1285  'attempt to get compression dimensions from axis "'//trim(axis%name)//&
1286  '" which is not compressed (does not have "compress" attributes).', fatal)
1287  if (axis%attributes(iatt)%type/=nf90_char) call error_mesg(tag, &
1288  'attempt to get compression dimensions from axis "'//trim(axis%name)//&
1289  '" but the axis attribute "compress" has incorrect type.', fatal)
1290 
1291  ! parse the "compress" attribute
1292  ! calculate the number of compression axes
1293  space = .true.; n=0
1294  do k = 1, len(axis%attributes(iatt)%catt)
1295  if (space.and.(axis%attributes(iatt)%catt(k:k)/=' ')) then
1296  n = n+1
1297  endif
1298  space = (axis%attributes(iatt)%catt(k:k)==' ')
1299  enddo
1300 
1301  allocate(r(n))
1302  ! make array of compression axes indices. Go from the last to the first to get the
1303  ! array in FORTRAN order: they are listed in "compress" attribute C order (fastest
1304  ! dimension last)
1305  k2 = 0
1306  do k = n, 1, -1
1307  do k1 = k2+1, len(axis%attributes(iatt)%catt)
1308  if (axis%attributes(iatt)%catt(k1:k1)/=' ') exit
1309  enddo
1310  do k2 = k1+1, len(axis%attributes(iatt)%catt)
1311  if (axis%attributes(iatt)%catt(k2:k2)==' ') exit
1312  enddo
1313  r(k) = get_axis_num(axis%attributes(iatt)%catt(k1:k2),axis_sets(axis%set))
1314  if (r(k)<=0) call error_mesg(tag, &
1315  'compression dimension "'//trim(axis%attributes(iatt)%catt(k1:k2))//&
1316  '" not found among the axes of set "'//trim(axis_sets(axis%set))//'".', fatal)
1317  enddo
1318  end associate ! axis
1319  end subroutine get_compressed_axes_ids
1320 END MODULE diag_axis_mod
1321 !> @}
1322 ! close documentation grouping
subroutine, public get_compressed_axes_ids(id, r)
given an index of compressed-by-gathering axis, return an array of axes used in compression....
Definition: diag_axis.F90:1262
integer(i4_kind) function, public axis_compatible_check(id, varname)
Checks if the axes are compatible.
Definition: diag_axis.F90:745
integer function, public get_axis_length(id)
Return the length of the axis.
Definition: diag_axis.F90:612
subroutine, public get_diag_axis(id, name, units, long_name, cart_name, direction, edges, Domain, DomainU, array_data, num_attributes, attributes, domain_position)
Return information about the axis with index ID.
Definition: diag_axis.F90:468
subroutine, public get_diag_axis_cart(id, cart_name)
Return the axis cartesian.
Definition: diag_axis.F90:566
subroutine valid_id_check(id, routine_name)
Check to see if the given axis id is a valid id. If the axis id is invalid, call a FATAL error....
Definition: diag_axis.F90:885
integer function, public get_axis_num(axis_name, set_name)
Returns index into axis table corresponding to a given axis name.
Definition: diag_axis.F90:847
subroutine attribute_init_axis(out_axis, err_msg)
Allocates memory in out_file for the attributes. Will FATAL if err_msg is not included in the subrout...
Definition: diag_axis.F90:1123
integer function, public get_axis_global_length(id)
Return the global length of the axis.
Definition: diag_axis.F90:650
type(domain2d) function, public get_domain2d(ids)
Return the 2D domain for the axis IDs given.
Definition: diag_axis.F90:699
integer function get_axis_set_num(set_name)
Returns index in axis set table corresponding to a given axis set name.
Definition: diag_axis.F90:869
type(domainug) function, public get_domainug(id)
Retrun the 1D domain for the axis ID given.
Definition: diag_axis.F90:731
subroutine, public get_diag_axis_data(id, axis_data)
Return the axis data.
Definition: diag_axis.F90:575
character(len=128) function, public get_axis_reqfld(id)
Return the required field names for the axis.
Definition: diag_axis.F90:641
subroutine diag_axis_attribute_init(diag_axis_id, name, type, cval, ival, rval)
Definition: diag_axis.F90:901
integer function, public get_tile_count(ids)
Return the tile count for the axis.
Definition: diag_axis.F90:659
subroutine, public get_diag_axis_name(id, axis_name)
Return the short name of the axis.
Definition: diag_axis.F90:589
subroutine prepend_attribute_axis(out_axis, att_name, prepend_value, err_msg)
Prepends the attribute value to an already existing attribute. If the attribute isn't yet defined,...
Definition: diag_axis.F90:1152
character(len=128) function, public get_axis_aux(id)
Return the auxiliary name for the axis.
Definition: diag_axis.F90:632
integer(i4_kind), parameter, public diag_axis_2ddomain
For unstructured grid support.
Definition: diag_axis.F90:63
type(diag_axis_type), dimension(:), allocatable, save axes
global storage for all defined axes
Definition: diag_axis.F90:73
integer function, public diag_axis_init(name, array_data, units, cart_name, long_name, direction, set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count, domain_position)
Initialize the axis, and return the axis ID.
Definition: diag_axis.F90:111
logical function, public axis_is_compressed(id)
given an axis, returns TRUE if the axis uses compression-by-gathering: that is, if this is an axis fo...
Definition: diag_axis.F90:1242
subroutine, public get_axes_shift(ids, ishift, jshift)
Return the value of the shift for the axis IDs given.
Definition: diag_axis.F90:823
integer, dimension(:), allocatable num_subaxes
counter of number of axes defined
Definition: diag_axis.F90:67
integer(i4_kind), parameter, public diag_axis_ugdomain
For unstructured grid support.
Definition: diag_axis.F90:64
subroutine, public get_diag_axis_domain_name(id, name)
Return the name of the axis' domain.
Definition: diag_axis.F90:602
type(domain1d) function, public get_domain1d(id)
Retrun the 1D domain for the axis ID given.
Definition: diag_axis.F90:686
integer function, public diag_subaxes_init(axis, subdata, start_indx, end_indx, domain_2d)
Create a subaxis on a parent axis.
Definition: diag_axis.F90:383
character(len=128), dimension(:), allocatable, save axis_sets
storage for axis set names
Definition: diag_axis.F90:70
Add an arbitrary attribute and value to the diagnostic axis.
Definition: diag_axis.F90:89
integer max_axis_attributes
Maximum number of user definable attributes per axis.
Definition: diag_data.F90:385
logical use_modern_diag
Namelist flag to use the modernized diag_manager code.
Definition: diag_data.F90:391
integer max_axes
Maximum number of independent axes.
Definition: diag_data.F90:361
Attribute type for diagnostic fields.
Definition: diag_data.F90:157
Type to hold the diagnostic axis description.
Definition: diag_data.F90:306
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:758
logical function, public fms_error_handler(routine, message, err_msg)
Facilitates the control of fatal error conditions.
Definition: fms.F90:525
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:498
subroutine mpp_get_domain_components(domain, x, y, tile_count)
Retrieve 1D components of 2D decomposition.
character(len=name_length) function mpp_get_domain_name(domain)
Set user stack size.
These routines retrieve the axis specifications associated with the compute domains....
These routines retrieve the axis specifications associated with the global domains....
One dimensional domain used to manage shared data access between pes.
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.