ESMF_CFIOSdfFileOpen Subroutine

public subroutine ESMF_CFIOSdfFileOpen(cfio, fmode, rc, expid, cyclic)

ESMF_CFIOSdfFileOpen – open a CFIO file, and get CFIO meta data into a cfio Object.

Arguments

Type IntentOptional Attributes Name
type(ESMF_CFIO), intent(inout) :: cfio

a CFIO object

integer, intent(in) :: fmode

0 for READ-WRITE non-zero for READ-ONLY

integer, intent(out), optional :: rc

Error return code: 0 all is well -1 invalid count -2 type mismatch -12 error determining default precision -10 ngatts is incompatible with file -11 character string not long enough -19 unable to identify coordinate variable -36 error from NF90_PUT_ATT (global attribute) -39 error from ncopn (file open) -40 error from NF90_INQ_VARID -41 error from NF90_INQ_DIMID (lat or lon) -42 error from NF90_INQ_DIMID (lev) -43 error from NF90_INQ_VARID (time variable) -47 error from NF90_INQ_DIMID (time) -48 error from NF90_INQUIRE -51 error from NF90_GET_ATT (global attribute) -52 error from NF90_INQUIRE_VARIABLE -53 error from NF90_GET_ATT -57 error from NF90_INQ_ATTNAME -58 error from NF90_INQUIRE_ATTRIBUTE

character(len=*), intent(in), optional :: expid

Experiment ID

logical, intent(in), optional :: cyclic

cyclic input file


Calls

proc~~esmf_cfiosdffileopen~~CallsGraph proc~esmf_cfiosdffileopen ESMF_CFIOSdfFileOpen nf90_get_att nf90_get_att proc~esmf_cfiosdffileopen->nf90_get_att nf90_get_var nf90_get_var proc~esmf_cfiosdffileopen->nf90_get_var nf90_inq_attname nf90_inq_attname proc~esmf_cfiosdffileopen->nf90_inq_attname nf90_inq_varid nf90_inq_varid proc~esmf_cfiosdffileopen->nf90_inq_varid nf90_inquire nf90_inquire proc~esmf_cfiosdffileopen->nf90_inquire nf90_inquire_attribute nf90_inquire_attribute proc~esmf_cfiosdffileopen->nf90_inquire_attribute nf90_inquire_dimension nf90_inquire_dimension proc~esmf_cfiosdffileopen->nf90_inquire_dimension nf90_inquire_variable nf90_inquire_variable proc~esmf_cfiosdffileopen->nf90_inquire_variable proc~cfio_attinquire CFIO_AttInquire proc~esmf_cfiosdffileopen->proc~cfio_attinquire proc~cfio_diminquire CFIO_DimInquire proc~esmf_cfiosdffileopen->proc~cfio_diminquire proc~cfio_getattnames CFIO_GetAttNames proc~esmf_cfiosdffileopen->proc~cfio_getattnames proc~cfio_getcharatt CFIO_GetCharAtt proc~esmf_cfiosdffileopen->proc~cfio_getcharatt proc~cfio_getintatt CFIO_GetIntAtt proc~esmf_cfiosdffileopen->proc~cfio_getintatt proc~cfio_getrealatt CFIO_GetRealAtt proc~esmf_cfiosdffileopen->proc~cfio_getrealatt proc~cfio_open CFIO_Open proc~esmf_cfiosdffileopen->proc~cfio_open proc~err err proc~esmf_cfiosdffileopen->proc~err proc~esmf_cfiovarinfocreate ESMF_CFIOVarInfoCreate proc~esmf_cfiosdffileopen->proc~esmf_cfiovarinfocreate proc~getbegdatetime GetBegDateTime proc~esmf_cfiosdffileopen->proc~getbegdatetime proc~identifydim IdentifyDim proc~esmf_cfiosdffileopen->proc~identifydim proc~strtemplate_ strTemplate_ proc~esmf_cfiosdffileopen->proc~strtemplate_ proc~cfio_attinquire->nf90_inquire_attribute proc~cfio_attinquire->proc~err proc~cfio_diminquire->nf90_get_att proc~cfio_diminquire->nf90_inq_varid proc~cfio_diminquire->nf90_inquire proc~cfio_diminquire->nf90_inquire_dimension proc~cfio_diminquire->nf90_inquire_variable proc~cfio_diminquire->proc~err proc~cfio_diminquire->proc~identifydim proc~cfio_getattnames->nf90_inq_attname proc~cfio_getattnames->nf90_inquire proc~cfio_getattnames->proc~err proc~cfio_getcharatt->nf90_get_att proc~cfio_getcharatt->nf90_inquire_attribute proc~cfio_getcharatt->proc~err proc~cfio_getintatt->nf90_get_att proc~cfio_getintatt->nf90_inquire_attribute proc~cfio_getintatt->proc~err proc~cfio_getrealatt->nf90_get_att proc~cfio_getrealatt->nf90_inquire_attribute proc~cfio_getrealatt->proc~err proc~cfio_open->proc~err ncopn ncopn proc~cfio_open->ncopn proc~esmf_cfiogridcreate ESMF_CFIOGridCreate proc~esmf_cfiovarinfocreate->proc~esmf_cfiogridcreate proc~getbegdatetime->nf90_get_att proc~getbegdatetime->nf90_get_var proc~getbegdatetime->nf90_inq_varid proc~getbegdatetime->nf90_inquire proc~getbegdatetime->nf90_inquire_dimension proc~getbegdatetime->nf90_inquire_variable proc~getbegdatetime->proc~err proc~getbegdatetime->proc~identifydim interface~getdate GetDate proc~getbegdatetime->interface~getdate proc~cfio_parseinttime CFIO_parseIntTime proc~getbegdatetime->proc~cfio_parseinttime proc~parsetimeunits ParseTimeUnits proc~getbegdatetime->proc~parsetimeunits proc~gx_ GX_ proc~strtemplate_->proc~gx_ proc~getdateint4 GetDateInt4 interface~getdate->proc~getdateint4 proc~getdateint8 GetDateInt8 interface~getdate->proc~getdateint8 proc~getdateint4->proc~getdateint8 proc~getdateint8->proc~cfio_parseinttime proc~caldat CALDAT proc~getdateint8->proc~caldat proc~julday julday proc~getdateint8->proc~julday

