1 !/****************************************************************************
4 !* National Center for Supercomputing Applications *
5 !* University of Illinois at Urbana-Champaign *
6 !* 605 E. Springfield, Champaign IL 61820 *
8 !* For conditions of distribution and use, see the accompanying *
11 !****************************************************************************/
13 module get_attrid_routine
16 module procedure get_attrid
21 subroutine get_attrid(DataHandle,Element,h5_attrid,Status,VAR)
24 use ext_phdf5_support_routines
25 USE HDF5 ! This module contains all necessary modules
27 include 'wrf_status_codes.h'
29 character*(*) ,intent(in) :: Element
30 integer ,intent(in) :: DataHandle
31 integer(hid_t) ,intent(out) :: h5_attrid
32 integer(hid_t) :: dset_id
33 integer ,intent(out) :: Status
34 character*(*) ,intent(in),optional :: VAR
35 integer(hid_t) :: hdf5err
36 type(wrf_phdf5_data_handle),pointer :: DH
38 character(Len = MaxTimeSLen) :: tname
39 character(Len = 256) :: tgroupname
40 integer(hid_t) :: tgroup_id
42 call GetDH(DataHandle,DH,Status)
43 if(Status /= WRF_NO_ERR) then
44 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
45 call wrf_debug ( WARN , msg)
50 call numtochar(1,tname)
51 tgroupname = 'TIME_STAMP_'//tname
52 call h5gopen_f(DH%GroupID,tgroupname,tgroup_id,hdf5err)
53 call h5dopen_f(tgroup_id,VAR,dset_id,hdf5err)
55 Status = WRF_HDF5_ERR_DATASET_OPEN
56 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
57 call wrf_debug ( WARN , msg)
61 call h5aopen_name_f(dset_id,Element,h5_attrid,hdf5err)
63 Status = WRF_HDF5_ERR_ATTRIBUTE_OPEN
64 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
65 call wrf_debug ( WARN , msg)
68 call h5dclose_f(dset_id,hdf5err)
69 call h5gclose_f(tgroup_id,hdf5err)
71 call h5aopen_name_f(DH%GroupID,Element,h5_attrid,hdf5err)
73 Status = WRF_HDF5_ERR_ATTRIBUTE_OPEN
74 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
75 call wrf_debug ( WARN , msg)
80 end subroutine get_attrid
81 end module get_attrid_routine
83 subroutine create_phdf5_objid(DataHandle,obj_id,routine_type,var,Status)
86 use ext_phdf5_support_routines
89 include 'wrf_status_codes.h'
92 integer ,intent(in) :: DataHandle
93 integer(hid_t) ,intent(out) :: obj_id
94 character*3 ,intent(in) :: routine_type
95 character*(*) ,intent(in) :: var
96 integer ,intent(out) :: Status
97 integer(hid_t) :: hdf5err
98 type(wrf_phdf5_data_handle),pointer :: DH
101 call GetDH(DataHandle,DH,Status)
102 if(Status /= WRF_NO_ERR) then
103 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,' ',ROUTINE_TYPE,', line', __LINE__
104 call wrf_debug ( WARN , msg)
108 if(routine_type == 'DOM') then
110 if(DH%FileStatus == WRF_FILE_OPENED_AND_COMMITTED) then
114 else if(routine_type == 'VAR') then
116 if(DH%FileStatus == WRF_FILE_OPENED_AND_COMMITTED) then
118 if(DH%VarNames(i) == var) then
119 obj_id = DH%dsetids(i)
120 write(*,*) "obj_id at attribute",obj_id
121 write(*,*) "DH%VarNames(i)",DH%VarNames(i)
128 Status = WRF_HDF5_ERR_DATA_ID_NOTFOUND
129 write(msg,*) 'CANNOT FIND DATASET ID of the attribute in',__FILE__,&
134 end subroutine create_phdf5_objid
137 subroutine create_phdf5_adtypeid(h5_atypeid,routine_datatype,Count,Status,DataHandle)
140 use ext_phdf5_support_routines
143 include 'wrf_status_codes.h'
146 integer(hid_t) ,intent(out) :: h5_atypeid
147 integer ,intent(in) :: routine_datatype
148 integer ,intent(in) :: Count
149 integer ,intent(out) :: Status
150 integer(hid_t) :: hdf5err
151 integer, intent(in), optional :: DataHandle
152 integer(size_t) :: count_size
154 type(wrf_phdf5_data_handle),pointer :: DH
156 if(routine_datatype == WRF_LOGICAL)then
157 call GetDH(DataHandle,DH,Status)
158 if(Status /= WRF_NO_ERR) then
159 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
160 call wrf_debug ( WARN , msg)
166 select case(routine_datatype)
168 h5_atypeid = H5T_NATIVE_REAL
170 h5_atypeid = H5T_NATIVE_DOUBLE
172 h5_atypeid = H5T_NATIVE_INTEGER
174 h5_atypeid = DH%EnumID
177 call h5tcopy_f(H5T_NATIVE_CHARACTER,h5_atypeid,hdf5err)
178 if(hdf5err.lt.0) then
179 Status = WRF_HDF5_ERR_DATATYPE
180 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
181 call wrf_debug ( WARN , msg)
186 call h5tset_size_f(h5_atypeid,count_size,hdf5err)
187 if(hdf5err.lt.0) then
188 Status = WRF_HDF5_ERR_DATATYPE
189 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
190 call wrf_debug ( WARN , msg)
194 call h5tset_strpad_f(h5_atypeid,H5T_STR_SPACEPAD_F,hdf5err)
195 if(hdf5err.lt.0) then
196 Status = WRF_HDF5_ERR_DATATYPE
197 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
198 call wrf_debug ( WARN , msg)
203 Status = WRF_HDF5_ERR_DATATYPE
210 end subroutine create_phdf5_adtypeid
212 subroutine create_phdf5_adspaceid(Count,str_flag,h5_aspaceid,Status)
217 include 'wrf_status_codes.h'
220 integer ,intent(in) :: Count
221 integer ,intent(in) :: str_flag
222 integer ,intent(out) :: Status
224 integer(hsize_t) , dimension(1) :: adims
225 integer(hid_t) :: hdf5err
226 integer(hid_t) ,intent(out) :: h5_aspaceid
229 ! if string, count is always 1
230 if(str_flag == 1) then
232 call h5screate_simple_f(arank,adims,h5_aspaceid,hdf5err)
233 if(hdf5err.lt.0) then
234 Status = WRF_HDF5_ERR_DATASPACE
235 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
236 call wrf_debug ( WARN , msg)
242 call h5screate_simple_f(arank,adims,h5_aspaceid,hdf5err)
243 if(hdf5err.lt.0) then
244 Status = WRF_HDF5_ERR_DATASPACE
245 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
246 call wrf_debug ( WARN , msg)
255 end subroutine create_phdf5_adspaceid
258 subroutine clean_phdf5_attrids(h5_attr_typeid,h5_space_typeid, &
259 h5_attrid,str_flag,Status)
264 include 'wrf_status_codes.h'
265 integer ,intent(out) :: Status
266 integer(hid_t) ,intent(in) :: h5_attr_typeid
267 integer(hid_t) ,intent(in) :: h5_space_typeid
268 integer(hid_t) ,intent(in) :: h5_attrid
269 integer ,intent(in) :: str_flag
272 if(str_flag == 1) then
273 call h5tclose_f(h5_attr_typeid,hdf5err)
274 if(hdf5err.lt.0) then
275 Status = WRF_HDF5_ERR_CLOSE_GENERAL
276 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
277 call wrf_debug ( WARN , msg)
282 call h5sclose_f(h5_space_typeid,hdf5err)
283 if(hdf5err.lt.0) then
284 Status = WRF_HDF5_ERR_CLOSE_GENERAL
285 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
286 call wrf_debug ( WARN , msg)
289 call h5aclose_f(h5_attrid,hdf5err)
290 if(hdf5err.lt.0) then
291 Status = WRF_HDF5_ERR_ATTRIBUTE_CLOSE
292 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
293 call wrf_debug ( WARN , msg)
301 end subroutine clean_phdf5_attrids
304 subroutine create_h5filetype(dtype_id,Status)
307 use ext_phdf5_support_routines
310 include 'wrf_status_codes.h'
312 integer(hid_t),intent(out) :: dtype_id
313 integer(hid_t) :: dtstr_id
314 integer(size_t) :: type_size
315 integer(size_t) :: type_sizes
316 integer(size_t) :: type_sizei
317 integer(size_t) :: offset
318 integer, intent(out) :: Status
319 integer(hid_t) :: hdf5err
321 call h5tcopy_f(H5T_NATIVE_CHARACTER,dtstr_id,hdf5err)
322 if(hdf5err.lt.0) then
323 Status = WRF_HDF5_ERR_DATATYPE
324 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
325 call wrf_debug ( WARN , msg)
330 call h5tset_size_f(dtstr_id,type_size,hdf5err)
331 if(hdf5err.lt.0) then
332 Status = WRF_HDF5_ERR_DATATYPE
333 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
334 call wrf_debug ( WARN , msg)
338 call h5tget_size_f(dtstr_id,type_sizes,hdf5err)
339 if(hdf5err.lt.0) then
340 Status = WRF_HDF5_ERR_DATATYPE
341 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
342 call wrf_debug ( WARN , msg)
346 call h5tget_size_f(H5T_NATIVE_INTEGER,type_sizei,hdf5err)
347 if(hdf5err.lt.0) then
348 Status = WRF_HDF5_ERR_DATATYPE
349 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
350 call wrf_debug ( WARN , msg)
354 type_size = type_sizes + 2*type_sizei
355 call h5tcreate_f(H5T_COMPOUND_F,type_size,dtype_id,hdf5err)
356 if(hdf5err.lt.0) then
357 Status = WRF_HDF5_ERR_DATATYPE
358 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
359 call wrf_debug ( WARN , msg)
365 call h5tinsert_f(dtype_id,"dim_name",offset,dtstr_id,hdf5err)
366 if(hdf5err.lt.0) then
367 Status = WRF_HDF5_ERR_DATATYPE
368 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
369 call wrf_debug ( WARN , msg)
373 offset = offset + type_sizes
374 call h5tinsert_f(dtype_id,"dim_length",offset,H5T_NATIVE_INTEGER,&
376 if(hdf5err.lt.0) then
377 Status = WRF_HDF5_ERR_DATATYPE
378 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
379 call wrf_debug ( WARN , msg)
383 offset = offset + type_sizei
385 call h5tinsert_f(dtype_id,"dim_unlimited",offset,H5T_NATIVE_INTEGER,&
387 if(hdf5err.lt.0) then
388 Status = WRF_HDF5_ERR_DATATYPE
389 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
390 call wrf_debug ( WARN , msg)
395 call h5tclose_f(dtstr_id,hdf5err)
396 if(hdf5err.lt.0) then
397 Status = WRF_HDF5_ERR_CLOSE_GENERAL
398 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
399 call wrf_debug ( WARN , msg)
405 end subroutine create_h5filetype
407 ! check whether two types are equal, attr_type and h5_attrid
408 subroutine check_type(DataHandle,attr_type,h5_attrid,Status)
411 use ext_phdf5_support_routines
412 USE HDF5 ! This module contains all necessary modules
414 include 'wrf_status_codes.h'
416 integer ,intent(in) :: DataHandle
417 integer(hid_t) ,intent(in) :: attr_type
418 integer(hid_t) ,intent(in) :: h5_attrid
419 integer ,intent(out) :: Status
420 integer(hid_t) :: h5_atypeid
421 integer(hid_t) :: h5_classid
422 integer(hid_t) :: h5_wrfclassid
425 type(wrf_phdf5_data_handle),pointer :: DH
427 call GetDH(DataHandle,DH,Status)
428 if(Status /= WRF_NO_ERR) then
429 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
430 call wrf_debug ( WARN , msg)
434 call h5aget_type_f(h5_attrid,h5_atypeid,hdf5err)
435 if(hdf5err.lt.0) then
436 Status = WRF_HDF5_ERR_DATATYPE
437 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
438 call wrf_debug ( WARN , msg)
442 call h5tget_class_f(h5_atypeid,h5_classid,hdf5err)
443 if(hdf5err.lt.0) then
444 Status = WRF_HDF5_ERR_DATATYPE
445 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
446 call wrf_debug ( WARN , msg)
450 call h5tget_class_f(attr_type,h5_wrfclassid,hdf5err)
451 if(hdf5err.lt.0) then
452 Status = WRF_HDF5_ERR_DATATYPE
453 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
454 call wrf_debug ( WARN , msg)
458 if((h5_classid==H5T_STRING_F).AND.&
459 (attr_type==H5T_NATIVE_CHARACTER)) then
462 if(h5_classid .NE. h5_wrfclassid) then
469 if(flag .EQV. .FALSE.) then
470 Status = WRF_HDF5_ERR_TYPE_MISMATCH
476 end subroutine check_type
479 subroutine retrieve_ti_info(DataHandle,h5_attrid,h5_atypeid,Count,OutCount,Status)
482 use ext_phdf5_support_routines
483 USE HDF5 ! This module contains all necessary modules
485 include 'wrf_status_codes.h'
487 integer ,intent(in) :: DataHandle
488 integer ,intent(in) :: h5_attrid
489 integer(hid_t) ,intent(out) :: h5_atypeid
490 integer ,intent(in) :: Count
491 integer ,intent(out) :: OutCount
492 integer ,intent(out) :: Status
493 integer(hid_t) :: h5_aspaceid
496 integer(hsize_t) :: npoints
498 type(wrf_phdf5_data_handle),pointer :: DH
501 call GetDH(DataHandle,DH,Status)
502 if(Status /= WRF_NO_ERR) then
503 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
504 call wrf_debug ( WARN , msg)
508 if(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
510 call h5aget_type_f(h5_attrid,h5_atypeid,hdf5err)
511 if(hdf5err.lt.0) then
512 Status = WRF_HDF5_ERR_DATATYPE
513 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
514 call wrf_debug ( WARN , msg)
518 call h5aget_space_f(h5_attrid,h5_aspaceid,hdf5err)
519 if(hdf5err.lt.0) then
520 Status = WRF_HDF5_ERR_DATASPACE
521 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
522 call wrf_debug ( WARN , msg)
526 call h5sget_simple_extent_ndims_f(h5_aspaceid,rank,hdf5err)
527 if(hdf5err.lt.0) then
528 Status = WRF_HDF5_ERR_DATASPACE
529 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
530 call wrf_debug ( WARN , msg)
535 ! The rank can be either 0 or 1
536 Status = WRF_HDF5_ERR_OTHERS
537 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
538 call wrf_debug ( WARN , msg)
542 call h5sget_simple_extent_npoints_f(h5_aspaceid,npoints,hdf5err)
543 if(hdf5err.lt.0) then
544 Status = WRF_HDF5_ERR_DATASPACE
545 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
546 call wrf_debug ( WARN , msg)
551 if(npoints > Count) then
553 Status = WRF_ERR_WARN_MORE_DATA_IN_FILE
559 end subroutine retrieve_ti_info
561 subroutine setup_wrtd_dataset(DataHandle,DataSetName,dtypeid,countmd,&
562 dsetid,dspace_id,fspace_id,tgroupid, &
566 use ext_phdf5_support_routines
567 USE HDF5 ! This module contains all necessary modules
569 include 'wrf_status_codes.h'
571 integer ,intent(in) :: DataHandle
572 character*(*) ,intent(in) :: DataSetName
573 integer(hid_t) ,intent(in) :: dtypeid
574 integer ,intent(in) :: countmd
575 integer ,intent(in) :: TimeIndex
577 integer(hid_t) ,intent(out) :: dsetid
578 integer(hid_t) ,intent(out) :: dspace_id
579 integer(hid_t) ,intent(out) :: fspace_id
580 integer(hid_t) ,intent(out) :: tgroupid
581 integer(hid_t) :: crp_list
582 integer ,intent(out) :: Status
584 integer(hsize_t) ,dimension(1) :: sizes
585 integer(hsize_t) ,dimension(1) :: chunk_dims
586 integer(hsize_t) ,dimension(1) :: dims
587 integer(hsize_t) ,dimension(1) :: hdf5_maxdims
588 integer(hsize_t) ,dimension(1) :: offset
589 integer(hsize_t) ,dimension(1) :: count
590 type(wrf_phdf5_data_handle),pointer :: DH
592 character(Len = MaxTimeSLen) :: tname
593 character(Len = 256) :: tgroupname
598 call GetDH(DataHandle,DH,Status)
599 if(Status /= WRF_NO_ERR) then
600 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
601 call wrf_debug ( WARN , msg)
606 chunk_dims(1) = countmd
616 hdf5_maxdims(1) = countmd
618 ! create the memory space id
619 call h5screate_simple_f(1,dims,dspace_id,hdf5err,dims)
620 if(hdf5err.lt.0) then
621 Status = WRF_HDF5_ERR_DATASPACE
622 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
623 call wrf_debug ( WARN , msg)
627 ! create file space(for parallel module, each dataset per time step)
628 call h5screate_simple_f(1,dims,fspace_id,hdf5err,hdf5_maxdims)
629 if(hdf5err.lt.0) then
630 Status = WRF_HDF5_ERR_DATASPACE
631 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
632 call wrf_debug ( WARN , msg)
636 ! obtain the absolute name of the group where the dataset is located
637 call numtochar(TimeIndex,tname)
638 tgroupname = 'TIME_STAMP_'//tname
639 call h5gopen_f(DH%GroupID,tgroupname,tgroupid,hdf5err)
642 call h5dcreate_f(tgroupid,DatasetName,H5T_NATIVE_REAL,fspace_id,&
644 if(hdf5err.lt.0) then
645 Status = WRF_HDF5_ERR_DATASET_CREATE
646 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
647 call wrf_debug ( WARN , msg)
648 write(*,*) "cannot create an HDF5 dataset "
653 end subroutine setup_wrtd_dataset
655 subroutine extend_wrtd_dataset(DataHandle,TimeIndex,countmd,dsetid,dspaceid,&
659 use ext_phdf5_support_routines
660 USE HDF5 ! This module contains all necessary modules
662 include 'wrf_status_codes.h'
664 integer ,intent(in) :: DataHandle
665 integer ,intent(in) :: countmd
666 integer ,intent(in) :: TimeIndex
668 integer(hid_t) ,intent(out) :: dsetid
669 integer(hid_t) ,intent(out) :: dspaceid
670 integer(hid_t) ,intent(out) :: fspaceid
671 integer ,intent(out) :: Status
673 integer(hsize_t) ,dimension(2) :: sizes
674 integer(hsize_t) ,dimension(2) :: chunk_dims
675 integer(hsize_t) ,dimension(2) :: dims
676 integer(hsize_t) ,dimension(2) :: hdf5_maxdims
677 integer(hsize_t) ,dimension(2) :: offset
678 integer(hsize_t) ,dimension(2) :: count
685 offset(2) = TimeIndex - 1
692 CALL h5dextend_f(dsetid,sizes,hdf5err)
693 if(hdf5err.lt.0) then
694 Status = WRF_HDF5_ERR_DATASET_GENERAL
695 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
696 call wrf_debug ( WARN , msg)
697 write(*,*) "cannot extend an HDF5 dataset in index ",TimeIndex
701 ! obtain file space id
702 CALL h5dget_space_f(dsetid,fspaceid,hdf5err)
703 if(hdf5err.lt.0) then
704 Status = WRF_HDF5_ERR_DATASPACE
705 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
706 call wrf_debug ( WARN , msg)
711 CALL h5sselect_hyperslab_f(fspaceid, H5S_SELECT_SET_F, &
712 offset, count, hdf5err)
713 if(hdf5err.lt.0) then
714 Status = WRF_HDF5_ERR_DATASET_GENERAL
715 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
716 call wrf_debug ( WARN , msg)
717 write(*,*) "cannot select hyperslab"
721 ! create the memory space id
722 call h5screate_simple_f(2,dims,dspaceid,hdf5err,dims)
723 if(hdf5err.lt.0) then
724 Status = WRF_HDF5_ERR_DATASPACE
725 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
726 call wrf_debug ( WARN , msg)
727 write(*,*) "cannot create HDF5 memory data space"
732 end subroutine extend_wrtd_dataset
734 subroutine setup_rdtd_dataset(DataHandle,DataSetName,mtypeid,TimeIndex,&
735 countmd,outcountmd,dset_id,memspaceid, &
736 dspace_id,tgroupid,Status)
739 use ext_phdf5_support_routines
740 USE HDF5 ! This module contains all necessary modules
742 include 'wrf_status_codes.h'
744 integer ,intent(in) :: DataHandle
745 character*(*) ,intent(in) :: DataSetName
746 integer ,intent(in) :: countmd
747 integer ,intent(out) :: outcountmd
748 integer ,intent(inout) :: mtypeid
749 integer ,intent(in) :: TimeIndex
751 integer(hid_t) ,intent(out) :: dset_id
752 integer(hid_t) ,intent(out) :: dspace_id
753 integer(hid_t) ,intent(out) :: memspaceid
754 integer(hid_t) ,intent(out) :: tgroupid
755 integer ,intent(out) :: Status
757 integer(hid_t) :: dtype_id
758 integer(hid_t) :: class_type
759 integer(hsize_t) ,dimension(1) :: sizes
760 integer(hsize_t) ,dimension(1) :: dims
761 integer(hsize_t) ,dimension(1) :: h5_dims
762 integer(hsize_t) ,dimension(1) :: hdf5_maxdims
763 integer(hsize_t) ,dimension(1) :: offset
764 integer(hsize_t) ,dimension(1) :: count
766 type(wrf_phdf5_data_handle),pointer :: DH
771 character(Len = MaxTimeSLen) :: tname
772 character(Len = 256) :: tgroupname
774 call GetDH(DataHandle,DH,Status)
775 if(Status /= WRF_NO_ERR) then
776 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
777 call wrf_debug ( WARN , msg)
781 ! obtain the absolute name of the group where the dataset is located
782 call numtochar(TimeIndex,tname)
783 tgroupname = 'TIME_STAMP_'//tname
784 call h5gopen_f(DH%GroupID,tgroupname,tgroupid,hdf5err)
786 ! Obtain HDF5 dataset id
787 call h5dopen_f(tgroupid,DataSetName,dset_id,hdf5err)
788 if(hdf5err.lt.0) then
789 Status = WRF_HDF5_ERR_DATASET_OPEN
790 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
791 call wrf_debug ( WARN , msg)
795 ! Obtain the datatype
796 call h5dget_type_f(dset_id,dtype_id,hdf5err)
797 if(hdf5err.lt.0) then
798 Status = WRF_HDF5_ERR_DATASET_GENERAL
799 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
800 call wrf_debug ( WARN , msg)
804 call h5tget_class_f(dtype_id,class_type,hdf5err)
805 if(hdf5err.lt.0) then
806 Status = WRF_HDF5_ERR_DATATYPE
807 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
808 call wrf_debug ( WARN , msg)
813 if(mtypeid == H5T_NATIVE_REAL .or. mtypeid == H5T_NATIVE_DOUBLE) then
814 if( class_type /= H5T_FLOAT_F) then
815 Status = WRF_HDF5_ERR_TYPE_MISMATCH
816 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
817 call wrf_debug ( WARN , msg)
820 else if(mtypeid ==H5T_NATIVE_CHARACTER) then
821 if(class_type /= H5T_STRING_F) then
822 Status = WRF_HDF5_ERR_TYPE_MISMATCH
823 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
824 call wrf_debug ( WARN , msg)
827 else if(mtypeid == H5T_NATIVE_INTEGER) then
828 if(class_type /= H5T_INTEGER_F) then
829 Status = WRF_HDF5_ERR_TYPE_MISMATCH
830 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
831 call wrf_debug ( WARN , msg)
834 else if(mtypeid == DH%EnumID) then
835 if(class_type /= H5T_ENUM_F) then
836 Status = WRF_HDF5_ERR_TYPE_MISMATCH
837 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
838 call wrf_debug ( WARN , msg)
841 call h5tequal_f(dtype_id,DH%EnumID,flag,hdf5err)
842 if(hdf5err.lt.0) then
843 Status = WRF_HDF5_ERR_DATASET_OPEN
844 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
845 call wrf_debug ( WARN , msg)
848 if(flag .EQV. .FALSE.) then
849 Status = WRF_HDF5_ERR_TYPE_MISMATCH
850 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
851 call wrf_debug ( WARN , msg)
855 Status = WRF_HDF5_ERR_BAD_DATA_TYPE
856 write(msg,*)'Fatal Non-WRF supported TYPE in ',__FILE__,', line',__LINE__
857 call wrf_debug(FATAL, msg)
861 if(mtypeid == H5T_NATIVE_CHARACTER) then
865 ! Obtain the dataspace
866 call h5dget_space_f(dset_id,dspace_id,hdf5err)
867 if(hdf5err.lt.0) then
868 Status = WRF_HDF5_ERR_DATASPACE
869 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
870 call wrf_debug ( WARN , msg)
874 ! Obtain the rank of the dimension
875 call h5sget_simple_extent_ndims_f(dspace_id,StoredDim,hdf5err)
876 if(hdf5err.lt.0) then
877 Status = WRF_HDF5_ERR_DATASPACE
878 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
879 call wrf_debug ( WARN , msg)
884 if(StoredDim /=1) then
885 Status = WRF_HDF5_ERR_DATASPACE
886 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
887 call wrf_debug ( WARN , msg)
892 call h5sget_simple_extent_dims_f(dspace_id,h5_dims,hdf5_maxdims,hdf5err)
893 if(hdf5err.lt.0) then
894 Status = WRF_HDF5_ERR_DATASPACE
895 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
896 call wrf_debug ( WARN , msg)
901 if(countmd <= 0) then
902 Status = WRF_HDF5_ERR_ZERO_LENGTH_READ
903 write(msg,*) 'Warning ZERO LENGTH READ in ',__FILE__,', line', __LINE__
904 call wrf_debug ( WARN , msg)
908 if(countmd .lt. h5_dims(1)) then
911 outcountmd = h5_dims(1)
917 call h5screate_simple_f(1,dims,memspaceid,hdf5err)
918 if(hdf5err.lt.0) then
919 Status = WRF_HDF5_ERR_DATASPACE
920 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
921 call wrf_debug ( WARN , msg)
926 count(1) = outcountmd
929 end subroutine setup_rdtd_dataset
931 subroutine make_strid(str_len,str_id,Status)
934 USE HDF5 ! This module contains all necessary modules
936 include 'wrf_status_codes.h'
938 integer ,intent(in) :: str_len
939 integer(hid_t),intent(out) :: str_id
940 integer ,intent(out) :: Status
941 integer(size_t) :: str_lensize
945 if(str_len <= 0) then
946 Status = WRF_HDF5_ERR_ATTRIBUTE_GENERAL
947 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
948 call wrf_debug ( WARN , msg)
952 call h5tcopy_f(H5T_NATIVE_CHARACTER,str_id,hdf5err)
953 if(hdf5err.lt.0) then
954 Status = WRF_HDF5_ERR_DATATYPE
955 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
956 call wrf_debug ( WARN , msg)
960 str_lensize = str_len
961 call h5tset_size_f(str_id,str_lensize,hdf5err)
962 if(hdf5err.lt.0) then
963 Status = WRF_HDF5_ERR_DATATYPE
964 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
965 call wrf_debug ( WARN , msg)
969 call h5tset_strpad_f(str_id,H5T_STR_SPACEPAD_F,hdf5err)
970 if(hdf5err.lt.0) then
971 Status = WRF_HDF5_ERR_DATATYPE
972 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
973 call wrf_debug ( WARN , msg)
977 end subroutine make_strid