GetBegDateTime
- Get begin date/time of file
This routine returns the begin date/begin time on file. For a native CFIO file, it simply returns the value of global attributes begin_date/begin_time. If these do no exist then it attempts to parse the COARDS compliant time unit.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer | :: | fid |
file ID |
|||
integer | :: | begDate |
beginning date |
|||
integer | :: | begTime |
beginning time |
|||
integer | :: | incSecs |
time increment in secs |
|||
integer | :: | rc |
error return code |
subroutine GetBegDateTime ( fid, begDate, begTime, incSecs, rc ) ! ! !USES: ! Implicit NONE ! ! !INPUT PARAMETERS: ! integer fid !! file ID ! ! !OUTPUT PARAMETERS: ! integer begDate !! beginning date integer begTime !! beginning time integer incSecs !! time increment in secs integer rc !! error return code ! !------------------------------------------------------------------------- integer i, timeId, hour, min, sec, corner(1), timInc integer year, month, day character(len=MAXCHR) timeUnits, dimUnits, stdName character(len=MAXCHR) varName, dimName integer type, nvDims, vdims(MAXVDIMS), nvAtts, dimSize integer nDims, nvars, ngatts, dimId ! Time conversion local variables real(kind=REAL32) :: rtime, rtime_array(1) real(kind=REAL64) :: dtime, dtime_array(1) integer(kind=INT16) :: itime, itime_array(1) integer(kind=INT32) :: ltime, ltime_array(1) integer t1, t2, tMult, newDate, newTime ! Start by determing the ID of the time coordinate variable ! --------------------------------------------------------- timeId = -1 rc = NF90_INQUIRE (fid, nDims, nvars, ngatts, dimId) if (err("GetBegDateTime: NF90_INQUIRE failed",rc,-48) .NE. 0)return do i=1,nDims rc = NF90_INQUIRE_DIMENSION (fid, i, dimName, dimSize) if (err("GetBegDateTime: can't get dim info",rc,-41) .NE. 0) return if (index(dimName,'station') .gt. 0) cycle if (trim(dimName) .eq. 'nv') cycle if ( index(dimName,'edges') .gt. 0 ) cycle if (dimName=='nf') cycle if (dimName=='orientationStrLen') cycle if (dimName=='ncontact') cycle if (trim(dimName) .eq. 'XCdim') cycle if (trim(dimName) .eq. 'YCdim') cycle rc = NF90_INQ_VARID (fid, dimName, dimId) if (err("GetBegDateTime: NF90_INQ_VARID failed",rc,-40) .NE. 0) return ! If it has the standard_name attribute, use that instead rc = NF90_GET_ATT(fid,dimId,'standard_name',stdName) if (rc .ne. 0) stdName = Trim(dimName) rc = NF90_GET_ATT(fid,dimId,'units',dimUnits) if (err("GetBegDateTime: could not get units for dimension",rc,-53)& .NE. 0) return if ( IdentifyDim (stdName, dimUnits) .eq. 3 ) then timeId = dimId timeUnits = dimUnits exit end if end do if ( timeId .lt. 0 ) then rc = -43 print *, "GetBegDateTime: could not find time coord" return end if ! Try assuming this file has been written with CFIO ! ------------------------------------------------- rc = NF90_GET_ATT(fid,timeId,'begin_date',begDate) if ( rc .eq. 0 ) then rc = NF90_GET_ATT(fid,timeId,'begin_time',begTime) end if ! Well, it must be a native CFIO file ! ----------------------------------- if ( rc .eq. 0 ) then rc = NF90_GET_ATT(fid,timeId,'time_increment',timInc) if (err("GetBegDateTime: missing time increment",rc,-44) .NE. 0) & return !ams write (strTmp,'(i6)') timinc !ams read (strTmp,'(3I2)') hour, min, sec call CFIO_parseIntTime ( timinc, hour, min, sec ) incSecs = hour*3600 + min*60 + sec !ams print *, 'begdate, begtime, incsecs: ',begdate, begtime, incsecs return ! all done. end if ! If could not find begin_date/begin_time attributes ! then this is not a native CFIO file. In this case ! attempt to parse the COARDS compliant time units ! -------------------------------------------------- !ams rc = NF90_GET_ATT(fid,timeId,'units',timeUnits) !ams if (err("GetBegDateTime: missing time.units",rc,-44) .NE. 0) return i = index(timeUnits,'since') if ( i .le. 0 ) then if (err("GetBegDateTime: invalid time units",1,-44) .NE. 0) return endif ! Call to ParseTimeUnits replaces an internal read, that made assumptions ! about the format of the Time Units string that were not always true. ! (RL: 10/2000) call ParseTimeUnits ( timeUnits, year, month, day, hour, min, sec, rc ) begDate = year*10000 + month*100 + day begTime = hour*10000 + min*100 + sec ! Determine time increment. ! ------------------------- rc = NF90_INQUIRE_VARIABLE (fid, timeID, varName, type, nvDims, vDims, & nvAtts) if (err("GetBegDateTime: error in time variable inquire",& rc,-52) .NE. 0) return if ( type .eq. NF90_FLOAT ) then corner(1) = 1 rc = NF90_GET_VAR(fid,timeID,rtime_array,corner,(/1/)) rtime = rtime_array(1) t1 = int(rtime) corner(1) = 2 rc = NF90_GET_VAR(fid,timeID,rtime_array,corner,(/1/)) rtime = rtime_array(1) t2 = int(rtime) else if ( type .eq. NF90_DOUBLE ) then corner(1) = 1 rc = NF90_GET_VAR(fid,timeID,dtime_array,corner,(/1/)) dtime = dtime_array(1) t1 = int(dtime) !ams print *, t1, dtime, rc corner(1) = 2 rc = NF90_GET_VAR(fid,timeID,dtime_array,corner,(/1/)) dtime = dtime_array(1) t2 = int(dtime) !ams print *, t2, dtime, rc else if ( type .eq. NF90_SHORT ) then corner(1) = 1 rc = NF90_GET_VAR(fid,timeID,itime_array,corner,(/1/)) itime = itime_array(1) t1 = itime corner(1) = 2 rc = NF90_GET_VAR(fid,timeID,itime_array,corner,(/1/)) itime = itime_array(1) t2 = itime else if ( type .eq. NF90_INT ) then corner(1) = 1 rc = NF90_GET_VAR(fid,timeID,ltime_array,corner,(/1/)) ltime = ltime_array(1) t1 = ltime corner(1) = 2 rc = NF90_GET_VAR(fid,timeID,ltime_array,corner,(/1/)) ltime = ltime_array(1) t2 = ltime else if (err("GetBegDateTime: invalid time data type",& 1,-44) .NE. 0) return endif ! Convert time increment to seconds if necessary ! ---------------------------------------------- incSecs = t2 - t1 if ( timeUnits(1:6) .eq. 'minute' ) then tMult = 60 else if ( timeUnits(1:4) .eq. 'hour' ) then tMult = 60 * 60 else if ( timeUnits(1:3) .eq. 'day' ) then tMult = 60 * 60 * 24 else if (err("GetBegDateTime: invalid time unit name",& 1,-44) .NE. 0) return endif incSecs = incSecs*tMult ! Combine the first time offset with the reference time to get the beginning ! date and time ! ----------------------------------------------------------------------------- Call GetDate(begDate,begTime,t1*tMult,newDate,newTime,rc) begDate=newDate begTime=newTime !ams print *, 'begdate, begtime, incsecs: ',begdate, begtime, incsecs incSecs = max ( 1, incSecs ) rc = 0 ! all done return end subroutine GetBegDateTime