1 !***************************************************************************
2 !* The HDF5 WRF IO module was written by the the HDF Group at NCSA, the *
3 !* National Center for Supercomputing Applications. *
5 !* National Center for Supercomputing Applications *
6 !* University of Illinois at Urbana-Champaign *
7 !* 605 E. Springfield, Champaign IL 61820 *
8 !* http://hdf.ncsa.uiuc.edu/ *
10 !* Copyright 2004 by the Board of Trustees, University of Illinois, *
12 !* Redistribution or use of this IO module, with or without modification, *
13 !* is permitted for any purpose, including commercial purposes. *
15 !* This software is an unsupported prototype. Use at your own risk. *
16 !* http://hdf.ncsa.uiuc.edu/apps/WRF-ROMS *
18 !* This work was funded by the MEAD expedition at the National Center *
19 !* for Supercomputing Applications, NCSA. For more information see: *
20 !* http://www.ncsa.uiuc.edu/expeditions/MEAD *
23 !****************************************************************************/
25 module get_attrid_routine
28 module procedure get_attrid
33 subroutine get_attrid(DataHandle,Element,h5_attrid,Status,VAR)
36 use ext_phdf5_support_routines
37 USE HDF5 ! This module contains all necessary modules
39 include 'wrf_status_codes.h'
41 character*(*) ,intent(in) :: Element
42 integer ,intent(in) :: DataHandle
43 integer(hid_t) ,intent(out) :: h5_attrid
44 integer(hid_t) :: dset_id
45 integer ,intent(out) :: Status
46 character*(*) ,intent(in),optional :: VAR
47 integer(hid_t) :: hdf5err
48 type(wrf_phdf5_data_handle),pointer :: DH
50 character(Len = MaxTimeSLen) :: tname
51 character(Len = 512) :: tgroupname
52 integer(hid_t) :: tgroup_id
54 call GetDH(DataHandle,DH,Status)
55 if(Status /= WRF_NO_ERR) then
56 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
57 call wrf_debug ( WARN , msg)
62 call numtochar(1,tname)
63 tgroupname = 'TIME_STAMP_'//tname
64 call h5gopen_f(DH%GroupID,tgroupname,tgroup_id,hdf5err)
66 Status = WRF_HDF5_ERR_GROUP
67 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
68 call wrf_debug ( WARN , msg)
71 call h5dopen_f(tgroup_id,VAR,dset_id,hdf5err)
73 Status = WRF_HDF5_ERR_DATASET_OPEN
74 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
75 call wrf_debug ( WARN , msg)
79 call h5aopen_name_f(dset_id,Element,h5_attrid,hdf5err)
81 Status = WRF_HDF5_ERR_ATTRIBUTE_OPEN
82 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
83 call wrf_debug ( WARN , msg)
86 call h5dclose_f(dset_id,hdf5err)
87 call h5gclose_f(tgroup_id,hdf5err)
89 call h5aopen_name_f(DH%GroupID,Element,h5_attrid,hdf5err)
90 write(*,*) "Element ",Element
92 Status = WRF_HDF5_ERR_ATTRIBUTE_OPEN
93 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
94 call wrf_debug ( WARN , msg)
99 end subroutine get_attrid
100 end module get_attrid_routine
102 subroutine create_phdf5_objid(DataHandle,obj_id,routine_type,var,Status)
105 use ext_phdf5_support_routines
108 include 'wrf_status_codes.h'
111 integer ,intent(in) :: DataHandle
112 integer(hid_t) ,intent(out) :: obj_id
113 character*3 ,intent(in) :: routine_type
114 character*(*) ,intent(in) :: var
115 integer ,intent(out) :: Status
116 integer(hid_t) :: hdf5err
117 type(wrf_phdf5_data_handle),pointer :: DH
120 call GetDH(DataHandle,DH,Status)
121 if(Status /= WRF_NO_ERR) then
122 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,' ',ROUTINE_TYPE,', line', __LINE__
123 call wrf_debug ( WARN , msg)
127 if(routine_type == 'DOM') then
129 if(DH%FileStatus == WRF_FILE_OPENED_AND_COMMITTED) then
133 else if(routine_type == 'VAR') then
135 if(DH%FileStatus == WRF_FILE_OPENED_AND_COMMITTED) then
137 if(DH%VarNames(i) == var) then
138 obj_id = DH%dsetids(i)
145 Status = WRF_HDF5_ERR_DATA_ID_NOTFOUND
146 write(msg,*) 'CANNOT FIND DATASET ID of the attribute in',__FILE__,&
151 end subroutine create_phdf5_objid
154 subroutine create_phdf5_adtypeid(h5_atypeid,routine_datatype,Count,Status,DataHandle)
157 use ext_phdf5_support_routines
160 include 'wrf_status_codes.h'
163 integer(hid_t) ,intent(out) :: h5_atypeid
164 integer ,intent(in) :: routine_datatype
165 integer ,intent(in) :: Count
166 integer ,intent(out) :: Status
167 integer(hid_t) :: hdf5err
168 integer, intent(in), optional :: DataHandle
169 integer(size_t) :: count_size
171 type(wrf_phdf5_data_handle),pointer :: DH
173 if(routine_datatype == WRF_LOGICAL)then
174 call GetDH(DataHandle,DH,Status)
175 if(Status /= WRF_NO_ERR) then
176 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
177 call wrf_debug ( WARN , msg)
183 select case(routine_datatype)
185 h5_atypeid = H5T_NATIVE_REAL
187 h5_atypeid = H5T_NATIVE_DOUBLE
189 h5_atypeid = H5T_NATIVE_INTEGER
191 h5_atypeid = DH%EnumID
194 call h5tcopy_f(H5T_NATIVE_CHARACTER,h5_atypeid,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 call h5tset_size_f(h5_atypeid,count_size,hdf5err)
204 if(hdf5err.lt.0) then
205 Status = WRF_HDF5_ERR_DATATYPE
206 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
207 call wrf_debug ( WARN , msg)
211 call h5tset_strpad_f(h5_atypeid,H5T_STR_SPACEPAD_F,hdf5err)
212 if(hdf5err.lt.0) then
213 Status = WRF_HDF5_ERR_DATATYPE
214 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
215 call wrf_debug ( WARN , msg)
220 Status = WRF_HDF5_ERR_DATATYPE
227 end subroutine create_phdf5_adtypeid
229 subroutine create_phdf5_adspaceid(Count,str_flag,h5_aspaceid,Status)
234 include 'wrf_status_codes.h'
237 integer ,intent(in) :: Count
238 integer ,intent(in) :: str_flag
239 integer ,intent(out) :: Status
241 integer(hsize_t) , dimension(1) :: adims
242 integer(hid_t) :: hdf5err
243 integer(hid_t) ,intent(out) :: h5_aspaceid
246 ! if string, count is always 1
247 if(str_flag == 1) then
249 call h5screate_simple_f(arank,adims,h5_aspaceid,hdf5err)
250 if(hdf5err.lt.0) then
251 Status = WRF_HDF5_ERR_DATASPACE
252 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
253 call wrf_debug ( WARN , msg)
259 call h5screate_simple_f(arank,adims,h5_aspaceid,hdf5err)
260 if(hdf5err.lt.0) then
261 Status = WRF_HDF5_ERR_DATASPACE
262 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
263 call wrf_debug ( WARN , msg)
272 end subroutine create_phdf5_adspaceid
275 subroutine clean_phdf5_attrids(h5_attr_typeid,h5_space_typeid, &
276 h5_attrid,str_flag,Status)
281 include 'wrf_status_codes.h'
282 integer ,intent(out) :: Status
283 integer(hid_t) ,intent(in) :: h5_attr_typeid
284 integer(hid_t) ,intent(in) :: h5_space_typeid
285 integer(hid_t) ,intent(in) :: h5_attrid
286 integer ,intent(in) :: str_flag
289 if(str_flag == 1) then
290 call h5tclose_f(h5_attr_typeid,hdf5err)
291 if(hdf5err.lt.0) then
292 Status = WRF_HDF5_ERR_CLOSE_GENERAL
293 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
294 call wrf_debug ( WARN , msg)
299 call h5sclose_f(h5_space_typeid,hdf5err)
300 if(hdf5err.lt.0) then
301 Status = WRF_HDF5_ERR_CLOSE_GENERAL
302 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
303 call wrf_debug ( WARN , msg)
306 call h5aclose_f(h5_attrid,hdf5err)
307 if(hdf5err.lt.0) then
308 Status = WRF_HDF5_ERR_ATTRIBUTE_CLOSE
309 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
310 call wrf_debug ( WARN , msg)
318 end subroutine clean_phdf5_attrids
321 subroutine create_h5filetype(dtype_id,Status)
324 use ext_phdf5_support_routines
327 include 'wrf_status_codes.h'
329 integer(hid_t),intent(out) :: dtype_id
330 integer(hid_t) :: dtstr_id
331 integer(size_t) :: type_size
332 integer(size_t) :: type_sizes
333 integer(size_t) :: type_sizei
334 integer(size_t) :: offset
335 integer, intent(out) :: Status
336 integer(hid_t) :: hdf5err
338 call h5tcopy_f(H5T_NATIVE_CHARACTER,dtstr_id,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)
348 call h5tset_size_f(dtstr_id,type_size,hdf5err)
349 if(hdf5err.lt.0) then
350 Status = WRF_HDF5_ERR_DATATYPE
351 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
352 call wrf_debug ( WARN , msg)
356 call h5tget_size_f(dtstr_id,type_sizes,hdf5err)
357 if(hdf5err.lt.0) then
358 Status = WRF_HDF5_ERR_DATATYPE
359 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
360 call wrf_debug ( WARN , msg)
364 call h5tget_size_f(H5T_NATIVE_INTEGER,type_sizei,hdf5err)
365 if(hdf5err.lt.0) then
366 Status = WRF_HDF5_ERR_DATATYPE
367 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
368 call wrf_debug ( WARN , msg)
372 type_size = type_sizes + 2*type_sizei
374 call h5tcreate_f(H5T_COMPOUND_F,type_size,dtype_id,hdf5err)
375 if(hdf5err.lt.0) then
376 Status = WRF_HDF5_ERR_DATATYPE
377 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
378 call wrf_debug ( WARN , msg)
385 call h5tinsert_f(dtype_id,"dim_name",offset,dtstr_id,hdf5err)
386 if(hdf5err.lt.0) then
387 Status = WRF_HDF5_ERR_DATATYPE
388 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
389 call wrf_debug ( WARN , msg)
393 offset = offset + type_sizes
394 call h5tinsert_f(dtype_id,"dim_length",offset,H5T_NATIVE_INTEGER,&
396 if(hdf5err.lt.0) then
397 Status = WRF_HDF5_ERR_DATATYPE
398 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
399 call wrf_debug ( WARN , msg)
403 offset = offset + type_sizei
405 call h5tinsert_f(dtype_id,"dim_unlimited",offset,H5T_NATIVE_INTEGER,&
407 if(hdf5err.lt.0) then
408 Status = WRF_HDF5_ERR_DATATYPE
409 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
410 call wrf_debug ( WARN , msg)
415 call h5tclose_f(dtstr_id,hdf5err)
416 if(hdf5err.lt.0) then
417 Status = WRF_HDF5_ERR_CLOSE_GENERAL
418 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
419 call wrf_debug ( WARN , msg)
425 end subroutine create_h5filetype
427 ! check whether two types are equal, attr_type and h5_attrid
428 subroutine check_type(DataHandle,attr_type,h5_attrid,Status)
431 use ext_phdf5_support_routines
432 USE HDF5 ! This module contains all necessary modules
434 include 'wrf_status_codes.h'
436 integer ,intent(in) :: DataHandle
437 integer(hid_t) ,intent(in) :: attr_type
438 integer(hid_t) ,intent(in) :: h5_attrid
439 integer ,intent(out) :: Status
440 integer(hid_t) :: h5_atypeid
441 integer(hid_t) :: h5_classid
442 integer(hid_t) :: h5_wrfclassid
445 type(wrf_phdf5_data_handle),pointer :: DH
447 call GetDH(DataHandle,DH,Status)
448 if(Status /= WRF_NO_ERR) then
449 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
450 call wrf_debug ( WARN , msg)
454 call h5aget_type_f(h5_attrid,h5_atypeid,hdf5err)
455 if(hdf5err.lt.0) then
456 Status = WRF_HDF5_ERR_DATATYPE
457 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
458 call wrf_debug ( WARN , msg)
462 call h5tget_class_f(h5_atypeid,h5_classid,hdf5err)
463 if(hdf5err.lt.0) then
464 Status = WRF_HDF5_ERR_DATATYPE
465 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
466 call wrf_debug ( WARN , msg)
470 call h5tget_class_f(attr_type,h5_wrfclassid,hdf5err)
471 if(hdf5err.lt.0) then
472 Status = WRF_HDF5_ERR_DATATYPE
473 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
474 call wrf_debug ( WARN , msg)
478 if((h5_classid==H5T_STRING_F).AND.&
479 (attr_type==H5T_NATIVE_CHARACTER)) then
482 if(h5_classid .NE. h5_wrfclassid) then
489 if(flag .EQV. .FALSE.) then
490 Status = WRF_HDF5_ERR_TYPE_MISMATCH
496 end subroutine check_type
499 subroutine retrieve_ti_info(DataHandle,h5_attrid,h5_atypeid,Count,OutCount,Status)
502 use ext_phdf5_support_routines
503 USE HDF5 ! This module contains all necessary modules
505 include 'wrf_status_codes.h'
507 integer ,intent(in) :: DataHandle
508 integer ,intent(in) :: h5_attrid
509 integer(hid_t) ,intent(out) :: h5_atypeid
510 integer ,intent(in) :: Count
511 integer ,intent(out) :: OutCount
512 integer ,intent(out) :: Status
513 integer(hid_t) :: h5_aspaceid
517 integer(hsize_t) :: npoints
519 type(wrf_phdf5_data_handle),pointer :: DH
522 call GetDH(DataHandle,DH,Status)
523 if(Status /= WRF_NO_ERR) then
524 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
525 call wrf_debug ( WARN , msg)
529 if(DH%FileStatus == WRF_FILE_OPENED_FOR_READ) then
531 call h5aget_type_f(h5_attrid,h5_atypeid,hdf5err)
532 if(hdf5err.lt.0) then
533 Status = WRF_HDF5_ERR_DATATYPE
534 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
535 call wrf_debug ( WARN , msg)
539 call h5aget_space_f(h5_attrid,h5_aspaceid,hdf5err)
540 if(hdf5err.lt.0) then
541 Status = WRF_HDF5_ERR_DATASPACE
542 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
543 call wrf_debug ( WARN , msg)
547 call h5sget_simple_extent_ndims_f(h5_aspaceid,rank,hdf5err)
548 if(hdf5err.lt.0) then
549 Status = WRF_HDF5_ERR_DATASPACE
550 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
551 call wrf_debug ( WARN , msg)
556 ! The rank can be either 0 or 1
557 Status = WRF_HDF5_ERR_OTHERS
558 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
559 call wrf_debug ( WARN , msg)
563 call h5sget_simple_extent_npoints_f(h5_aspaceid,npoints,hdf5err)
564 if(hdf5err.lt.0) then
565 Status = WRF_HDF5_ERR_DATASPACE
566 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
567 call wrf_debug ( WARN , msg)
572 call h5tget_class_f(h5_atypeid,typeclass,hdf5err)
573 if(hdf5err.lt.0) then
574 Status = WRF_HDF5_ERR_DATATYPE
575 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
576 call wrf_debug ( WARN , msg)
579 if((npoints > Count).and.(typeclass.ne.H5T_STRING_F)) then
581 Status = WRF_HDF5_ERR_MORE_DATA_IN_FILE
587 end subroutine retrieve_ti_info
589 subroutine setup_wrtd_dataset(DataHandle,DataSetName,dtypeid,countmd,&
590 dsetid,dspace_id,fspace_id,tgroupid, &
594 use ext_phdf5_support_routines
595 USE HDF5 ! This module contains all necessary modules
597 include 'wrf_status_codes.h'
599 integer ,intent(in) :: DataHandle
600 character*(*) ,intent(in) :: DataSetName
601 integer(hid_t) ,intent(in) :: dtypeid
602 integer ,intent(in) :: countmd
603 integer ,intent(in) :: TimeIndex
605 integer(hid_t) ,intent(out) :: dsetid
606 integer(hid_t) ,intent(out) :: dspace_id
607 integer(hid_t) ,intent(out) :: fspace_id
608 integer(hid_t) ,intent(out) :: tgroupid
609 integer(hid_t) :: crp_list
610 integer ,intent(out) :: Status
612 integer(hsize_t) ,dimension(1) :: sizes
613 integer(hsize_t) ,dimension(1) :: chunk_dims
614 integer(hsize_t) ,dimension(1) :: dims
615 integer(hsize_t) ,dimension(1) :: hdf5_maxdims
616 integer(hsize_t) ,dimension(1) :: offset
617 integer(hsize_t) ,dimension(1) :: count
618 type(wrf_phdf5_data_handle),pointer :: DH
620 character(Len = MaxTimeSLen) :: tname
621 character(Len = 512) :: tgroupname
626 call GetDH(DataHandle,DH,Status)
627 if(Status /= WRF_NO_ERR) then
628 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
629 call wrf_debug ( WARN , msg)
634 chunk_dims(1) = countmd
644 hdf5_maxdims(1) = countmd
646 ! create the memory space id
647 call h5screate_simple_f(1,dims,dspace_id,hdf5err,dims)
648 if(hdf5err.lt.0) then
649 Status = WRF_HDF5_ERR_DATASPACE
650 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
651 call wrf_debug ( WARN , msg)
655 ! create file space(for parallel module, each dataset per time step)
656 call h5screate_simple_f(1,dims,fspace_id,hdf5err,hdf5_maxdims)
657 if(hdf5err.lt.0) then
658 Status = WRF_HDF5_ERR_DATASPACE
659 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
660 call wrf_debug ( WARN , msg)
664 ! obtain the absolute name of the group where the dataset is located
665 call numtochar(TimeIndex,tname)
666 tgroupname = 'TIME_STAMP_'//tname
667 if(DH%TgroupIDs(TimeIndex) /= -1) then
668 tgroupid = DH%TgroupIDs(TimeIndex)
669 ! call h5gopen_f(DH%GroupID,tgroupname,tgroupid,hdf5err)
671 call h5gcreate_f(DH%GroupID,tgroupname,tgroupid,hdf5err)
672 if(hdf5err.lt.0) then
673 Status = WRF_HDF5_ERR_GROUP
674 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
675 call wrf_debug(WARN,msg)
678 DH%TgroupIDs(TimeIndex) = tgroupid
682 call h5dcreate_f(tgroupid,DatasetName,dtypeid,fspace_id,&
684 if(hdf5err.lt.0) then
685 Status = WRF_HDF5_ERR_DATASET_CREATE
686 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
687 call wrf_debug ( WARN , msg)
692 end subroutine setup_wrtd_dataset
694 subroutine extend_wrtd_dataset(DataHandle,TimeIndex,countmd,dsetid,dspaceid,&
698 use ext_phdf5_support_routines
699 USE HDF5 ! This module contains all necessary modules
701 include 'wrf_status_codes.h'
703 integer ,intent(in) :: DataHandle
704 integer ,intent(in) :: countmd
705 integer ,intent(in) :: TimeIndex
707 integer(hid_t) ,intent(out) :: dsetid
708 integer(hid_t) ,intent(out) :: dspaceid
709 integer(hid_t) ,intent(out) :: fspaceid
710 integer ,intent(out) :: Status
712 integer(hsize_t) ,dimension(2) :: sizes
713 integer(hsize_t) ,dimension(2) :: chunk_dims
714 integer(hsize_t) ,dimension(2) :: dims
715 integer(hsize_t) ,dimension(2) :: hdf5_maxdims
716 integer(hsize_t) ,dimension(2) :: offset
717 integer(hsize_t) ,dimension(2) :: count
724 offset(2) = TimeIndex - 1
731 CALL h5dextend_f(dsetid,sizes,hdf5err)
732 if(hdf5err.lt.0) then
733 Status = WRF_HDF5_ERR_DATASET_GENERAL
734 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
735 call wrf_debug ( WARN , msg)
739 ! obtain file space id
740 CALL h5dget_space_f(dsetid,fspaceid,hdf5err)
741 if(hdf5err.lt.0) then
742 Status = WRF_HDF5_ERR_DATASPACE
743 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
744 call wrf_debug ( WARN , msg)
749 CALL h5sselect_hyperslab_f(fspaceid, H5S_SELECT_SET_F, &
750 offset, count, hdf5err)
751 if(hdf5err.lt.0) then
752 Status = WRF_HDF5_ERR_DATASET_GENERAL
753 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
754 call wrf_debug ( WARN , msg)
758 ! create the memory space id
759 call h5screate_simple_f(2,dims,dspaceid,hdf5err,dims)
760 if(hdf5err.lt.0) then
761 Status = WRF_HDF5_ERR_DATASPACE
762 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
763 call wrf_debug ( WARN , msg)
768 end subroutine extend_wrtd_dataset
770 subroutine setup_rdtd_dataset(DataHandle,DataSetName,mtypeid,TimeIndex,&
771 countmd,outcountmd,dset_id,memspaceid, &
772 dspace_id,tgroupid,Status)
775 use ext_phdf5_support_routines
776 USE HDF5 ! This module contains all necessary modules
778 include 'wrf_status_codes.h'
780 integer ,intent(in) :: DataHandle
781 character*(*) ,intent(in) :: DataSetName
782 integer ,intent(in) :: countmd
783 integer ,intent(out) :: outcountmd
784 integer ,intent(inout) :: mtypeid
785 integer ,intent(in) :: TimeIndex
787 integer(hid_t) ,intent(out) :: dset_id
788 integer(hid_t) ,intent(out) :: dspace_id
789 integer(hid_t) ,intent(out) :: memspaceid
790 integer(hid_t) ,intent(out) :: tgroupid
791 integer ,intent(out) :: Status
793 integer(hid_t) :: dtype_id
794 integer(hid_t) :: class_type
795 integer(hsize_t) ,dimension(1) :: sizes
796 integer(hsize_t) ,dimension(1) :: dims
797 integer(hsize_t) ,dimension(1) :: h5_dims
798 integer(hsize_t) ,dimension(1) :: hdf5_maxdims
799 integer(hsize_t) ,dimension(1) :: offset
800 integer(hsize_t) ,dimension(1) :: count
802 type(wrf_phdf5_data_handle),pointer :: DH
807 character(Len = MaxTimeSLen) :: tname
808 character(Len = 512) :: tgroupname
810 call GetDH(DataHandle,DH,Status)
811 if(Status /= WRF_NO_ERR) then
812 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
813 call wrf_debug ( WARN , msg)
817 ! obtain the absolute name of the group where the dataset is located
818 call numtochar(TimeIndex,tname)
819 tgroupname = 'TIME_STAMP_'//tname
820 call h5gopen_f(DH%GroupID,tgroupname,tgroupid,hdf5err)
822 ! Obtain HDF5 dataset id
823 call h5dopen_f(tgroupid,DataSetName,dset_id,hdf5err)
824 if(hdf5err.lt.0) then
825 Status = WRF_HDF5_ERR_DATASET_OPEN
826 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
827 call wrf_debug ( WARN , msg)
831 ! Obtain the datatype
832 call h5dget_type_f(dset_id,dtype_id,hdf5err)
833 if(hdf5err.lt.0) then
834 Status = WRF_HDF5_ERR_DATASET_GENERAL
835 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
836 call wrf_debug ( WARN , msg)
840 call h5tget_class_f(dtype_id,class_type,hdf5err)
841 if(hdf5err.lt.0) then
842 Status = WRF_HDF5_ERR_DATATYPE
843 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
844 call wrf_debug ( WARN , msg)
849 if(mtypeid == H5T_NATIVE_REAL .or. mtypeid == H5T_NATIVE_DOUBLE) then
850 if( class_type /= H5T_FLOAT_F) then
851 Status = WRF_HDF5_ERR_TYPE_MISMATCH
852 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
853 call wrf_debug ( WARN , msg)
856 else if(mtypeid ==H5T_NATIVE_CHARACTER) then
857 if(class_type /= H5T_STRING_F) then
858 Status = WRF_HDF5_ERR_TYPE_MISMATCH
859 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
860 call wrf_debug ( WARN , msg)
863 else if(mtypeid == H5T_NATIVE_INTEGER) then
864 if(class_type /= H5T_INTEGER_F) then
865 Status = WRF_HDF5_ERR_TYPE_MISMATCH
866 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
867 call wrf_debug ( WARN , msg)
870 else if(mtypeid == DH%EnumID) then
871 if(class_type /= H5T_ENUM_F) then
872 Status = WRF_HDF5_ERR_TYPE_MISMATCH
873 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
874 call wrf_debug ( WARN , msg)
877 call h5tequal_f(dtype_id,DH%EnumID,flag,hdf5err)
878 if(hdf5err.lt.0) then
879 Status = WRF_HDF5_ERR_DATASET_OPEN
880 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
881 call wrf_debug ( WARN , msg)
884 if(flag .EQV. .FALSE.) then
885 Status = WRF_HDF5_ERR_TYPE_MISMATCH
886 write(msg,*) 'Warning TYPE MISMATCH in ',__FILE__,', line', __LINE__
887 call wrf_debug ( WARN , msg)
891 Status = WRF_HDF5_ERR_BAD_DATA_TYPE
892 write(msg,*)'Fatal Non-WRF supported TYPE in ',__FILE__,', line',__LINE__
893 call wrf_debug(FATAL, msg)
897 if(mtypeid == H5T_NATIVE_CHARACTER) then
901 ! Obtain the dataspace
902 call h5dget_space_f(dset_id,dspace_id,hdf5err)
903 if(hdf5err.lt.0) then
904 Status = WRF_HDF5_ERR_DATASPACE
905 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
906 call wrf_debug ( WARN , msg)
910 ! Obtain the rank of the dimension
911 call h5sget_simple_extent_ndims_f(dspace_id,StoredDim,hdf5err)
912 if(hdf5err.lt.0) then
913 Status = WRF_HDF5_ERR_DATASPACE
914 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
915 call wrf_debug ( WARN , msg)
920 if(StoredDim /=1) then
921 Status = WRF_HDF5_ERR_DATASPACE
922 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
923 call wrf_debug ( WARN , msg)
928 call h5sget_simple_extent_dims_f(dspace_id,h5_dims,hdf5_maxdims,hdf5err)
929 if(hdf5err.lt.0) then
930 Status = WRF_HDF5_ERR_DATASPACE
931 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
932 call wrf_debug ( WARN , msg)
937 if(countmd <= 0) then
938 Status = WRF_HDF5_ERR_ZERO_LENGTH_READ
939 write(msg,*) 'Warning ZERO LENGTH READ in ',__FILE__,', line', __LINE__
940 call wrf_debug ( WARN , msg)
944 if(countmd .lt. h5_dims(1)) then
947 outcountmd = h5_dims(1)
953 call h5screate_simple_f(1,dims,memspaceid,hdf5err)
954 if(hdf5err.lt.0) then
955 Status = WRF_HDF5_ERR_DATASPACE
956 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
957 call wrf_debug ( WARN , msg)
962 count(1) = outcountmd
965 end subroutine setup_rdtd_dataset
967 subroutine make_strid(str_len,str_id,Status)
970 USE HDF5 ! This module contains all necessary modules
972 include 'wrf_status_codes.h'
974 integer ,intent(in) :: str_len
975 integer(hid_t),intent(out) :: str_id
976 integer ,intent(out) :: Status
977 integer(size_t) :: str_lensize
981 if(str_len <= 0) then
982 Status = WRF_HDF5_ERR_ATTRIBUTE_GENERAL
983 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
984 call wrf_debug ( WARN , msg)
988 call h5tcopy_f(H5T_NATIVE_CHARACTER,str_id,hdf5err)
989 if(hdf5err.lt.0) then
990 Status = WRF_HDF5_ERR_DATATYPE
991 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
992 call wrf_debug ( WARN , msg)
996 str_lensize = str_len
997 call h5tset_size_f(str_id,str_lensize,hdf5err)
998 if(hdf5err.lt.0) then
999 Status = WRF_HDF5_ERR_DATATYPE
1000 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1001 call wrf_debug ( WARN , msg)
1005 call h5tset_strpad_f(str_id,H5T_STR_SPACEPAD_F,hdf5err)
1006 if(hdf5err.lt.0) then
1007 Status = WRF_HDF5_ERR_DATATYPE
1008 write(msg,*) 'Warning Status = ',Status,' in ',__FILE__,', line', __LINE__
1009 call wrf_debug ( WARN , msg)
1013 end subroutine make_strid