36 & north, east, center, &
43 use fms_diag_object_mod,
only:fms_diag_object
44 USE netcdf,
ONLY: nf90_int, nf90_float, nf90_char
59 #include<file_version.h>
62 integer(I4_KIND),
parameter,
public :: DIAG_AXIS_NODOMAIN = 0
68 INTEGER :: num_def_axes = 0
70 CHARACTER(len=128),
DIMENSION(:),
ALLOCATABLE,
SAVE ::
axis_sets
71 INTEGER :: num_axis_sets = 0
74 LOGICAL :: module_is_initialized = .false.
90 MODULE PROCEDURE diag_axis_add_attribute_scalar_r
91 MODULE PROCEDURE diag_axis_add_attribute_scalar_i
92 MODULE PROCEDURE diag_axis_add_attribute_scalar_c
93 MODULE PROCEDURE diag_axis_add_attribute_r1d
94 MODULE PROCEDURE diag_axis_add_attribute_i1d
109 INTEGER FUNCTION diag_axis_init(name, array_data, units, cart_name, long_name, direction,&
110 & set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count, domain_position )
111 CHARACTER(len=*),
INTENT(in) :: name
112 CLASS(*),
DIMENSION(:),
INTENT(in) :: array_data
113 CHARACTER(len=*),
INTENT(in) :: units
114 CHARACTER(len=*),
INTENT(in) :: cart_name
115 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: long_name
116 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: set_name
117 INTEGER,
INTENT(in),
OPTIONAL :: direction
118 INTEGER,
INTENT(in),
OPTIONAL :: edges
119 TYPE(
domain1d),
INTENT(in),
OPTIONAL :: domain
120 TYPE(
domain2d),
INTENT(in),
OPTIONAL :: domain2
121 TYPE(
domainug),
INTENT(in),
OPTIONAL :: domainu
122 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: aux
123 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: req
124 INTEGER,
INTENT(in),
OPTIONAL :: tile_count
125 INTEGER,
INTENT(in),
OPTIONAL :: domain_position
128 TYPE(
domain1d) :: domain_x, domain_y
129 INTEGER :: ierr, axlen
130 INTEGER :: i, set, tile
131 INTEGER :: isc, iec, isg, ieg
132 CHARACTER(len=128) :: emsg
134 IF ( .NOT.module_is_initialized )
THEN
141 diag_axis_init = fms_diag_object%fms_diag_axis_init(name, array_data, units, cart_name,
size(array_data(:)), &
142 & long_name=long_name, direction=direction, set_name=set_name, edges=edges, domain=domain, domain2=domain2, &
143 & domainu=domainu, aux=aux, req=req, tile_count=tile_count, domain_position=domain_position)
146 IF (
PRESENT(tile_count))
THEN
161 IF (
PRESENT(set_name) )
THEN
165 num_axis_sets = num_axis_sets + 1
166 IF ( num_axis_sets > max_num_axis_sets )
THEN
167 WRITE (emsg, fmt=
'("num_axis_sets (",I2,") exceeds max_num_axis_sets (",I2,"). ")')&
168 & num_axis_sets, max_num_axis_sets
173 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
174 & trim(emsg)//
' Increase max_num_axis_sets via diag_manager_nml.', fatal)
186 DO i = 1, num_def_axes
187 IF ( trim(name) ==
axes(i)%name )
THEN
188 IF ( trim(name) ==
'Stations' .OR. trim(name) ==
'Levels')
THEN
191 ELSE IF ( set ==
axes(i)%set )
THEN
192 IF ( trim(lowercase(name)) ==
'time' .OR.&
193 & trim(lowercase(cart_name)) ==
't' .OR.&
194 & trim(lowercase(name)) ==
'nv' .OR.&
195 & trim(lowercase(cart_name)) ==
'n' )
THEN
198 ELSE IF ( (lowercase(cart_name) /=
'x' .AND. lowercase(cart_name) /=
'y')&
199 & .OR. tile /=
axes(i)%tile_count)
THEN
201 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
202 &
'axis_name '//trim(name)//
' and axis_set already exist.', fatal)
209 num_def_axes = num_def_axes + 1
212 &
'max_axes exceeded, increase via diag_manager_nml', fatal)
216 IF ( trim(uppercase(cart_name)) ==
'X' .OR.&
217 & trim(uppercase(cart_name)) ==
'Y' .OR.&
218 & trim(uppercase(cart_name)) ==
'Z' .OR.&
219 & trim(uppercase(cart_name)) ==
'T' .OR.&
220 & trim(uppercase(cart_name)) ==
'U' .OR.&
221 & trim(uppercase(cart_name)) ==
'N' )
THEN
225 CALL error_mesg(
'diag_axis_mod::diag_axis_init',
'Invalid cart_name name. '//trim(uppercase(cart_name)), fatal)
232 axlen =
SIZE(array_data(:))
238 SELECT TYPE (array_data)
239 TYPE IS (real(kind=r4_kind))
241 TYPE IS (real(kind=r8_kind))
244 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
245 &
'The axis data is not one of the supported types of real(kind=4) or real(kind=8)', fatal)
257 IF (
PRESENT(long_name) )
THEN
263 IF (
PRESENT(aux) )
THEN
269 IF (
PRESENT(req) )
THEN
274 IF (
PRESENT(domain_position) )
THEN
275 if (domain_position == north .or. domain_position == east .or. domain_position == center)
then
278 CALL error_mesg(
'diag_axis_mod::diag_axis_init',
"Position must be NORTH, EAST, or CENTER" ,&
286 IF (
PRESENT(direction) )
THEN
287 IF ( abs(direction) /= 1 .AND. direction /= 0 )&
289 &
CALL error_mesg(
'diag_axis_mod::diag_axis_init',
'direction must be 0, +1 or -1', fatal)
296 IF (
present(domainu) .AND. (
PRESENT(domain2) .OR.
PRESENT(domain)) )
THEN
298 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
299 &
'Presence of DomainU and another Domain at the same time is prohibited', fatal)
301 ELSE IF (
PRESENT(domain2) .AND.
PRESENT(domain))
THEN
303 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
304 &
'Presence of both Domain and Domain2 at the same time is prohibited', fatal)
305 ELSE IF (
PRESENT(domain2) .OR.
PRESENT(domain))
THEN
308 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
309 &
'A Structured Domain must not be present for an axis which is not in the X or Y direction', fatal)
312 CALL error_mesg(
'diag_axis_mod::diag_axis_init',&
313 &
'In the unstructured domain, the axis cart_name must be U', fatal)
318 IF (
PRESENT(domain2) )
THEN
324 ELSE IF (
PRESENT(domain))
THEN
329 ELSE IF (
present(domainu))
THEN
350 IF (
PRESENT(edges) )
THEN
351 IF ( edges > 0 .AND. edges < num_def_axes )
THEN
358 WRITE (emsg,
'("Edges axis does not match axis (code ",I1,").")') ierr
359 CALL error_mesg(
'diag_axis_mod::diag_axis_init', emsg, fatal)
364 CALL error_mesg(
'diag_axis_mod::diag_axis_init',
'Edges axis is not defined', fatal)
369 module_is_initialized = .true.
383 INTEGER,
INTENT(in) :: axis
384 REAL,
DIMENSION(:),
INTENT(in) :: subdata
385 INTEGER,
INTENT(in) :: start_indx
386 INTEGER,
INTENT(in) :: end_indx
387 TYPE(
domain2d),
INTENT(in),
OPTIONAL :: domain_2d
389 INTEGER :: i, nsub_axis, direction
390 INTEGER :: xbegin, xend, ybegin, yend
391 INTEGER :: ad_xbegin, ad_xend, ad_ybegin, ad_yend
392 CHARACTER(len=128) :: name, nsub_name
393 CHARACTER(len=128) :: units
394 CHARACTER(len=128) :: cart_name
395 CHARACTER(len=128) :: long_name
396 CHARACTER(len=128) :: emsg
397 LOGICAL :: subaxis_set, hasdomain
401 subaxis_set = .false.
403 IF (
PRESENT(domain_2d) )
THEN
410 IF ( start_indx ==
axes(axis)%start(i) .AND. end_indx ==
axes(axis)%end(i) )
THEN
411 IF ( hasdomain )
THEN
413 IF ( .NOT.((xbegin == ad_xbegin .AND. xend == ad_xend) .AND.&
414 & (ybegin == ad_ybegin .AND. yend == ad_yend)) )
THEN
420 name = trim(
axes(axis)%subaxis_name(nsub_axis))
425 IF ( nsub_axis == 0 )
THEN
429 WRITE (emsg,
'("max_subaxes (value ",I4,") is too small. Consider increasing max_subaxes.")') max_subaxes
430 CALL error_mesg(
'diag_axis_mod::diag_subaxes_init', emsg, fatal)
433 axes(axis)%start(nsub_axis) = start_indx
434 axes(axis)%end(nsub_axis) = end_indx
435 if ( hasdomain )
axes(axis)%subaxis_domain2(nsub_axis) = domain_2d
441 IF (
axes(axis)%set > 0 )
THEN
449 WRITE (nsub_name,
'(I2.2)') nsub_axis
450 name = trim(
axes(axis)%name)//
'_sub'//trim(nsub_name)
451 axes(axis)%subaxis_name(nsub_axis) = name
452 long_name = trim(
axes(axis)%long_name)
453 units = trim(
axes(axis)%units)
454 cart_name = trim(
axes(axis)%cart_name)
455 direction =
axes(axis)%direction
456 IF (
axes(axis)%set > 0)
THEN
458 & set_name=trim(
axis_sets(
axes(axis)%set)), direction=direction, domain2=domain_2d)
461 & direction=direction, domain2=domain_2d)
467 & direction, edges, Domain, DomainU, array_data, num_attributes, attributes, domain_position)
468 CHARACTER(len=*),
INTENT(out) :: name, units, long_name, cart_name
469 INTEGER,
INTENT(in) :: id
470 TYPE(
domain1d),
INTENT(out) :: domain
471 TYPE(
domainug),
INTENT(out) :: domainu
472 INTEGER,
INTENT(out) :: direction
474 INTEGER,
INTENT(out) :: edges
475 CLASS(*),
DIMENSION(:),
INTENT(out) :: array_data
476 INTEGER,
INTENT(out),
OPTIONAL :: num_attributes
477 TYPE(
diag_atttype),
ALLOCATABLE,
DIMENSION(:),
INTENT(out),
OPTIONAL :: attributes
478 INTEGER,
INTENT(out),
OPTIONAL :: domain_position
480 INTEGER :: i, j, istat
484 units =
axes(id)%units
485 long_name =
axes(id)%long_name
486 cart_name =
axes(id)%cart_name
487 direction =
axes(id)%direction
488 edges =
axes(id)%edges
489 domain =
axes(id)%Domain
490 domainu =
axes(id)%DomainUG
491 if (
present(domain_position)) domain_position =
axes(id)%domain_position
492 IF (
axes(id)%length >
SIZE(array_data(:)) )
THEN
494 CALL error_mesg(
'diag_axis_mod::get_diag_axis',
'array data is too small', fatal)
496 SELECT TYPE (array_data)
497 TYPE IS (real(kind=r4_kind))
498 array_data(1:
axes(id)%length) = real(
axes(id)%diag_type_data(1:
axes(id)%length), kind=r4_kind)
499 TYPE IS (real(kind=r8_kind))
500 array_data(1:
axes(id)%length) =
axes(id)%diag_type_data(1:
axes(id)%length)
502 CALL error_mesg(
'diag_axis_mod::get_diag_axis',&
503 &
'The axis data is not one of the supported types of real(kind=4) or real(kind=8)', fatal)
506 IF (
PRESENT(num_attributes) )
THEN
507 num_attributes =
axes(id)%num_attributes
509 IF (
PRESENT(attributes) )
THEN
510 IF (
allocated(
axes(id)%attributes) )
THEN
511 IF (
ALLOCATED(attributes) )
THEN
513 IF (
axes(id)%num_attributes .GT.
SIZE(attributes(:)) )
THEN
514 CALL error_mesg(
'diag_axis_mod::get_diag_axis',
'array attribute is too small', fatal)
518 ALLOCATE(attributes(
axes(id)%num_attributes), stat=istat)
519 IF ( istat .NE. 0 )
THEN
520 CALL error_mesg(
'diag_axis_mod::get_diag_axis',
'Unable to allocate memory for attribute', fatal)
523 DO i=1,
axes(id)%num_attributes
525 IF (
allocated(attributes(i)%fatt) )
THEN
526 DEALLOCATE(attributes(i)%fatt)
528 IF (
allocated(attributes(i)%iatt) )
THEN
529 DEALLOCATE(attributes(i)%iatt)
533 attributes(i)%type =
axes(id)%attributes(i)%type
534 attributes(i)%len =
axes(id)%attributes(i)%len
535 attributes(i)%name =
axes(id)%attributes(i)%name
536 attributes(i)%catt =
axes(id)%attributes(i)%catt
538 IF (
allocated(
axes(id)%attributes(i)%fatt) )
THEN
539 ALLOCATE(attributes(i)%fatt(
SIZE(
axes(id)%attributes(i)%fatt(:))), stat=istat)
540 IF ( istat .NE. 0 )
THEN
541 CALL error_mesg(
'diag_axis_mod::get_diag_axis', &
542 &
'Unable to allocate memory for attribute%fatt', fatal)
544 DO j=1,
SIZE(attributes(i)%fatt(:))
545 attributes(i)%fatt(j) =
axes(id)%attributes(i)%fatt(j)
549 IF (
allocated(
axes(id)%attributes(i)%iatt) )
THEN
550 ALLOCATE(attributes(i)%iatt(
SIZE(
axes(id)%attributes(i)%iatt(:))), stat=istat)
551 IF ( istat .NE. 0 )
THEN
552 CALL error_mesg(
'diag_axis_mod::get_diag_axis', &
553 &
'Unable to allocate memory for attribute%iatt', fatal)
555 DO j=1,
SIZE(attributes(i)%iatt(:))
556 attributes(i)%iatt(j) =
axes(id)%attributes(i)%iatt(j)
566 INTEGER,
INTENT(in) :: id
567 CHARACTER(len=*),
INTENT(out) :: cart_name
570 cart_name =
axes(id)%cart_name
575 INTEGER,
INTENT(in) :: id
576 REAL,
DIMENSION(:),
INTENT(out) :: axis_data
579 IF (
axes(id)%length >
SIZE(axis_data(:)))
THEN
581 CALL error_mesg(
'diag_axis_mod::get_diag_axis_data',
'array data is too small', fatal)
583 axis_data(1:
axes(id)%length) =
axes(id)%diag_type_data
589 INTEGER ,
INTENT(in) :: id
590 CHARACTER(len=*),
INTENT(out) :: axis_name
593 axis_name = fms_diag_object%fms_get_axis_name_from_id(id)
596 axis_name =
axes(id)%name
602 INTEGER,
INTENT(in) :: id
603 CHARACTER(len=*),
INTENT(out) :: name
612 INTEGER,
INTENT(in) :: id
619 IF (
axes(id)%Domain .NE. null_domain1d )
THEN
632 INTEGER,
INTENT(in) :: id
641 INTEGER,
INTENT(in) :: id
650 INTEGER,
INTENT(in) :: id
659 INTEGER,
DIMENSION(:),
INTENT(in) :: ids
662 INTEGER :: i, id, flag
664 IF (
SIZE(ids(:)) < 1 )
THEN
666 CALL error_mesg(
'diag_axis_mod::get_tile_count',
'input argument has incorrect size', fatal)
670 DO i = 1,
SIZE(ids(:))
673 IF (
axes(id)%cart_name ==
'X' .OR. &
674 axes(id)%cart_name ==
'Y' ) flag = flag + 1
676 IF ( flag == 2 )
THEN
686 INTEGER,
INTENT(in) :: id
689 IF (
axes(id)%Domain .NE. null_domain1d)
THEN
699 INTEGER,
DIMENSION(:),
INTENT(in) :: ids
702 INTEGER :: i, id, flag
704 IF (
SIZE(ids(:)) < 1 )
THEN
706 CALL error_mesg(
'diag_axis_mod::get_domain2d',
'input argument has incorrect size', fatal)
709 if (use_modern_diag)
then
716 DO i = 1,
SIZE(ids(:))
719 IF (
axes(id)%cart_name ==
'X' .OR.
axes(id)%cart_name ==
'Y' ) flag = flag + 1
721 IF ( flag == 2 )
THEN
731 INTEGER,
INTENT(in) :: id
734 IF (
axes(id)%DomainUG .NE. null_domainug)
THEN
747 integer,
dimension(:),
intent(in) :: id
748 character(*),
intent(in),
optional :: varname
749 integer(I4_KIND) :: domain_type
757 logical :: uses_domain2d
758 logical :: uses_domainug
763 uses_domain2d = .false.
764 uses_domainug = .false.
770 "axis_compatible_check")
771 if (
axes(id(n))%cart_name .eq.
"X" .or. &
772 axes(id(n))%cart_name .eq.
"Y")
then
774 elseif (
axes(id(n))%cart_name .eq.
"U")
then
777 if (
axes(id(n))%Domain2 .ne. null_domain2d)
then
778 uses_domain2d = .true.
779 elseif (
axes(id(n))%DomainUG .ne. null_domainug)
then
780 uses_domainug = .true.
783 if (ug .and. xory)
then
784 if (
present(varname))
then
785 call error_mesg(
"axis_compatible_check", &
786 "Can not use an unstructured grid with a "// &
787 "horizontal cartesian coordinate for the field " &
791 call error_mesg(
"axis_compatible_check", &
792 "Can not use an unstructured grid with a horizontal "// &
793 "cartesian coordinate", &
797 if (uses_domain2d .and. uses_domainug)
then
798 if (
present(varname))
then
799 call error_mesg(
"axis_compatible_check", &
800 "Can not use an unstructured grid with a"// &
801 "structured grid for the field "//trim(varname), &
804 call error_mesg(
"axis_compatible_check", &
805 "Can not use an unstructured grid with a"// &
806 "structured grid.", &
810 if (uses_domain2d)
then
812 elseif (uses_domainug)
then
815 domain_type = diag_axis_nodomain
823 INTEGER,
DIMENSION(:),
INTENT(in) :: ids
824 INTEGER,
INTENT(out) :: ishift
825 INTEGER,
INTENT(out) :: jshift
832 DO i = 1,
SIZE(ids(:))
835 SELECT CASE (
axes(id)%cart_name)
837 ishift =
axes(id)%shift
839 jshift =
axes(id)%shift
847 CHARACTER(len=*),
INTENT(in) :: axis_name
848 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: set_name
852 IF (
PRESENT(set_name) )
THEN
858 DO n = 1, num_def_axes
859 IF ( trim(axis_name) == trim(
axes(n)%name) .AND.
axes(n)%set == set )
THEN
869 CHARACTER(len=*),
INTENT(in) :: set_name
874 DO iset = 1, num_axis_sets
885 INTEGER,
INTENT(in) :: id
886 CHARACTER(len=*),
INTENT(in) :: routine_name
888 CHARACTER(len=5) :: emsg
890 IF ( id < 1 .OR. id > num_def_axes)
THEN
894 WRITE (emsg,
'(I2)') id
895 CALL error_mesg(
'diag_axis_mod::'//trim(routine_name),&
896 &
'Illegal value for axis_id used (value '//trim(emsg)//
').', fatal)
901 INTEGER,
INTENT(in) :: diag_axis_id
902 CHARACTER(len=*) :: name
903 INTEGER,
INTENT(in) ::
type
904 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: cval
905 INTEGER,
DIMENSION(:),
INTENT(in),
OPTIONAL :: ival
906 REAL,
DIMENSION(:),
INTENT(in),
OPTIONAL :: rval
908 INTEGER :: istat, length, i, this_attribute
909 CHARACTER(len=1024) :: err_msg
911 IF ( .NOT.first_send_data_call )
THEN
917 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',
'Attempting to add attribute "'&
918 &//trim(name)//
'" to axis after first send_data call. Too late.', fatal)
922 IF ( diag_axis_id .LE. 0 )
THEN
924 ELSE IF ( diag_axis_id .GT. num_def_axes )
THEN
929 WRITE(err_msg,
'(I5)') diag_axis_id
930 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',
'Attempting to add attribute "'&
931 &//trim(name)//
'" to axis ID "'//trim(err_msg)//
'", however ID unknown.', fatal)
939 DO i=1,
axes(diag_axis_id)%num_attributes
940 IF ( trim(
axes(diag_axis_id)%attributes(i)%name) .EQ. trim(name) )
THEN
946 IF ( this_attribute.NE.0 .AND. (type.EQ.nf90_int .OR. type.EQ.nf90_float) )
THEN
951 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',&
952 &
'Attribute "'//trim(name)//
'" already defined for axis "'&
953 &//trim(
axes(diag_axis_id)%name)//
'". Contact the developers.', fatal)
954 ELSE IF ( this_attribute.NE.0 .AND. type.EQ.nf90_char .AND. debug_diag_manager )
THEN
959 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',&
960 &
'Attribute "'//trim(name)//
'" already defined for axis"'&
961 &//trim(
axes(diag_axis_id)%name)//
'". Prepending.', note)
965 this_attribute =
axes(diag_axis_id)%num_attributes + 1
967 IF ( this_attribute .GT. max_axis_attributes )
THEN
972 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',&
973 &
'Number of attributes exceeds max_axis_attributes for attribute "'&
974 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)&
975 & //
'". Increase diag_manager_nml:max_axis_attributes.',&
978 axes(diag_axis_id)%num_attributes = this_attribute
980 axes(diag_axis_id)%attributes(this_attribute)%name = name
981 axes(diag_axis_id)%attributes(this_attribute)%type =
type
983 axes(diag_axis_id)%attributes(this_attribute)%catt =
''
989 IF ( .NOT.
PRESENT(ival) )
THEN
994 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',&
995 &
'Attribute type claims INTEGER, but ival not present for attribute "'&
996 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)&
997 & //
'". Contact the developers.', fatal)
1001 ALLOCATE(
axes(diag_axis_id)%attributes(this_attribute)%iatt(length), stat=istat)
1002 IF ( istat.NE.0 )
THEN
1006 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',
'Unable to allocate iatt for attribute "'&
1007 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)//
'"', fatal)
1010 axes(diag_axis_id)%attributes(this_attribute)%len = length
1011 axes(diag_axis_id)%attributes(this_attribute)%iatt = ival
1013 IF ( .NOT.
PRESENT(rval) )
THEN
1018 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',&
1019 &
'Attribute type claims REAL, but rval not present for attribute "'&
1020 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)&
1021 & //
'". Contact the developers.', fatal)
1025 ALLOCATE(
axes(diag_axis_id)%attributes(this_attribute)%fatt(length), stat=istat)
1026 IF ( istat.NE.0 )
THEN
1030 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',
'Unable to allocate fatt for attribute "'&
1031 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)&
1035 axes(diag_axis_id)%attributes(this_attribute)%len = length
1036 axes(diag_axis_id)%attributes(this_attribute)%fatt = rval
1038 IF ( .NOT.
PRESENT(cval) )
THEN
1043 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',&
1044 &
'Attribute type claims CHARACTER, but cval not present for attribute "'&
1045 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)&
1046 & //
'". Contact the developers.', fatal)
1054 CALL error_mesg(
'diag_manager_mod::diag_axis_add_attribute',
'Unknown attribute type for attribute "'&
1055 & //trim(name)//
'" for axis "'//trim(
axes(diag_axis_id)%name)&
1056 & //
'". Contact the developers.', fatal)
1061 SUBROUTINE diag_axis_add_attribute_scalar_r(diag_axis_id, att_name, att_value)
1062 INTEGER,
INTENT(in) :: diag_axis_id
1063 CHARACTER(len=*),
INTENT(in) :: att_name
1064 REAL,
INTENT(in) :: att_value
1066 if (use_modern_diag)
then
1067 call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, (/ att_value /))
1069 CALL diag_axis_add_attribute_r1d(diag_axis_id, att_name, (/ att_value /))
1071 END SUBROUTINE diag_axis_add_attribute_scalar_r
1073 SUBROUTINE diag_axis_add_attribute_scalar_i(diag_axis_id, att_name, att_value)
1074 INTEGER,
INTENT(in) :: diag_axis_id
1075 CHARACTER(len=*),
INTENT(in) :: att_name
1076 INTEGER,
INTENT(in) :: att_value
1078 if (use_modern_diag)
then
1079 call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, (/ att_value /))
1081 CALL diag_axis_add_attribute_i1d(diag_axis_id, att_name, (/ att_value /))
1083 END SUBROUTINE diag_axis_add_attribute_scalar_i
1085 SUBROUTINE diag_axis_add_attribute_scalar_c(diag_axis_id, att_name, att_value)
1086 INTEGER,
INTENT(in) :: diag_axis_id
1087 CHARACTER(len=*),
INTENT(in) :: att_name
1088 CHARACTER(len=*),
INTENT(in) :: att_value
1090 if (use_modern_diag)
then
1091 call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, (/ att_value /))
1095 END SUBROUTINE diag_axis_add_attribute_scalar_c
1097 SUBROUTINE diag_axis_add_attribute_r1d(diag_axis_id, att_name, att_value)
1098 INTEGER,
INTENT(in) :: diag_axis_id
1099 CHARACTER(len=*),
INTENT(in) :: att_name
1100 REAL,
DIMENSION(:),
INTENT(in) :: att_value
1102 if (use_modern_diag)
then
1103 call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, att_value)
1107 END SUBROUTINE diag_axis_add_attribute_r1d
1109 SUBROUTINE diag_axis_add_attribute_i1d(diag_axis_id, att_name, att_value)
1110 INTEGER,
INTENT(in) :: diag_axis_id
1111 CHARACTER(len=*),
INTENT(in) :: att_name
1112 INTEGER,
DIMENSION(:),
INTENT(in) :: att_value
1113 if (use_modern_diag)
then
1114 call fms_diag_object%fms_diag_axis_add_attribute(diag_axis_id, att_name, att_value)
1118 END SUBROUTINE diag_axis_add_attribute_i1d
1123 TYPE(diag_axis_type),
INTENT(inout) :: out_axis
1124 CHARACTER(LEN=*),
INTENT(out),
OPTIONAL :: err_msg
1129 IF (
PRESENT(err_msg) ) err_msg =
''
1132 IF ( .NOT.
allocated(out_axis%attributes) )
THEN
1133 ALLOCATE(out_axis%attributes(max_axis_attributes), stat=istat)
1134 IF ( istat.NE.0 )
THEN
1138 IF ( fms_error_handler(
'diag_util_mod::attribute_init_axis', &
1139 &
'Unable to allocate memory for diag axis attributes', err_msg) )
THEN
1144 out_axis%num_attributes = 0
1152 TYPE(diag_axis_type),
INTENT(inout) :: out_axis
1153 CHARACTER(len=*),
INTENT(in) :: att_name
1154 CHARACTER(len=*),
INTENT(in) :: prepend_value
1155 CHARACTER(len=*),
INTENT(out) ,
OPTIONAL :: err_msg
1157 INTEGER :: length, i, this_attribute
1158 CHARACTER(len=512) :: err_msg_local
1162 IF (
PRESENT(err_msg) ) err_msg =
''
1166 IF ( trim(err_msg_local) .NE.
'' )
THEN
1167 IF ( fms_error_handler(
'diag_util_mod::prepend_attribute_axis', trim(err_msg_local), err_msg) )
THEN
1174 DO i=1, out_axis%num_attributes
1175 IF ( trim(out_axis%attributes(i)%name) .EQ. trim(att_name) )
THEN
1181 IF ( this_attribute > 0 )
THEN
1182 IF ( out_axis%attributes(this_attribute)%type .NE. nf90_char )
THEN
1186 IF ( fms_error_handler(
'diag_util_mod::prepend_attribute_axis',&
1187 &
'Attribute "'//trim(att_name)//
'" is not a character attribute.',&
1195 this_attribute = out_axis%num_attributes + 1
1196 IF ( this_attribute .GT. max_axis_attributes )
THEN
1201 IF ( fms_error_handler(
'diag_util_mod::prepend_attribute_axis',&
1202 &
'Number of attributes exceeds max_axis_attributes for attribute "'&
1203 &//trim(att_name)//
'". Increase diag_manager_nml:max_axis_attributes.',&
1208 out_axis%num_attributes = this_attribute
1210 out_axis%attributes(this_attribute)%name = att_name
1211 out_axis%attributes(this_attribute)%type = nf90_char
1213 out_axis%attributes(this_attribute)%catt =
''
1218 IF ( index(trim(out_axis%attributes(this_attribute)%catt), trim(prepend_value)).EQ.0 )
THEN
1220 length = len_trim(trim(prepend_value)//
" "//trim(out_axis%attributes(this_attribute)%catt))
1221 IF ( length.GT.len(out_axis%attributes(this_attribute)%catt) )
THEN
1225 IF ( fms_error_handler(
'diag_util_mod::prepend_attribute_file',&
1226 &
'Prepend length for attribute "'//trim(att_name)//
'" is longer than allowed.',&
1232 out_axis%attributes(this_attribute)%catt =&
1233 & trim(prepend_value)//
' '//trim(out_axis%attributes(this_attribute)%catt)
1234 out_axis%attributes(this_attribute)%len = length
1242 integer,
intent(in) :: id
1249 if (.not.
allocated(
axes(id)%attributes))
return
1250 do i = 1,
axes(id)%num_attributes
1251 if (trim(
axes(id)%attributes(i)%name)==
'compress')
then
1262 integer,
intent(in) :: id
1263 integer,
intent(out),
allocatable :: r(:)
1265 integer iatt, k, k1, k2, n
1268 character(*),
parameter :: tag =
'get_compressed_axes_ids'
1272 associate(axis=>
axes(id))
1273 if (.not.
allocated(axis%attributes))
call error_mesg(tag, &
1274 'attempt to get compression dimensions from axis "'//trim(axis%name)// &
1275 &
'" which is not compressed (does not have any attributes)', fatal)
1278 do k = 1,axis%num_attributes
1279 if (trim(axis%attributes(k)%name)==
'compress')
then
1284 if (iatt == 0)
call error_mesg(tag, &
1285 'attempt to get compression dimensions from axis "'//trim(axis%name)//&
1286 '" which is not compressed (does not have "compress" attributes).', fatal)
1287 if (axis%attributes(iatt)%type/=nf90_char)
call error_mesg(tag, &
1288 'attempt to get compression dimensions from axis "'//trim(axis%name)//&
1289 '" but the axis attribute "compress" has incorrect type.', fatal)
1294 do k = 1, len(axis%attributes(iatt)%catt)
1295 if (space.and.(axis%attributes(iatt)%catt(k:k)/=
' '))
then
1298 space = (axis%attributes(iatt)%catt(k:k)==
' ')
1307 do k1 = k2+1, len(axis%attributes(iatt)%catt)
1308 if (axis%attributes(iatt)%catt(k1:k1)/=
' ')
exit
1310 do k2 = k1+1, len(axis%attributes(iatt)%catt)
1311 if (axis%attributes(iatt)%catt(k2:k2)==
' ')
exit
1314 if (r(k)<=0)
call error_mesg(tag, &
1315 'compression dimension "'//trim(axis%attributes(iatt)%catt(k1:k2))//&
1316 '" not found among the axes of set "'//trim(
axis_sets(axis%set))//
'".', fatal)
1320 END MODULE diag_axis_mod
subroutine, public get_compressed_axes_ids(id, r)
given an index of compressed-by-gathering axis, return an array of axes used in compression....
integer(i4_kind) function, public axis_compatible_check(id, varname)
Checks if the axes are compatible.
integer function, public get_axis_length(id)
Return the length of the axis.
subroutine, public get_diag_axis(id, name, units, long_name, cart_name, direction, edges, Domain, DomainU, array_data, num_attributes, attributes, domain_position)
Return information about the axis with index ID.
subroutine, public get_diag_axis_cart(id, cart_name)
Return the axis cartesian.
subroutine valid_id_check(id, routine_name)
Check to see if the given axis id is a valid id. If the axis id is invalid, call a FATAL error....
integer function, public get_axis_num(axis_name, set_name)
Returns index into axis table corresponding to a given axis name.
subroutine attribute_init_axis(out_axis, err_msg)
Allocates memory in out_file for the attributes. Will FATAL if err_msg is not included in the subrout...
integer function, public get_axis_global_length(id)
Return the global length of the axis.
type(domain2d) function, public get_domain2d(ids)
Return the 2D domain for the axis IDs given.
integer function get_axis_set_num(set_name)
Returns index in axis set table corresponding to a given axis set name.
type(domainug) function, public get_domainug(id)
Retrun the 1D domain for the axis ID given.
subroutine, public get_diag_axis_data(id, axis_data)
Return the axis data.
character(len=128) function, public get_axis_reqfld(id)
Return the required field names for the axis.
subroutine diag_axis_attribute_init(diag_axis_id, name, type, cval, ival, rval)
integer function, public get_tile_count(ids)
Return the tile count for the axis.
subroutine, public get_diag_axis_name(id, axis_name)
Return the short name of the axis.
subroutine prepend_attribute_axis(out_axis, att_name, prepend_value, err_msg)
Prepends the attribute value to an already existing attribute. If the attribute isn't yet defined,...
character(len=128) function, public get_axis_aux(id)
Return the auxiliary name for the axis.
integer(i4_kind), parameter, public diag_axis_2ddomain
For unstructured grid support.
type(diag_axis_type), dimension(:), allocatable, save axes
global storage for all defined axes
integer function, public diag_axis_init(name, array_data, units, cart_name, long_name, direction, set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count, domain_position)
Initialize the axis, and return the axis ID.
logical function, public axis_is_compressed(id)
given an axis, returns TRUE if the axis uses compression-by-gathering: that is, if this is an axis fo...
subroutine, public get_axes_shift(ids, ishift, jshift)
Return the value of the shift for the axis IDs given.
integer, dimension(:), allocatable num_subaxes
counter of number of axes defined
integer(i4_kind), parameter, public diag_axis_ugdomain
For unstructured grid support.
subroutine, public get_diag_axis_domain_name(id, name)
Return the name of the axis' domain.
type(domain1d) function, public get_domain1d(id)
Retrun the 1D domain for the axis ID given.
integer function, public diag_subaxes_init(axis, subdata, start_indx, end_indx, domain_2d)
Create a subaxis on a parent axis.
character(len=128), dimension(:), allocatable, save axis_sets
storage for axis set names
Add an arbitrary attribute and value to the diagnostic axis.
integer max_axis_attributes
Maximum number of user definable attributes per axis.
logical use_modern_diag
Namelist flag to use the modernized diag_manager code.
integer max_axes
Maximum number of independent axes.
Attribute type for diagnostic fields.
Type to hold the diagnostic axis description.
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
logical function, public fms_error_handler(routine, message, err_msg)
Facilitates the control of fatal error conditions.
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
subroutine mpp_get_domain_components(domain, x, y, tile_count)
Retrieve 1D components of 2D decomposition.
character(len=name_length) function mpp_get_domain_name(domain)
Set user stack size.
These routines retrieve the axis specifications associated with the compute domains....
These routines retrieve the axis specifications associated with the global domains....
One dimensional domain used to manage shared data access between pes.
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.