Called by

proc~~esmf_cfiosdffileopen~~CalledByGraph proc~esmf_cfiosdffileopen ESMF_CFIOSdfFileOpen proc~esmf_cfiofileopen ESMF_CFIOFileOpen proc~esmf_cfiofileopen->proc~esmf_cfiosdffileopen none~find~4 CFIOCollection%find none~find~4->proc~esmf_cfiofileopen proc~mapl_cfioopenwrite MAPL_CFIOOpenWrite proc~mapl_cfioopenwrite->proc~esmf_cfiofileopen program~test~11 test program~test~11->proc~esmf_cfiofileopen program~test~2 test program~test~2->proc~esmf_cfiofileopen program~test~4 test program~test~4->proc~esmf_cfiofileopen program~test~6 test program~test~6->proc~esmf_cfiofileopen proc~mapl_cfiocreatefromfile MAPL_CFIOCreateFromFile proc~mapl_cfiocreatefromfile->none~find~4 proc~mapl_cfioreadbundleread MAPL_CFIOReadBundleRead proc~mapl_cfioreadbundleread->none~find~4

Source Code

      subroutine ESMF_CFIOSdfFileOpen (cfio, fmode, rc, expid, cyclic)
!
! !INPUT PARAMETERS:
!
      integer, intent(in) :: fmode              !! 0 for READ-WRITE
                                                !! non-zero for READ-ONLY
      character(len=*), intent(in), OPTIONAL :: expid   !! Experiment ID
      logical, intent(in), OPTIONAL :: cyclic           !! cyclic input file
!
! !OUTPUT PARAMETERS:
!
      integer, intent(out), OPTIONAL :: rc      !! Error return code:
                                                !! 0   all is well
                         !! -1   invalid count
                         !! -2   type mismatch
                         !! -12  error determining default precision
                         !! -10  ngatts is incompatible with file
                         !! -11  character string not long enough
                         !! -19  unable to identify coordinate variable
                         !! -36  error from NF90_PUT_ATT (global attribute)
                         !! -39  error from ncopn (file open)
                         !! -40  error from NF90_INQ_VARID
                         !! -41  error from NF90_INQ_DIMID (lat or lon)
                         !! -42  error from NF90_INQ_DIMID (lev)
                         !! -43  error from NF90_INQ_VARID (time variable)
                         !! -47  error from NF90_INQ_DIMID (time)
                         !! -48  error from NF90_INQUIRE
                         !! -51  error from NF90_GET_ATT (global attribute)
                         !! -52  error from NF90_INQUIRE_VARIABLE
                         !! -53  error from NF90_GET_ATT
                         !! -57  error from NF90_INQ_ATTNAME
                         !! -58  error from NF90_INQUIRE_ATTRIBUTE

!
! !INPUT/OUTPUT PARAMETERS:
!
      type(ESMF_CFIO), intent(inout) :: cfio    !! a CFIO object
!
!------------------------------------------------------------------------------
      integer :: ngatts, lm, i, ii, iv
      real(kind=REAL32) :: amiss
      real(kind=REAL32) :: vRange32(2)
      real(kind=REAL32), pointer :: lon(:), lat(:)
      real(kind=REAL32), pointer :: lev(:) => null()
      real(kind=REAL64), pointer :: lon_64(:), lat_64(:), lev_64(:)
      integer :: coXType = NF90_FLOAT
      integer :: coYType = NF90_FLOAT
      integer :: coZType = NF90_FLOAT
      character(len=MVARLEN) :: vAttName
      character(len=MLEN), pointer :: attNames(:)
      integer :: iCnt, rCnt, cCnt
      integer :: iMaxLen, rMaxLen, cMaxLen
      integer :: type, count, rtcode
      integer :: varId
      integer :: datatype         ! variable type
      integer :: vtype            ! variable type
      integer :: nvDims           ! number of dimensions
      integer :: vDims(MAXVDIMS)  ! variable shape
      integer :: nvatts           ! number of attributes
      real(kind=REAL32), pointer :: rtmp(:)
      integer, pointer :: itmp(:)
      logical :: new_grid
      integer :: nDims, allVars, recdim
      integer :: im, jm, km
      integer :: hour, minute, seconds
      integer :: fid, nVars, dimSize(8), myIndex
      character(len=MVARLEN) :: dimName(8), dimUnits(8), vnameTemp
      character(len=MVARLEN) :: nameAk, nameBk, namePtop
      character(len=MVARLEN) :: levunits
      integer :: loc1, loc2
      integer :: akid, bkid, ptopid
      integer :: icount
      integer :: vdir
      real(kind=REAL32), pointer :: ak(:), bk(:)
      real(kind=REAL32) :: ptop
      real(kind=REAL32) :: ptop_array(1)
      real(kind=REAL32) :: scale, offset
      character, pointer ::  globalAtt(:)
      character(len=MLEN) :: fNameTmp     ! file name
      character(len=MVARLEN) :: fversion
      logical :: cs_found
      integer :: nf

      fNameTmp = ''
