29 module fms_diag_axis_object_mod
34 use platform_mod,
only: r8_kind, r4_kind, i4_kind, i8_kind
41 use fms2_io_mod,
only: fmsnetcdffile_t, fmsnetcdfdomainfile_t, fmsnetcdfunstructureddomainfile_t, &
67 procedure :: get_ntiles
88 INTEGER ,
private :: axis_id
106 class(fmsDiagAxis_type),
allocatable :: axis
112 CHARACTER(len=:),
ALLOCATABLE ,
private :: subaxis_name
113 INTEGER ,
private :: starting_index
115 INTEGER ,
private :: ending_index
117 INTEGER ,
private :: parent_axis_id
118 INTEGER ,
private :: compute_idx(2)
119 INTEGER,
allocatable,
private :: global_idx(:)
120 real(kind=r4_kind),
allocatable,
private :: zbounds(:)
133 INTEGER ,
private :: ndiurnal_samples
134 CHARACTER(len=:),
ALLOCATABLE,
private :: axis_name
135 CHARACTER(len=:),
ALLOCATABLE,
private :: long_name
136 CHARACTER(len=:),
ALLOCATABLE,
private :: units
137 INTEGER ,
private :: edges_id
138 CHARACTER(len=:),
ALLOCATABLE,
private :: edges_name
139 CLASS(*),
ALLOCATABLE,
private :: diurnal_data(:)
149 CHARACTER(len=:),
ALLOCATABLE,
private :: axis_name
150 CHARACTER(len=:),
ALLOCATABLE,
private :: units
151 CHARACTER(len=:),
ALLOCATABLE,
private :: long_name
152 CHARACTER(len=1) ,
private :: cart_name
153 CLASS(*),
ALLOCATABLE,
private :: axis_data(:)
154 CHARACTER(len=:),
ALLOCATABLE,
private :: type_of_data
156 integer,
ALLOCATABLE,
private :: subaxis(:)
157 integer ,
private :: nsubaxis
159 INTEGER ,
private :: type_of_domain
161 INTEGER ,
private :: length
162 INTEGER ,
private :: direction
163 INTEGER,
ALLOCATABLE,
private :: edges_id
165 CHARACTER(len=:),
ALLOCATABLE,
private :: edges_name
167 CHARACTER(len=:),
ALLOCATABLE,
private :: aux
169 CHARACTER(len=128) ,
private :: req
170 INTEGER ,
private :: tile_count
172 INTEGER ,
private :: num_attributes
173 INTEGER ,
private :: domain_position
174 integer,
allocatable ,
private :: structured_ids(:)
176 CHARACTER(len=:),
ALLOCATABLE,
private :: set_name
206 & set_name, Domain, Domain2, DomainU, aux, req, tile_count, domain_position, axis_length )
208 CHARACTER(len=*),
INTENT(in) :: axis_name
209 class(*),
INTENT(in) :: axis_data(:)
210 CHARACTER(len=*),
INTENT(in) :: units
211 CHARACTER(len=1),
INTENT(in) :: cart_name
212 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: long_name
213 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: set_name
214 INTEGER,
INTENT(in),
OPTIONAL :: direction
215 TYPE(
domain1d),
INTENT(in),
OPTIONAL :: Domain
216 TYPE(
domain2d),
INTENT(in),
OPTIONAL :: Domain2
217 TYPE(
domainug),
INTENT(in),
OPTIONAL :: DomainU
218 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: aux
220 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: req
221 INTEGER,
INTENT(in),
OPTIONAL :: tile_count
222 INTEGER,
INTENT(in),
OPTIONAL :: domain_position
223 integer,
intent(in),
optional :: axis_length
225 this%axis_name = trim(axis_name)
226 this%units = trim(units)
227 this%cart_name = uppercase(cart_name)
230 if (
present(long_name)) this%long_name = trim(long_name)
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"
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"
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.")
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))
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))
260 call this%axis_domain%set(domain2=domain2)
262 else if (
present(domainu))
then
264 call this%axis_domain%set(domainu=domainu)
269 if (
present(tile_count)) this%tile_count = tile_count
271 this%domain_position = center
272 if (
present(domain_position)) this%domain_position = domain_position
276 if (
present(direction)) this%direction = direction
279 if (
present(aux)) this%aux = trim(aux)
280 if (
present(req)) this%req = trim(req)
282 if (
present(set_name)) this%set_name = trim(set_name)
284 if (max_subaxes .gt. 0)
then
285 allocate(this%subaxis(max_subaxes))
286 this%subaxis = diag_null
290 this%num_attributes = 0
296 character(len=*),
intent(in) :: att_name
297 class(*),
intent(in) :: att_value(:)
301 if (.not.
allocated(this%attributes)) &
304 this%num_attributes = this%num_attributes + 1
306 j = this%num_attributes
307 call this%attributes(j)%add(att_name, att_value)
313 class(fmsnetcdffile_t),
INTENT(INOUT) :: fms2io_fileobj
314 logical,
INTENT(IN) :: edges_in_file
320 character(len=:),
ALLOCATABLE :: axis_edges_name
321 character(len=:),
pointer :: axis_name
322 integer :: axis_length
326 integer :: type_of_domain
327 logical :: is_subaxis
328 logical :: needs_domain_decomposition
330 integer :: domain_decomposition(4)
333 needs_domain_decomposition = .false.
337 axis_name => this%axis_name
338 axis_length = this%length
340 type_of_domain = this%type_of_domain
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
352 if (
present(parent_axis))
then
353 select type(parent_axis)
355 diag_axis => parent_axis
360 call this%write_diurnal_metadata(fms2io_fileobj)
365 select type (fms2io_fileobj)
368 type is (fmsnetcdffile_t)
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)
376 type is (fmsnetcdfdomainfile_t)
377 select case (type_of_domain)
382 call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
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/))
388 type is (fmsnetcdfunstructureddomainfile_t)
389 select case (type_of_domain)
393 call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
398 call register_field(fms2io_fileobj, axis_name, diag_axis%type_of_data, (/axis_name/))
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))
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)
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))
414 select case (diag_axis%direction)
416 call register_variable_attribute(fms2io_fileobj, axis_name,
"positive",
"up", str_len=2)
418 call register_variable_attribute(fms2io_fileobj, axis_name,
"positive",
"down", str_len=4)
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))
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)))
434 call register_variable_attribute(fms2io_fileobj, axis_name, diag_axis%attributes(i)%att_name, att_value)
444 class(fmsnetcdffile_t),
INTENT(INOUT) :: fms2io_fileobj
449 integer :: global_io_index(2)
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)))
455 i = this%starting_index
456 j = this%ending_index
458 if (
present(parent_axis))
then
459 select type(parent_axis)
461 call write_data(fms2io_fileobj, this%subaxis_name, parent_axis%axis_data(i:j))
465 call write_data(fms2io_fileobj, this%axis_name, this%diurnal_data)
473 integer,
intent(inout) :: naxis
475 integer,
intent(in) :: n_diurnal_samples
477 logical,
intent(in) :: is_edges
480 CHARACTER(32) :: axis_name
481 CHARACTER(32) :: long_name
482 CHARACTER(32) :: edges_name
483 CHARACTER(128) :: units
484 real(kind=r8_kind),
allocatable :: diurnal_data(:)
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
498 do i = 1, n_diurnal_samples
499 diurnal_data(i+1) = 24.0* real(i)/n_diurnal_samples
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))
506 do i = 1, n_diurnal_samples
507 diurnal_data(i) = 24.0*(real(i)-0.5)/n_diurnal_samples
509 WRITE (edges_name,
'(a,i2.2)')
'time_of_day_edges_', n_diurnal_samples
514 11
FORMAT(a,
' since ',i4.4,
'-',i2.2,
'-',i2.2,
' ',i2.2,
':',i2.2,
':',i2.2)
517 select type (diurnal_axis => diag_axis(naxis)%axis)
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
527 WRITE (edges_name,
'(a,i2.2)')
'time_of_day_edges_', n_diurnal_samples
528 diurnal_axis%edges_name = trim(edges_name)
547 integer,
intent(in) :: axis_ids(2)
551 allocate(this%structured_ids(2))
552 this%structured_ids = axis_ids
566 rslt = this%structured_ids
579 if (
allocated(this%edges_id))
get_edges_id = this%edges_id
586 integer,
intent(out) :: global_io_index(2)
587 logical,
intent(in) :: use_collective_writes
589 type(
domain2d),
pointer :: io_domain
591 global_io_index(1) = 1
592 global_io_index(2) = this%length
594 if (
allocated(this%axis_domain))
then
595 select type(domain => this%axis_domain)
597 if (use_collective_writes)
then
598 io_domain => domain%domain2
603 if (this%cart_name .eq.
"X")
then
605 position=this%domain_position)
606 elseif (this%cart_name .eq.
"Y")
then
608 position=this%domain_position)
622 if (
allocated(this%axis_domain))
then
623 axis_length = this%axis_domain%length(this%cart_name, this%domain_position, this%length)
639 if (
allocated(this%aux)) rslt = trim(this%aux) .ne.
""
650 if (
allocated(this%set_name)) rslt = trim(this%set_name) .ne.
""
658 integer,
optional,
intent(inout) :: x_or_y
662 select case (trim(this%cart_name))
671 if (
present(x_or_y)) x_or_y = diag_null
680 integer,
intent(out) :: dim_size
681 integer,
intent(out) :: layout
684 integer :: layout_xy(2)
686 select type (domain => this%axis_domain)
691 if (this%cart_name .eq.
"X")
then
693 layout = layout_xy(1)
694 else if (this%cart_name .eq.
"Y")
then
696 layout = layout_xy(2)
706 character(len=:),
allocatable :: rslt
716 character(len=:),
allocatable :: rslt
724 integer,
intent(in) :: axis_id
726 this%axis_id = axis_id
733 CHARACTER(len=*),
intent(in) :: edges_name
734 integer,
intent(in) :: edges_id
740 this%edges_name = edges_name
741 this%edges_id = edges_id
746 subroutine get_indices(this, compute_idx, corners_indices, starting_index, ending_index, need_to_define_axis)
748 integer,
intent(in) :: compute_idx(:)
749 class(*),
intent(in) :: corners_indices(:)
750 integer,
intent(out) :: starting_index
752 integer,
intent(out) :: ending_index
754 logical,
intent(out) :: need_to_define_axis
757 integer :: subregion_start
758 integer :: subregion_end
763 select type (corners_indices)
764 type is (
integer(kind=i4_kind))
765 subregion_start = minval(corners_indices)
766 subregion_end = maxval(corners_indices)
770 need_to_define_axis = .false.
771 starting_index = diag_null
772 ending_index = diag_null
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
778 need_to_define_axis = .true.
779 if (compute_idx(1) >= subregion_start .and. compute_idx(2) >= subregion_end)
then
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
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
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
793 starting_index = subregion_start
794 ending_index = subregion_end
797 if (this%domain_position .ne. center)
then
798 if (subregion_end - subregion_start + 1 .eq. 1)
then
800 if (ending_index .eq. compute_idx(2)) need_to_define_axis = .false.
802 if (ending_index - starting_index + 1 .eq. 1)
then
804 if (starting_index .eq. compute_idx(2) .or. ending_index .eq. compute_idx(1)) &
805 need_to_define_axis = .false.
815 integer,
intent(inout) :: compute_idx(:)
816 logical,
intent(out) :: need_to_define_axis
817 integer,
optional,
intent(in) :: tile_number
820 need_to_define_axis = .false.
821 compute_idx = diag_null
823 if (.not.
allocated(this%axis_domain))
then
825 if (this%cart_name .eq.
"X" .or. this%cart_name .eq.
"Y")
then
827 compute_idx(2) =
size(this%axis_data)
828 need_to_define_axis = .true.
833 select type(domain => this%axis_domain)
835 if (
present(tile_number))
then
839 need_to_define_axis = .false.
845 select case (this%cart_name)
848 & position=this%domain_position)
849 need_to_define_axis = .true.
852 & position=this%domain_position)
853 need_to_define_axis = .true.
861 subroutine fill_subaxis(this, starting_index, ending_index, axis_id, parent_id, parent_axis_name, compute_idx, &
862 global_idx, zbounds, nz_subaxis)
864 integer ,
intent(in) :: starting_index
865 integer ,
intent(in) :: ending_index
866 integer ,
intent(in) :: axis_id
867 integer ,
intent(in) :: parent_id
868 character(len=*) ,
intent(in) :: parent_axis_name
869 integer ,
intent(in) :: compute_idx(2)
871 integer,
optional,
intent(in) :: global_idx(2)
873 real(kind=r4_kind),
optional,
intent(in) :: zbounds(2)
874 integer,
optional,
intent(in) :: nz_subaxis
878 character(len=2) :: nsubaxis_char
881 if (
present(nz_subaxis)) nsubaxis = nz_subaxis
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
891 if (
present(zbounds))
then
893 allocate(this%zbounds(2))
894 this%zbounds = zbounds
897 if (
present(global_idx))
then
899 allocate(this%global_idx(2))
900 this%global_idx = global_idx
911 res = this%ending_index - this%starting_index + 1
919 indx = this%starting_index
927 indx = this%ending_index
935 indx = this%compute_idx
942 real(kind=r4_kind),
intent(in) :: zbounds(2)
945 is_same = zbounds(1) .eq. this%zbounds(1) .and. zbounds(2) .eq. this%zbounds(2)
964 function get_length(this, cart_axis, domain_position, global_length) &
967 character(len=*),
INTENT(IN) :: cart_axis
968 integer,
INTENT(IN) :: domain_position
969 integer,
INTENT(IN) :: global_length
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)
979 length = global_length
988 TYPE(
domain1d),
INTENT(in),
OPTIONAL :: Domain
989 TYPE(
domain2d),
INTENT(in),
OPTIONAL :: Domain2
990 TYPE(
domainug),
INTENT(in),
OPTIONAL :: DomainU
996 this%Domain2 = domain2
998 this%DomainUG = domainu
1007 if (
allocated(axis_array))
call mpp_error(fatal,
"The diag_axis containers is already allocated")
1019 if (
allocated(axis_array))
deallocate(axis_array)
1032 logical,
intent(in),
optional :: is_regional
1034 character(len=:),
allocatable :: axis_name
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"
1045 axis_name = this%subaxis_name
1056 if (this%cart_name .eq.
"Z")
is_z_axis = .true.
1062 character(len=*),
intent(in) :: cart_name
1064 select case (cart_name)
1065 case (
"X",
"Y",
"Z",
"T",
"U",
"N")
1067 call mpp_error(fatal,
"diag_axit_init: Invalid cart_name: "//cart_name//&
1068 "The acceptable values are X, Y, Z, T, U, N.")
1074 integer,
INTENT(IN) :: domain_position
1076 select case (domain_position)
1077 case (center, north, east)
1079 call mpp_error(fatal,
"diag_axit_init: Invalid domain_positon. &
1080 &The acceptable values are NORTH, EAST, CENTER")
1086 integer,
INTENT(IN) :: direction
1088 select case(direction)
1091 call mpp_error(fatal,
"diag_axit_init: Invalid direction. &
1092 &The acceptable values are-1 0 1")
1099 integer,
INTENT(IN) :: axis_id(:)
1100 integer,
INTENT(OUT) :: domain_type
1102 character(len=*),
INTENT(IN) :: var_name
1110 do i = 1,
size(axis_id)
1112 select type (axis => diag_axis(j)%axis)
1115 if (domain_type .ne. axis%type_of_domain)
then
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
1126 call mpp_error(fatal,
"The variable:"//trim(var_name)//
" has axis that are not in the same domain")
1137 integer,
intent(inout) :: naxis
1139 integer,
intent(in) :: is_x_or_y
1141 logical,
intent(out) :: write_on_this_pe
1143 integer :: compute_idx(2)
1144 integer :: global_idx(2)
1145 integer :: starting_index
1146 integer :: ending_index
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
1153 call parent_axis%get_indices(compute_idx, subregion%corners(:,is_x_or_y), starting_index, ending_index, &
1156 if (.not. write_on_this_pe)
return
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))
1165 call define_new_axis(diag_axis, parent_axis, naxis, parent_axis%axis_id, &
1166 starting_index, ending_index, compute_idx, global_idx)
1173 integer,
INTENT(in) :: axis_ids(:)
1174 integer,
intent(inout) :: naxis
1176 logical,
intent(in) :: is_cube_sphere
1177 logical,
intent(out) :: write_on_this_pe
1182 integer :: lat_indices(2)
1183 integer :: lon_indices(2)
1184 integer :: compute_idx(2)
1185 integer :: starting_index(2)
1187 integer :: ending_index(2)
1188 logical :: need_to_define_axis(2)
1190 integer :: parent_axis_ids(2)
1191 logical :: is_x_y_axis
1192 integer :: compute_idx_2(2, 2)
1193 integer :: global_idx (2, 2)
1195 write_on_this_pe = .false.
1196 need_to_define_axis = .true.
1197 parent_axis_ids = diag_null
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))
1210 if_is_cube_sphere:
if (is_cube_sphere)
then
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)
1218 call parent_axis%get_compute_domain(compute_idx, is_x_y_axis)
1221 if (.not. is_x_y_axis) cycle
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
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)
1245 call parent_axis%get_compute_domain(compute_idx, is_x_y_axis)
1248 if (.not. is_x_y_axis) cycle
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)
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)
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
1281 enddo loop_over_axis_ids2
1282 endif if_is_cube_sphere
1285 if (any(.not. need_to_define_axis ))
return
1288 write_on_this_pe = .true.
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)
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,:))
1303 starting_index, ending_index, compute_idx, global_idx, new_axis_id, zbounds, &
1308 integer,
intent(inout) :: naxis
1310 integer,
intent(in) :: parent_id
1311 integer,
intent(in) :: starting_index
1312 integer,
intent(in) :: ending_index
1313 integer,
intent(in) :: compute_idx(2)
1315 integer,
optional,
intent(in) :: global_idx(2)
1317 integer,
optional,
intent(out) :: new_axis_id
1318 real(kind=r4_kind),
optional,
intent(in) :: zbounds(2)
1319 integer,
optional,
intent(in) :: nz_subaxis
1325 parent_axis%nsubaxis = parent_axis%nsubaxis + 1
1326 parent_axis%subaxis(parent_axis%nsubaxis) = naxis
1330 diag_axis(naxis)%axis%axis_id = naxis
1331 if (
present(new_axis_id)) new_axis_id = naxis
1333 select type (sub_axis => diag_axis(naxis)%axis)
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)
1343 result(parent_axis_id)
1346 integer :: parent_axis_id
1350 parent_axis_id = diag_null
1352 parent_axis_id = this%parent_axis_id
1354 parent_axis_id = diag_null
1365 integer :: sub_axis_id
1367 sub_axis_id = this%axis_id
1370 if (this%cart_name .ne.
"Z") sub_axis_id = this%subaxis(this%nsubaxis)
1379 class(*),
intent(in) :: compress_att(:)
1380 character(len=120) :: axis_names(2)
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 =
""
1398 character(len=*),
intent(in) :: axis_name
1399 integer,
intent(in) :: naxis
1400 character(len=*),
intent(in) :: set_name
1407 select type(axis => diag_axis(i)%axis)
1409 if (trim(axis%axis_name) .eq. trim(axis_name))
then
1410 if (trim(axis%set_name) .eq. trim(set_name))
then
1423 result(n_diurnal_samples)
1426 integer :: n_diurnal_samples
1428 n_diurnal_samples = this%ndiurnal_samples
1434 class(fmsnetcdffile_t),
intent(inout) :: fms2io_fileobj
1436 call register_axis(fms2io_fileobj, this%axis_name,
size(this%diurnal_data))
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))
1450 real(kind=r4_kind),
intent(in) :: zbounds(2)
1451 integer,
intent(inout) :: var_axis_ids(:)
1453 integer,
intent(inout) :: naxis
1455 integer,
intent(inout) :: file_axis_id(:)
1456 integer,
intent(inout) :: nfile_axis
1458 integer,
intent(inout) :: nz_subaxis
1460 character(len=*),
intent(inout) :: error_mseg
1463 class(*),
pointer :: zaxis_data(:)
1464 integer :: subaxis_indices(2)
1467 integer :: subaxis_id
1468 integer :: parent_axis_id
1469 integer :: zaxis_index
1472 parent_axis_id = diag_null
1473 zaxis_index = diag_null
1476 do i = 1,
size(var_axis_ids)
1477 select type (parent_axis => diag_axis(var_axis_ids(i))%axis)
1479 if (parent_axis%cart_name .eq.
"Z")
then
1480 parent_axis_id = var_axis_ids(i)
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))
1491 do i = 1, nfile_axis
1492 select type (axis => diag_axis(file_axis_id(i))%axis)
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)
1502 select type (axis => diag_axis(parent_axis_id)%axis)
1504 zaxis_data => axis%axis_data
1508 select type(zaxis_data)
1509 type is (real(kind=r4_kind))
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))
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
1530 integer,
intent(in) :: axis_id
1531 integer,
intent(in) :: parent_axis_id
1537 select type(axis => diag_axis(axis_id)%axis)
1539 if (axis%parent_axis_id .eq. parent_axis_id) rslt = .true.
1546 character(len=*),
intent(inout) :: dim_name
1547 integer,
intent(in) :: parent_axis_id
1548 integer,
intent(in) :: file_axis_id(:)
1555 do i = 1,
size(file_axis_id)
1556 id = file_axis_id(i)
1557 select type (axis_ptr => diag_axis(id)%axis)
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
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()))
1572 end module fms_diag_axis_object_mod
integer, parameter direction_down
The axis points down if positive.
integer function get_base_minute()
gets the module variable base_minute
integer function get_base_year()
gets the module variable base_year
integer function get_base_hour()
gets the module variable base_hour
integer, parameter no_domain
Use the FmsNetcdfFile_t fileobj.
integer max_axis_attributes
Maximum number of user definable attributes per axis.
character(len=6) pack_size_str
Pack size as a string to be used in fms2_io register call set to "double" or "float".
integer max_axes
Maximum number of independent axes.
integer, parameter is_x_axis
integer indicating that it is a x axis
integer, parameter is_y_axis
integer indicating that it is a y axis
integer function get_base_day()
gets the module variable base_day
integer, parameter ug_domain
Use the FmsNetcdfUnstructuredDomainFile_t fileobj.
integer, parameter direction_up
The axis points up if positive.
integer function get_base_month()
gets the module variable base_month
integer function get_base_second()
gets the module variable base_second
integer, parameter two_d_domain
Use the FmsNetcdfDomainFile_t fileobj.
Attribute type for diagnostic fields.
Type to hold the attributes of the field/axis/file.
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.
Add a dimension to a given file.
Defines a new field within the given file Example usage:
Write data to a defined field within a file Example usage:
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.
integer function mpp_pe()
Returns processor ID.
Type to hold the 1d domain.
Type to hold the 2d domain.
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 diurnal axis.
Type to hold the diagnostic axis description.
Type to hold the subaxis.
type to hold the info a diag_field
type to hold the sub region information about a file