FMS 2025.01-dev
Flexible Modeling System
Loading...
Searching...
No Matches
diag_output.F90
1!***********************************************************************
2!* GNU Lesser General Public License
3!*
4!* This file is part of the GFDL Flexible Modeling System (FMS).
5!*
6!* FMS is free software: you can redistribute it and/or modify it under
7!* the terms of the GNU Lesser General Public License as published by
8!* the Free Software Foundation, either version 3 of the License, or (at
9!* your option) any later version.
10!*
11!* FMS is distributed in the hope that it will be useful, but WITHOUT
12!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14!* for more details.
15!*
16!* You should have received a copy of the GNU Lesser General Public
17!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18!***********************************************************************
19!> @defgroup diag_output_mod diag_output_mod
20!> @ingroup diag_manager
21!! @brief diag_output_mod is an integral part of
22!! diag_manager_mod. Its function is to write axis-meta-data,
23!! field-meta-data and field data.
24!! @author Seth Underwood
25
26!> @addtogroup diag_output_mod
27!> @{
28MODULE diag_output_mod
29
30use platform_mod
31use,intrinsic :: iso_fortran_env, only: real128
32use,intrinsic :: iso_c_binding, only: c_double,c_float,c_int64_t, &
33 c_int32_t,c_int16_t,c_intptr_t
34 USE mpp_domains_mod, ONLY: domain1d, domain2d, mpp_define_domains, mpp_get_pelist,&
35 & mpp_get_global_domain, mpp_get_compute_domains, null_domain1d, null_domain2d,&
36 & domainug, null_domainug, center, east, north, mpp_get_compute_domain,&
37 & OPERATOR(.NE.), mpp_get_layout, OPERATOR(.EQ.), mpp_get_io_domain, &
39 USE mpp_mod, ONLY: mpp_npes, mpp_pe, mpp_root_pe, mpp_get_current_pelist
40 USE diag_axis_mod, ONLY: diag_axis_init, get_diag_axis, get_axis_length,&
44 USE time_manager_mod, ONLY: get_calendar_type, valid_calendar_types
45 USE fms_mod, ONLY: error_mesg, mpp_pe, write_version_number, fms_error_handler, fatal, note
46
47 USE netcdf, ONLY: nf90_int, nf90_float, nf90_char
48
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
52 use mpp_mod, only: mpp_gather
53 use mpp_mod, only: uppercase,lowercase
54 use fms2_io_mod
55
56 IMPLICIT NONE
57
58 PRIVATE
62 TYPE(diag_global_att_type), SAVE :: diag_global_att
63
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 /)
70
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
75
76 LOGICAL :: module_is_initialized = .false.
77
78 ! Include variable "version" to be written to log file.
79 character(len=*), parameter :: version = '2020.03'
80 !> @}
81
82!> @addtogroup diag_output_mod
83!> @{
84CONTAINS
85
86 !> @brief Opens the output file.
87 SUBROUTINE diag_output_init (file_name, file_title, file_unit,&
88 & domain, domainU, fileobj, fileobjU, fileobjND, fnum_domain, &
89 & attributes)
90 CHARACTER(len=*), INTENT(in) :: file_name !< Output file name
91 CHARACTER(len=*), INTENT(in) :: file_title !< Descriptive title for the file
92 INTEGER , INTENT(out) :: file_unit !< File unit number assigned to the output file.
93 !! Needed for subsuquent calls to
94 !! diag_output_mod
95 TYPE(domain2d) , INTENT(in) :: domain !< Domain associated with file, if domain decomposed
96 TYPE(domainug) , INTENT(in) :: domainu !< The unstructure domain
97 type(fmsnetcdfdomainfile_t),intent(inout),target :: fileobj !< Domain decomposed fileobj
98 type(fmsnetcdfunstructureddomainfile_t),intent(inout),target :: fileobju !< Unstructured domain fileobj
99 type(fmsnetcdffile_t),intent(inout),target :: fileobjnd !< Non domain decomposed fileobj
100 character(*),intent(out) :: fnum_domain !< String indicating the type of fileobj was used:
101 !! "2d" domain decomposed
102 !! "ug" unstrucuted domain decomposed
103 !! "nd" no domain
104 TYPE(diag_atttype), INTENT(in), OPTIONAL :: attributes(:) !< Array of global attributes to be written to file
105
106 class(fmsnetcdffile_t), pointer :: fileob => null()
107 integer :: i !< For looping through number of attributes
108 TYPE(diag_global_att_type) :: gatt
109 integer, allocatable, dimension(:) :: current_pelist
110 integer :: mype !< The pe you are on
111 character(len=9) :: mype_string !< a string to store the pe
112 character(len=FMS_FILE_LEN) :: filename_tile !< Filename with the tile number included
113 !! It is needed for subregional diagnostics
114
115 !---- initialize mpp_io ----
116 IF ( .NOT.module_is_initialized ) THEN
117 module_is_initialized = .true.
118 CALL write_version_number("DIAG_OUTPUT_MOD", version)
119 END IF
120
121!> Checks to make sure that only domain2D or domainUG is used. If both are not null, then FATAL
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)
125
126 !---- open output file (return file_unit id) -----
127 IF ( domain .NE. null_domain2d ) THEN
128 !> Check if there is an io_domain
129 iF ( associated(mpp_get_io_domain(domain)) ) then
130 fileob => fileobj
131 if (.not.check_if_open(fileob)) call open_check(open_file(fileobj, trim(file_name)//".nc", "overwrite", &
132 domain, is_restart=.false.))
133 fnum_domain = "2d" ! 2d domain
134 file_unit = 2
135 elSE !< No io domain, so every core is going to write its own file.
136 fileob => fileobjnd
137 mype = mpp_pe()
138 write(mype_string,'(I0.4)') mype
139 !! Add the tile number to the subregional file
140 !! This is needed for the combiner to work correctly
141 call get_mosaic_tile_file(file_name, filename_tile, .true., domain)
142 filename_tile = trim(filename_tile)//"."//trim(mype_string)
143
144 if (.not.check_if_open(fileob)) then
145 call open_check(open_file(fileobjnd, trim(filename_tile), "overwrite", &
146 is_restart=.false.))
147 !< For regional subaxis add the NumFilesInSet attribute, which is added by fms2_io for (other)
148 !< domains with sufficient decomposition info. Note mppnccombine will work with an entry of zero.
149 call register_global_attribute(fileobjnd, "NumFilesInSet", 0)
150 endif
151 fnum_domain = "nd" ! no domain
152 if (file_unit < 0) file_unit = 10
153 endiF
154 ELSE IF (domainu .NE. null_domainug) THEN
155 fileob => fileobju
156 if (.not.check_if_open(fileob)) call open_check(open_file(fileobju, trim(file_name)//".nc", "overwrite", &
157 domainu, is_restart=.false.))
158 fnum_domain = "ug" ! unstructured grid
159 file_unit=3
160 ELSE
161 fileob => fileobjnd
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.))
167 endif
168 fnum_domain = "nd" ! no domain
169 if (file_unit < 0) file_unit = 10
170 deallocate(current_pelist)
171 END IF
172
173 !---- write global attributes ----
174 IF ( file_title(1:1) /= ' ' ) THEN
175 call register_global_attribute(fileob, 'title', trim(file_title), str_len=len_trim(file_title))
176 END IF
177
178 IF ( PRESENT(attributes) ) THEN
179 DO i=1, SIZE(attributes)
180 SELECT CASE (attributes(i)%type)
181 CASE (nf90_int)
182 call register_global_attribute(fileob, trim(attributes(i)%name), attributes(i)%iatt)
183 CASE (nf90_float)
184
185 call register_global_attribute(fileob, trim(attributes(i)%name), attributes(i)%fatt)
186 CASE (nf90_char)
187
188 call register_global_attribute(fileob, trim(attributes(i)%name), attributes(i)%catt, &
189 & str_len=len_trim(attributes(i)%catt))
190 CASE default
191 ! <ERROR STATUS="FATAL">
192 ! Unknown attribute type for attribute <name> to module/input_field <module_name>/<field_name>.
193 ! Contact the developers.
194 ! </ERROR>
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)
197 END SELECT
198 END DO
199 END IF
200 !---- write grid type (mosaic or regular)
201 CALL get_diag_global_att(gatt)
202
203 call register_global_attribute(fileob, 'grid_type', trim(gatt%grid_type), str_len=len_trim(gatt%grid_type))
204
205 call register_global_attribute(fileob, 'grid_tile', trim(gatt%tile_name), str_len=len_trim(gatt%tile_name))
206
207 END SUBROUTINE diag_output_init
208
209 !> @brief Write the axis meta data to file.
210 SUBROUTINE write_axis_meta_data(file_unit, axes, fileob, time_ops, time_axis_registered)
211 INTEGER, INTENT(in) :: file_unit !< File unit number
212 INTEGER, INTENT(in) :: axes(:) !< Array of axis ID's, including the time axis
213 class(fmsnetcdffile_t) , intent(inout) :: fileob !< FMS2_io fileobj
214 LOGICAL, INTENT(in), OPTIONAL :: time_ops !< .TRUE. if this file contains any min, max, time_rms, or time_average
215 logical, intent(inout) , optional :: time_axis_registered !< .TRUE. if the time axis was already
216 !! written to the file
217
218 TYPE(domain1d) :: domain
219 TYPE(domainug) :: domainu
220
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(:)
226 integer :: axis_pos
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
231 INTEGER :: gend !< End index of global io_domain
232 LOGICAL :: time_ops1
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 !< Start and end of global and compute domains
238 integer :: clength !< Length of compute domain
239 character(len=32) :: type_str !< Str indicating the type of the axis data
240
241 ! Make sure err_msg is initialized
242 err_msg = ''
243 IF ( PRESENT(time_ops) ) THEN
244 time_ops1 = time_ops
245 ELSE
246 time_ops1 = .false.
247 END IF
248 if (present(time_axis_registered)) then
249 is_time_axis_registered = time_axis_registered
250 else
251 is_time_axis_registered = .false.
252 endif
253 !---- save the current file_unit ----
254 IF ( num_axis_in_file == 0 ) current_file_unit = file_unit
255
256 !---- dummy checks ----
257 num = SIZE(axes(:))
258 ! <ERROR STATUS="FATAL">number of axes < 1 </ERROR>
259 IF ( num < 1 ) CALL error_mesg('write_axis_meta_data', 'number of axes < 1.', fatal)
260
261 ! <ERROR STATUS="FATAL">writing meta data out-of-order to different files.</ERROR>
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)
264
265 IF (pack_size .eq. 1) then
266 type_str = "double"
267 ELSE IF (pack_size .eq. 2) then
268 type_str = "float"
269 ENDIF
270
271 !---- check all axes ----
272 !---- write axis meta data for new axes ----
273 DO i = 1, num
274 id_axis = axes(i)
275 index = get_axis_index( id_axis )
276
277 !---- skip axes already written -----
278 IF ( index > 0 ) cycle
279
280 !---- create new axistype (then point to) -----
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.
284 length = get_axis_global_length(id_axis)
285 ALLOCATE(axis_data(length))
286
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)
290
291 IF ( domain .NE. null_domain1d ) THEN
292 select type (fileob)
293 type is (fmsnetcdffile_t)
294 !> If the axis is domain decomposed and the type is FmsNetcdfFile_t, this is regional diagnostic
295 !! So treat it as any other dimension
296 call mpp_get_global_domain(domain, begin=gstart, end=gend) !< Get the global indicies
297 call mpp_get_compute_domain(domain, begin=cstart, end=cend, size=clength) !< Get the compute indicies
298 iend = cend - gstart + 1 !< Get the array indicies for the axis data
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/) )
302
303 !> For regional subaxis add the "domain_decomposition" attribute, which is added
304 !> fms2_io for (other) domains with sufficient decomposition info.
305 call register_variable_attribute(fileob, axis_name, "domain_decomposition", &
306 (/gstart, gend, cstart, cend/))
307 type is (fmsnetcdfdomainfile_t)
308 !> If the axis is domain decomposed and the type is FmsNetcdfDomainFile_t,
309 !! this is a domain decomposed dimension
310 !! so register it as one
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/) )
314 end select
315
316 ELSE IF ( domainu .NE. null_domainug) THEN
317 select type(fileob)
318 type is (fmsnetcdfunstructureddomainfile_t)
319 !> If the axis is in unstructured domain and the type is FmsNetcdfUnstructuredDomainFile_t,
320 !! this is an unstrucutred axis
321 !! so register it as one
322 call register_axis(fileob, axis_name )
323 end select
324 call register_field(fileob, axis_name, type_str, (/axis_name/) )
325 istart = lbound(axis_data,1)
326 iend = ubound(axis_data,1)
327 ELSE
328 !> If the axis is not in a domain, register it as a normal dimension
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)
333 ENDIF !! IF ( Domain .NE. null_domain1d )
334
335 if (length <= 0) then
336 !> @note Check if the time variable is registered. It's possible that is_time_axis_registered
337 !! is set to true if using
338 !! time-templated files because they aren't closed when done writing. An alternative
339 !! to this set up would be to put
340 !! variable_exists into the if statement with an .or. so that it gets registered.
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
347 endif !! if (.not. is_time_axis_registered)
348 endif !! if (length <= 0)
349
350 !> Add the attributes
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))
357
358 if (length > 0 ) then
359 !> If not a time axis, add the positive attribute and write the data
360 select case (axis_direction)
361 case (1)
362 call register_variable_attribute(fileob, axis_name, "positive", "up", str_len=len_trim("up"))
363 case (-1)
364 call register_variable_attribute(fileob, axis_name, "positive", "down", str_len=len_trim("down"))
365 end select
366 call write_data(fileob, axis_name, axis_data(istart:iend) )
367 endif
368
369 !> Write additional axis attributes, from diag_axis_add_attribute calls
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)
373 END IF
374
375 !> Write additional attribute (calendar_type) for time axis ----
376 !> @note calendar attribute is compliant with CF convention
377 !! http://www.cgd.ucar.edu/cms/eaton/netcdf/CF-current.htm#cal
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()
382
383
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'))
393 END IF
394 call set_fileobj_time_name(fileob, axis_name)
395 ELSE
396 time_axis_flag(num_axis_in_file) = .false.
397 END IF
398
399 DEALLOCATE(axis_data)
400
401 !> Deallocate attributes
402 IF ( ALLOCATED(attributes) ) THEN
403 DO j=1, num_attributes
404 IF ( allocated(attributes(j)%fatt ) ) THEN
405 DEALLOCATE(attributes(j)%fatt)
406 END IF
407 IF ( allocated(attributes(j)%iatt ) ) THEN
408 DEALLOCATE(attributes(j)%iatt)
409 END IF
410 END DO
411 DEALLOCATE(attributes)
412 END IF
413
414 !------------- write axis containing edge information ---------------
415
416 ! --- this axis has no edges -----
417 IF ( axis_edges <= 0 ) cycle
418
419 ! --- was this axis edge previously defined? ---
420 id_axis_current = id_axis
421 axis_name_current = axis_name
422 id_axis = axis_edges
423 edges_index = get_axis_index(id_axis)
424 IF ( edges_index > 0 ) cycle
425
426 ! ---- get data for axis edges ----
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)
431
432 ! ---- write edges attribute to original axis ----
433 call register_variable_attribute(fileob, axis_name_current, "edges",trim(axis_name), &
434 & str_len=len_trim(axis_name))
435 ! ---- add edges index to axis list ----
436 ! ---- assume this is not a time axis ----
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.
441
442!> Add edges axis with fms2_io
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)
453 case (1)
454 call register_variable_attribute(fileob, axis_name, "positive", "up", str_len=len_trim("up"))
455 case (-1)
456 call register_variable_attribute(fileob, axis_name, "positive", "down", str_len=len_trim("down"))
457 end select
458 call write_data(fileob, axis_name, axis_data)
459 endif !! if (.not.variable_exists(fileob, axis_name))
460
461 DEALLOCATE (axis_data)
462 END DO
463 END SUBROUTINE write_axis_meta_data
464
465 !> @brief Write the field meta data to file.
466 !! @return diag_fieldtype Field
467 !! @details The meta data for the field is written to the file indicated by file_unit
468 FUNCTION write_field_meta_data ( file_unit, name, axes, units, long_name, range, pack, mval,&
469 & avg_name, time_method, standard_name, interp_method, attributes, num_attributes, &
470 & use_UGdomain, fileob) result ( Field )
471 INTEGER, INTENT(in) :: file_unit !< Output file unit number
472 INTEGER, INTENT(in) :: axes(:) !< Array of axis IDs
473 CHARACTER(len=*), INTENT(in) :: name !< Field name
474 CHARACTER(len=*), INTENT(in) :: units !< Field units
475 CHARACTER(len=*), INTENT(in) :: long_name !< Field's long name
476 REAL, OPTIONAL, INTENT(in) :: range(2) !< Valid range (min, max). If min > max, the range will be ignored
477 REAL, OPTIONAL, INTENT(in) :: mval !< Missing value, must be within valid range
478 INTEGER, OPTIONAL, INTENT(in) :: pack !< Packing flag. Only valid when range specified. Valid values:
479 !! Flag | Size
480 !! --- | ---
481 !! 1 | 64bit
482 !! 2 | 32bit
483 !! 4 | 16bit
484 !! 8 | 8bit
485 CHARACTER(len=*), OPTIONAL, INTENT(in) :: avg_name !< Name of variable containing time averaging info
486 CHARACTER(len=*), OPTIONAL, INTENT(in) :: time_method !< Name of transformation applied to the time-varying data,
487 !! i.e. "avg", "min", "max"
488 CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard name of field
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
493class(fmsnetcdffile_t), intent(inout) :: fileob
494
495 logical :: is_time_bounds !< Flag indicating if the variable 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
502
503character(len=128),dimension(size(axes)) :: axis_names
504 REAL :: scale, add
505 INTEGER :: i, indexx, num, ipack, np
506 LOGICAL :: use_range
507 INTEGER :: axis_indices(size(axes))
508 logical :: use_ugdomain_local
509 !---- Initialize err_msg to bank ----
510 err_msg = ''
511
512 !---- dummy checks ----
513 coord_present = .false.
514 IF( PRESENT(standard_name) ) THEN
515 standard_name2 = standard_name
516 ELSE
517 standard_name2 = 'none'
518 END IF
519
520 use_ugdomain_local = .false.
521 if(present(use_ugdomain)) use_ugdomain_local = use_ugdomain
522
523 num = SIZE(axes(:))
524 ! <ERROR STATUS="FATAL">number of axes < 1</ERROR>
525 IF ( num < 1 ) CALL error_mesg ( 'write_meta_data', 'number of axes < 1', fatal)
526 ! <ERROR STATUS="FATAL">writing meta data out-of-order to different files</ERROR>
527 IF ( file_unit /= current_file_unit ) CALL error_mesg ( 'write_meta_data', &
528 & 'writing meta data out-of-order to different files', fatal)
529
530 IF (trim(name) .eq. "time_bnds") then
531 is_time_bounds = .true.
532 ELSE
533 is_time_bounds = .false.
534 ENDIF
535
536 !---- check all axes for this field ----
537 !---- set up indexing to axistypes ----
538 DO i = 1, num
539 indexx = get_axis_index(axes(i))
540 !---- point to existing axistype -----
541 IF ( indexx > 0 ) THEN
542 axis_indices(i) = indexx
543 ELSE
544 ! <ERROR STATUS="FATAL">axis data not written for field</ERROR>
545 CALL error_mesg ('write_field_meta_data',&
546 & 'axis data not written for field '//trim(name), fatal)
547 END IF
548 !Get the axes names
549 call get_diag_axis_name(axes(i),axis_names(i))
550 END DO
551
552 ! Create coordinate attribute
553 IF ( num >= 2 .OR. (num==1 .and. use_ugdomain_local) ) THEN
554 coord_att = ' '
555 DO i = 1, num
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))
560 ELSE
561 coord_att = trim(coord_att)// ' '//trim(aux_axes(i))
562 ENDIF
563 coord_present = .true.
564 END IF
565 END DO
566 END IF
567
568 !--------------------- write field meta data ---------------------------
569
570 !---- select packing? ----
571 !(packing option only valid with range option)
572 IF ( PRESENT(pack) ) THEN
573 ipack = pack
574 ELSE
575 ipack = 2
576 END IF
577
578 !---- check range ----
579 use_range = .false.
580 add = 0.0
581 scale = 1.0
582 IF ( PRESENT(range) ) THEN
583 IF ( range(2) > range(1) ) THEN
584 use_range = .true.
585 !---- set packing parameters ----
586 IF ( ipack > 2 ) THEN
587 np = ipack/4
588 add = 0.5*(range(1)+range(2))
589 scale = (range(2)-range(1)) / real(max_range(2,np)-max_range(1,np))
590 END IF
591 END IF
592 END IF
593
594 !---- select packing? ----
595 IF ( PRESENT(mval) ) THEN
596 field%miss = mval
597 field%miss_present = .true.
598 IF ( ipack > 2 ) THEN
599 np = ipack/4
600 field%miss_pack = real(missval(np))*scale+add
601 field%miss_pack_present = .true.
602 ELSE
603 field%miss_pack = mval
604 field%miss_pack_present = .false.
605 END IF
606 ELSE
607 field%miss_present = .false.
608 field%miss_pack_present = .false.
609 END IF
610
611 !< Save the fieldname in the diag_fieldtype, so it can be used later
612 field%fieldname = name
613
614 if (.not. variable_exists(fileob,name)) then
615 ! ipack Valid values:
616 ! 1 = 64bit </LI>
617 ! 2 = 32bit </LI>
618 ! 4 = 16bit </LI>
619 ! 8 = 8bit </LI>
620 select case (ipack)
621 case (1)
622 call register_field(fileob,name,"double",axis_names)
623 !< Don't write the _FillValue, missing_value if the variable is
624 !time_bounds to be cf compliant
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))
629 ELSE
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))
632 ENDIF
633 IF ( use_range ) then
634 call register_variable_attribute(fileob,name,"valid_range", real(range,8))
635 ENDIF
636 endif !< if (.not. is_time_bounds)
637 case (2) !default
638 call register_field(fileob,name,"float",axis_names)
639 !< Don't write the _FillValue, missing_value if the variable is
640 !time_bounds to be cf compliant
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))
645 ELSE
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))
648 ENDIF
649 IF ( use_range ) then
650 call register_variable_attribute(fileob,name,"valid_range", real(range,4))
651 ENDIF
652 endif !< if (.not. is_time_bounds)
653 case default
654 CALL error_mesg('diag_output_mod::write_field_meta_data',&
655 &"Pack values must be 1 or 2. Contact the developers.", fatal)
656 end select
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)))
663 ENDIF
664 endif
665 !---- write user defined attributes -----
666 IF ( PRESENT(num_attributes) ) THEN
667 IF ( PRESENT(attributes) ) THEN
668 IF ( num_attributes .GT. 0 .AND. allocated(attributes) ) THEN
669 CALL write_attribute_meta(file_unit, num_attributes, attributes, time_method, err_msg, &
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)
674 END IF
675 ELSE
676 ! Catch some bad cases
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)
685 END IF
686 END IF
687 ELSE
688 ! More edge error cases
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)
692 END IF
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)
697 END IF
698
699 !---- write additional attribute for time averaging -----
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'))
705 END IF
706 END IF
707
708 ! write coordinates attribute for CF compliance
709 IF ( coord_present ) then
710 call register_variable_attribute(fileob,name,'coordinates',trim(coord_att), str_len=len_trim(coord_att))
711 ENDIF
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))
715 ENDIF
716 !---- write attribute for interp_method ----
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))
720 END IF
721
722 !---- get axis domain ----
723 field%Domain = get_domain2d( axes )
724 field%tile_count = get_tile_count( axes )
725 field%DomainU = get_domainug( axes(1) )
726
727 END FUNCTION write_field_meta_data
728
729 !> \brief Write out attribute meta data to file
730 !!
731 !! Write out the attribute meta data to file, for field and axes
732 SUBROUTINE write_attribute_meta(file_unit, num_attributes, attributes, time_method, err_msg, varname, fileob)
733 INTEGER, INTENT(in) :: file_unit !< File unit number
734 INTEGER, INTENT(in) :: num_attributes !< Number of attributes to write
735 TYPE(diag_atttype), DIMENSION(:), INTENT(in) :: attributes !< Array of attributes
736 CHARACTER(len=*), INTENT(in), OPTIONAL :: time_method !< To include in cell_methods attribute if present
737 CHARACTER(len=*), INTENT(out), OPTIONAL :: err_msg !< Return error message
738 CHARACTER(len=*), INTENT(IN), OPTIONAL :: varname !< The name of the variable
739 class(fmsnetcdffile_t), intent(inout) :: fileob !< FMS2_io fileobj
740
741 INTEGER :: i, att_len
742 CHARACTER(len=1280) :: att_str
743
744 ! Clear err_msg if present
745 IF ( PRESENT(err_msg) ) err_msg = ''
746
747 DO i = 1, num_attributes
748 SELECT CASE (attributes(i)%type)
749 CASE (nf90_int)
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
754 RETURN
755 END IF
756 END IF
757 if (present(varname))call register_variable_attribute(fileob, varname,trim(attributes(i)%name), &
758 & attributes(i)%iatt)
759 CASE (nf90_float)
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
764 RETURN
765 END IF
766 END IF
767 if (present(varname))call register_variable_attribute(fileob, varname,trim(attributes(i)%name), &
768 & real(attributes(i)%fatt,4) )
769 CASE (nf90_char)
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
773 ! Append ",time: time_method" if time_method present
774 att_str = attributes(i)%catt(1:attributes(i)%len)//' time: '//time_method
775 att_len = len_trim(att_str)
776 END IF
777 if (present(varname))&
778 call register_variable_attribute(fileob, varname,trim(attributes(i)%name) , att_str(1:att_len), &
779 & str_len=att_len)
780
781 CASE default
782 IF ( fms_error_handler('diag_output_mod::write_attribute_meta', 'Invalid type for attribute '&
783 &//trim(attributes(i)%name)//'.', err_msg) ) THEN
784 RETURN
785 END IF
786 END SELECT
787 END DO
788 END SUBROUTINE write_attribute_meta
789
790 !> @brief Writes axis data to file.
791 !! @details Writes axis data to file. This subroutine is to be called once per file
792 !! after all <TT>write_meta_data</TT> calls, and before the first
793 !! <TT>diag_field_out</TT> call.
794 SUBROUTINE done_meta_data(file_unit)
795 INTEGER, INTENT(in) :: file_unit !< Output file unit number
796
797 !---- write data for all non-time axes ----
798 num_axis_in_file = 0
799 END SUBROUTINE done_meta_data
800
801 !> \brief Writes diagnostic data out using fms2_io routine.
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 !< Variable name
805 REAL , INTENT(inout) :: buffer(:,:,:,:) !< Buffer containing the variable data
806 logical, intent(in) :: static !< Flag indicating if a variable is static
807 integer, intent(in) :: file_num !< Index in the fileobj* types array
808 type(fmsnetcdfunstructureddomainfile_t), intent(inout) :: fileobju(:) !< Array of non domain decomposed fileobj
809 type(fmsnetcdfdomainfile_t), intent(inout) :: fileobj(:) !< Array of domain decomposed fileobj
810 type(fmsnetcdffile_t), intent(inout) :: fileobjnd(:) !< Array of unstructured domain fileobj
811 character(len=2), intent(in) :: fnum_for_domain !< String indicating the type of domain
812 !! "2d" domain decomposed
813 !! "ug" unstructured domain decomposed
814 !! "nd" no domain
815 INTEGER, OPTIONAL, INTENT(in) :: time_in !< Time index
816
817 integer :: time !< Time index
818 real,allocatable :: local_buffer(:,:,:,:) !< Buffer containing the data will be sent to fms2io
819
820!> Set up the time. Static field and default time is 0
821 if ( static ) then
822 time = 0
823 elseif (present(time_in)) then
824 time = time_in
825 else
826 time = 0
827 endif
828
829 if (size(buffer,3) .eq. 1 .and. size(buffer,2) .eq. 1) then
830 !> If the variable is 1D, switch the buffer so that n_diurnal_samples is
831 !! the second dimension (nx, n_diurnal_samples, 1, 1)
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
835 !> If the variable is 2D, switch the n_diurnal_samples and nz dimension, so local_buffer has
836 !! dimension (nx, ny, n_diurnal_samples, 1).
837 allocate(local_buffer(size(buffer,1),size(buffer,2),size(buffer,4),size(buffer,3)))
838 local_buffer(:,:,:,1) = buffer(:,:,1,:)
839 else
840 allocate(local_buffer(size(buffer,1),size(buffer,2),size(buffer,3),size(buffer,4)))
841 local_buffer = buffer(:,:,:,:)
842 endif
843
844 !> Figure out which file object to write output to
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 )
848 endif
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)
852 endif
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)
856 endif
857 else
858 call error_mesg("diag_field_write","fnum_for_domain must be '2d', 'nd', or 'ug'",fatal)
859 endif
860
861 deallocate(local_buffer)
862 end subroutine diag_field_write
863
864!> \brief Writes the time data to the history file
865 subroutine diag_write_time (fileob,rtime_value,time_index,time_name)
866 class(fmsnetcdffile_t), intent(inout) :: fileob !< fms2_io file object
867 real, intent(in) :: rtime_value !< The value of time to be written
868 integer, intent(in) :: time_index !< The index of the time variable
869 character(len=*), intent(in), optional :: time_name !< The name of the time variable
870 character(len=:),allocatable :: name_time !< The name of the time variable
871
872!> Get the name of the time variable
873 if (present(time_name)) then
874 allocate(character(len=len(time_name)) :: name_time)
875 name_time = time_name
876 else
877 allocate(character(len=4) :: name_time)
878 name_time = "time"
879 endif
880!> Write the time data
881 call write_data (fileob, trim(name_time), rtime_value, unlim_dim_level=time_index)
882!> Cleanup
883 if (allocated(name_time)) deallocate(name_time)
884 end subroutine diag_write_time
885
886 !> @brief Return the axis index number.
887 !! @return Integer index
888 FUNCTION get_axis_index(num) RESULT ( index )
889 INTEGER, INTENT(in) :: num
890
891 INTEGER :: index
892 INTEGER :: i
893
894 !---- get the array index for this axis type ----
895 !---- set up pointers to axistypes ----
896 !---- write axis meta data for new axes ----
897 index = 0
898 DO i = 1, num_axis_in_file
899 IF ( num == axis_in_file(i) ) THEN
900 index = i
901 EXIT
902 END IF
903 END DO
904 END FUNCTION get_axis_index
905
906 !> @brief Return the global attribute type.
907 SUBROUTINE get_diag_global_att(gAtt)
908 TYPE(diag_global_att_type), INTENT(out) :: gatt
909
910 gatt=diag_global_att
911 END SUBROUTINE get_diag_global_att
912
913 !> @brief Set the global attribute type.
914 SUBROUTINE set_diag_global_att(component, gridType, tileName)
915 CHARACTER(len=*),INTENT(in) :: component, gridtype, tilename
916
917 ! The following two lines are set to remove compile time warnings
918 ! about 'only used once'.
919 CHARACTER(len=64) :: component_tmp
920 component_tmp = component
921 ! Don't know how to set these for specific component
922 ! Want to be able to say
923 ! if(output_file has component) then
924 diag_global_att%grid_type = gridtype
925 diag_global_att%tile_name = tilename
926 ! endif
927 END SUBROUTINE set_diag_global_att
928
929 !> @brief Flushes the file into disk
930 subroutine diag_flush(file_num, fileobjU, fileobj, fileobjND, fnum_for_domain)
931 integer, intent(in) :: file_num !< Index in the fileobj* types array
932 type(fmsnetcdfunstructureddomainfile_t),intent(inout) :: fileobju(:) !< Array of non domain decomposed fileobj
933 type(fmsnetcdfdomainfile_t), intent(inout) :: fileobj(:) !< Array of domain decomposed fileobj
934 type(fmsnetcdffile_t), intent(inout) :: fileobjnd(:) !< Array of unstructured domain fileobj
935 character(len=2), intent(in) :: fnum_for_domain !< String indicating the type of domain
936 !! "2d" domain decomposed
937 !! "ug" unstructured domain decomposed
938 !! "nd" no 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))
945 else
946 call error_mesg("diag_field_write","No file object is associated with this file number",fatal)
947 endif
948 end subroutine diag_flush
949END MODULE diag_output_mod
950!> @}
951! close documentation grouping
type(domain2d) function, public get_domain2d(ids)
Return the 2D domain for the axis IDs given.
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.
type(domainug) function, public get_domainug(id)
Retrun the 1D domain for the axis ID given.
character(len=128) function, public get_axis_aux(id)
Return the auxiliary name for 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.
type(diag_axis_type), dimension(:), allocatable, save axes
global storage for all defined axes
Definition diag_axis.F90:73
integer function, public get_tile_count(ids)
Return the tile count for the axis.
integer function, public get_axis_length(id)
Return the length of the axis.
integer function, public get_axis_global_length(id)
Return the global length of the axis.
subroutine, public get_diag_axis_name(id, axis_name)
Return the short name of the axis.
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.
Diagnostic field type.
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 done_meta_data(file_unit)
Writes axis data to file.
subroutine write_attribute_meta(file_unit, num_attributes, attributes, time_method, err_msg, varname, fileob)
Write out attribute meta data to file.
subroutine, public set_diag_global_att(component, gridtype, tilename)
Set the global attribute type.
subroutine, public write_axis_meta_data(file_unit, axes, fileob, time_ops, time_axis_registered)
Write the axis meta data to file.
subroutine, public diag_flush(file_num, fileobju, fileobj, fileobjnd, fnum_for_domain)
Flushes the file into disk.
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.
subroutine, public diag_write_time(fileob, rtime_value, time_index, time_name)
Writes the time data to the history file.
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 get_diag_global_att(gatt)
Return the global attribute type.
integer function get_axis_index(num)
Return the axis index number.
Opens a given netcdf or domain file.
Definition fms2_io.F90:122
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.
Gather data sent from pelist onto the root pe Wrapper for MPI_gather, can be used with and without in...
Definition mpp.F90:698
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.