!     checking file name template
      if (present(expid)) cfio%expid = expid
      if (present(cyclic)) cfio%isCyclic = cyclic
      if (present(expid) .and. cfio%date .gt. 0 .and. cfio%begTime .ge. 0) then
         call strTemplate_(fNameTmp,cfio%fName,xid=expid,nymd=cfio%date, &
                           nhms=cfio%begTime, stat=rtcode)
      else
         if (cfio%date .gt. 0 .and. cfio%begTime .ge. 0) then
            call strTemplate_(fNameTmp,cfio%fName,nymd=cfio%date, &
                              nhms=cfio%begTime, stat=rtcode)
         else
            if (present(expid)) then
               call strTemplate_(fNameTmp,cfio%fName,xid=expid, stat=rtcode)
            end if
         end if
      end if
      if (trim(fNameTmp) .ne. trim(cfio%fName) .and. len(trim(fNameTmp)) .gt. 0) then
         cfio%fNameTmplt = cfio%fName
         cfio%fName = fNameTmp
      end if

!     open a cfio file
      call CFIO_Open ( cfio%fName, fmode, cfio%fid, rtcode )
      if (err("problem in CFIO_Open",rtcode,rtcode) .lt. 0 ) then
         if ( present(rc) ) rc = rtcode
         return
      end if
      cfio%isOpen = .true.
      if (fmode == 0) then
         rc = 0
         return
      endif
      fid =cfio%fid

!     get grid information and global meta data

      call CFIO_DimInquire (cfio%fid, im, jm, km, lm, &
                            cfio%mVars, ngatts, vdir=vdir, rc=rtcode)
      if (err("CFIO_DimInquire failed",rtcode,rtcode) .lt. 0) then
         if ( present(rc) ) rc = rtcode
         return
      end if
      cfio%tSteps = lm
      cfio%vDir = vdir

      rtcode = NF90_INQUIRE (cfio%fid,nDims,allVars,ngatts,recdim)
      if (err("FileOpen: NF90_INQUIRE failed",rtcode,-48) .NE. 0) then
         if ( present(rc) ) rc = rtcode
         return
      end if

      allocate(cfio%varObjs(cfio%mVars))
      do i=1,cfio%mVars
         cfio%varObjs(i)=ESMF_CFIOVarInfoCreate(rc=rc)
         cfio%varObjs(i)%timAve = .false.
      end do
      nVars = 0
      cfio%mGrids = 0
      cs_found = .false.
      do i=1,allVars
        rtcode = NF90_INQUIRE_VARIABLE(fid,i,vnameTemp,vtype,nvDims,vDims,nvAtts)
        if (err("Inquire: variable inquire error",rtcode,-52) .NE. 0) then
           if ( present(rc) ) rc = rtcode
           return
        end if
        if (nvDims .EQ. 1 .and. (index(vnameTemp, 'lon') .gt. 0 .or.  &
            vnameTemp == 'Xdim' .or. &
            index(vnameTemp, 'XDim:EOSGRID') .gt. 0) ) then
           coXType = vtype
           cfio%mGrids = cfio%mGrids + 1
        end if
        if (nvDims .EQ. 1 .and. (index(vnameTemp, 'lat') .gt. 0 .or.  &
            index(vnameTemp, 'YDim:EOSGRID') .gt. 0) ) then
           coYType = vtype
        end if
        if (nvDims .EQ. 1 .and. (index(vnameTemp, 'lev') .gt. 0 .or.  &
            index(vnameTemp, 'Height:EOSGRID') .gt. 0) ) then
           coZType = vtype
        end if
        if(vnameTemp == 'nf') then
           cs_found = .true.
           cycle
        endif

      end do
      if (cs_found) then
         fversion="2.91"
         read(fversion, *) cfio%formatVersion
      end if

      do i=1,allVars
        rtcode = NF90_INQUIRE_VARIABLE(fid,i,vnameTemp,vtype,nvDims,vDims,nvAtts)
        if (err("Inquire: variable inquire error",rtcode,-52) .NE. 0) then
           if ( present(rc) ) rc = rtcode
           return
        end if
        if (trim(vnameTemp) .eq. 'time_bnds') then
           cfio%varObjs(nVars)%timAve = .true.
           cycle
        end if
        if (cs_found) then
           if(vnameTemp == 'ncontact') cycle
           if(vnameTemp == 'lons') cycle
           if(vnameTemp == 'lats') cycle
           if(vnameTemp == 'contacts') cycle
           if(vnameTemp == 'orientation') cycle
           if(vnameTemp == 'anchor') cycle
           if(vnameTemp == 'corner_lons') cycle
           if(vnameTemp == 'corner_lats') cycle
        end if
        if (nvDims .EQ. 1) cycle
        nVars = nVars + 1
        cfio%varObjs(nVars)%vName = trim(vnameTemp)
        cfio%varObjs(nVars)%grid%km = 0
