FMS  2026.01
Flexible Modeling System
fms_diag_axis_object.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 
19 !> @defgroup fms_diag_axis_object_mod fms_diag_axis_object_mod
20 !> @ingroup diag_manager
21 !! @brief fms_diag_axis_object_mod stores the diag axis object, a diag domain
22 !! object, and a subaxis object.
23 
24 !> @file
25 !> @brief File for @ref diag_axis_object_mod
26 
27 !> @addtogroup fms_diag_axis_object_mod
28 !> @{
29 module fms_diag_axis_object_mod
30 #ifdef use_yaml
31  use mpp_domains_mod, only: domain1d, domain2d, domainug, mpp_get_compute_domain, center, &
32  & mpp_get_global_domain, north, east, mpp_get_tile_id, &
34  use platform_mod, only: r8_kind, r4_kind, i4_kind, i8_kind
35  use diag_data_mod, only: diag_atttype, max_axes, no_domain, two_d_domain, ug_domain, &
37  diag_null, index_gridtype, latlon_gridtype, pack_size_str, &
40  use mpp_mod, only: fatal, mpp_error, uppercase, mpp_pe, mpp_root_pe, stdout
41  use fms2_io_mod, only: fmsnetcdffile_t, fmsnetcdfdomainfile_t, fmsnetcdfunstructureddomainfile_t, &
42  & register_axis, register_field, register_variable_attribute, write_data
43  use fms_diag_yaml_mod, only: subregion_type, diag_yaml, max_subaxes, diagyamlfilesvar_type
44  use diag_grid_mod, only: get_local_indices_cubesphere => get_local_indexes
45  use axis_utils2_mod, only: nearest_index
46  implicit none
47 
48  PRIVATE
49 
56 
57  !> @}
58 
59  !> @brief Type to hold the domain info for an axis
60  !! This type was created to avoid having to send in "Domain", "Domain2", "DomainUG" as arguments into subroutines
61  !! and instead only 1 class(diagDomain_t) argument can be send
62  !> @ingroup diag_axis_object_mod
64  contains
65  procedure :: set => set_axis_domain
66  procedure :: length => get_length
67  procedure :: get_ntiles
68  end type diagdomain_t
69 
70  !> @brief Type to hold the 1d domain
71  type, extends(diagdomain_t) :: diagdomain1d_t
72  type(domain1d) :: domain !< 1d Domain of the axis
73  end type
74 
75  !> @brief Type to hold the 2d domain
76  type, extends(diagdomain_t) :: diagdomain2d_t
77  type(domain2d) :: domain2 !< 2d Domain of an "X" or "Y" axis
78  end type
79 
80  !> @brief Type to hold the unstructured domain
81  type, extends(diagdomain_t) :: diagdomainug_t
82  type(domainug) :: domainug !< Domain of "U" axis
83  end type
84 
85  !> @brief Type to hold the diagnostic axis description.
86  !> @ingroup diag_axis_object_mod
88  INTEGER , private :: axis_id !< ID of the axis
89 
90  contains
91  procedure :: get_parent_axis_id
92  procedure :: get_subaxes_id
93  procedure :: get_axis_name
94  procedure :: is_z_axis
95  procedure :: write_axis_metadata
96  procedure :: write_axis_data
97  procedure :: add_structured_axis_ids
98  procedure :: get_structured_axis
99  procedure :: is_unstructured_grid
100  procedure :: get_edges_id
101  END TYPE fmsdiagaxis_type
102 
103  !> @brief Type to hold the diag_axis (either subaxis or a full axis)
104  !> @ingroup diag_axis_object_mod
106  class(fmsDiagAxis_type), allocatable :: axis
107  end type
108 
109  !> @brief Type to hold the subaxis
110  !> @ingroup diag_axis_object_mod
112  CHARACTER(len=:), ALLOCATABLE , private :: subaxis_name !< Name of the subaxis
113  INTEGER , private :: starting_index !< Starting index of the subaxis relative to the
114  !! parent axis
115  INTEGER , private :: ending_index !< Ending index of the subaxis relative to the
116  !! parent axis
117  INTEGER , private :: parent_axis_id !< Id of the parent_axis
118  INTEGER , private :: compute_idx(2) !< Starting and ending index of the compute domain
119  INTEGER, allocatable, private :: global_idx(:) !< Starting and ending index of the global domain
120  real(kind=r4_kind), allocatable, private :: zbounds(:) !< Bounds of the Z axis
121  contains
122  procedure :: fill_subaxis
123  procedure :: axis_length
124  procedure :: get_starting_index
125  procedure :: get_ending_index
126  procedure :: get_compute_indices
127  procedure :: is_same_zbounds
128  END TYPE fmsdiagsubaxis_type
129 
130  !> @brief Type to hold the diurnal axis
131  !> @ingroup diag_axis_object_mod
133  INTEGER , private :: ndiurnal_samples !< The number of diurnal samples
134  CHARACTER(len=:), ALLOCATABLE, private :: axis_name !< The diurnal axis name
135  CHARACTER(len=:), ALLOCATABLE, private :: long_name !< The longname of the diurnal axis
136  CHARACTER(len=:), ALLOCATABLE, private :: units !< The units
137  INTEGER , private :: edges_id !< The id of the diurnal edges
138  CHARACTER(len=:), ALLOCATABLE, private :: edges_name !< The name of the edges axis
139  CLASS(*), ALLOCATABLE, private :: diurnal_data(:) !< The diurnal data
140 
141  contains
142  procedure :: get_diurnal_axis_samples
143  procedure :: write_diurnal_metadata
144  END TYPE fmsdiagdiurnalaxis_type
145 
146  !> @brief Type to hold the diagnostic axis description.
147  !> @ingroup diag_axis_object_mod
149  CHARACTER(len=:), ALLOCATABLE, private :: axis_name !< Name of the axis
150  CHARACTER(len=:), ALLOCATABLE, private :: units !< Units of the axis
151  CHARACTER(len=:), ALLOCATABLE, private :: long_name !< Long_name attribute of the axis
152  CHARACTER(len=1) , private :: cart_name !< Cartesian name "X", "Y", "Z", "T", "U", "N"
153  CLASS(*), ALLOCATABLE, private :: axis_data(:) !< Data of the axis
154  CHARACTER(len=:), ALLOCATABLE, private :: type_of_data !< The type of the axis_data ("float" or "double")
155  !< TO DO this can be a dlinked to avoid having limits
156  integer, ALLOCATABLE, private :: subaxis(:) !< Array of subaxis
157  integer , private :: nsubaxis !< Number of subaxis
158  class(diagdomain_t),ALLOCATABLE, private :: axis_domain !< Domain
159  INTEGER , private :: type_of_domain !< The type of domain ("NO_DOMAIN", "TWO_D_DOMAIN",
160  !! or "UG_DOMAIN")
161  INTEGER , private :: length !< Global axis length
162  INTEGER , private :: direction !< Direction of the axis 0, 1, -1
163  INTEGER, ALLOCATABLE, private :: edges_id !< Axis ID for the edges axis
164  !! This axis will be written to the file
165  CHARACTER(len=:), ALLOCATABLE, private :: edges_name !< Name for the previously defined "edges axis"
166  !! This will be written as an attribute
167  CHARACTER(len=:), ALLOCATABLE, private :: aux !< Auxiliary name, can only be <TT>geolon_t</TT>
168  !! or <TT>geolat_t</TT>
169  CHARACTER(len=128) , private :: req !< Required field names.
170  INTEGER , private :: tile_count !< The number of tiles
171  TYPE(fmsdiagattribute_type),allocatable , private :: attributes(:) !< Array to hold user definable attributes
172  INTEGER , private :: num_attributes !< Number of defined attibutes
173  INTEGER , private :: domain_position !< The position in the doman (NORTH, EAST or CENTER)
174  integer, allocatable , private :: structured_ids(:) !< If the axis is in the unstructured grid,
175  !! this is the axis ids of the structured axis
176  CHARACTER(len=:), ALLOCATABLE, private :: set_name !< Name of the axis set. This is to distinguish
177  !! two axis with the same name
178 
179  contains
180 
181  PROCEDURE :: add_axis_attribute
182  PROCEDURE :: register => register_diag_axis_obj
183  PROCEDURE :: axis_length => get_axis_length
184  PROCEDURE :: set_edges
185  PROCEDURE :: set_axis_id
186  PROCEDURE :: get_compute_domain
187  PROCEDURE :: get_indices
188  PROCEDURE :: get_global_io_domain
189  PROCEDURE :: get_aux
190  PROCEDURE :: has_aux
191  PROCEDURE :: get_set_name
192  PROCEDURE :: has_set_name
193  PROCEDURE :: is_x_or_y_axis
194  PROCEDURE :: get_dim_size_layout
195  ! TO DO:
196  ! Get/has/is subroutines as needed
197  END TYPE fmsdiagfullaxis_type
198 
199  !> @addtogroup fms_diag_yaml_mod
200  !> @{
201  contains
202 
203  !!!!!!!!!!!!!!!!! DIAG AXIS PROCEDURES !!!!!!!!!!!!!!!!!
204  !> @brief Initialize the axis
205  subroutine register_diag_axis_obj(this, axis_name, axis_data, units, cart_name, long_name, direction,&
206  & set_name, Domain, Domain2, DomainU, aux, req, tile_count, domain_position, axis_length )
207  class(fmsdiagfullaxis_type),INTENT(inout):: this !< Diag_axis obj
208  CHARACTER(len=*), INTENT(in) :: axis_name !< Name of the axis
209  class(*), INTENT(in) :: axis_data(:) !< Array of coordinate values
210  CHARACTER(len=*), INTENT(in) :: units !< Units for the axis
211  CHARACTER(len=1), INTENT(in) :: cart_name !< Cartesian axis ("X", "Y", "Z", "T", "U", "N")
212  CHARACTER(len=*), INTENT(in), OPTIONAL :: long_name !< Long name for the axis.
213  CHARACTER(len=*), INTENT(in), OPTIONAL :: set_name !< Name of the parent axis, if it is a subaxis
214  INTEGER, INTENT(in), OPTIONAL :: direction !< Indicates the direction of the axis
215  TYPE(domain1d), INTENT(in), OPTIONAL :: Domain !< 1D domain
216  TYPE(domain2d), INTENT(in), OPTIONAL :: Domain2 !< 2D domain
217  TYPE(domainug), INTENT(in), OPTIONAL :: DomainU !< Unstructured domain
218  CHARACTER(len=*), INTENT(in), OPTIONAL :: aux !< Auxiliary name, can only be <TT>geolon_t</TT>
219  !! or <TT>geolat_t</TT>
220  CHARACTER(len=*), INTENT(in), OPTIONAL :: req !< Required field names.
221  INTEGER, INTENT(in), OPTIONAL :: tile_count !< Number of tiles
222  INTEGER, INTENT(in), OPTIONAL :: domain_position !< Domain position, "NORTH" or "EAST"
223  integer, intent(in), optional :: axis_length !< The length of the axis size(axis_data(:))
224 
225  this%axis_name = trim(axis_name)
226  this%units = trim(units)
227  this%cart_name = uppercase(cart_name)
228  call check_if_valid_cart_name(this%cart_name)
229 
230  if (present(long_name)) this%long_name = trim(long_name)
231 
232  select type (axis_data)
233  type is (real(kind=r8_kind))
234  allocate(real(kind=r8_kind) :: this%axis_data(axis_length))
235  this%axis_data = axis_data
236  this%length = axis_length
237  this%type_of_data = "double" !< This is what fms2_io expects in the register_field call
238  type is (real(kind=r4_kind))
239  allocate(real(kind=r4_kind) :: this%axis_data(axis_length))
240  this%axis_data = axis_data
241  this%length = axis_length
242  this%type_of_data = "float" !< This is what fms2_io expects in the register_field call
243  class default
244  call mpp_error(fatal, "The axis_data in your diag_axis_init call is not a supported type. &
245  & Currently only r4 and r8 data is supported.")
246  end select
247 
248  this%type_of_domain = no_domain
249  if (present(domain)) then
250  if (present(domain2) .or. present(domainu)) call mpp_error(fatal, &
251  "The presence of Domain with any other domain type is prohibited. "//&
252  "Check you diag_axis_init call for axis_name:"//trim(axis_name))
253  allocate(diagdomain1d_t :: this%axis_domain)
254  call this%axis_domain%set(domain=domain)
255  else if (present(domain2)) then
256  if (present(domainu)) call mpp_error(fatal, &
257  "The presence of Domain2 with any other domain type is prohibited. "//&
258  "Check you diag_axis_init call for axis_name:"//trim(axis_name))
259  allocate(diagdomain2d_t :: this%axis_domain)
260  call this%axis_domain%set(domain2=domain2)
261  this%type_of_domain = two_d_domain
262  else if (present(domainu)) then
263  allocate(diagdomainug_t :: this%axis_domain)
264  call this%axis_domain%set(domainu=domainu)
265  this%type_of_domain = ug_domain
266  endif
267 
268  this%tile_count = 1
269  if (present(tile_count)) this%tile_count = tile_count
270 
271  this%domain_position = center
272  if (present(domain_position)) this%domain_position = domain_position
273  call check_if_valid_domain_position(this%domain_position)
274 
275  this%direction = 0
276  if (present(direction)) this%direction = direction
277  call check_if_valid_direction(this%direction)
278 
279  if (present(aux)) this%aux = trim(aux)
280  if (present(req)) this%req = trim(req)
281  this%set_name = ""
282  if (present(set_name)) this%set_name = trim(set_name)
283 
284  if (max_subaxes .gt. 0) then
285  allocate(this%subaxis(max_subaxes))
286  this%subaxis = diag_null
287  endif
288 
289  this%nsubaxis = 0
290  this%num_attributes = 0
291  end subroutine register_diag_axis_obj
292 
293  !> @brief Add an attribute to an axis
294  subroutine add_axis_attribute(this, att_name, att_value)
295  class(fmsdiagfullaxis_type),INTENT(INOUT) :: this !< diag_axis obj
296  character(len=*), intent(in) :: att_name !< Name of the attribute
297  class(*), intent(in) :: att_value(:) !< The attribute value to add
298 
299  integer :: j !< obj%num_attributes (for less typing)
300 
301  if (.not. allocated(this%attributes)) &
302  allocate(this%attributes(max_axis_attributes))
303 
304  this%num_attributes = this%num_attributes + 1
305 
306  j = this%num_attributes
307  call this%attributes(j)%add(att_name, att_value)
308  end subroutine add_axis_attribute
309 
310  !> @brief Write the axis meta data to an open fileobj
311  subroutine write_axis_metadata(this, fms2io_fileobj, edges_in_file, parent_axis)
312  class(fmsdiagaxis_type), target, INTENT(IN) :: this !< diag_axis obj
313  class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write the data to
314  logical, INTENT(IN) :: edges_in_file !< .True. if the edges to this axis are
315  !! already in the file
316  class(fmsdiagaxis_type), OPTIONAL, target, INTENT(IN) :: parent_axis !< If the axis is a subaxis, axis object
317  !! for the parent axis (this will be used
318  !! to get some of the metadata info)
319 
320  character(len=:), ALLOCATABLE :: axis_edges_name !< Name of the edges, if it exist
321  character(len=:), pointer :: axis_name !< Name of the axis
322  integer :: axis_length !< Size of the axis
323  integer :: i !< For do loops
324  type(fmsdiagfullaxis_type), pointer :: diag_axis !< Local pointer to the diag_axis
325 
326  integer :: type_of_domain !< The type of domain the current axis is in
327  logical :: is_subaxis !< .true. if the axis is a subaxis
328  logical :: needs_domain_decomposition !< .True. if the axis needs the domain decomposition attribute
329  !! (i.e for "X" and "Y" subaxis)
330  integer :: domain_decomposition(4) !< indices of the global (1:2) and compute (3:4) domain for a "X" and "Y" subaxis
331 
332  is_subaxis = .false.
333  needs_domain_decomposition = .false.
334 
335  select type(this)
336  type is (fmsdiagfullaxis_type)
337  axis_name => this%axis_name
338  axis_length = this%length
339  diag_axis => this
340  type_of_domain = this%type_of_domain
341  type is (fmsdiagsubaxis_type)
342  is_subaxis = .true.
343  axis_name => this%subaxis_name
344  axis_length = this%ending_index - this%starting_index + 1
345  if (allocated(this%global_idx)) then
346  needs_domain_decomposition = .true.
347  domain_decomposition(1:2) = this%global_idx
348  domain_decomposition(3) = this%starting_index
349  domain_decomposition(4) = this%ending_index
350  endif
351  !< Get all the other information from the parent axis (i.e the cart_name, units, etc)
352  if (present(parent_axis)) then
353  select type(parent_axis)
354  type is (fmsdiagfullaxis_type)
355  diag_axis => parent_axis
356  end select
357  endif
358  type_of_domain = no_domain !< All subaxes are treated as non-domain decomposed (each rank writes it own file)
359  type is (fmsdiagdiurnalaxis_type)
360  call this%write_diurnal_metadata(fms2io_fileobj)
361  return
362  end select
363 
364  !< Add the axis as a dimension in the netcdf file based on the type of axis_domain and the fileobj type
365  select type (fms2io_fileobj)
366  !< The register_field calls need to be inside the select type block so that it can go inside the correct
367  !! register_field interface
368  type is (fmsnetcdffile_t)
369  !< Here the axis is not domain decomposed (i.e z_axis)
370  call register_axis(fms2io_fileobj, axis_name, axis_length)
371  call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
372  if (needs_domain_decomposition) then
373  call register_variable_attribute(fms2io_fileobj, axis_name, "domain_decomposition", &
374  domain_decomposition)
375  endif
376  type is (fmsnetcdfdomainfile_t)
377  select case (type_of_domain)
378  case (no_domain)
379  !< Here the fms2io_fileobj is domain decomposed, but the axis is not
380  !! Domain decomposed fileobjs can have axis that are not domain decomposed (i.e "Z" axis)
381  call register_axis(fms2io_fileobj, axis_name, axis_length)
382  call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
383  case (two_d_domain)
384  !< Here the axis is domain decomposed
385  call register_axis(fms2io_fileobj, axis_name, diag_axis%cart_name, domain_position=diag_axis%domain_position)
386  call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
387  end select
388  type is (fmsnetcdfunstructureddomainfile_t)
389  select case (type_of_domain)
390  case (ug_domain)
391  !< Here the axis is in a unstructured domain
392  call register_axis(fms2io_fileobj, axis_name)
393  call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
394  case default
395  !< Here the fms2io_fileobj is in the unstructured domain, but the axis is not
396  !< Unstructured domain fileobjs can have axis that are not domain decomposed (i.e "Z" axis)
397  call register_axis(fms2io_fileobj, axis_name, axis_length)
398  call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
399  end select
400  end select
401 
402  !< Write its metadata
403  if(allocated(diag_axis%long_name)) &
404  call register_variable_attribute(fms2io_fileobj, axis_name, "long_name", diag_axis%long_name, &
405  str_len=len_trim(diag_axis%long_name))
406 
407  if (diag_axis%cart_name .NE. "N") &
408  call register_variable_attribute(fms2io_fileobj, axis_name, "axis", diag_axis%cart_name, str_len=1)
409 
410  if (trim(diag_axis%units) .NE. "none") &
411  call register_variable_attribute(fms2io_fileobj, axis_name, "units", diag_axis%units, &
412  str_len=len_trim(diag_axis%units))
413 
414  select case (diag_axis%direction)
415  case (direction_up)
416  call register_variable_attribute(fms2io_fileobj, axis_name, "positive", "up", str_len=2)
417  case (direction_down)
418  call register_variable_attribute(fms2io_fileobj, axis_name, "positive", "down", str_len=4)
419  end select
420 
421  !< Ignore the edges attribute, if the edges are already in the file or if it is subaxis
422  if (.not. edges_in_file .and. allocated(diag_axis%edges_name) .and. .not. is_subaxis) then
423  call register_variable_attribute(fms2io_fileobj, axis_name, "edges", diag_axis%edges_name, &
424  str_len=len_trim(diag_axis%edges_name))
425  endif
426 
427  if(allocated(diag_axis%attributes)) then
428  do i = 1, diag_axis%num_attributes
429  select type (att_value => diag_axis%attributes(i)%att_value)
430  type is (character(len=*))
431  call register_variable_attribute(fms2io_fileobj, axis_name, diag_axis%attributes(i)%att_name, &
432  trim(att_value(1)), str_len=len_trim(att_value(1)))
433  class default
434  call register_variable_attribute(fms2io_fileobj, axis_name, diag_axis%attributes(i)%att_name, att_value)
435  end select
436  enddo
437  endif
438 
439  end subroutine write_axis_metadata
440 
441  !> @brief Write the axis data to an open fms2io_fileobj
442  subroutine write_axis_data(this, fms2io_fileobj, parent_axis)
443  class(fmsdiagaxis_type), target, INTENT(IN) :: this !< diag_axis obj
444  class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write the data to
445  class(fmsdiagaxis_type), OPTIONAL, target, INTENT(IN) :: parent_axis !< The parent axis if this is a subaxis
446 
447  integer :: i !< Starting index of a sub_axis
448  integer :: j !< Ending index of a sub_axis
449  integer :: global_io_index(2)!< Global io domain starting and ending index
450  select type(this)
451  type is (fmsdiagfullaxis_type)
452  call this%get_global_io_domain(global_io_index, fms2io_fileobj%is_file_using_netcdf_mpi())
453  call write_data(fms2io_fileobj, this%axis_name, this%axis_data(global_io_index(1):global_io_index(2)))
454  type is (fmsdiagsubaxis_type)
455  i = this%starting_index
456  j = this%ending_index
457 
458  if (present(parent_axis)) then
459  select type(parent_axis)
460  type is (fmsdiagfullaxis_type)
461  call write_data(fms2io_fileobj, this%subaxis_name, parent_axis%axis_data(i:j))
462  end select
463  endif
464  type is (fmsdiagdiurnalaxis_type)
465  call write_data(fms2io_fileobj, this%axis_name, this%diurnal_data)
466  end select
467  end subroutine write_axis_data
468 
469 
470  !> @brief Defined a new diurnal axis
471  subroutine define_diurnal_axis(diag_axis, naxis, n_diurnal_samples, is_edges)
472  class(fmsdiagaxiscontainer_type), target, intent(inout) :: diag_axis(:) !< Array of axis containers
473  integer, intent(inout) :: naxis !< Number of axis that have
474  !! been defined
475  integer, intent(in) :: n_diurnal_samples !< The number of diurnal samples
476  !! for the curent axis
477  logical, intent(in) :: is_edges !< Flag indicating if this is
478  !! an edge axis
479 
480  CHARACTER(32) :: axis_name !< name of the axis
481  CHARACTER(32) :: long_name !< long name of the axis
482  CHARACTER(32) :: edges_name !< name of the axis edge
483  CHARACTER(128) :: units !< units of the axis
484  real(kind=r8_kind), allocatable :: diurnal_data(:) !< Data for the axis
485  integer :: edges_id !< Id of the axis edge
486  integer :: i !< For do loops
487 
488  naxis = naxis + 1
489 
490  axis_name = ''
491  edges_name = ''
492  if (is_edges) then
493  WRITE (axis_name,'(a,i2.2)') 'time_of_day_edges_', n_diurnal_samples
494  long_name = "time of day edges"
495  allocate(diurnal_data(n_diurnal_samples + 1))
496  diurnal_data(1) = 0.0
497  edges_id = diag_null
498  do i = 1, n_diurnal_samples
499  diurnal_data(i+1) = 24.0* real(i)/n_diurnal_samples
500  enddo
501  else
502  WRITE (axis_name,'(a,i2.2)') 'time_of_day_', n_diurnal_samples
503  long_name = "time of day"
504  allocate(diurnal_data(n_diurnal_samples))
505  edges_id = naxis -1 !< The diurnal edges is the last defined axis
506  do i = 1, n_diurnal_samples
507  diurnal_data(i) = 24.0*(real(i)-0.5)/n_diurnal_samples
508  enddo
509  WRITE (edges_name,'(a,i2.2)') 'time_of_day_edges_', n_diurnal_samples
510  endif
511 
512  WRITE (units,11) 'hours', get_base_year(), get_base_month(), &
514 11 FORMAT(a,' since ',i4.4,'-',i2.2,'-',i2.2,' ',i2.2,':',i2.2,':',i2.2)
515 
516  allocate(fmsdiagdiurnalaxis_type :: diag_axis(naxis)%axis)
517  select type (diurnal_axis => diag_axis(naxis)%axis)
518  type is (fmsdiagdiurnalaxis_type)
519  diurnal_axis%axis_id = naxis
520  diurnal_axis%ndiurnal_samples = n_diurnal_samples
521  diurnal_axis%axis_name = trim(axis_name)
522  diurnal_axis%long_name = trim(long_name)
523  diurnal_axis%units = trim(units)
524  diurnal_axis%diurnal_data = diurnal_data
525  diurnal_axis%edges_id = edges_id
526  if (is_edges) &
527  WRITE (edges_name,'(a,i2.2)') 'time_of_day_edges_', n_diurnal_samples
528  diurnal_axis%edges_name = trim(edges_name)
529  end select
530  end subroutine define_diurnal_axis
531 
532  !< @brief Determine if the axis is in the unstructured grid
533  !! @return .True. if the axis is in unstructured grid
534  pure logical function is_unstructured_grid(this)
535  class(fmsdiagaxis_type), target, INTENT(in) :: this !< diag_axis obj
536 
537  is_unstructured_grid = .false.
538  select type (this)
539  type is (fmsdiagfullaxis_type)
540  is_unstructured_grid = trim(this%cart_name) .eq. "U"
541  end select
542  end function is_unstructured_grid
543 
544  !< @brief Adds the structured axis ids to the axis object
545  subroutine add_structured_axis_ids(this, axis_ids)
546  class(fmsdiagaxis_type), target, INTENT(inout) :: this !< diag_axis obj
547  integer, intent(in) :: axis_ids(2) !< axis ids to add to the axis object
548 
549  select type (this)
550  type is (fmsdiagfullaxis_type)
551  allocate(this%structured_ids(2))
552  this%structured_ids = axis_ids
553  end select
554  end subroutine add_structured_axis_ids
555 
556  !< @brief Get the structured axis ids from the axis object
557  !! @return the structured axis ids
558  pure function get_structured_axis(this) &
559  result(rslt)
560  class(fmsdiagaxis_type), target, INTENT(in) :: this !< diag_axis obj
561  integer :: rslt(2)
562 
563  rslt = diag_null
564  select type (this)
565  type is (fmsdiagfullaxis_type)
566  rslt = this%structured_ids
567  end select
568  end function get_structured_axis
569 
570 
571  !< @brief Get the edges_id of an axis_object
572  !! @return The edges_id of an axis object
573  pure integer function get_edges_id(this)
574  class(fmsdiagaxis_type), INTENT(in) :: this !< diag_axis obj
575 
576  get_edges_id = diag_null
577  select type (this)
578  type is (fmsdiagfullaxis_type)
579  if (allocated(this%edges_id)) get_edges_id = this%edges_id
580  end select
581  end function
582 
583  !> @brief Get the starting and ending indices of the global io domain of the axis
584  subroutine get_global_io_domain(this, global_io_index, use_collective_writes)
585  class(fmsdiagfullaxis_type), target, intent(in) :: this !< diag_axis obj
586  integer, intent(out) :: global_io_index(2) !< Global io domain starting and ending index
587  logical, intent(in) :: use_collective_writes !< .True. if using collective writes
588 
589  type(domain2d), pointer :: io_domain !< pointer to the io domain
590 
591  global_io_index(1) = 1
592  global_io_index(2) = this%length
593 
594  if (allocated(this%axis_domain)) then
595  select type(domain => this%axis_domain)
596  type is (diagdomain2d_t)
597  if (use_collective_writes) then
598  io_domain => domain%domain2
599  else
600  io_domain => mpp_get_io_domain(domain%domain2)
601  endif
602 
603  if (this%cart_name .eq. "X") then
604  call mpp_get_global_domain(io_domain, xbegin=global_io_index(1), xend=global_io_index(2), &
605  position=this%domain_position)
606  elseif (this%cart_name .eq. "Y") then
607  call mpp_get_global_domain(io_domain, ybegin=global_io_index(1), yend=global_io_index(2), &
608  position=this%domain_position)
609  endif
610  end select
611  endif
612  end subroutine get_global_io_domain
613 
614  !> @brief Get the length of the axis
615  !> @return axis length
616  function get_axis_length(this) &
617  result(axis_length)
618  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
619  integer :: axis_length
620 
621  !< If the axis is domain decomposed axis_length will be set to the length for the current PE:
622  if (allocated(this%axis_domain)) then
623  axis_length = this%axis_domain%length(this%cart_name, this%domain_position, this%length)
624  else
625  axis_length = this%length
626  endif
627 
628  end function
629 
630 
631  !> @brief Determine if an axis object has an auxiliary name
632  !! @return .true. if an axis object has an auxiliary name
633  pure function has_aux(this) &
634  result(rslt)
635  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
636  logical :: rslt
637 
638  rslt = .false.
639  if (allocated(this%aux)) rslt = trim(this%aux) .ne. ""
640  end function has_aux
641 
642  !> @brief Determine if an axis object has a set_name
643  !! @return .true. if an axis object has a set_name
644  pure function has_set_name(this) &
645  result(rslt)
646  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
647  logical :: rslt
648 
649  rslt = .false.
650  if (allocated(this%set_name)) rslt = trim(this%set_name) .ne. ""
651  end function has_set_name
652 
653  !> @brief Determine if an axis object is an x or y axis
654  !! @return .true. if an axis object is an x or y axis, optionally return a flag indicating which it is
655  function is_x_or_y_axis(this, x_or_y) &
656  result(rslt)
657  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
658  integer, optional, intent(inout) :: x_or_y !< returns is_x_axis if it is a x axis
659  !! is_y_axis if it is a y axis
660  logical :: rslt
661 
662  select case (trim(this%cart_name))
663  case ("X")
664  if (present(x_or_y)) x_or_y = is_x_axis
665  rslt = .true.
666  case ("Y")
667  if (present(x_or_y)) x_or_y = is_y_axis
668  rslt = .true.
669  case default
670  rslt = .false.
671  if (present(x_or_y)) x_or_y = diag_null
672  end select
673  end function is_x_or_y_axis
674 
675  !< @brief Get the global size of the axis, and the layout
676  !! It is assumed that this function is only called on "X" and "Y" axes
677  !! using the `is_x_or_y_axis` function from above
678  subroutine get_dim_size_layout(this, dim_size, layout)
679  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
680  integer, intent(out) :: dim_size !< Size of the dimension
681  integer, intent(out) :: layout !< Layout of the dimension
682 
683  integer :: nx, ny
684  integer :: layout_xy(2)
685 
686  select type (domain => this%axis_domain)
687  type is (diagdomain2d_t)
688  call mpp_get_global_domain(domain%Domain2, xsize=nx, ysize=ny)
689  call mpp_get_layout(domain%Domain2, layout_xy)
690 
691  if (this%cart_name .eq. "X") then
692  dim_size = nx
693  layout = layout_xy(1)
694  else if (this%cart_name .eq. "Y") then
695  dim_size = ny
696  layout = layout_xy(2)
697  endif
698  end select
699  end subroutine get_dim_size_layout
700 
701  !> @brief Get the set name of an axis object
702  !! @return the set name of an axis object
703  pure function get_set_name(this) &
704  result(rslt)
705  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
706  character(len=:), allocatable :: rslt
707 
708  rslt = this%set_name
709  end function get_set_name
710 
711  !> @brief Get the auxiliary name of an axis object
712  !! @return the auxiliary name of an axis object
713  pure function get_aux(this) &
714  result(rslt)
715  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
716  character(len=:), allocatable :: rslt
717 
718  rslt = this%aux
719  end function get_aux
720 
721  !> @brief Set the axis_id
722  subroutine set_axis_id(this, axis_id)
723  class(fmsdiagfullaxis_type), intent(inout) :: this !< diag_axis obj
724  integer, intent(in) :: axis_id !< Axis_id
725 
726  this%axis_id = axis_id
727 
728  end subroutine set_axis_id
729 
730  !> @brief Set the name and ids of the edges
731  subroutine set_edges(this, edges_name, edges_id)
732  class(fmsdiagfullaxis_type), intent(inout) :: this !< diag_axis obj
733  CHARACTER(len=*), intent(in) :: edges_name !< Name of the edges
734  integer, intent(in) :: edges_id !< Axis id of the edges
735 
736  !< Saving the name and the id of the edges axis because it will make it easier to use
737  !! downstream (i.e you need the edges name to write the attribute to the current axis,
738  !! and you need the edges id to add to the diag file object so that you can write the edges
739  !! to the file)
740  this%edges_name = edges_name
741  this%edges_id = edges_id
742  end subroutine set_edges
743 
744  !> @brief Determine if the subRegion is in the current PE.
745  !! If it is, determine the starting and ending indices of the current PE that belong to the subRegion
746  subroutine get_indices(this, compute_idx, corners_indices, starting_index, ending_index, need_to_define_axis)
747  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
748  integer, intent(in) :: compute_idx(:) !< Current PE's compute domain
749  class(*), intent(in) :: corners_indices(:) !< The indices of the corners of the subRegion
750  integer, intent(out) :: starting_index !< Starting index of the subRegion
751  !! for the current PE
752  integer, intent(out) :: ending_index !< Ending index of the subRegion
753  !! for the current PE
754  logical, intent(out) :: need_to_define_axis !< .true. if it is needed to define
755  !! an axis
756 
757  integer :: subregion_start !< Starting index of the subRegion
758  integer :: subregion_end !< Ending index of the subRegion
759 
760  !< Get the rectangular coordinates of the subRegion
761  !! If the subRegion is not rectangular, the points outside of the subRegion will be masked
762  !! out later
763  select type (corners_indices)
764  type is (integer(kind=i4_kind))
765  subregion_start = minval(corners_indices)
766  subregion_end = maxval(corners_indices)
767  end select
768 
769  !< Initiliaze the output
770  need_to_define_axis = .false.
771  starting_index = diag_null
772  ending_index = diag_null
773 
774  !< If the compute domain of the current PE is outisde of the range of sub_axis, return
775  if (compute_idx(1) < subregion_start .and. compute_idx(2) < subregion_start) return
776  if (compute_idx(1) > subregion_end .and. compute_idx(2) > subregion_end) return
777 
778  need_to_define_axis = .true.
779  if (compute_idx(1) >= subregion_start .and. compute_idx(2) >= subregion_end) then
780  !< In this case all the point of the current PE are inside the range of the sub_axis
781  starting_index = compute_idx(1)
782  ending_index = subregion_end
783  else if (compute_idx(1) >= subregion_start .and. compute_idx(2) <= subregion_end) then
784  !< In this case all the points of the current PE are valid up to the end point
785  starting_index = compute_idx(1)
786  ending_index = compute_idx(2)
787  else if (compute_idx(1) <= subregion_start .and. compute_idx(2) <= subregion_end) then
788  !< In this case all the points of the current PE are valid starting with t subregion_start
789  starting_index = subregion_start
790  ending_index = compute_idx(2)
791  else if (compute_idx(1) <= subregion_start .and. compute_idx(2) >= subregion_end) then
792  !< In this case only the points in the current PE ar valid
793  starting_index = subregion_start
794  ending_index = subregion_end
795  endif
796 
797  if (this%domain_position .ne. center) then
798  if (subregion_end - subregion_start + 1 .eq. 1) then
799  !< If your subregion consitsts of just 1 one, only include 1 PE
800  if (ending_index .eq. compute_idx(2)) need_to_define_axis = .false.
801  else
802  if (ending_index - starting_index + 1 .eq. 1) then
803  !< If the PEs section is only 1, only include 1 PE
804  if (starting_index .eq. compute_idx(2) .or. ending_index .eq. compute_idx(1)) &
805  need_to_define_axis = .false.
806  endif
807  endif
808  endif
809 
810  end subroutine get_indices
811 
812  !< Get the compute domain of the axis
813  subroutine get_compute_domain(this, compute_idx, need_to_define_axis, tile_number)
814  class(fmsdiagfullaxis_type), intent(in) :: this !< diag_axis obj
815  integer, intent(inout) :: compute_idx(:) !< Compute domain of the axis
816  logical, intent(out) :: need_to_define_axis !< .true. if it needed to define the axis
817  integer, optional, intent(in) :: tile_number !< The tile number of the axis
818 
819  !< Initialize the output
820  need_to_define_axis = .false.
821  compute_idx = diag_null
822 
823  if (.not. allocated(this%axis_domain)) then
824  !< If the axis is not domain decomposed, use the whole axis as the compute domain
825  if (this%cart_name .eq. "X" .or. this%cart_name .eq. "Y") then
826  compute_idx(1) = 1
827  compute_idx(2) = size(this%axis_data)
828  need_to_define_axis = .true.
829  endif
830  return
831  endif
832 
833  select type(domain => this%axis_domain)
834  type is (diagdomain2d_t)
835  if (present(tile_number)) then
836  !< If the tile number is present and the current PE is not on the tile, then there is no need
837  !! to define the axis
838  if (any(mpp_get_tile_id(domain%Domain2) .ne. tile_number)) then
839  need_to_define_axis = .false.
840  return
841  endif
842  endif
843 
844  !< Get the compute domain for the current PE if it is an "X" or "Y" axis
845  select case (this%cart_name)
846  case ("X")
847  call mpp_get_compute_domain(domain%Domain2, xbegin=compute_idx(1), xend=compute_idx(2), &
848  & position=this%domain_position)
849  need_to_define_axis = .true.
850  case ("Y")
851  call mpp_get_compute_domain(domain%Domain2, ybegin=compute_idx(1), yend=compute_idx(2), &
852  & position=this%domain_position)
853  need_to_define_axis = .true.
854  end select
855  end select
856 
857  end subroutine get_compute_domain
858 
859  !!!!!!!!!!!!!!!!!! SUB AXIS PROCEDURES !!!!!!!!!!!!!!!!!
860  !> @brief Fills in the information needed to define a subaxis
861  subroutine fill_subaxis(this, starting_index, ending_index, axis_id, parent_id, parent_axis_name, compute_idx, &
862  global_idx, zbounds, nz_subaxis)
863  class(fmsdiagsubaxis_type) , INTENT(INOUT) :: this !< diag_sub_axis obj
864  integer , intent(in) :: starting_index !< Starting index of the subRegion for the PE
865  integer , intent(in) :: ending_index !< Ending index of the subRegion for the PE
866  integer , intent(in) :: axis_id !< Axis id to assign to the subaxis
867  integer , intent(in) :: parent_id !< The id of the parent axis the subaxis belongs to
868  character(len=*) , intent(in) :: parent_axis_name !< Name of the parent_axis
869  integer , intent(in) :: compute_idx(2) !< Starting and ending index of
870  !! the axis's compute domain
871  integer, optional, intent(in) :: global_idx(2) !< Starting and ending index of
872  !! the axis's compute domain
873  real(kind=r4_kind), optional, intent(in) :: zbounds(2) !< Bounds of the z-axis
874  integer, optional, intent(in) :: nz_subaxis !< The number of z subaxis that have been defined
875  !! in the file
876 
877  integer :: nsubaxis !< The subaxis number in the axis name subXX
878  character(len=2) :: nsubaxis_char !< nsubaxis converted to a string
879 
880  nsubaxis = 1
881  if (present(nz_subaxis)) nsubaxis = nz_subaxis
882 
883  this%axis_id = axis_id
884  this%starting_index = starting_index
885  this%ending_index = ending_index
886  this%parent_axis_id = parent_id
887  write(nsubaxis_char, '(i2.2)') nsubaxis
888  this%subaxis_name = trim(parent_axis_name)//"_sub"//nsubaxis_char
889  this%compute_idx = compute_idx
890 
891  if (present(zbounds)) then
892  ! This is needed to avoid duplicating z sub axis!
893  allocate(this%zbounds(2))
894  this%zbounds = zbounds
895  endif
896 
897  if (present(global_idx)) then
898  ! This is needed for the "domain_decomposition" attribute which is needed for the combiner
899  allocate(this%global_idx(2))
900  this%global_idx = global_idx
901  endif
902  end subroutine fill_subaxis
903 
904  !> @brief Get the axis length of a subaxis
905  !> @return the axis length
906  function axis_length(this) &
907  result(res)
908  class(fmsdiagsubaxis_type) , INTENT(IN) :: this !< diag_sub_axis obj
909  integer :: res
910 
911  res = this%ending_index - this%starting_index + 1
912  end function
913 
914  !> @brief Accesses its member starting_index
915  !! @return a copy of the starting_index
916  function get_starting_index(this) result(indx)
917  class(fmsdiagsubaxis_type), intent(in) :: this !< diag_sub_axis object
918  integer :: indx !< Result to return
919  indx = this%starting_index
920  end function get_starting_index
921 
922  !> @brief Accesses its member ending_index
923  !! @return a copy of the ending_index
924  function get_ending_index(this) result(indx)
925  class(fmsdiagsubaxis_type), intent(in) :: this !< diag_sub_axis object
926  integer :: indx !< Result to return
927  indx = this%ending_index
928  end function get_ending_index
929 
930  !> @brief Accesses its member compute_indices
931  !! @return a copy of the ending_index
932  function get_compute_indices(this) result(indx)
933  class(fmsdiagsubaxis_type), intent(in) :: this !< diag_sub_axis object
934  integer :: indx(2) !< Result to return
935  indx = this%compute_idx
936  end function get_compute_indices
937 
938  !> @brief Determines if the zbounds passed in are the same as those in the file
939  !! @return .True. if the zbounds are the same
940  function is_same_zbounds(this, zbounds) result(is_same)
941  class(fmsdiagsubaxis_type), intent(in) :: this !< diag_sub_axis object
942  real(kind=r4_kind), intent(in) :: zbounds(2) !< zbounds to compare with
943  logical :: is_same
944 
945  is_same = zbounds(1) .eq. this%zbounds(1) .and. zbounds(2) .eq. this%zbounds(2)
946  end function
947 
948  !> @brief Get the ntiles in a domain
949  !> @return the number of tiles in a domain
950  function get_ntiles(this) &
951  result(ntiles)
952  class(diagdomain_t), INTENT(IN) :: this !< diag_axis obj
953 
954  integer :: ntiles
955 
956  select type (this)
957  type is (diagdomain2d_t)
958  ntiles = mpp_get_ntile_count(this%domain2)
959  end select
960  end function get_ntiles
961 
962  !> @brief Get the length of a 2D domain
963  !> @return Length of the 2D domain
964  function get_length(this, cart_axis, domain_position, global_length) &
965  result(length)
966  class(diagdomain_t), INTENT(IN) :: this !< diag_axis obj
967  character(len=*), INTENT(IN) :: cart_axis !< cart_axis of the axis
968  integer, INTENT(IN) :: domain_position !< Domain position (CENTER, NORTH, EAST)
969  integer, INTENT(IN) :: global_length !< global_length of the axis
970 
971  integer :: length
972 
973  select type (this)
974  type is(diagdomain2d_t)
975  if (trim(cart_axis) == "X") call mpp_get_compute_domain(this%Domain2, xsize=length, position=domain_position)
976  if (trim(cart_axis) == "Y") call mpp_get_compute_domain(this%Domain2, ysize=length, position=domain_position)
977  class default
978  !< If domain is 1D or UG, just set it to the global length
979  length = global_length
980  end select
981  end function get_length
982 
983  !!!!!!!!!!!!!!!!! FMS_DOMAIN PROCEDURES !!!!!!!!!!!!!!!!!
984 
985  !> @brief Set the axis domain
986  subroutine set_axis_domain(this, Domain, Domain2, DomainU)
987  class(diagdomain_t) :: this !< fms_domain obj
988  TYPE(domain1d), INTENT(in), OPTIONAL :: Domain !< 1d domain
989  TYPE(domain2d), INTENT(in), OPTIONAL :: Domain2 !< 2d domain
990  TYPE(domainug), INTENT(in), OPTIONAL :: DomainU !< Unstructured domain
991 
992  select type(this)
993  type is (diagdomain1d_t)
994  this%Domain = domain
995  type is (diagdomain2d_t)
996  this%Domain2 = domain2
997  type is (diagdomainug_t)
998  this%DomainUG = domainu
999  end select
1000  end subroutine set_axis_domain
1001 
1002  !< @brief Allocates the array of axis/subaxis objects
1003  !! @return true if there the aray of axis/subaxis objects is allocated
1004  logical function fms_diag_axis_object_init(axis_array)
1005  class(fmsdiagaxiscontainer_type) , allocatable, intent(inout) :: axis_array(:) !< Array of diag_axis
1006 
1007  if (allocated(axis_array)) call mpp_error(fatal, "The diag_axis containers is already allocated")
1008  allocate(axis_array(max_axes))
1009  !axis_array%axis_id = DIAG_NULL
1010 
1011  fms_diag_axis_object_init = .true.
1012  end function fms_diag_axis_object_init
1013 
1014  !< @brief Deallocates the array of axis/subaxis objects
1015  !! @return false if the aray of axis/subaxis objects was allocated
1016  logical function fms_diag_axis_object_end(axis_array)
1017  class(fmsdiagaxiscontainer_type) , allocatable, intent(inout) :: axis_array(:) !< Array of diag_axis
1018 
1019  if (allocated(axis_array)) deallocate(axis_array)
1020  fms_diag_axis_object_end = .false.
1021 
1022  end function fms_diag_axis_object_end
1023 
1024  !< @brief Determine the axis name of an axis_object
1025  !! @return The name of the axis
1026  !! @note This function may be called from the field object (i.e. to determine the dimension names for io),
1027  !! The field object only contains the parent axis ids, because the subregion is defined in a per file basis,
1028  !! so the is_regional flag is needed so that the correct axis name can be used
1029  pure function get_axis_name(this, is_regional) &
1030  result(axis_name)
1031  class(fmsdiagaxis_type), intent(in) :: this !< Axis object
1032  logical, intent(in), optional :: is_regional !< Flag indicating if the axis is regional
1033 
1034  character(len=:), allocatable :: axis_name
1035 
1036  select type (this)
1037  type is (fmsdiagfullaxis_type)
1038  axis_name = this%axis_name
1039  if (present(is_regional)) then
1040  if (is_regional) then
1041  if (this%cart_name .eq. "X" .or. this%cart_name .eq. "Y") axis_name = axis_name//"_sub01"
1042  endif
1043  endif
1044  type is (fmsdiagsubaxis_type)
1045  axis_name = this%subaxis_name
1046  end select
1047  end function get_axis_name
1048 
1049  !< @brief Determine if the axis is a Z axis by looking at the cartesian name
1050  !! @return .True. if the axis is a Z axis
1051  pure logical function is_z_axis(this)
1052  class(fmsdiagaxis_type), intent(in) :: this !< Axis object
1053  is_z_axis = .false.
1054  select type (this)
1055  type is (fmsdiagfullaxis_type)
1056  if (this%cart_name .eq. "Z") is_z_axis = .true.
1057  end select
1058  end function
1059 
1060  !> @brief Check if a cart_name is valid and crashes if it isn't
1061  subroutine check_if_valid_cart_name(cart_name)
1062  character(len=*), intent(in) :: cart_name
1063 
1064  select case (cart_name)
1065  case ("X", "Y", "Z", "T", "U", "N")
1066  case default
1067  call mpp_error(fatal, "diag_axit_init: Invalid cart_name: "//cart_name//&
1068  "The acceptable values are X, Y, Z, T, U, N.")
1069  end select
1070  end subroutine check_if_valid_cart_name
1071 
1072  !> @brief Check if a domain_position is valid and crashes if it isn't
1073  subroutine check_if_valid_domain_position(domain_position)
1074  integer, INTENT(IN) :: domain_position
1075 
1076  select case (domain_position)
1077  case (center, north, east)
1078  case default
1079  call mpp_error(fatal, "diag_axit_init: Invalid domain_positon. &
1080  &The acceptable values are NORTH, EAST, CENTER")
1081  end select
1082  end subroutine check_if_valid_domain_position
1083 
1084  !> @brief Check if a direction is valid and crashes if it isn't
1085  subroutine check_if_valid_direction(direction)
1086  integer, INTENT(IN) :: direction
1087 
1088  select case(direction)
1089  case(-1, 0, 1)
1090  case default
1091  call mpp_error(fatal, "diag_axit_init: Invalid direction. &
1092  &The acceptable values are-1 0 1")
1093  end select
1094  end subroutine check_if_valid_direction
1095 
1096  !> @brief Loop through a variable's axis_id to determine and return the domain type and domain to use
1097  subroutine get_domain_and_domain_type(diag_axis, axis_id, domain_type, domain, var_name)
1098  class(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Array of diag_axis
1099  integer, INTENT(IN) :: axis_id(:) !< Array of axis ids
1100  integer, INTENT(OUT) :: domain_type !< fileobj_type to use
1101  CLASS(diagdomain_t), POINTER, INTENT(OUT) :: domain !< Domain
1102  character(len=*), INTENT(IN) :: var_name !< Name of the variable (for error messages)
1103 
1104  integer :: i !< For do loops
1105  integer :: j !< axis_id(i) (for less typing)
1106 
1107  domain_type = no_domain
1108  domain => null()
1109 
1110  do i = 1, size(axis_id)
1111  j = axis_id(i)
1112  select type (axis => diag_axis(j)%axis)
1113  type is (fmsdiagfullaxis_type)
1114  !< Check that all the axis are in the same domain
1115  if (domain_type .ne. axis%type_of_domain) then
1116  !< If they are different domains, one of them can be NO_DOMAIN
1117  !! i.e a variable can have axis that are domain decomposed (x,y) and an axis that isn't (z)
1118  if (domain_type .eq. no_domain .or. axis%type_of_domain .eq. no_domain ) then
1119  !< Update the domain_type and domain, if needed
1120  if ((axis%type_of_domain .eq. two_d_domain .and. size(axis_id) > 1) &
1121  & .or. axis%type_of_domain .eq. ug_domain) then
1122  domain_type = axis%type_of_domain
1123  domain => axis%axis_domain
1124  endif
1125  else
1126  call mpp_error(fatal, "The variable:"//trim(var_name)//" has axis that are not in the same domain")
1127  endif
1128  endif
1129  end select
1130  enddo
1131  end subroutine get_domain_and_domain_type
1132 
1133  !> @brief Fill in the subaxis object for a subRegion defined by index
1134  subroutine define_new_subaxis_index(parent_axis, subRegion, diag_axis, naxis, is_x_or_y, write_on_this_pe)
1135  class(fmsdiagaxiscontainer_type), target, intent(inout) :: diag_axis(:) !< Diag_axis object
1136  type(fmsdiagfullaxis_type), intent(inout) :: parent_axis !< axis object of the parent
1137  integer, intent(inout) :: naxis !< Number of axis registered
1138  type(subregion_type), intent(in) :: subregion !< SubRegion definition from the yaml
1139  integer, intent(in) :: is_x_or_y !< Flag indicating if it is
1140  !! a x or y axis
1141  logical, intent(out) :: write_on_this_pe !< .true. if the subregion
1142  !! is on this PE
1143  integer :: compute_idx(2) !< Indices of the compute domain
1144  integer :: global_idx(2) !< Indices of the "global" domain
1145  integer :: starting_index !< starting index of the subregion
1146  integer :: ending_index !< ending index of the subregion
1147 
1148  call parent_axis%get_compute_domain(compute_idx, write_on_this_pe, tile_number=subregion%tile)
1149  if (.not. write_on_this_pe) return
1150 
1151  !< Determine if the PE's compute domain is inside the subRegion
1152  !! If it is get the starting and ending indices for that PE
1153  call parent_axis%get_indices(compute_idx, subregion%corners(:,is_x_or_y), starting_index, ending_index, &
1154  write_on_this_pe)
1155 
1156  if (.not. write_on_this_pe) return
1157 
1158  select type(corners=> subregion%corners)
1159  type is (integer(kind=i4_kind))
1160  global_idx(1) = minval(corners(:,is_x_or_y))
1161  global_idx(2) = maxval(corners(:,is_x_or_y))
1162  end select
1163 
1164  !< If it made it to this point, the current PE is in the subRegion!
1165  call define_new_axis(diag_axis, parent_axis, naxis, parent_axis%axis_id, &
1166  starting_index, ending_index, compute_idx, global_idx)
1167 
1168  end subroutine define_new_subaxis_index
1169 
1170  !> @brief Fill in the subaxis object for a subRegion defined by lat lon
1171  subroutine define_new_subaxis_latlon(diag_axis, axis_ids, naxis, subRegion, is_cube_sphere, write_on_this_pe)
1172  class(fmsdiagaxiscontainer_type), target, intent(inout) :: diag_axis(:) !< Diag_axis object
1173  integer, INTENT(in) :: axis_ids(:) !< Array of axes_ids
1174  integer, intent(inout) :: naxis !< Number of axis registered
1175  type(subregion_type), intent(in) :: subregion !< SubRegion definition from the yaml
1176  logical, intent(in) :: is_cube_sphere !< .true. if this is a cubesphere
1177  logical, intent(out) :: write_on_this_pe !< .true. if the subregion
1178  !! is on this PE
1179 
1180  real :: lat(2) !< Starting and ending lattiude of the subRegion
1181  real :: lon(2) !< Starting and ending longitude or the subRegion
1182  integer :: lat_indices(2) !< Starting and ending latitude indices of the subRegion
1183  integer :: lon_indices(2) !< Starting and ending longitude indices of the subRegion
1184  integer :: compute_idx(2) !< Compute domain of the current axis
1185  integer :: starting_index(2) !< Starting index of the subRegion for the current PE for the "x" and "y"
1186  !! direction
1187  integer :: ending_index(2) !< Ending index of the subRegion for the current PE for the "x" and "y" direction
1188  logical :: need_to_define_axis(2) !< .true. if it is needed to define the subaxis for the "x" and "y" direction
1189  integer :: i !< For do loops
1190  integer :: parent_axis_ids(2) !< The axis id of the parent axis for the "x" and "y" direction
1191  logical :: is_x_y_axis !< .true. if the axis is x or y
1192  integer :: compute_idx_2(2, 2) !< Starting and ending indices of the compute domain for the "x" and "y" direction
1193  integer :: global_idx (2, 2) !< Starting and ending indices of the global domain for the "x" and "y" direction
1194 
1195  write_on_this_pe = .false.
1196  need_to_define_axis = .true.
1197  parent_axis_ids = diag_null
1198 
1199  !< Get the rectangular coordinates of the subRegion
1200  !! If the subRegion is not rectangular, the points outside of the subRegion will be masked
1201  !! out later
1202  select type (corners => subregion%corners)
1203  type is (real(kind=r4_kind))
1204  lon(1) = minval(corners(:,1))
1205  lon(2) = maxval(corners(:,1))
1206  lat(1) = minval(corners(:,2))
1207  lat(2) = maxval(corners(:,2))
1208  end select
1209 
1210  if_is_cube_sphere: if (is_cube_sphere) then
1211  !< Get the starting and ending indices of the subregion in the cubesphere relative to the global domain
1212  call get_local_indices_cubesphere(lat(1), lat(2), lon(1), lon(2),&
1213  & lon_indices(1), lon_indices(2), lat_indices(1), lat_indices(2))
1214  loop_over_axis_ids: do i = 1, size(axis_ids)
1215  select_axis_type: select type (parent_axis => diag_axis(axis_ids(i))%axis)
1216  type is (fmsdiagfullaxis_type)
1217  !< Get the PEs compute domain
1218  call parent_axis%get_compute_domain(compute_idx, is_x_y_axis)
1219 
1220  !< If this is not a "X" or "Y" axis go to the next axis
1221  if (.not. is_x_y_axis) cycle
1222 
1223  !< Determine if the PE's compute domain is inside the subRegion
1224  !! If it is get the starting and ending indices for that PE
1225  if (parent_axis%cart_name .eq. "X") then
1226  call parent_axis%get_indices(compute_idx, lon_indices, starting_index(1), ending_index(1), &
1227  need_to_define_axis(1))
1228  parent_axis_ids(1) = axis_ids(i)
1229  compute_idx_2(1,:) = compute_idx
1230  global_idx(1,:) = lon_indices
1231  else if (parent_axis%cart_name .eq. "Y") then
1232  call parent_axis%get_indices(compute_idx, lat_indices, starting_index(2), ending_index(2), &
1233  need_to_define_axis(2))
1234  parent_axis_ids(2) = axis_ids(i)
1235  compute_idx_2(2,:) = compute_idx
1236  global_idx(2,:) = lat_indices
1237  endif
1238  end select select_axis_type
1239  enddo loop_over_axis_ids
1240  else if_is_cube_sphere
1241  loop_over_axis_ids2: do i = 1, size(axis_ids)
1242  select type (parent_axis => diag_axis(axis_ids(i))%axis)
1243  type is (fmsdiagfullaxis_type)
1244  !< Get the PEs compute domain
1245  call parent_axis%get_compute_domain(compute_idx, is_x_y_axis)
1246 
1247  !< If this is not a "X" or "Y" axis go to the next axis
1248  if (.not. is_x_y_axis) cycle
1249 
1250  !< Get the starting and ending indices of the subregion relative to the global grid
1251  if (parent_axis%cart_name .eq. "X") then
1252  select type(adata=>parent_axis%axis_data)
1253  type is (real(kind=r8_kind))
1254  lon_indices(1) = nearest_index(real(lon(1), kind=r8_kind), adata)
1255  lon_indices(2) = nearest_index(real(lon(2), kind=r8_kind), adata)
1256  type is (real(kind=r4_kind))
1257  lon_indices(1) = nearest_index(real(lon(1), kind=r4_kind), adata)
1258  lon_indices(2) = nearest_index(real(lon(2), kind=r4_kind), adata)
1259  end select
1260  call parent_axis%get_indices(compute_idx, lon_indices, starting_index(1), ending_index(1), &
1261  need_to_define_axis(1))
1262  parent_axis_ids(1) = axis_ids(i)
1263  compute_idx_2(1,:) = compute_idx
1264  global_idx(1,:) = lon_indices
1265  else if (parent_axis%cart_name .eq. "Y") then
1266  select type(adata=>parent_axis%axis_data)
1267  type is (real(kind=r8_kind))
1268  lat_indices(1) = nearest_index(real(lat(1), kind=r8_kind), adata)
1269  lat_indices(2) = nearest_index(real(lat(2), kind=r8_kind), adata)
1270  type is (real(kind=r4_kind))
1271  lat_indices(1) = nearest_index(real(lat(1), kind=r4_kind), adata)
1272  lat_indices(2) = nearest_index(real(lat(2), kind=r4_kind), adata)
1273  end select
1274  call parent_axis%get_indices(compute_idx, lat_indices, starting_index(2), ending_index(2), &
1275  need_to_define_axis(2))
1276  parent_axis_ids(2) = axis_ids(i)
1277  compute_idx_2(2,:) = compute_idx
1278  global_idx(2,:) = lat_indices
1279  endif
1280  end select
1281  enddo loop_over_axis_ids2
1282  endif if_is_cube_sphere
1283 
1284  !< If the PE's compute is not inside the subRegion move to the next axis
1285  if (any(.not. need_to_define_axis )) return
1286 
1287  !< If it made it to this point, the current PE is in the subRegion!
1288  write_on_this_pe = .true.
1289 
1290  do i = 1, size(parent_axis_ids)
1291  if (parent_axis_ids(i) .eq. diag_null) cycle
1292  select type (parent_axis => diag_axis(parent_axis_ids(i))%axis)
1293  type is (fmsdiagfullaxis_type)
1294  call define_new_axis(diag_axis, parent_axis, naxis, parent_axis_ids(i), &
1295  starting_index(i), ending_index(i), compute_idx_2(i,:), global_idx(i,:))
1296  end select
1297  enddo
1298 
1299  end subroutine define_new_subaxis_latlon
1300 
1301  !> @brief Creates a new subaxis and fills it will all the information it needs
1302  subroutine define_new_axis(diag_axis, parent_axis, naxis, parent_id, &
1303  starting_index, ending_index, compute_idx, global_idx, new_axis_id, zbounds, &
1304  nz_subaxis)
1305 
1306  class(fmsdiagaxiscontainer_type), target, intent(inout) :: diag_axis(:) !< Diag_axis object
1307  class(fmsdiagfullaxis_type), intent(inout) :: parent_axis !< The parent axis
1308  integer, intent(inout) :: naxis !< The number of axis that
1309  !! have been defined
1310  integer, intent(in) :: parent_id !< Id of the parent axis
1311  integer, intent(in) :: starting_index !< PE's Starting index
1312  integer, intent(in) :: ending_index !< PE's Ending index
1313  integer, intent(in) :: compute_idx(2) !< Starting and ending index of
1314  !! the axis's compute domain
1315  integer, optional, intent(in) :: global_idx(2) !< Starting and ending index of
1316  !! the axis's global domain
1317  integer, optional, intent(out) :: new_axis_id !< Axis id of the axis this is creating
1318  real(kind=r4_kind), optional, intent(in) :: zbounds(2) !< Bounds of the Z axis
1319  integer, optional, intent(in) :: nz_subaxis !< The number of z subaxis that have
1320  !! been defined in the file
1321 
1322  naxis = naxis + 1 !< This is the axis id of the new axis!
1323 
1324  !< Add the axis_id of the new subaxis to the parent axis
1325  parent_axis%nsubaxis = parent_axis%nsubaxis + 1
1326  parent_axis%subaxis(parent_axis%nsubaxis) = naxis
1327 
1328  !< Allocate the new axis as a subaxis and fill it
1329  allocate(fmsdiagsubaxis_type :: diag_axis(naxis)%axis)
1330  diag_axis(naxis)%axis%axis_id = naxis
1331  if (present(new_axis_id)) new_axis_id = naxis
1332 
1333  select type (sub_axis => diag_axis(naxis)%axis)
1334  type is (fmsdiagsubaxis_type)
1335  call sub_axis%fill_subaxis(starting_index, ending_index, naxis, parent_id, &
1336  parent_axis%axis_name, compute_idx, global_idx=global_idx, zbounds=zbounds, nz_subaxis=nz_subaxis)
1337  end select
1338  end subroutine define_new_axis
1339 
1340  !< @brief Determine the parent_axis_id of a subaxis
1341  !! @return parent_axis_id if it is a subaxis and diag_null if is not a subaxis
1342  pure function get_parent_axis_id(this) &
1343  result(parent_axis_id)
1344 
1345  class(fmsdiagaxis_type), intent(in) :: this !< Axis Object
1346  integer :: parent_axis_id
1347 
1348  select type (this)
1349  type is (fmsdiagfullaxis_type)
1350  parent_axis_id = diag_null
1351  type is (fmsdiagsubaxis_type)
1352  parent_axis_id = this%parent_axis_id
1353  type is (fmsdiagdiurnalaxis_type)
1354  parent_axis_id = diag_null
1355  end select
1356 
1357  end function
1358 
1359  !< @brief Determine the most recent subaxis id in a diag_axis object
1360  !! @return the most recent subaxis id in a diag_axis object
1361  pure function get_subaxes_id(this) &
1362  result(sub_axis_id)
1363 
1364  class(fmsdiagaxis_type), intent(in) :: this !< Axis Object
1365  integer :: sub_axis_id
1366 
1367  sub_axis_id = this%axis_id
1368  select type (this)
1369  type is (fmsdiagfullaxis_type)
1370  if (this%cart_name .ne. "Z") sub_axis_id = this%subaxis(this%nsubaxis)
1371  end select
1372 
1373  end function
1374 
1375  !< @brief Parses the "compress" attribute to get the names of the two axis
1376  !! @return the names of the structured axis
1377  pure function parse_compress_att(compress_att) &
1378  result(axis_names)
1379  class(*), intent(in) :: compress_att(:) !< The compress attribute to parse
1380  character(len=120) :: axis_names(2)
1381 
1382  integer :: ios !< Errorcode after parsing the compress attribute
1383 
1384  select type (compress_att)
1385  type is (character(len=*))
1386  read(compress_att(1),*, iostat=ios) axis_names
1387  if (ios .ne. 0) axis_names = ""
1388  class default
1389  axis_names = ""
1390  end select
1391  end function parse_compress_att
1392 
1393  !< @brief Determine the axis id of a axis
1394  !! @return Axis id
1395  pure function get_axis_id_from_name(axis_name, diag_axis, naxis, set_name) &
1396  result(axis_id)
1397  class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Array of axis object
1398  character(len=*), intent(in) :: axis_name !< Name of the axis
1399  integer, intent(in) :: naxis !< Number of axis that have been registered
1400  character(len=*), intent(in) :: set_name !< Name of the axis set
1401  integer :: axis_id
1402 
1403  integer :: i !< For do loops
1404 
1405  axis_id = diag_null
1406  do i = 1, naxis
1407  select type(axis => diag_axis(i)%axis)
1408  type is (fmsdiagfullaxis_type)
1409  if (trim(axis%axis_name) .eq. trim(axis_name)) then
1410  if (trim(axis%set_name) .eq. trim(set_name)) then
1411  axis_id = i
1412  return
1413  endif
1414  endif
1415  end select
1416  enddo
1417 
1418  end function get_axis_id_from_name
1419 
1420  !< @brief Get the number of diurnal samples for a diurnal axis
1421  !! @return The number of diurnal samples
1422  pure function get_diurnal_axis_samples(this) &
1423  result(n_diurnal_samples)
1424 
1425  class(fmsdiagdiurnalaxis_type), intent(in) :: this !< Axis Object
1426  integer :: n_diurnal_samples
1427 
1428  n_diurnal_samples = this%ndiurnal_samples
1429  end function get_diurnal_axis_samples
1430 
1431  !< @brief Writes out the metadata for a diurnal axis
1432  subroutine write_diurnal_metadata(this, fms2io_fileobj)
1433  class(fmsdiagdiurnalaxis_type), intent(in) :: this !< Diurnal axis Object
1434  class(fmsnetcdffile_t), intent(inout) :: fms2io_fileobj !< Fms2_io fileobj to write the data to
1435 
1436  call register_axis(fms2io_fileobj, this%axis_name, size(this%diurnal_data))
1437  call register_field(fms2io_fileobj, this%axis_name, pack_size_str, (/trim(this%axis_name)/))
1438  call register_variable_attribute(fms2io_fileobj, this%axis_name, "units", &
1439  &trim(this%units), str_len=len_trim(this%units))
1440  call register_variable_attribute(fms2io_fileobj, this%axis_name, "long_name", &
1441  &trim(this%long_name), str_len=len_trim(this%long_name))
1442  if (this%edges_id .ne. diag_null) &
1443  call register_variable_attribute(fms2io_fileobj, this%axis_name, "edges", &
1444  &trim(this%edges_name), str_len=len_trim(this%edges_name))
1445  end subroutine write_diurnal_metadata
1446 
1447  !> @brief Creates a new z subaxis to use
1448  subroutine create_new_z_subaxis(zbounds, var_axis_ids, diag_axis, naxis, file_axis_id, nfile_axis, nz_subaxis, &
1449  error_mseg)
1450  real(kind=r4_kind), intent(in) :: zbounds(2) !< Bounds of the Z axis
1451  integer, intent(inout) :: var_axis_ids(:) !< The variable's axis_ids
1452  class(fmsdiagaxiscontainer_type), target, intent(inout) :: diag_axis(:) !< Array of diag_axis objects
1453  integer, intent(inout) :: naxis !< Number of axis that have been
1454  !! registered
1455  integer, intent(inout) :: file_axis_id(:) !< The file's axis_ids
1456  integer, intent(inout) :: nfile_axis !< Number of axis that have been
1457  !! defined in file
1458  integer, intent(inout) :: nz_subaxis !< The number of z subaxis currently
1459  !! defined in the file
1460  character(len=*), intent(inout) :: error_mseg !! Message to include in error message
1461  !! if there is an error
1462 
1463  class(*), pointer :: zaxis_data(:) !< The data of the full zaxis
1464  integer :: subaxis_indices(2) !< The starting and ending indices of the subaxis relative to the full
1465  !! axis
1466  integer :: i !< For do loops
1467  integer :: subaxis_id !< The id of the new z subaxis
1468  integer :: parent_axis_id !< Id of parent axis id
1469  integer :: zaxis_index !< Index of the z axis (i.e 3 if the variable is x,y,z)
1470  type(fmsdiagfullaxis_type), pointer :: parent_axis !< Pointer to the parent axis
1471 
1472  parent_axis_id = diag_null
1473  zaxis_index = diag_null
1474 
1475  !< Determine which axis is the z axis:
1476  do i = 1, size(var_axis_ids)
1477  select type (parent_axis => diag_axis(var_axis_ids(i))%axis)
1478  type is (fmsdiagfullaxis_type)
1479  if (parent_axis%cart_name .eq. "Z") then
1480  parent_axis_id = var_axis_ids(i)
1481  zaxis_index = i
1482  endif
1483  end select
1484  enddo
1485 
1486  if (parent_axis_id .eq. diag_null) then
1487  call mpp_error(fatal, "create_new_z_subaxis:: unable to find the zaxis for "//trim(error_mseg))
1488  endif
1489 
1490  !< Determine if the axis was already created
1491  do i = 1, nfile_axis
1492  select type (axis => diag_axis(file_axis_id(i))%axis)
1493  type is (fmsdiagsubaxis_type)
1494  if (axis%parent_axis_id .ne. parent_axis_id) cycle
1495  if (axis%zbounds(1) .eq. zbounds(1) .and. axis%zbounds(2) .eq. zbounds(2)) then
1496  var_axis_ids(zaxis_index) = file_axis_id(i)
1497  return
1498  endif
1499  end select
1500  enddo
1501 
1502  select type (axis => diag_axis(parent_axis_id)%axis)
1503  type is (fmsdiagfullaxis_type)
1504  zaxis_data => axis%axis_data
1505  parent_axis => axis
1506  end select
1507 
1508  select type(zaxis_data)
1509  type is (real(kind=r4_kind))
1510  !TODO need to include the conversion to "real" because nearest_index doesn't take r4s and r8s
1511  subaxis_indices(1) = nearest_index(real(zbounds(1)), real(zaxis_data))
1512  subaxis_indices(2) = nearest_index(real(zbounds(2)), real(zaxis_data))
1513  type is (real(kind=r8_kind))
1514  subaxis_indices(1) = nearest_index(real(zbounds(1)), real(zaxis_data))
1515  subaxis_indices(2) = nearest_index(real(zbounds(2)), real(zaxis_data))
1516  end select
1517 
1518  nz_subaxis = nz_subaxis + 1
1519  call define_new_axis(diag_axis, parent_axis, naxis, parent_axis%axis_id, &
1520  &subaxis_indices(1), subaxis_indices(2), (/lbound(zaxis_data,1), ubound(zaxis_data,1)/), &
1521  &new_axis_id=subaxis_id, zbounds=zbounds, nz_subaxis=nz_subaxis)
1522  var_axis_ids(zaxis_index) = subaxis_id
1523 
1524  end subroutine
1525 
1526  !> @brief Determine if the diag_axis(parent_axis_id) is the parent of diag_axis(axis_id)
1527  !! @return .True. if diag_axis(parent_axis_id) is the parent of diag_axis(axis_id)
1528  function is_parent_axis(axis_id, parent_axis_id, diag_axis) &
1529  result(rslt)
1530  integer, intent(in) :: axis_id !< Axis id to check
1531  integer, intent(in) :: parent_axis_id !< Axis id of the parent to check
1532  class(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Array of diag_axis objects
1533 
1534  logical :: rslt
1535 
1536  rslt = .false.
1537  select type(axis => diag_axis(axis_id)%axis)
1538  type is (fmsdiagsubaxis_type)
1539  if (axis%parent_axis_id .eq. parent_axis_id) rslt = .true.
1540  end select
1541  end function is_parent_axis
1542 
1543  !> @brief Determine the name of the z subaxis by matching the parent axis id and the zbounds
1544  !! in the diag table yaml
1545  subroutine find_z_sub_axis_name(dim_name, parent_axis_id, file_axis_id, field_yaml, diag_axis)
1546  character(len=*), intent(inout) :: dim_name !< Name of z subaxis
1547  integer, intent(in) :: parent_axis_id !< Axis id of the parent
1548  integer, intent(in) :: file_axis_id(:) !< Axis ids of the file
1549  type(diagyamlfilesvar_type), intent(in) :: field_yaml !< Field info from diag_table yaml
1550  class(fmsdiagaxiscontainer_type),intent(in) :: diag_axis(:) !< Array of axis objections
1551 
1552  integer :: id
1553  integer :: i
1554 
1555  do i = 1, size(file_axis_id)
1556  id = file_axis_id(i)
1557  select type (axis_ptr => diag_axis(id)%axis)
1558  type is (fmsdiagsubaxis_type)
1559  if (axis_ptr%parent_axis_id .eq. parent_axis_id) then
1560  if (axis_ptr%is_same_zbounds(field_yaml%get_var_zbounds())) then
1561  dim_name = axis_ptr%subaxis_name
1562  return
1563  endif
1564  endif
1565  end select
1566  enddo
1567  call mpp_error(fatal, "Unable to determine the z subaxis name for field "//&
1568  trim(field_yaml%get_var_varname())//" in file: "//&
1569  trim(field_yaml%get_var_fname()))
1570  end subroutine
1571 #endif
1572 end module fms_diag_axis_object_mod
1573 !> @}
1574 ! close documentation grouping
integer, parameter direction_down
The axis points down if positive.
Definition: diag_data.F90:106
integer function get_base_minute()
gets the module variable base_minute
Definition: diag_data.F90:551
integer function get_base_year()
gets the module variable base_year
Definition: diag_data.F90:519
integer function get_base_hour()
gets the module variable base_hour
Definition: diag_data.F90:543
integer, parameter no_domain
Use the FmsNetcdfFile_t fileobj.
Definition: diag_data.F90:100
integer max_axis_attributes
Maximum number of user definable attributes per axis.
Definition: diag_data.F90:385
character(len=6) pack_size_str
Pack size as a string to be used in fms2_io register call set to "double" or "float".
Definition: diag_data.F90:409
integer max_axes
Maximum number of independent axes.
Definition: diag_data.F90:361
integer, parameter is_x_axis
integer indicating that it is a x axis
Definition: diag_data.F90:130
integer, parameter is_y_axis
integer indicating that it is a y axis
Definition: diag_data.F90:131
integer function get_base_day()
gets the module variable base_day
Definition: diag_data.F90:535
integer, parameter ug_domain
Use the FmsNetcdfUnstructuredDomainFile_t fileobj.
Definition: diag_data.F90:102
integer, parameter direction_up
The axis points up if positive.
Definition: diag_data.F90:105
integer function get_base_month()
gets the module variable base_month
Definition: diag_data.F90:527
integer function get_base_second()
gets the module variable base_second
Definition: diag_data.F90:559
integer, parameter two_d_domain
Use the FmsNetcdfDomainFile_t fileobj.
Definition: diag_data.F90:101
Attribute type for diagnostic fields.
Definition: diag_data.F90:157
Type to hold the attributes of the field/axis/file.
Definition: diag_data.F90:334
subroutine, public get_local_indexes(latStart, latEnd, lonStart, lonEnd, istart, iend, jstart, jend)
Find the local start and local end indexes on the local PE for regional output.
Definition: diag_grid.F90:391
Add a dimension to a given file.
Definition: fms2_io.F90:189
Defines a new field within the given file Example usage:
Definition: fms2_io.F90:211
Write data to a defined field within a file Example usage:
Definition: fms2_io.F90:261
subroutine, public define_new_subaxis_latlon(diag_axis, axis_ids, naxis, subRegion, is_cube_sphere, write_on_this_pe)
Fill in the subaxis object for a subRegion defined by lat lon.
subroutine get_dim_size_layout(this, dim_size, layout)
pure character(len=:) function, allocatable get_axis_name(this, is_regional)
integer function get_length(this, cart_axis, domain_position, global_length)
Get the length of a 2D domain.
integer function get_ending_index(this)
Accesses its member ending_index.
logical function, public fms_diag_axis_object_end(axis_array)
subroutine check_if_valid_domain_position(domain_position)
Check if a domain_position is valid and crashes if it isn't.
logical function is_x_or_y_axis(this, x_or_y)
Determine if an axis object is an x or y axis.
pure logical function is_unstructured_grid(this)
subroutine set_axis_id(this, axis_id)
Set the axis_id.
pure logical function is_z_axis(this)
integer function axis_length(this)
Get the axis length of a subaxis.
pure integer function get_diurnal_axis_samples(this)
integer function get_ntiles(this)
Get the ntiles in a domain.
subroutine, public define_new_subaxis_index(parent_axis, subRegion, diag_axis, naxis, is_x_or_y, write_on_this_pe)
Fill in the subaxis object for a subRegion defined by index.
subroutine, public find_z_sub_axis_name(dim_name, parent_axis_id, file_axis_id, field_yaml, diag_axis)
Determine the name of the z subaxis by matching the parent axis id and the zbounds in the diag table ...
pure integer function, public get_axis_id_from_name(axis_name, diag_axis, naxis, set_name)
subroutine get_compute_domain(this, compute_idx, need_to_define_axis, tile_number)
subroutine, public create_new_z_subaxis(zbounds, var_axis_ids, diag_axis, naxis, file_axis_id, nfile_axis, nz_subaxis, error_mseg)
Creates a new z subaxis to use.
subroutine, public define_diurnal_axis(diag_axis, naxis, n_diurnal_samples, is_edges)
Defined a new diurnal axis.
pure integer function, dimension(2) get_structured_axis(this)
pure logical function has_aux(this)
Determine if an axis object has an auxiliary name.
logical function is_same_zbounds(this, zbounds)
Determines if the zbounds passed in are the same as those in the file.
subroutine add_structured_axis_ids(this, axis_ids)
subroutine check_if_valid_cart_name(cart_name)
Check if a cart_name is valid and crashes if it isn't.
subroutine write_axis_metadata(this, fms2io_fileobj, edges_in_file, parent_axis)
Write the axis meta data to an open fileobj.
subroutine fill_subaxis(this, starting_index, ending_index, axis_id, parent_id, parent_axis_name, compute_idx, global_idx, zbounds, nz_subaxis)
Fills in the information needed to define a subaxis.
subroutine get_global_io_domain(this, global_io_index, use_collective_writes)
Get the starting and ending indices of the global io domain of the axis.
logical function, public fms_diag_axis_object_init(axis_array)
pure integer function get_subaxes_id(this)
integer function get_starting_index(this)
Accesses its member starting_index.
subroutine get_indices(this, compute_idx, corners_indices, starting_index, ending_index, need_to_define_axis)
Determine if the subRegion is in the current PE. If it is, determine the starting and ending indices ...
pure integer function get_parent_axis_id(this)
subroutine check_if_valid_direction(direction)
Check if a direction is valid and crashes if it isn't.
subroutine register_diag_axis_obj(this, axis_name, axis_data, units, cart_name, long_name, direction, set_name, Domain, Domain2, DomainU, aux, req, tile_count, domain_position, axis_length)
Initialize the axis.
subroutine write_axis_data(this, fms2io_fileobj, parent_axis)
Write the axis data to an open fms2io_fileobj.
subroutine, public define_new_axis(diag_axis, parent_axis, naxis, parent_id, starting_index, ending_index, compute_idx, global_idx, new_axis_id, zbounds, nz_subaxis)
Creates a new subaxis and fills it will all the information it needs.
integer function get_axis_length(this)
Get the length of the axis.
pure character(len=:) function, allocatable get_set_name(this)
Get the set name of an axis object.
pure character(len=120) function, dimension(2), public parse_compress_att(compress_att)
subroutine set_edges(this, edges_name, edges_id)
Set the name and ids of the edges.
logical function, public is_parent_axis(axis_id, parent_axis_id, diag_axis)
Determine if the diag_axis(parent_axis_id) is the parent of diag_axis(axis_id)
subroutine add_axis_attribute(this, att_name, att_value)
Add an attribute to an axis.
subroutine write_diurnal_metadata(this, fms2io_fileobj)
pure logical function has_set_name(this)
Determine if an axis object has a set_name.
pure integer function get_edges_id(this)
integer function, dimension(2) get_compute_indices(this)
Accesses its member compute_indices.
subroutine, public get_domain_and_domain_type(diag_axis, axis_id, domain_type, domain, var_name)
Loop through a variable's axis_id to determine and return the domain type and domain to use.
subroutine set_axis_domain(this, Domain, Domain2, DomainU)
Set the axis domain.
pure character(len=:) function, allocatable get_aux(this)
Get the auxiliary name of an axis object.
integer function, dimension(size(domain%tile_id(:))) mpp_get_tile_id(domain)
Returns the tile_id on current pe.
integer function mpp_get_ntile_count(domain)
Returns number of tiles in mosaic.
type(domain2d) function, pointer mpp_get_io_domain(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....
Retrieve layout associated with a domain decomposition The 1D version of this call returns the number...
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.
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:42
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381
Type to hold the domain info for an axis This type was created to avoid having to send in "Domain",...
Type to hold the unstructured domain.
Type to hold the diagnostic axis description.
Type to hold the diag_axis (either subaxis or a full axis)
Type to hold the diagnostic axis description.
type to hold the info a diag_field
type to hold the sub region information about a file