UGridReadCoords Subroutine

public subroutine UGridReadCoords(filename, meshvar, meshloc, lonarray, latarray, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
character(len=*), intent(in) :: meshvar
type(ESMF_MeshLoc) :: meshloc
real(kind=ESMF_KIND_R8), TARGET :: lonarray(:)
real(kind=ESMF_KIND_R8), TARGET :: latarray(:)
integer :: rc

Source Code

  subroutine UGridReadCoords(filename, meshvar, meshloc, lonarray, latarray, rc)

  character(len=*), intent(in) :: filename
  character(len=*), intent(in) :: meshvar
  type(ESMF_MeshLoc)           :: meshloc
  real(ESMF_KIND_R8), TARGET  :: lonarray(:)
  real(ESMF_KIND_R8), TARGET  :: latarray(:)
  integer                      :: rc

    character(len=128) :: errmsg
    integer, parameter :: nf90_noerror = 0
    integer :: ncStatus
    integer ::  gridid, varid, tempids(1), varids(2), meshid, len
    character(len=128) :: attvalue, coordsname, varnames(2)
    integer :: ndims, dimids(2)
    integer :: i, nvars, pos
    real(ESMF_KIND_R8), parameter :: d2r = 3.141592653589793238/180

#ifdef ESMF_NETCDF
    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, meshvar, meshid)
    errmsg = 'Dummy variable '//trim(meshvar)//' does not exist'
    if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg, &
            rc)) return

    if (meshloc .eq. ESMF_MESHLOC_ELEMENT) then
        ncStatus=nf90_get_att(gridid, meshid, 'face_coordinates', coordsname) 
        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', coordsname)
        errmsg = 'node_coordinates attribute does not exist '//trim(filename)
        if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
    endif
    pos = INDEX(coordsname, ' ')
    varnames(1)=coordsname(1:pos-1)
    varnames(2)=coordsname(pos+1:)
      
    do i=1,2
       ncStatus = nf90_inq_varid(gridid, varnames(i), varid)
       errmsg = 'Coordinate variable '//trim(varnames(i))//' does not exist'
       if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg, &
            rc)) return
       ! Check it is a longitude or a latitude
       ncStatus = nf90_inquire_attribute(gridid, varid, 'units', len=len)
       errmsg = 'Attribute units of '//trim(varnames(i))//' does not exist'
       if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg, &
            rc)) return
       ncStatus=nf90_get_att(gridid, varid, 'units', attvalue)
       errmsg = 'Attribute units of '//trim(varnames(i))//' does not exist'
       if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg, &
            rc)) return
        if (attvalue(len:len) .eq. achar(0)) len = len-1
        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
             ncStatus = nf90_get_var(gridid, varid,lonarray)
             errmsg = 'Read variable failed: '//trim(varnames(i))
             if (CDFCheckError (ncStatus, &
                 ESMF_METHOD, &
                 ESMF_SRCLINE,&
                 errmsg, &
                 rc)) return
        else
            ncStatus = nf90_get_var(gridid, varid, latarray)
            errmsg = 'Read variable failed: '//trim(varnames(i))
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
        endif
    enddo  
    latarray = latarray*d2r
    lonarray = lonarray*d2r
    rc = ESMF_SUCCESS
    return
#else
    call ESMF_LogSetError(rcToCheck=ESMF_RC_LIB_NOT_PRESENT, &
      msg="- ESMF_NETCDF not defined when lib was compiled", &
      ESMF_CONTEXT, rcToReturn=rc)
    return
#endif
  end subroutine UGridReadCoords