!        cfio%varObjs(nVars)%grid%km = 1
        cfio%varObjs(nVars)%grid%stnGrid = .false.
        do iv = 1, nvDims
           rtcode = NF90_INQUIRE_DIMENSION(fid, vDims(iv), dimName(iv), dimSize(iv))
           if (err("problem in NF90_INQUIRE_DIMENSION",rtcode,-41) .NE. 0) then
              if ( present(rc) ) rc = rtcode
              return
           end if
           if (index(dimName(iv),'station') .gt. 0) then
              cfio%varObjs(nVars)%grid%im = dimSize(iv)
              cfio%varObjs(nVars)%grid%jm = dimSize(iv)
              cfio%varObjs(nVars)%grid%stnGrid = .true.
              cycle
           end if
           if (cs_found) then
              if (dimName(iv)=='nf') cycle
              if (dimName(iv)=='orientationStrLen') cycle
              if (dimName(iv)=='ncontact') cycle
              if (dimName(iv)=='XCdim') cycle
              if (dimName(iv)=='YCdim') cycle
           end if
           rtcode = NF90_INQ_VARID(fid,dimName(iv),varId)
           dimUnits(iv) = ' '
           rtcode = NF90_GET_ATT(fid,varId,'units',dimUnits(iv))
           if (err("problem in NF90_GET_ATT",rtcode,-53) .NE. 0) then
              if ( present(rc) ) rc = rtcode
              return
           end if
           myIndex = IdentifyDim (dimName(iv), dimUnits(iv))
           if (myIndex .EQ. 0) then
              cfio%varObjs(nVars)%grid%im = dimSize(iv)
              if (.not. associated(cfio%varObjs(nVars)%grid%lon)) then
                 allocate(cfio%varObjs(nVars)%grid%lon(dimSize(iv)))
              end if
              allocate(lon(dimSize(iv)))
              if ( coXType .eq. NF90_FLOAT ) then
                 rtcode =  NF90_GET_VAR (fid, varId, lon, (/1/), (/dimSize(iv)/))
              else
                 allocate(lon_64(dimSize(iv)))
                 rtcode = NF90_GET_VAR (fid, varId, lon_64, (/1/), (/dimSize(iv)/))
                 lon =lon_64
                 deallocate(lon_64)
              end if
              if (err("problem in NF90_GET_VAR",rtcode,-53) .NE. 0) then
                 if ( present(rc) ) rc = rtcode
                 return
              end if
              cfio%varObjs(nVars)%grid%lon = lon
              deallocate(lon)
           end if
           nf = 1
           if (cs_found) then
              nf = 6
           end if
           if (myIndex .EQ. 1) then
              cfio%varObjs(nVars)%grid%jm = dimSize(iv)*nf
              if (.not. associated(cfio%varObjs(nVars)%grid%lat)) then
                 allocate(cfio%varObjs(nVars)%grid%lat(dimSize(iv)))
              end if
              allocate(lat(dimSize(iv)))
              if ( coYType .eq. NF90_FLOAT ) then
                 rtcode = NF90_GET_VAR(fid, varId, lat, (/1/), (/dimSize(iv)/))
              else
                 allocate(lat_64(dimSize(iv)))
                 rtcode = NF90_GET_VAR(fid, varId, lat_64, (/1/), (/dimSize(iv)/))
                 lat = lat_64
                 deallocate(lat_64)
              end if
