CheckVarInfo Subroutine

public subroutine CheckVarInfo(filename, varname, varexist, filetype, meshVar, attstr, ndims, dims, dimids, haveMask, missingval, locstr, vartype, coordnames, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
character(len=*), intent(in) :: varname
logical, intent(out) :: varexist
type(ESMF_FileFormat_Flag), intent(in) :: filetype
character(len=*), intent(in) :: meshVar
character(len=*), intent(out) :: attstr
integer, intent(out) :: ndims
integer, intent(out) :: dims(:)
integer, intent(out) :: dimids(:)
logical, intent(out) :: haveMask
real(kind=ESMF_KIND_R8) :: missingval
character(len=*), intent(inout), optional :: locstr
integer, intent(out), optional :: vartype
character(len=*), intent(in), optional :: coordnames(:)
integer, intent(out), optional :: rc

Source Code

  subroutine CheckVarInfo(filename, varname, varexist, filetype, meshVar, attstr, &
   ndims, dims, dimids, haveMask, missingval, locstr, vartype, coordnames, rc)

! ARGUMENTS:
    character(len=*),  intent(in) :: filename
    character(len=*),  intent(in) :: varname
    logical, intent(out)          :: varexist
    type(ESMF_Fileformat_Flag), intent(in) :: filetype
    character(len=*), intent(in)  :: meshVar
    character(len=*),  intent(out):: attstr  ! the coordinate variable names
    integer, intent(out) :: ndims
    integer, intent(out) :: dims(:)
    integer, intent(out) :: dimids(:)
    logical, intent(out) :: haveMask
    real(ESMF_KIND_R8) :: missingval
    character(len=*),  intent(inout), optional :: locstr
    integer, intent(out), optional :: vartype
    character(len=*), intent(in), optional :: coordnames(:)
    integer, intent(out), optional :: rc

    logical :: foundlon, foundlat
    integer :: ncStatus
    integer ::  gridid, varid, tempids(1), varids(2), meshid, len
    character(len=128) :: attvalue, locallocstr, varnames(2)
    integer :: i, nvars, pos
    character(len=128) :: errmsg
    integer, parameter :: nf90_noerror = 0

    varexist = .true.
#ifdef ESMF_NETCDF
    rc = ESMF_FAILURE
 
    ! check if varname exist for GRIDSPEC and UGRID
    ! varname could be a list of variables separated by comma, need to check all of 
    ! them
    ncStatus = nf90_open (path=filename, mode=nf90_nowrite, ncid=gridid)
      errmsg = 'Fail to open '//trim(filename)
      if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg, &
            rc)) return
      ncStatus = nf90_inq_varid(gridid, varname, varid)
      if (ncStatus /= nf90_noerror) then
         varexist = .false.
         !print *, varname, ' does not exist in ',trim(filename)
      endif
      if (varexist) then
         if (filetype == ESMF_FILEFORMAT_UGRID) then
           ! get location attribute
           ncStatus = nf90_inquire_attribute(gridid, varid, "location", len=len)
           if (ncStatus == nf90_noerror) then
              ncStatus = nf90_get_att(gridid, varid, "location", locallocstr)
              errmsg = 'Fail to get attribute location '//trim(varname)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD, &
                 ESMF_SRCLINE,&
                 errmsg, &
                 rc)) return
              if (present(locstr)) then
                 locstr = locallocstr(1:4)
              endif
           else
              call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
                msg = " location attribute not defined for the variable", &
                ESMF_CONTEXT, rcToReturn=rc)
              return
           endif
           attstr = MeshVar
         else
           !  GRIDSPEC or MOSAIC fileformat
           ncStatus = nf90_inquire_attribute(gridid, varid, "coordinates", len=len)
           if (ncStatus == nf90_noerror) then
              ncStatus = nf90_get_att(gridid, varid, "coordinates", attstr)
              errmsg = 'fail to get attribute coordinates for '//trim(varname)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD, &
                 ESMF_SRCLINE,&
                 errmsg, &
                 rc)) return
              if (present(locStr)) then
                locStr = "face"
              endif
           else
             call ESMF_LogSetError(ESMF_FAILURE, &
                 msg=trim(filename)//' is not a GRIDSPEC, MOSAIC or a UGRID file', &
                 ESMF_CONTEXT, rcToReturn=rc)
             return
           endif
         endif

         ! Get missing value attribute
         ncStatus = nf90_get_att(gridid, varid, "_FillValue", missingval)
         if (ncStatus /= nf90_noerror) then
            ncStatus = nf90_get_att(gridid, varid, "missing_value", missingval)
            if (ncStatus == nf90_noerror) then
                   haveMask = .TRUE.
            else
                   haveMask = .FALSE.
            endif
         else
            haveMask = .TRUE.
         endif

         if (.not. haveMask) missingval = 0.0

         ! get the dimension info for the variable
         if (present(vartype)) then
           ncStatus = nf90_inquire_variable(gridid, varid, xtype=vartype, ndims=ndims, &
               dimids=dimids)
         else
           ncStatus = nf90_inquire_variable(gridid, varid, ndims=ndims, dimids=dimids)
         endif
         errmsg = 'nf90_inquire_variable failed '//trim(varname)
         if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
         do i=1, ndims
           ncStatus = nf90_inquire_dimension(gridid, dimids(i), len=dims(i))
           errmsg = 'nf90_inquire_dimension failed '//trim(filename)
           if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
         enddo
      else ! variable does not exist      
         if (fileType == ESMF_FILEFORMAT_UGRID) then
           if (present(locStr)) then
              ncStatus = nf90_inq_varid(gridid, meshVar,meshid)
              errmsg = 'mesh_topology dummy varible does not exist '//trim(filename)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD, &
                 ESMF_SRCLINE,&
                 errmsg, &
                 rc)) return
              if (locStr .eq. 'face') then
                 ncStatus=nf90_get_att(gridid, meshid, 'face_coordinates', locallocstr) 
                 errmsg = 'face_coordinates attribute does not exist '//trim(filename)
                 if (CDFCheckError (ncStatus, &
                    ESMF_METHOD, &
                    ESMF_SRCLINE,&
                    errmsg, &
                    rc)) return
              else ! node, default for non-conservative
                 ncStatus=nf90_get_att(gridid, meshid, 'node_coordinates', locallocstr)
                 errmsg = 'node_coordinates attribute does not exist '//trim(filename)
                 if (CDFCheckError (ncStatus, &
                    ESMF_METHOD, &
                    ESMF_SRCLINE,&
                    errmsg, &
                    rc)) return
              endif
           else
              ! signal an error if the variable does not exist and locStr is not present
              call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
                 msg = " locStr is not present when the variable is not defined", &
                 ESMF_CONTEXT, rcToReturn=rc)
              return
           endif 
        
           ! split the locallocstr to find the coordinate var name
           pos = INDEX(locallocstr, ' ')
           varnames(1)=locallocstr(1:pos-1)
           varnames(2)=locallocstr(pos+1:)
           ncStatus=nf90_inq_varid(gridid, varnames(1), varids(1))
           errmsg = trim(varnames(1))//' does not exist '
           if (CDFCheckError (ncStatus, &
                  ESMF_METHOD, &
                  ESMF_SRCLINE,&
                  errmsg, &
                  rc)) return
           ncStatus=nf90_inq_varid(gridid, varnames(2), varids(2))
           errmsg = trim(varnames(2))//' does not exist '
           if (CDFCheckError (ncStatus, &
                  ESMF_METHOD, &
                  ESMF_SRCLINE,&
                  errmsg, &
                  rc)) return
           attstr = MeshVar
         elseif (fileType == ESMF_FILEFORMAT_GRIDSPEC .or. &
            fileType == ESMF_FILEFORMAT_MOSAIC) then
           !GRIDSPEC, find the coordinate variables using units
           !Check if the optional coordinate argument exist
           if (present(coordnames)) then
             ! check if the coordinate variables exist or not
             ncStatus=nf90_inq_varid(gridid, coordnames(1), varids(1))
             errmsg = trim(coordnames(1))//' does not exist '
             if (CDFCheckError (ncStatus, &
                  ESMF_METHOD, &
                  ESMF_SRCLINE,&
                  errmsg, &
                  rc)) return
             ncStatus=nf90_inq_varid(gridid, coordnames(2), varids(2))
             errmsg = trim(coordnames(2))//' does not exist '
             if (CDFCheckError (ncStatus, &
                  ESMF_METHOD, &
                  ESMF_SRCLINE,&
                  errmsg, &
                  rc)) return
             ! check the unit attribute and make sure it is a latitude var
             ncStatus=nf90_get_att(gridid, varids(1), 'units', attvalue)
             if (ncStatus == nf90_noerror) then
               if (attvalue(len:len) .eq. achar(0)) len = len-1
                if (.not. (attvalue(1:len) .eq. "degrees_east" .or. &
                    attvalue(1:len) .eq. "degree_east" .or. &
                    attvalue(1:len) .eq. "degree_E" .or. &
                    attvalue(1:len) .eq. "degrees_E" .or. &
                    attvalue(1:len) .eq. "degreeE" .or. &
                    attvalue(1:len) .eq. "degreesE"))  then
                  call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                    msg = " not a valid latitude variable - units attribute incorrect", &
                    ESMF_CONTEXT, rcToReturn=rc)
                  return
               endif
              else
                  call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                    msg = " not a valid latitude variable - units attribute does not exist", &
                    ESMF_CONTEXT, rcToReturn=rc)
                  return
              endif
              foundlon = .true.
              foundlat = .true.
              ! check the unit attribute and make sure it is a longitude var
              ncStatus=nf90_get_att(gridid, varids(2), 'units', attvalue)
              if (ncStatus == nf90_noerror) then
                if (attvalue(len:len) .eq. achar(0)) len = len-1
               if (.not. (attvalue(1:len) .eq. "degrees_north" .or. &
                   attvalue(1:len) .eq. "degree_north" .or. &
                   attvalue(1:len) .eq. "degree_N" .or. &
                   attvalue(1:len) .eq. "degrees_N" .or. &
                   attvalue(1:len) .eq. "degreeN" .or. &
                   attvalue(1:len) .eq. "degreesN"))  then
                  call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                    msg = " not a valid longitude variable - units attribute incorrect", &
                    ESMF_CONTEXT, rcToReturn=rc)
                  return
                endif
              else
                  call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                    msg = " not a valid longitude variable - units attribute does not exist", &
                    ESMF_CONTEXT, rcToReturn=rc)
                  return
              endif
              attstr = trim(coordnames(1))//' '//trim(coordnames(2))
           else  ! did not specify the optional coordnames, find them
              ncStatus = nf90_inquire(gridid, nVariables=nvars)
              errmsg = 'nf90_inquire failed '//trim(filename)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD, &
                 ESMF_SRCLINE,&
                 errmsg, &
                 rc)) return
              foundlon = .false.
              foundlat = .false.
              do i=1,nvars
                ncStatus = nf90_inquire_attribute(gridid, i, "units", len=len)
                if (ncStatus /= nf90_noerror) cycle
                ncStatus=nf90_get_att(gridid, i, 'units', attvalue)
                if (ncStatus /= nf90_noerror) then
                  print '("NetCDF error: ", A)', trim(nf90_strerror(ncStatus))
                  return
                endif
                if (attvalue(len:len) .eq. achar(0)) len = len-1
                if (len >= 6 .and. (attvalue(1:6) .eq. "degree")) then
                  if (attvalue(1:len) .eq. "degrees_east" .or. &
                    attvalue(1:len) .eq. "degree_east" .or. &
                    attvalue(1:len) .eq. "degree_E" .or. &
                    attvalue(1:len) .eq. "degrees_E" .or. &
                    attvalue(1:len) .eq. "degreeE" .or. &
                    attvalue(1:len) .eq. "degreesE")  then
                    if (foundlon) then
                      call ESMF_LogSetError(ESMF_FAILURE, &
                        msg="- Duplicate longitude variables defined", &
                        ESMF_CONTEXT, rcToReturn=rc)
                      return
                    else
                      ncStatus = nf90_inquire_variable(gridid,i,name=varnames(1))
                      varids(1)=i
                      foundlon = .true.
                    endif
                  else if (attvalue(1:len) .eq. "degrees_north" .or. &
                    attvalue(1:len) .eq. "degree_north" .or. &
                    attvalue(1:len) .eq. "degree_N" .or. &
                    attvalue(1:len) .eq. "degrees_N" .or. &
                    attvalue(1:len) .eq. "degreeN" .or. &
                    attvalue(1:len) .eq. "degreesN")  then
                    if (foundlat) then
                      call ESMF_LogSetError(ESMF_FAILURE, &
                         msg="- Duplicate latitude variables defined", &
                         ESMF_CONTEXT, rcToReturn=rc)
                      return
                    else
                      ncStatus = nf90_inquire_variable(gridid,i,name=varnames(2))
                      varids(2)=i
                      foundlat = .true.
                    endif
                  endif
                endif
              enddo ! loop through all the variables
              if (foundlon .and. foundlat) then
                 attstr = trim(varnames(1))//' '//trim(varnames(2))
              else
                 call ESMF_LogSetError(ESMF_FAILURE, &
                    msg="- Did not find the coordinate variables", &
                    ESMF_CONTEXT, rcToReturn=rc)
                 return
              endif
            endif !coordname is not present
            if (present(locStr)) then
               locStr = "face"
            endif
         endif ! not UGRID
         ! find the lon, lat dimension
         ncStatus=nf90_inquire_variable(gridid, varids(1), ndims=ndims, dimids=dimids)
         do i=1, ndims
              ncStatus = nf90_inquire_dimension(gridid, dimids(i), len=dims(i))
              errmsg = 'nf90_inquire_dimension '//trim(filename)
              if (CDFCheckError (ncStatus, &
                     ESMF_METHOD, &
                     ESMF_SRCLINE,&
                     errmsg, &
                     rc)) return
         enddo
         if (filetype == ESMF_FILEFORMAT_GRIDSPEC .and. ndims==1) then
             ! find the dimension of the latitude coordinate
             ncStatus=nf90_inquire_variable(gridid, varids(2), ndims=ndims, dimids=tempids)
             dimids(2)=tempids(1)
             ncStatus = nf90_inquire_dimension(gridid, dimids(2), len=dims(2))
              errmsg = 'nf90_inquire_dimension '//trim(filename)
              if (CDFCheckError (ncStatus, &
                     ESMF_METHOD, &
                     ESMF_SRCLINE,&
                     errmsg, &
                     rc)) return
             ndims=2
         endif
     endif ! not varExist
     ncStatus = nf90_close(gridid)
     if (CDFCheckError (ncStatus, &
              ESMF_METHOD, &
              ESMF_SRCLINE,&
              errmsg, &
              rc)) return
     rc = ESMF_SUCCESS
     return
#else
    call ESMF_LogSetError(ESMF_RC_LIB_NOT_PRESENT, &
                 msg="- ESMF_NETCDF not defined when lib was compiled", &
                 ESMF_CONTEXT, rcToReturn=rc)
    return
#endif

  end subroutine checkVarInfo