ESMF_ScripGetVar Subroutine

public subroutine ESMF_ScripGetVar(filename, grid_center_lon, grid_center_lat, grid_imask, grid_corner_lon, grid_corner_lat, grid_area, convertToDeg, start, count, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
real(kind=ESMF_KIND_R8), intent(inout), optional :: grid_center_lon(:)
real(kind=ESMF_KIND_R8), intent(inout), optional :: grid_center_lat(:)
integer, intent(inout), optional :: grid_imask(:)
real(kind=ESMF_KIND_R8), intent(inout), optional :: grid_corner_lon(:,:)
real(kind=ESMF_KIND_R8), intent(inout), optional :: grid_corner_lat(:,:)
real(kind=ESMF_KIND_R8), intent(inout), optional :: grid_area(:)
logical, intent(in), optional :: convertToDeg
integer, intent(in), optional :: start
integer, intent(in), optional :: count
integer, intent(out), optional :: rc

Source Code

  subroutine ESMF_ScripGetVar(filename, grid_center_lon, grid_center_lat, &
                        grid_imask, grid_corner_lon, grid_corner_lat, &
                        grid_area, convertToDeg, start, count, rc)
!
! !ARGUMENTS:
    character(len=*), intent(in)  :: filename
    real(ESMF_KIND_R8), intent(inout),optional:: grid_center_lon(:)
    real(ESMF_KIND_R8), intent(inout),optional:: grid_center_lat(:)
    integer, intent(inout), optional:: grid_imask(:)
    real(ESMF_KIND_R8), intent(inout),optional:: grid_corner_lon(:,:)
    real(ESMF_KIND_R8), intent(inout),optional:: grid_corner_lat(:,:)
    real(ESMF_KIND_R8), intent(inout),optional:: grid_area(:)
    logical, intent(in), optional:: convertToDeg
    integer, intent(in), optional:: start, count
    integer, intent(out), optional:: rc

    integer:: ncStatus
    integer:: ncid, VarId
    integer:: len
    character(len=128) :: units
    character(len=256) :: errmsg
    integer:: start1(1), count1(1), start2(2), count2(2)
    integer:: totalcells, grid_corners
    logical :: convertToDegLocal
    integer:: localrc
    integer, parameter :: nf90_noerror = 0

#ifdef ESMF_NETCDF
    convertToDegLocal = .false.
    if (present(convertToDeg)) convertToDegLocal = convertToDeg

    call ESMF_ScripInq(trim(filename), grid_size=totalcells, &
           grid_corners=grid_corners, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
           ESMF_CONTEXT, rcToReturn=rc)) return

    ncStatus = nf90_open (path=trim(filename), mode=nf90_nowrite, ncid=ncid)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD, &
      ESMF_SRCLINE,&
      trim(filename), &
      rc)) return

    if (present(start)) then
        start1(1)=start
        start2(2)=start
        start2(1)=1
    else
        start1(1)=1
        start2(:)=1
    endif
    if (present(count)) then
        count1(1)=count
        count2(2)=count
        count2(1)=grid_corners
    else
        count1(1)=totalcells
        count2(2)=totalcells
        count2(1)=grid_corners
    endif

    ! Read in grid_center_lon and grid_center_lat
    if (present(grid_center_lon)) then
      ncStatus = nf90_inq_varid (ncid, "grid_center_lon", VarId)
      errmsg = "variable grid_center_lon in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ncStatus = nf90_get_var (ncid, VarId, grid_center_lon, start1, count1)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return

      ! get the attribute 'units'
      ncStatus = nf90_inquire_attribute(ncid, VarId, "units", len=len)
      errmsg = "Attribute units for grid_center_lon in "//trim(filename)
      if (CDFCheckError (ncStatus, &
          ESMF_METHOD, &
          ESMF_SRCLINE,&
          errmsg,&
          rc)) return
      ncStatus = nf90_get_att(ncid, VarId, "units", units)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! if units is not "degrees" or "radians" return errors
      if (units(len:len) .eq. achar(0)) len = len-1
      units = ESMF_UtilStringLowerCase(units(1:len))
      if (units(1:len) .ne. 'degrees' .and. units(1:len) .ne. 'radians') then
          call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- units attribute is not degrees or radians", &
                 ESMF_CONTEXT, rcToReturn=rc)
          return
      endif
      ! if units is "radians", convert it to degrees
      if (convertToDegLocal) then
         if (units(1:len) .eq. "radians") then
            grid_center_lon(:) = &
                   grid_center_lon(:)*ESMF_COORDSYS_RAD2DEG
         endif
      endif     
    endif
    if (present(grid_center_lat)) then
      ncStatus = nf90_inq_varid (ncid, "grid_center_lat", VarId)
      errmsg = "variable grid_center_lat in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ncStatus = nf90_get_var (ncid, VarId, grid_center_lat, start1, count1)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! get the attribute 'units'
      ncStatus = nf90_inquire_attribute(ncid, VarId, "units", len=len)
      errmsg = "Attribute units for grid_center_lat in "//trim(filename)
      if (CDFCheckError (ncStatus, &
          ESMF_METHOD, &
          ESMF_SRCLINE,&
          errmsg,&
          rc)) return
      ncStatus = nf90_get_att(ncid, VarId, "units", units)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! if units is not "degrees" or "radians" return errors
     if (units(len:len) .eq. achar(0)) len = len-1
      units = ESMF_UtilStringLowerCase(units(1:len))
      if (units(1:len) .ne. 'degrees' .and. units(1:len) .ne. 'radians') then
          call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- units attribute is not degrees or radians", &
                 ESMF_CONTEXT, rcToReturn=rc)
          return
      endif
      ! if units is "radians", convert it to degree
      if (convertToDegLocal) then
         if (units(1:len) .eq. "radians") then
            grid_center_lat(:) = &
                 grid_center_lat(:)*ESMF_COORDSYS_RAD2DEG
         endif
      endif
    endif

    if (present(grid_imask)) then
      ncStatus = nf90_inq_varid (ncid, "grid_imask", VarId)
      errmsg = "variable grid_imask is not defined in "//trim(filename)
      if (ncStatus /= nf90_noerror) then
         print *, 'Warning:', errmsg
         grid_imask = 1
      else
         ncStatus = nf90_get_var (ncid, VarId, grid_imask, start1, count1)
         if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg,&
            rc)) return
      endif
    endif

    ! Read in grid_corner_lon and grid_corner_lat
    if (present(grid_corner_lon)) then
      ncStatus = nf90_inq_varid (ncid, "grid_corner_lon", VarId)
      errmsg = "variable grid_corner_lon in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ncStatus = nf90_get_var (ncid, VarId, grid_corner_lon, start2, count2)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! get the attribute 'units'
      ncStatus = nf90_inquire_attribute(ncid, VarId, "units", len=len)
      errmsg = "attribute units for grid_center_lon in "//trim(filename)
      if (CDFCheckError (ncStatus, &
          ESMF_METHOD, &
          ESMF_SRCLINE,&
          errmsg,&
          rc)) return
      ncStatus = nf90_get_att(ncid, VarId, "units", units)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! if units is not "degrees" or "radians" return errors
     if (units(len:len) .eq. achar(0)) len = len-1
      units = ESMF_UtilStringLowerCase(units(1:len))
      if (units(1:len) .ne. 'degrees' .and. units(1:len) .ne. 'radians') then
          call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- units attribute is not degrees or radians", &
                 ESMF_CONTEXT, rcToReturn=rc)
          return
      endif
      ! if units is "radians", convert it to degree
      if (convertToDegLocal) then
         if (units(1:len) .eq. "radians") then
            grid_corner_lon(:,:) = &
               grid_corner_lon(:,:)*ESMF_COORDSYS_RAD2DEG
         endif
      endif     
    endif
    if (present(grid_corner_lat)) then
      ncStatus = nf90_inq_varid (ncid, "grid_corner_lat", VarId)
      errmsg = "variable grid_corner_lat in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ncStatus = nf90_get_var (ncid, VarId, grid_corner_lat, start2, count2)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! get the attribute 'units'
      ncStatus = nf90_inquire_attribute(ncid, VarId, "units", len=len)
      errmsg = "Attribute units for grid_center_lon in "//trim(filename)
      if (CDFCheckError (ncStatus, &
          ESMF_METHOD, &
          ESMF_SRCLINE,&
          errmsg,&
          rc)) return
      ncStatus = nf90_get_att(ncid, VarId, "units", units)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ! if units is not "degrees" or "radians" return errors
     if (units(len:len) .eq. achar(0)) len = len-1
      units = ESMF_UtilStringLowerCase(units(1:len))
      if (units(1:len) .ne. 'degrees' .and. units(1:len) .ne. 'radians') then
          call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- units attribute is not degrees or radians", &
                 ESMF_CONTEXT, rcToReturn=rc)
          return
      endif
      ! if units is "radians", convert it to degree
      if (convertToDegLocal) then
         if (units(1:len) .eq. "radians") then
            grid_corner_lat(:,:) = &
                grid_corner_lat(:,:)*ESMF_COORDSYS_RAD2DEG
         endif
      endif
    endif

    ! Read in grid_area and grid_corner_lat
    if (present(grid_area)) then
      ncStatus = nf90_inq_varid (ncid, "grid_area", VarId)
      errmsg = "variable grid_area in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
      ncStatus = nf90_get_var (ncid, VarId, grid_area, start1, count1)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg,&
        rc)) return
    endif

    ncStatus = nf90_close(ncid)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD, &
      ESMF_SRCLINE,&
      trim(filename),&
      rc)) return

    if(present(rc)) rc = ESMF_SUCCESS
    return
#else
    if (ESMF_LogFoundError(ESMF_RC_LIB_NOT_PRESENT, &
                msg="- ESMF_NETCDF not defined when lib was compiled", &
                ESMF_CONTEXT, rcToReturn=rc)) return
#endif

end subroutine ESMF_ScripGetVar