!print *, "vDims(iv) varId: ", vDims(iv), varId
!print *, "dimName dimUnits: ", trim(dimName(iv)), trim(dimUnits(iv))
              if (err("problem in NF90_GET_VAR",rtcode,-51) .NE. 0) then
                 if ( present(rc) ) rc = rtcode
                 return
              end if
              cfio%varObjs(nVars)%grid%lat = lat
              deallocate(lat)
           end if
           if (myIndex .EQ. 2) then
              cfio%varObjs(nVars)%grid%km = dimSize(iv)
              rtcode = NF90_GET_ATT(fid,varId,'standard_name',cfio%varObjs(nVars)%grid%standardName)
              if (rtcode /= 0) cfio%varObjs(nVars)%grid%standardName="pressure"
              if ( index(cfio%varObjs(nVars)%grid%standardName,        &
                   'atmosphere_sigma_coordinate') .gt. 0  .or.         &
                   index(cfio%varObjs(nVars)%grid%standardName,        &
                   'atmosphere_hybrid_sigma_pressure_coordinate' )     &
                   .gt.  0 ) then

                 rtcode = NF90_GET_ATT(fid,varId,'formula_term',cfio%varObjs(nVars)%grid%formulaTerm)
                 if ( index(cfio%varObjs(nVars)%grid%standardName,     &
                   'atmosphere_sigma_coordinate') .gt. 0 ) then
                    loc1 = index(cfio%varObjs(nVars)%grid%formulaTerm,'ptop:')
                    icount = loc1  + 5
                    do icount = loc1+5, len(cfio%varObjs(nVars)%grid%formulaTerm)
                      if (cfio%varObjs(nVars)%grid%formulaTerm(icount:icount) &
                           .ne. ' ') exit
                    end do
                    namePtop=trim(cfio%varObjs(nVars)%grid%formulaTerm    &
                         (icount:len(cfio%varObjs(nVars)%grid%formulaTerm)))
                    rtcode = NF90_INQ_VARID(cfio%fid,trim(namePtop),ptopid)
                    if (rtcode .ne. 0) print *, "problem in getting ptopid in NF90_INQ_VARID"
                    if (rtcode .eq. 0) then
                       rtcode = NF90_GET_VAR(cfio%fid,ptopid,ptop_array,(/1/), (/1/))
                       ptop=ptop_array(1)
                    end if
                    if (rtcode .eq. 0) cfio%varObjs(nVars)%grid%ptop = ptop
                 end if
              end if
              if (index(cfio%varObjs(nVars)%grid%standardName,             &
                        'atmosphere_hybrid_sigma_pressure_coordinate')     &
                        .gt. 0)  then
                 loc1 = index(cfio%varObjs(nVars)%grid%formulaTerm,'a:')
                 loc2 = index(cfio%varObjs(nVars)%grid%formulaTerm,'b:')
                 icount = 0
                 do icount = loc1+2, loc2
                   if (cfio%varObjs(nVars)%grid%formulaTerm(icount:icount) &
                        .ne. ' ') exit
                 end do
                 nameAk=trim(cfio%varObjs(nVars)%grid%formulaTerm          &
                             (icount:loc2-1))
                 loc1 = index(cfio%varObjs(nVars)%grid%formulaTerm,'b:')
                 loc2 = index(cfio%varObjs(nVars)%grid%formulaTerm,'ps:')
                 do icount = loc1+2, loc2
                   if (cfio%varObjs(nVars)%grid%formulaTerm(icount:icount) &
                        .ne. ' ') exit
                 end do
                 nameBk=trim(cfio%varObjs(nVars)%grid%formulaTerm          &
                             (icount:loc2-1))
                 loc1 = index(cfio%varObjs(nVars)%grid%formulaTerm,'p0:')
                 icount = loc1  + 4
                 namePtop=trim(cfio%varObjs(nVars)%grid%formulaTerm        &
                         (icount:len(cfio%varObjs(nVars)%grid%formulaTerm)))

                 rtcode = NF90_INQ_VARID(cfio%fid, trim(nameAk), akid)
                 if (rtcode .ne. 0) print *, "problem in getting akid in NF90_INQ_VARID"

                 allocate(cfio%varObjs(nVars)%grid%ak                      &
                          (cfio%varObjs(nVars)%grid%km+1),                 &
                          ak(cfio%varObjs(nVars)%grid%km+1))
                 rtcode = NF90_GET_VAR(cfio%fid,akid,ak,(/1/),(/cfio%varObjs(nVars)%grid%km+1/))
                 if (rtcode .ne. 0) print *, "problem in getting ak in NF90_GET_VAR"
                 cfio%varObjs(nVars)%grid%ak = ak
                 deallocate(ak)
                 rtcode = NF90_INQ_VARID(cfio%fid, trim(nameBk), bkid)
                 if (rtcode .ne. 0) print *, "problem in getting bkid in NF90_INQ_VARID"
                 allocate(cfio%varObjs(nVars)%grid%bk                      &
                          (cfio%varObjs(nVars)%grid%km+1),                 &
                          bk(cfio%varObjs(nVars)%grid%km+1))
                 rtcode = NF90_GET_VAR(cfio%fid,bkid,bk,(/1/),(/cfio%varObjs(nVars)%grid%km+1/))
                 if (rtcode .ne. 0) print *, "problem in getting bk in NF90_GET_VAR"
                 cfio%varObjs(nVars)%grid%bk = bk
                 deallocate(bk)

                 rtcode = NF90_INQ_VARID(cfio%fid, trim(namePtop), ptopid)
                 if (rtcode .ne. 0) print *, "problem in getting ptopid in NF90_INQ_VARID"
                 rtcode = NF90_GET_VAR(cfio%fid,ptopid,ptop_array,(/1/), (/1/))
                 ptop = ptop_array(1)
                 if (rtcode .ne. 0) print *, "problem in getting ptop in NF90_GET_VAR"
                 cfio%varObjs(nVars)%grid%ptop = ptop
             end if
              rtcode = NF90_GET_ATT(fid,varId,'coordinate',cfio%varObjs(nVars)%grid%coordinate)
              if (rtcode .ne. 0) cfio%varObjs(nVars)%grid%coordinate = "pressure"
              cfio%varObjs(nVars)%grid%levUnits = trim(dimUnits(iv))

              allocate(cfio%varObjs(nVars)%grid%lev(dimSize(iv)))
              if (.not.associated(lev))  allocate(lev(dimSize(iv)))

              if ( coZType .eq. NF90_FLOAT ) then
                 rtcode = NF90_GET_VAR(fid, varId, lev, (/1/), (/dimSize(iv)/))
!print *, "Lev from CFIO SDFFileOpen: ", lev
              else
                 allocate(lev_64(dimSize(iv)))
                 rtcode = NF90_GET_VAR(fid, varId, lev_64, (/1/), (/dimSize(iv)/))
                 lev =lev_64
                 deallocate(lev_64)
              end if
              cfio%varObjs(nVars)%grid%lev = 0.0
              cfio%varObjs(nVars)%grid%lev = lev
