subroutine ESMF_UGridInqVarLoc (ncid, VarId, varname,location, rc)
integer, intent(in) :: ncid
integer, intent(in) :: VarId
character(len=*), intent(in) :: varname
integer, intent(out) :: location ! 1 for node, 2 for face
integer :: rc
integer :: meshId
integer :: ncStatus
character(len=256) :: errmsg
character(len=80) :: attstr, attstr1, locationStr
integer :: len,len1
integer, parameter :: nf90_noerror=0
#ifdef ESMF_NETCDF
ncStatus = nf90_get_att(ncid, VarId, "location", locationStr)
if (ncStatus /= nf90_noerror) then
! location attribute does not exist, check coordinates attribute
! check if the coordinate attribute exist or not
ncStatus = nf90_inquire_attribute(ncid, VarId, "coordinates", len=len)
errmsg ="No location or coordinates attributes defined for "//varname
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
ncStatus = nf90_get_att(ncid, VarId, "coordinates", attstr)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
! check if it matches with the node_coordinates or face_coordinates defined
! on the topology variable
ncStatus = nf90_inquire_attribute(ncid, VarId, "mesh", len=len1)
errmsg ="No mesh attribute defined for "//varname
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
ncStatus = nf90_get_att(ncid, VarId, "mesh", attstr1)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
ncStatus = nf90_inq_varid (ncid, attstr1(1:len1), meshId)
errmsg = "Dummy Variable "//attstr1(1:len1)//" does not exist"
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
errmsg = "Attribute node_coordinates in variable "//attstr1(1:len1)
ncStatus = nf90_inquire_attribute(ncid, meshId, "node_coordinates", len=len1)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
ncStatus = nf90_get_att (ncid, meshId, "node_coordinates", attstr1)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
if (attstr1(1:len1) .eq. attstr(1:len)) then
location = 1
else
! check if the coordinates match with face_coordinates
ncStatus = nf90_inquire_attribute(ncid, meshId, "face_coordinates", len=len1)
errmsg = "Attribute face_coordinates"
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
ncStatus = nf90_get_att (ncid, meshId, "face_coordinates", attstr1)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
if (attstr1(1:len1) .eq. attstr(1:len)) then
location = 2
else
errmsg = "- coordinates attribute does not match with face or node coordinates"
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg=errmsg, ESMF_CONTEXT, rcToReturn=rc)
return
endif
endif
else
ncStatus = nf90_inquire_attribute(ncid, VarId, "location", len=len)
if (locationStr(1:len) .eq. 'node') then
location = 1
elseif (locationStr(1:len) .eq. 'face') then
location = 2
else
errmsg = "- location attribute is not recognizable: "//locationStr(1:len)
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg=errmsg, ESMF_CONTEXT, rcToReturn=rc)
return
endif
endif
rc = ESMF_SUCCESS
#else
call ESMF_LogSetError(rcToCheck=ESMF_RC_LIB_NOT_PRESENT, &
msg="- ESMF_NETCDF not defined when lib was compiled", &
ESMF_CONTEXT, rcToReturn=rc)
#endif
return
end subroutine ESMF_UGridInqVarLoc