28 MODULE diag_output_mod
31 use,
intrinsic :: iso_fortran_env, only: real128
32 use,
intrinsic :: iso_c_binding, only: c_double,c_float,c_int64_t, &
33 c_int32_t,c_int16_t,c_intptr_t
39 USE mpp_mod,
ONLY:
mpp_npes,
mpp_pe, mpp_root_pe, mpp_get_current_pelist
47 USE netcdf,
ONLY: nf90_int, nf90_float, nf90_char
49 use mpp_domains_mod,
only: mpp_get_ug_io_domain
50 use mpp_domains_mod,
only: mpp_get_ug_domain_npes
51 use mpp_domains_mod,
only: mpp_get_ug_domain_pelist
53 use mpp_mod,
only: uppercase,lowercase
62 TYPE(diag_global_att_type),
SAVE :: diag_global_att
64 INTEGER,
PARAMETER :: NETCDF1 = 1
65 INTEGER,
PARAMETER :: mxch = 128
66 INTEGER,
PARAMETER :: mxchl = 256
67 INTEGER :: current_file_unit = -1
68 INTEGER,
DIMENSION(2,2) :: max_range = reshape((/ -32767, 32767, -127, 127 /),(/2,2/))
69 INTEGER,
DIMENSION(2) :: missval = (/ -32768, -128 /)
71 INTEGER,
PARAMETER :: max_axis_num = 20
72 INTEGER :: num_axis_in_file = 0
73 INTEGER,
DIMENSION(max_axis_num) :: axis_in_file
74 LOGICAL,
DIMENSION(max_axis_num) :: time_axis_flag, edge_axis_flag
76 LOGICAL :: module_is_initialized = .false.
79 character(len=*),
parameter :: version =
'2020.03'
88 & domain, domainU, fileobj, fileobjU, fileobjND, fnum_domain, &
90 CHARACTER(len=*),
INTENT(in) :: file_name
91 CHARACTER(len=*),
INTENT(in) :: file_title
92 INTEGER ,
INTENT(out) :: file_unit
95 TYPE(
domain2d) ,
INTENT(in) :: domain
96 TYPE(
domainug) ,
INTENT(in) :: domainu
97 type(fmsnetcdfdomainfile_t),
intent(inout),
target :: fileobj
98 type(fmsnetcdfunstructureddomainfile_t),
intent(inout),
target :: fileobju
99 type(fmsnetcdffile_t),
intent(inout),
target :: fileobjnd
100 character(*),
intent(out) :: fnum_domain
104 TYPE(
diag_atttype),
INTENT(in),
OPTIONAL :: attributes(:)
106 class(fmsnetcdffile_t),
pointer :: fileob => null()
109 integer,
allocatable,
dimension(:) :: current_pelist
111 character(len=9) :: mype_string
112 character(len=FMS_FILE_LEN) :: filename_tile
116 IF ( .NOT.module_is_initialized )
THEN
117 module_is_initialized = .true.
122 if (domain .NE. null_domain2d .AND. domainu .NE. null_domainug)&
123 &
CALL error_mesg(
'diag_output_init',
"Domain2D and DomainUG can not be used at the same time in "//&
124 & trim(file_name), fatal)
127 IF ( domain .NE. null_domain2d )
THEN
131 if (.not.check_if_open(fileob))
call open_check(
open_file(fileobj, trim(file_name)//
".nc",
"overwrite", &
132 domain, is_restart=.false.))
138 write(mype_string,
'(I0.4)') mype
141 call get_mosaic_tile_file(file_name, filename_tile, .true., domain)
142 filename_tile = trim(filename_tile)//
"."//trim(mype_string)
144 if (.not.check_if_open(fileob))
then
145 call open_check(
open_file(fileobjnd, trim(filename_tile),
"overwrite", &
149 call register_global_attribute(fileobjnd,
"NumFilesInSet", 0)
152 if (file_unit < 0) file_unit = 10
154 ELSE IF (domainu .NE. null_domainug)
THEN
156 if (.not.check_if_open(fileob))
call open_check(
open_file(fileobju, trim(file_name)//
".nc",
"overwrite", &
157 domainu, is_restart=.false.))
162 allocate(current_pelist(
mpp_npes()))
163 call mpp_get_current_pelist(current_pelist)
164 if (.not.check_if_open(fileob))
then
165 call open_check(
open_file(fileobjnd, trim(file_name)//
".nc",
"overwrite", &
166 pelist=current_pelist, is_restart=.false.))
169 if (file_unit < 0) file_unit = 10
170 deallocate(current_pelist)
174 IF ( file_title(1:1) /=
' ' )
THEN
175 call register_global_attribute(fileob,
'title', trim(file_title), str_len=len_trim(file_title))
178 IF (
PRESENT(attributes) )
THEN
179 DO i=1,
SIZE(attributes)
180 SELECT CASE (attributes(i)%type)
182 call register_global_attribute(fileob, trim(attributes(i)%name), attributes(i)%iatt)
185 call register_global_attribute(fileob, trim(attributes(i)%name), attributes(i)%fatt)
188 call register_global_attribute(fileob, trim(attributes(i)%name), attributes(i)%catt, &
189 & str_len=len_trim(attributes(i)%catt))
195 CALL error_mesg(
'diag_output_mod::diag_output_init',
'Unknown attribute type for global attribute "'&
196 &//trim(attributes(i)%name)//
'" in file "'//trim(file_name)//
'". Contact the developers.', fatal)
203 call register_global_attribute(fileob,
'grid_type', trim(gatt%grid_type), str_len=len_trim(gatt%grid_type))
205 call register_global_attribute(fileob,
'grid_tile', trim(gatt%tile_name), str_len=len_trim(gatt%tile_name))
211 INTEGER,
INTENT(in) :: file_unit
212 INTEGER,
INTENT(in) :: axes(:)
213 class(fmsnetcdffile_t) ,
intent(inout) :: fileob
214 LOGICAL,
INTENT(in),
OPTIONAL :: time_ops
215 logical,
intent(inout) ,
optional :: time_axis_registered
221 CHARACTER(len=mxch) :: axis_name, axis_units, axis_name_current
222 CHARACTER(len=mxchl) :: axis_long_name
223 CHARACTER(len=1) :: axis_cart_name
224 INTEGER :: axis_direction, axis_edges
225 REAL,
ALLOCATABLE :: axis_data(:)
227 INTEGER :: num_attributes
228 TYPE(
diag_atttype),
DIMENSION(:),
ALLOCATABLE :: attributes
229 INTEGER :: calendar, id_axis, id_time_axis
230 INTEGER :: i, j, index, num, length, edges_index
233 CHARACTER(len=2048) :: err_msg
234 integer :: id_axis_current
235 logical :: is_time_axis_registered
236 integer :: istart, iend
237 integer :: gstart, cstart, cend
239 character(len=32) :: type_str
243 IF (
PRESENT(time_ops) )
THEN
248 if (
present(time_axis_registered))
then
249 is_time_axis_registered = time_axis_registered
251 is_time_axis_registered = .false.
254 IF ( num_axis_in_file == 0 ) current_file_unit = file_unit
259 IF ( num < 1 )
CALL error_mesg(
'write_axis_meta_data',
'number of axes < 1.', fatal)
262 IF ( file_unit /= current_file_unit )
CALL error_mesg(
'write_axis_meta_data',&
263 &
'writing meta data out-of-order to different files.', fatal)
278 IF ( index > 0 ) cycle
281 num_axis_in_file = num_axis_in_file + 1
282 axis_in_file(num_axis_in_file) = id_axis
283 edge_axis_flag(num_axis_in_file) = .false.
285 ALLOCATE(axis_data(length))
287 CALL get_diag_axis(id_axis, axis_name, axis_units, axis_long_name,&
288 & axis_cart_name, axis_direction, axis_edges, domain, domainu, axis_data,&
289 & num_attributes, attributes, domain_position=axis_pos)
291 IF ( domain .NE. null_domain1d )
THEN
293 type is (fmsnetcdffile_t)
298 iend = cend - gstart + 1
299 istart = cstart - gstart + 1
300 call register_axis(fileob, axis_name, dimension_length=clength)
301 call register_field(fileob, axis_name, type_str, (/axis_name/) )
305 call register_variable_attribute(fileob, axis_name,
"domain_decomposition", &
306 (/gstart, gend, cstart, cend/))
307 type is (fmsnetcdfdomainfile_t)
311 call register_axis(fileob, axis_name, lowercase(trim(axis_cart_name)), domain_position=axis_pos )
312 call get_global_io_domain_indices(fileob, trim(axis_name), istart, iend)
313 call register_field(fileob, axis_name, type_str, (/axis_name/) )
316 ELSE IF ( domainu .NE. null_domainug)
THEN
318 type is (fmsnetcdfunstructureddomainfile_t)
322 call register_axis(fileob, axis_name )
324 call register_field(fileob, axis_name, type_str, (/axis_name/) )
325 istart = lbound(axis_data,1)
326 iend = ubound(axis_data,1)
329 call register_axis(fileob, axis_name, dimension_length=
size(axis_data))
330 call register_field(fileob, axis_name, type_str, (/axis_name/) )
331 istart = lbound(axis_data,1)
332 iend = ubound(axis_data,1)
335 if (length <= 0)
then
341 is_time_axis_registered = variable_exists(fileob,trim(axis_name),.true.)
342 if (.not. is_time_axis_registered)
then
343 call register_axis(fileob, trim(axis_name), unlimited )
344 call register_field(fileob, axis_name, type_str, (/axis_name/) )
345 is_time_axis_registered = .true.
346 if (
present(time_axis_registered)) time_axis_registered = is_time_axis_registered
351 if(trim(axis_units) .ne.
"none")
call register_variable_attribute(fileob, axis_name,
"units", &
352 & trim(axis_units), str_len=len_trim(axis_units))
353 call register_variable_attribute(fileob, axis_name,
"long_name", trim(axis_long_name), &
354 & str_len=len_trim(axis_long_name))
355 if(trim(axis_cart_name).ne.
"N")
call register_variable_attribute(fileob, axis_name,
"axis", &
356 & trim(axis_cart_name), str_len=len_trim(axis_cart_name))
358 if (length > 0 )
then
360 select case (axis_direction)
362 call register_variable_attribute(fileob, axis_name,
"positive",
"up", str_len=len_trim(
"up"))
364 call register_variable_attribute(fileob, axis_name,
"positive",
"down", str_len=len_trim(
"down"))
366 call write_data(fileob, axis_name, axis_data(istart:iend) )
370 CALL write_attribute_meta(file_unit, num_attributes, attributes, err_msg, varname=axis_name, fileob=fileob)
371 IF ( len_trim(err_msg) .GT. 0 )
THEN
372 CALL error_mesg(
'diag_output_mod::write_axis_meta_data', trim(err_msg), fatal)
378 IF ( axis_cart_name ==
'T' )
THEN
379 time_axis_flag(num_axis_in_file) = .true.
380 id_time_axis = num_axis_in_file
381 calendar = get_calendar_type()
384 call register_variable_attribute(fileob, axis_name,
"calendar_type", &
385 uppercase(trim(valid_calendar_types(calendar))), &
386 & str_len=len_trim(valid_calendar_types(calendar)) )
387 call register_variable_attribute(fileob, axis_name,
"calendar", &
388 lowercase(trim(valid_calendar_types(calendar))), &
389 & str_len=len_trim(valid_calendar_types(calendar)) )
390 IF ( time_ops1 )
THEN
391 call register_variable_attribute(fileob, axis_name,
'bounds', trim(axis_name)//
'_bnds', &
392 & str_len=len_trim(trim(axis_name)//
'_bnds'))
394 call set_fileobj_time_name(fileob, axis_name)
396 time_axis_flag(num_axis_in_file) = .false.
399 DEALLOCATE(axis_data)
402 IF (
ALLOCATED(attributes) )
THEN
403 DO j=1, num_attributes
404 IF (
allocated(attributes(j)%fatt ) )
THEN
405 DEALLOCATE(attributes(j)%fatt)
407 IF (
allocated(attributes(j)%iatt ) )
THEN
408 DEALLOCATE(attributes(j)%iatt)
411 DEALLOCATE(attributes)
417 IF ( axis_edges <= 0 ) cycle
420 id_axis_current = id_axis
421 axis_name_current = axis_name
424 IF ( edges_index > 0 ) cycle
427 length = get_axis_global_length( id_axis )
428 ALLOCATE(axis_data(length))
429 CALL get_diag_axis(id_axis, axis_name, axis_units, axis_long_name, axis_cart_name,&
430 & axis_direction, axis_edges, domain, domainu, axis_data)
433 call register_variable_attribute(fileob, axis_name_current,
"edges",trim(axis_name), &
434 & str_len=len_trim(axis_name))
437 num_axis_in_file = num_axis_in_file + 1
438 axis_in_file(num_axis_in_file) = id_axis
439 edge_axis_flag(num_axis_in_file) = .true.
440 time_axis_flag(num_axis_in_file) = .false.
443 if (.not.variable_exists(fileob, axis_name))
then
444 call register_axis(fileob, axis_name,
size(axis_data) )
445 call register_field(fileob, axis_name, type_str, (/axis_name/) )
446 if(trim(axis_units) .ne.
"none")
call register_variable_attribute(fileob, axis_name,
"units", &
447 & trim(axis_units), str_len=len_trim(axis_units))
448 call register_variable_attribute(fileob, axis_name,
"long_name", trim(axis_long_name), &
449 & str_len=len_trim(axis_long_name))
450 if(trim(axis_cart_name).ne.
"N")
call register_variable_attribute(fileob, axis_name,
"axis", &
451 & trim(axis_cart_name), str_len=len_trim(axis_cart_name))
452 select case (axis_direction)
454 call register_variable_attribute(fileob, axis_name,
"positive",
"up", str_len=len_trim(
"up"))
456 call register_variable_attribute(fileob, axis_name,
"positive",
"down", str_len=len_trim(
"down"))
458 call write_data(fileob, axis_name, axis_data)
461 DEALLOCATE (axis_data)
469 & avg_name, time_method, standard_name, interp_method, attributes, num_attributes, &
470 & use_UGdomain, fileob)
result ( Field )
471 INTEGER,
INTENT(in) :: file_unit
472 INTEGER,
INTENT(in) :: axes(:)
473 CHARACTER(len=*),
INTENT(in) :: name
474 CHARACTER(len=*),
INTENT(in) :: units
475 CHARACTER(len=*),
INTENT(in) :: long_name
476 REAL,
OPTIONAL,
INTENT(in) :: range(2)
477 REAL,
OPTIONAL,
INTENT(in) :: mval
478 INTEGER,
OPTIONAL,
INTENT(in) :: pack
485 CHARACTER(len=*),
OPTIONAL,
INTENT(in) :: avg_name
486 CHARACTER(len=*),
OPTIONAL,
INTENT(in) :: time_method
488 CHARACTER(len=*),
OPTIONAL,
INTENT(in) :: standard_name
489 CHARACTER(len=*),
OPTIONAL,
INTENT(in) :: interp_method
490 TYPE(diag_atttype),
DIMENSION(:),
allocatable,
OPTIONAL,
INTENT(in) :: attributes
491 INTEGER,
OPTIONAL,
INTENT(in) :: num_attributes
492 LOGICAL,
OPTIONAL,
INTENT(in) :: use_ugdomain
493 class(fmsnetcdffile_t),
intent(inout) :: fileob
495 logical :: is_time_bounds
496 CHARACTER(len=256) :: standard_name2
497 TYPE(diag_fieldtype) :: field
498 LOGICAL :: coord_present
499 CHARACTER(len=128) :: aux_axes(size(axes))
500 CHARACTER(len=160) :: coord_att
501 CHARACTER(len=1024) :: err_msg
503 character(len=128),
dimension(size(axes)) :: axis_names
505 INTEGER :: i, indexx, num, ipack, np
507 INTEGER :: axis_indices(size(axes))
508 logical :: use_ugdomain_local
513 coord_present = .false.
514 IF(
PRESENT(standard_name) )
THEN
515 standard_name2 = standard_name
517 standard_name2 =
'none'
520 use_ugdomain_local = .false.
521 if(
present(use_ugdomain)) use_ugdomain_local = use_ugdomain
525 IF ( num < 1 )
CALL error_mesg (
'write_meta_data',
'number of axes < 1', fatal)
527 IF ( file_unit /= current_file_unit )
CALL error_mesg (
'write_meta_data', &
528 &
'writing meta data out-of-order to different files', fatal)
530 IF (trim(name) .eq.
"time_bnds")
then
531 is_time_bounds = .true.
533 is_time_bounds = .false.
541 IF ( indexx > 0 )
THEN
542 axis_indices(i) = indexx
545 CALL error_mesg (
'write_field_meta_data',&
546 &
'axis data not written for field '//trim(name), fatal)
549 call get_diag_axis_name(axes(i),axis_names(i))
553 IF ( num >= 2 .OR. (num==1 .and. use_ugdomain_local) )
THEN
556 aux_axes(i) = get_axis_aux(axes(i))
557 IF( trim(aux_axes(i)) /=
'none' )
THEN
558 IF(len_trim(coord_att) == 0)
THEN
559 coord_att = trim(aux_axes(i))
561 coord_att = trim(coord_att)//
' '//trim(aux_axes(i))
563 coord_present = .true.
572 IF (
PRESENT(pack) )
THEN
582 IF (
PRESENT(range) )
THEN
583 IF ( range(2) > range(1) )
THEN
586 IF ( ipack > 2 )
THEN
588 add = 0.5*(range(1)+range(2))
589 scale = (range(2)-range(1)) / real(max_range(2,np)-max_range(1,np))
595 IF (
PRESENT(mval) )
THEN
597 field%miss_present = .true.
598 IF ( ipack > 2 )
THEN
600 field%miss_pack = real(missval(np))*scale+add
601 field%miss_pack_present = .true.
603 field%miss_pack = mval
604 field%miss_pack_present = .false.
607 field%miss_present = .false.
608 field%miss_pack_present = .false.
612 field%fieldname = name
614 if (.not. variable_exists(fileob,name))
then
622 call register_field(fileob,name,
"double",axis_names)
625 if (.not. is_time_bounds)
then
626 IF ( field%miss_present )
THEN
627 call register_variable_attribute(fileob,name,
"_FillValue",real(field%miss_pack,8))
628 call register_variable_attribute(fileob,name,
"missing_value",real(field%miss_pack,8))
630 call register_variable_attribute(fileob,name,
"_FillValue",real(cmor_missing_value,8))
631 call register_variable_attribute(fileob,name,
"missing_value",real(cmor_missing_value,8))
633 IF ( use_range )
then
634 call register_variable_attribute(fileob,name,
"valid_range", real(range,8))
638 call register_field(fileob,name,
"float",axis_names)
641 if (.not. is_time_bounds)
then
642 IF ( field%miss_present )
THEN
643 call register_variable_attribute(fileob,name,
"_FillValue",real(field%miss_pack,4))
644 call register_variable_attribute(fileob,name,
"missing_value",real(field%miss_pack,4))
646 call register_variable_attribute(fileob,name,
"_FillValue",real(cmor_missing_value,4))
647 call register_variable_attribute(fileob,name,
"missing_value",real(cmor_missing_value,4))
649 IF ( use_range )
then
650 call register_variable_attribute(fileob,name,
"valid_range", real(range,4))
654 CALL error_mesg(
'diag_output_mod::write_field_meta_data',&
655 &
"Pack values must be 1 or 2. Contact the developers.", fatal)
657 if (trim(units) .ne.
"none")
call register_variable_attribute(fileob,name,
"units",trim(units), &
658 & str_len=len_trim(units))
659 call register_variable_attribute(fileob,name,
"long_name",long_name, str_len=len_trim(long_name))
660 IF (
present(time_method) )
then
661 call register_variable_attribute(fileob,name,
'cell_methods',
'time: '//trim(time_method), &
662 & str_len=len_trim(
'time: '//trim(time_method)))
666 IF (
PRESENT(num_attributes) )
THEN
667 IF (
PRESENT(attributes) )
THEN
668 IF ( num_attributes .GT. 0 .AND.
allocated(attributes) )
THEN
670 & fileob=fileob, varname=name)
671 IF ( len_trim(err_msg) .GT. 0 )
THEN
672 CALL error_mesg(
'diag_output_mod::write_field_meta_data',&
673 & trim(err_msg)//
" Contact the developers.", fatal)
677 IF ( num_attributes .GT. 0 .AND. .NOT.
allocated(attributes) )
THEN
678 CALL error_mesg(
'diag_output_mod::write_field_meta_data',&
679 &
'num_attributes > 0 but attributes is not allocated for attribute '&
680 &//trim(attributes(i)%name)//
' for field '//trim(name)//
'. Contact the developers.', fatal)
681 ELSE IF ( num_attributes .EQ. 0 .AND.
allocated(attributes) )
THEN
682 CALL error_mesg(
'diag_output_mod::write_field_meta_data',&
683 &
'num_attributes == 0 but attributes is allocated for attribute '&
684 &//trim(attributes(i)%name)//
' for field '//trim(name)//
'. Contact the developers.', fatal)
689 CALL error_mesg(
'diag_output_mod::write_field_meta_data',&
690 &
'num_attributes present but attributes missing for attribute '&
691 &//trim(attributes(i)%name)//
' for field '//trim(name)//
'. Contact the developers.', fatal)
693 ELSE IF (
PRESENT(attributes) )
THEN
694 CALL error_mesg(
'diag_output_mod::write_field_meta_data',&
695 &
'attributes present but num_attributes missing for attribute '&
696 &//trim(attributes(i)%name)//
' for field '//trim(name)//
'. Contact the developers.', fatal)
700 IF (
PRESENT(avg_name) )
THEN
701 IF ( avg_name(1:1) /=
' ' )
THEN
702 call register_variable_attribute(fileob,name,
'time_avg_info',&
703 & trim(avg_name)//
'_T1,'//trim(avg_name)//
'_T2,'//trim(avg_name)//
'_DT', &
704 & str_len=len_trim(trim(avg_name)//
'_T1,'//trim(avg_name)//
'_T2,'//trim(avg_name)//
'_DT'))
709 IF ( coord_present )
then
710 call register_variable_attribute(fileob,name,
'coordinates',trim(coord_att), str_len=len_trim(coord_att))
712 IF ( trim(standard_name2) /=
'none' )
then
713 call register_variable_attribute(fileob,name,
'standard_name',trim(standard_name2), &
714 & str_len=len_trim(standard_name2))
717 IF(
PRESENT(interp_method) )
THEN
718 call register_variable_attribute(fileob,name,
'interp_method', trim(interp_method), &
719 & str_len=len_trim(interp_method))
723 field%Domain = get_domain2d( axes )
724 field%tile_count = get_tile_count( axes )
725 field%DomainU = get_domainug( axes(1) )
733 INTEGER,
INTENT(in) :: file_unit
734 INTEGER,
INTENT(in) :: num_attributes
735 TYPE(diag_atttype),
DIMENSION(:),
INTENT(in) :: attributes
736 CHARACTER(len=*),
INTENT(in),
OPTIONAL :: time_method
737 CHARACTER(len=*),
INTENT(out),
OPTIONAL :: err_msg
738 CHARACTER(len=*),
INTENT(IN),
OPTIONAL :: varname
739 class(fmsnetcdffile_t),
intent(inout) :: fileob
741 INTEGER :: i, att_len
742 CHARACTER(len=1280) :: att_str
745 IF (
PRESENT(err_msg) ) err_msg =
''
747 DO i = 1, num_attributes
748 SELECT CASE (attributes(i)%type)
750 IF ( .NOT.
allocated(attributes(i)%iatt) )
THEN
751 IF ( fms_error_handler(
'diag_output_mod::write_attribute_meta',&
752 &
'Integer attribute type indicated, but array not allocated for attribute '&
753 &//trim(attributes(i)%name)//
'.', err_msg) )
THEN
757 if (
present(varname))
call register_variable_attribute(fileob, varname,trim(attributes(i)%name), &
758 & attributes(i)%iatt)
760 IF ( .NOT.
allocated(attributes(i)%fatt) )
THEN
761 IF ( fms_error_handler(
'diag_output_mod::write_attribute_meta',&
762 &
'Real attribute type indicated, but array not allocated for attribute '&
763 &//trim(attributes(i)%name)//
'.', err_msg) )
THEN
767 if (
present(varname))
call register_variable_attribute(fileob, varname,trim(attributes(i)%name), &
768 & real(attributes(i)%fatt,4) )
770 att_str = attributes(i)%catt
771 att_len = attributes(i)%len
772 IF ( trim(attributes(i)%name).EQ.
'cell_methods' .AND.
PRESENT(time_method) )
THEN
774 att_str = attributes(i)%catt(1:attributes(i)%len)//
' time: '//time_method
775 att_len = len_trim(att_str)
777 if (
present(varname))&
778 call register_variable_attribute(fileob, varname,trim(attributes(i)%name) , att_str(1:att_len), &
782 IF ( fms_error_handler(
'diag_output_mod::write_attribute_meta',
'Invalid type for attribute '&
783 &//trim(attributes(i)%name)//
'.', err_msg) )
THEN
795 INTEGER,
INTENT(in) :: file_unit
802 subroutine diag_field_write (varname, buffer, static, file_num, fileobjU, fileobj, fileobjND, &
803 & fnum_for_domain, time_in)
804 CHARACTER(len=*),
INTENT(in) :: varname
805 REAL ,
INTENT(inout) :: buffer(:,:,:,:)
806 logical,
intent(in) :: static
807 integer,
intent(in) :: file_num
808 type(fmsnetcdfunstructureddomainfile_t),
intent(inout) :: fileobju(:)
809 type(fmsnetcdfdomainfile_t),
intent(inout) :: fileobj(:)
810 type(fmsnetcdffile_t),
intent(inout) :: fileobjnd(:)
811 character(len=2),
intent(in) :: fnum_for_domain
815 INTEGER,
OPTIONAL,
INTENT(in) :: time_in
818 real,
allocatable :: local_buffer(:,:,:,:)
823 elseif (
present(time_in))
then
829 if (
size(buffer,3) .eq. 1 .and.
size(buffer,2) .eq. 1)
then
832 allocate(local_buffer(
size(buffer,1),
size(buffer,4),
size(buffer,2),
size(buffer,3)))
833 local_buffer(:,:,1,1) = buffer(:,1,1,:)
834 else if (
size(buffer,3) .eq. 1)
then
837 allocate(local_buffer(
size(buffer,1),
size(buffer,2),
size(buffer,4),
size(buffer,3)))
838 local_buffer(:,:,:,1) = buffer(:,:,1,:)
840 allocate(local_buffer(
size(buffer,1),
size(buffer,2),
size(buffer,3),
size(buffer,4)))
841 local_buffer = buffer(:,:,:,:)
845 if (fnum_for_domain ==
"2d" )
then
846 if (check_if_open(fileobj(file_num)))
then
847 call write_data (fileobj(file_num), trim(varname), local_buffer, unlim_dim_level=time )
849 elseif (fnum_for_domain ==
"nd")
then
850 if (check_if_open(fileobjnd(file_num)) )
then
851 call write_data (fileobjnd(file_num), trim(varname), local_buffer, unlim_dim_level=time)
853 elseif (fnum_for_domain ==
"ug")
then
854 if (check_if_open(fileobju(file_num)))
then
855 call write_data (fileobju(file_num), trim(varname), local_buffer, unlim_dim_level=time)
858 call error_mesg(
"diag_field_write",
"fnum_for_domain must be '2d', 'nd', or 'ug'",fatal)
861 deallocate(local_buffer)
866 class(fmsnetcdffile_t),
intent(inout) :: fileob
867 real,
intent(in) :: rtime_value
868 integer,
intent(in) :: time_index
869 character(len=*),
intent(in),
optional :: time_name
870 character(len=:),
allocatable :: name_time
873 if (
present(time_name))
then
874 allocate(
character(len=len(time_name)) :: name_time)
875 name_time = time_name
877 allocate(
character(len=4) :: name_time)
881 call write_data (fileob, trim(name_time), rtime_value, unlim_dim_level=time_index)
883 if (
allocated(name_time))
deallocate(name_time)
889 INTEGER,
INTENT(in) :: num
898 DO i = 1, num_axis_in_file
899 IF ( num == axis_in_file(i) )
THEN
908 TYPE(diag_global_att_type),
INTENT(out) :: gatt
915 CHARACTER(len=*),
INTENT(in) :: component, gridtype, tilename
919 CHARACTER(len=64) :: component_tmp
920 component_tmp = component
924 diag_global_att%grid_type = gridtype
925 diag_global_att%tile_name = tilename
930 subroutine diag_flush(file_num, fileobjU, fileobj, fileobjND, fnum_for_domain)
931 integer,
intent(in) :: file_num
932 type(fmsnetcdfunstructureddomainfile_t),
intent(inout) :: fileobju(:)
933 type(fmsnetcdfdomainfile_t),
intent(inout) :: fileobj(:)
934 type(fmsnetcdffile_t),
intent(inout) :: fileobjnd(:)
935 character(len=2),
intent(in) :: fnum_for_domain
939 if (fnum_for_domain ==
"2d" )
then
940 call flush_file (fileobj(file_num))
941 elseif (fnum_for_domain ==
"nd")
then
942 call flush_file (fileobjnd(file_num))
943 elseif (fnum_for_domain ==
"ug")
then
944 call flush_file (fileobju(file_num))
946 call error_mesg(
"diag_field_write",
"No file object is associated with this file number",fatal)
949 END MODULE diag_output_mod
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.
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.
type(domainug) function, public get_domainug(id)
Retrun the 1D domain for the axis ID given.
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.
character(len=128) function, public get_axis_aux(id)
Return the auxiliary name for the axis.
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.
type(domain1d) function, public get_domain1d(id)
Retrun the 1D domain for the axis ID given.
integer pack_size
1 for double and 2 for float
real(r8_kind), parameter cmor_missing_value
CMOR standard missing value.
Attribute type for diagnostic fields.
subroutine, public diag_field_write(varname, buffer, static, file_num, fileobjU, fileobj, fileobjND, fnum_for_domain, time_in)
Writes diagnostic data out using fms2_io routine.
subroutine, public diag_output_init(file_name, file_title, file_unit, domain, domainU, fileobj, fileobjU, fileobjND, fnum_domain, attributes)
Opens the output file.
subroutine, public get_diag_global_att(gAtt)
Return the global attribute type.
subroutine, public set_diag_global_att(component, gridType, tileName)
Set the global attribute type.
subroutine, public diag_flush(file_num, fileobjU, fileobj, fileobjND, fnum_for_domain)
Flushes the file into disk.
integer function get_axis_index(num)
Return the axis index number.
subroutine write_attribute_meta(file_unit, num_attributes, attributes, time_method, err_msg, varname, fileob)
Write out attribute meta data to file.
subroutine, public done_meta_data(file_unit)
Writes axis data to file.
subroutine, public diag_write_time(fileob, rtime_value, time_index, time_name)
Writes the time data to the history file.
subroutine, public write_axis_meta_data(file_unit, axes, fileob, time_ops, time_axis_registered)
Write the axis meta data to file.
type(diag_fieldtype) function, public write_field_meta_data(file_unit, name, axes, units, long_name, range, pack, mval, avg_name, time_method, standard_name, interp_method, attributes, num_attributes, use_UGdomain, fileob)
Write the field meta data to file.
Opens a given netcdf or domain file.
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....
type(domain2d) function, pointer mpp_get_io_domain(domain)
Set user stack size.
Set up a domain decomposition.
These routines retrieve the axis specifications associated with the compute domains....
Retrieve the entire array of compute domain extents associated with a decomposition.
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...
Retrieve list of PEs associated with a domain decomposition. The 1D version of this call returns an a...
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 mpp_npes()
Returns processor count for current pelist.
integer function mpp_pe()
Returns processor ID.
Gather data sent from pelist onto the root pe Wrapper for MPI_gather, can be used with and without in...
character(len=24) function, public valid_calendar_types(ncal, err_msg)
Returns a character string that describes the calendar type corresponding to the input integer.
integer function, public get_calendar_type()
Returns default calendar type for mapping from time to date.