!print *, "cfio%varObjs(nVars)%grid%lev from CFIO SDFFileOpen: ", cfio%varObjs(nVars)%grid%lev
              rtcode = NF90_GET_ATT(fid,varId,'units',levunits)
           end if
        end do
        rtcode = NF90_INQ_VARID (cfio%fid, cfio%varObjs(nVars)%vName, varId)
        if (rtcode .ne. 0) then
           print *, "problem in getting varId in NF90_INQ_VARID"
           if ( present(rc) ) rc = -40
           return
        end if
       rtcode = NF90_GET_ATT(fid,varId,'units',cfio%varObjs(nVars)%vunits)
        if (rtcode .ne. 0) then
           print *, "NF90_GET_ATT failed for units"
           if ( present(rc) ) rc = -53
           return
        end if
        cfio%varObjs(nVars)%vtitle = ' '
       rtcode = NF90_GET_ATT(fid,varId,'long_name',cfio%varObjs(nVars)%vtitle)
       rtcode = NF90_GET_ATT(fid,varId,'standard_name',cfio%varObjs(nVars)%standardName)
        if ( cfio%varObjs(nVars)%grid%km .gt. 0 ) then
            cfio%varObjs(nVars)%twoDimVar = .false.
        else
            cfio%varObjs(nVars)%twoDimVar = .true.
        end if
        rtcode = NF90_GET_ATT(fid,varId,'_FillValue',amiss)
        if (rtcode .NE. 0) then
           rtcode = NF90_GET_ATT(fid,varId,'missing_value',amiss)
        end if
        cfio%varObjs(nVars)%amiss = amiss
        rtcode = NF90_GET_ATT(fid,varId,'scale_factor',scale)
        if (rtcode .NE. 0) then
           cfio%varObjs(nVars)%scaleFactor = 1.0
        else
           cfio%varObjs(nVars)%scaleFactor = scale
        end if
        rtcode = NF90_GET_ATT(fid,varId,'add_offset',offset)
        if (rtcode .NE. 0) then
           cfio%varObjs(nVars)%addOffset = 0.0
        else
           cfio%varObjs(nVars)%addOffset = offset
        end if
        rtcode = NF90_GET_ATT(fid,varId,'vmin',vRange32(1))
        if (rtcode .NE. 0) then
          cfio%varObjs(nVars)%validRange(1) = cfio%varObjs(nVars)%amiss
        else
          cfio%varObjs(nVars)%validRange(1) = vRange32(1)
        endif
        rtcode = NF90_GET_ATT(fid,varId,'vmax',vRange32(2))
        if (rtcode .NE. 0) then
          cfio%varObjs(nVars)%validRange(2) = cfio%varObjs(nVars)%amiss
        else
          cfio%varObjs(nVars)%validRange(2) = vRange32(2)
        endif

      end do

      call GetBegDateTime(fid,cfio%date,cfio%begTime,cfio%timeInc,rtcode)
      if (rtcode .ne. 0) then
         print *, "GetBegDateTime failed to get data/time/timeInc"
         if ( present(rc) ) rc = rtcode
         return
      end if

      hour = cfio%timeInc/3600
      minute = (cfio%timeInc-(3600*hour))/60
      seconds = cfio%timeInc-(3600*hour) - (60 * minute)
      cfio%timeInc = hour*10000 + minute*100 + seconds

      allocate(attNames(ngatts))
      attNames = " "
      call CFIO_GetAttNames ( cfio%fid, ngatts, attNames, rtcode )
      if (err("CFIO_GetAttNames failed",rtcode,rtcode) .lt. 0) then
         if ( present(rc) ) rc = rtcode
         return
      end if

      iCnt = 0
      rCnt = 0
      cCnt = 0
      iMaxLen = 0
      rMaxLen = 0
      cMaxLen = 0

!     get how many int/real/char attributes in attNames
      do i =1, ngatts
         call CFIO_AttInquire (cfio%fid, attNames(i), type, count, rtcode)
         if (err("CFIO_AttInquire failed",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = rtcode
            return
         end if
         select case  (type)
            case ( 0 )
               iCnt = iCnt + 1
               if ( count .gt. iMaxLen ) iMaxLen = count
            case ( 1 )
               rCnt = rCnt + 1
               if ( count .gt. rMaxLen ) rMaxLen = count
            case ( 2 )
               cCnt = cCnt + 1
               if ( count .gt. cMaxLen ) cMaxLen = count
            case ( 3 )
               rCnt = rCnt + 1
               if ( count .gt. rMaxLen ) rMaxLen = count
            case ( 4 )
               iCnt = iCnt + 1
               if ( count .gt. iMaxLen ) iMaxLen = count
         end select
      end do

      cfio%nAttChar = cCnt
      cfio%nAttReal = rCnt
      cfio%nAttInt = iCnt

      allocate(cfio%attCharCnts(cCnt), cfio%attRealCnts(rCnt), &
               cfio%attIntCnts(iCnt))
      allocate(cfio%attCharNames(cCnt), cfio%attRealNames(rCnt), &
               cfio%attIntNames(iCnt))

      iCnt = 0
      rCnt = 0
      cCnt = 0
!     get attNames and count, then put them into a cfio obj
      do i =1, ngatts
         call CFIO_AttInquire (cfio%fid, attNames(i), type, count, rtcode)
         if (err("CFIO_AttInquire failed",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = rtcode
            return
         end if
         select case  (type)
            case ( 0 )
               iCnt = iCnt + 1
               cfio%attIntNames(iCnt) = attNames(i)
               cfio%attIntCnts(iCnt) = count
            case ( 1 )
               rCnt = rCnt + 1
               cfio%attRealNames(rCnt) = attNames(i)
               cfio%attRealCnts(rCnt) = count
            case ( 2 )
               cCnt = cCnt + 1
               cfio%attCharNames(cCnt) = attNames(i)
               cfio%attCharCnts(cCnt) = count
            case ( 3 )
               rCnt = rCnt + 1
               cfio%attRealNames(rCnt) = attNames(i)
               cfio%attRealCnts(rCnt) = count
            case ( 4 )
               iCnt = iCnt + 1
               cfio%attIntNames(iCnt) = attNames(i)
               cfio%attIntCnts(iCnt) = count
         end select
      end do

      deallocate(attNames)

      allocate(cfio%attReals(rCnt, rMaxLen), cfio%attInts(iCnt, iMaxLen),    &
               cfio%attChars(cCnt))
!     get global integer attributes
      do i = 1, iCnt
         call CFIO_GetIntAtt(cfio%fid,cfio%attIntNames(i),cfio%attIntCnts(i) &
                            , cfio%attInts(i,:), rtcode)
         if (err("CFIO_GetIntAtt failed",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = rtcode
            return
         end if
      end do

!     get global real attributes
      do i = 1, rCnt
         call CFIO_GetRealAtt(cfio%fid,cfio%attRealNames(i),               &
                              cfio%attRealCnts(i),                         &
                              cfio%attReals(i,:), rtcode)
         if (err("CFIO_GetRealAtt",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = rtcode
            return
         end if
      end do

!     get global char attributes
      do i = 1, cCnt
         allocate(globalAtt(cfio%attCharCnts(i)))
         call CFIO_GetCharAtt(cfio%fid,cfio%attCharNames(i),   &
                              cfio%attCharCnts(i),             &
                              globalAtt, rtcode)
         if (err("GetCharAtt",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = rtcode
            return
         end if
!        cfio%attChars(i) can only hold MLEN characters.
         do ii = 1, cfio%attCharCnts(i)
            cfio%attChars(i)(ii:ii) = globalAtt(ii)
            if (ii .ge. MLEN) then
               !bma commented out this warning as merra1 hdf4 files
               !have really long global attributes which prints too
               !many warnings
               !print *,"global attribute ",trim(cfio%attCharNames(i)), &
                       !" is longer than MLEN"
               exit
            end if
         end do
         cfio%attChars(i)(cfio%attCharCnts(i)+1:MLEN) = ' '

         if (index(cfio%attCharNames(i),'History') .gt. 0)  &
            cfio%History=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'Source') .gt. 0)  &
            cfio%source=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'Title') .gt. 0)  &
            cfio%title=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'Contact') .gt. 0)  &
            cfio%contact=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'Conventions') .gt. 0)  &
            cfio%convention=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'Institution') .gt. 0)  &
            cfio%institution=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'References') .gt. 0)  &
            cfio%references=cfio%attChars(i)
         if (index(cfio%attCharNames(i),'Comment') .gt. 0)  &
            cfio%comment=cfio%attChars(i)
         deallocate(globalAtt)
      end do


!     get variable meta data
      do i = 1, cfio%mVars
         rtcode = NF90_INQ_VARID (cfio%fid, cfio%varObjs(i)%vName, varId)
         if (err("NF90_INQ_VARID failed for vName",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = -40
            return
         end if
         rtcode = NF90_INQUIRE_VARIABLE(cfio%fid, varId, cfio%varObjs(i)%vName, datatype, &
                     nvdims, vdims, nvatts)
         if (err("NF90_INQUIRE_VARIABLE failed for vName",rtcode,rtcode) .lt. 0) then
            if ( present(rc) ) rc = -52
            return
         end if
         iCnt = 0
         rCnt = 0
         cCnt = 0
         iMaxLen = 0
         rMaxLen = 0
         cMaxLen = 0

!        get variable int/real/char attribute count
         do iv =1, nvatts
            rtcode = NF90_INQ_ATTNAME(cfio%fid,varId,iv, vAttName)
            if (err("NF90_INQ_ATTNAME failed for vName",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -57
               return
            end if
            rtcode = NF90_INQUIRE_ATTRIBUTE (cfio%fid,varId,vAttName,vtype,count)
            if (err("NF90_INQUIRE_ATTRIBUTE failed for vName",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -58
               return
            end if
            select case  (vtype)
               case ( NF90_SHORT )
                  iCnt = iCnt + 1
                  if ( count .gt. iMaxLen ) iMaxLen = count
               case ( NF90_FLOAT )
                  rCnt = rCnt + 1
                  if ( count .gt. rMaxLen ) rMaxLen = count
               case ( NF90_CHAR )
                  cCnt = cCnt + 1
                  if ( count .gt. cMaxLen ) cMaxLen = count
               case ( NF90_DOUBLE )
                  rCnt = rCnt + 1
                  if ( count .gt. rMaxLen ) rMaxLen = count
               case ( NF90_INT )
                  iCnt = iCnt + 1
                  if ( count .gt. iMaxLen ) iMaxLen = count
            end select
         end do

         cfio%varObjs(i)%nVarAttChar = cCnt
         cfio%varObjs(i)%nVarAttReal = rCnt
         cfio%varObjs(i)%nVarAttInt = iCnt

         allocate(cfio%varObjs(i)%attCharCnts(cCnt),  &
                  cfio%varObjs(i)%attRealCnts(rCnt),  &
                  cfio%varObjs(i)%attIntCnts(iCnt))
         allocate(cfio%varObjs(i)%attCharNames(cCnt), &
                  cfio%varObjs(i)%attRealNames(rCnt),&
                  cfio%varObjs(i)%attIntNames(iCnt))

         iCnt = 0
         rCnt = 0
         cCnt = 0
!        get variable int/real/char attribute names and counts
         do iv =1, nvatts
            rtcode = NF90_INQ_ATTNAME (cfio%fid, varId, iv, vAttName)
            if (err("NF90_INQ_ATTNAME failed for vName",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -57
               return
            end if
            rtcode = NF90_INQUIRE_ATTRIBUTE (cfio%fid,varId,vAttName,vtype,count)
            if (err("NF90_INQUIRE_ATTRIBUTE failed for vName",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -58
               return
            end if
            select case  (vtype)
               case ( NF90_SHORT )
                  iCnt = iCnt + 1
                  cfio%varObjs(i)%attIntNames(iCnt) = vAttName
                  cfio%varObjs(i)%attIntCnts(iCnt) = count
               case ( NF90_FLOAT )
                  rCnt = rCnt + 1
                  cfio%varObjs(i)%attRealNames(rCnt) = vAttName
                  cfio%varObjs(i)%attRealCnts(rCnt) = count
               case ( NF90_CHAR )
                  cCnt = cCnt + 1
                  cfio%varObjs(i)%attCharNames(cCnt) = vAttName
                  cfio%varObjs(i)%attCharCnts(cCnt) = count
               case ( NF90_DOUBLE )
                  rCnt = rCnt + 1
                  cfio%varObjs(i)%attRealNames(rCnt) = vAttName
                  cfio%varObjs(i)%attRealCnts(rCnt) = count
               case ( NF90_INT )
                  iCnt = iCnt + 1
                  cfio%varObjs(i)%attIntNames(iCnt) = vAttName
                  cfio%varObjs(i)%attIntCnts(iCnt) = count
            end select
         end do

         allocate(cfio%varObjs(i)%varAttReals(rCnt, rMaxLen), &
                  cfio%varObjs(i)%varAttInts(iCnt, iMaxLen),  &
                  cfio%varObjs(i)%varAttChars(cCnt))

!        get int variable attributes
         do ii = 1, iCnt
            allocate(itmp(cfio%varObjs(i)%attIntCnts(ii)))
            rtcode = NF90_GET_ATT(cfio%fid,varId,cfio%varObjs(i)%attIntNames(ii),itmp)
            if (err("NF90_GET_ATT failed for attIntNames",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -53
               return
            end if
            cfio%varObjs(i)%varAttInts(ii,1:cfio%varObjs(i)%attIntCnts(ii))&
                       = itmp
            deallocate(itmp)
         end do

!        get real variable attributes
         do ii = 1, rCnt
            allocate(rtmp(cfio%varObjs(i)%attRealCnts(ii)))
            rtcode = NF90_GET_ATT(cfio%fid,varId,cfio%varObjs(i)%attRealNames(ii),rtmp)
            if (err("NF90_GET_ATT failed for attRealNames",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -53
               return
            end if
            cfio%varObjs(i)%varAttReals(ii,1:cfio%varObjs(i)%attRealCnts(ii))&
                       = rtmp
            deallocate(rtmp)
         end do

!        get char variable attributes
         do ii = 1, cCnt
            rtcode = NF90_GET_ATT(cfio%fid,varId,cfio%varObjs(i)%attCharNames(ii),cfio%varObjs(i)%varAttChars(ii))
            if (err("NF90_GET_ATT failed for attCharNames",rtcode,rtcode) .lt. 0) then
               if ( present(rc) ) rc = -53
               return
            end if
            cfio%varObjs(i)%varAttChars(ii)  &
                 (cfio%varObjs(i)%attCharCnts(ii)+1:MLEN) = ' '
         end do

      end do

!     set grids objects in a CFIO object
      allocate( cfio%grids(cfio%mGrids), stat = rtcode)
      cfio%grids(1) = cfio%varObjs(1)%grid
      if (associated(lev)) then
         cfio%grids(1)%levunits = levUnits
         allocate(cfio%grids(1)%lev(size(lev)),stat = rtcode)
         cfio%grids(1)%lev = lev
         deallocate(lev)
      end if
      if ( cfio%mGrids .eq. 1 .and. cfio%varObjs(1)%grid%km .eq. 0) &
         cfio%grids(1)%km = km

      if ( cfio%mGrids .gt. 1 ) then
        do i = 2, cfio%mGrids
           iCnt = 1
           do iv = 2, cfio%mVars
              new_grid = .true.
              iCnt = iCnt + 1
              do ii = 2, i
                if (cfio%varObjs(iv)%grid%im .eq. cfio%grids(ii-1)%im .and.  &
                  cfio%varObjs(iv)%grid%jm .eq. cfio%grids(ii-1)%jm .and.  &
                  cfio%varObjs(iv)%grid%km .eq. cfio%grids(ii-1)%km ) then
                  new_grid = .false.
                end if
              end do
              if ( new_grid ) exit
           end do
           cfio%grids(i) = cfio%varObjs(iCnt)%grid
        end do
      end if

      rtcode = 0
      if ( present(rc) ) rc = rtcode

      end subroutine ESMF_CFIOSdfFileOpen