subroutine ESMF_GridspecReadStaggerR8(filename, nx, ny, lon, lat, staggerLoc, start, count, rc)
! !ARGUMENTS:
character(len=*), intent(in) :: filename
integer, intent(in) :: nx, ny
real(ESMF_KIND_R8), TARGET :: lon(:,:)
real(ESMF_KIND_R8), TARGET :: lat(:,:)
type(ESMF_StaggerLoc) :: staggerLoc
integer, optional, intent(in) :: start(2)
integer, optional, intent(in) :: count(2)
integer, optional, intent(out) :: rc
!EOPI
integer :: ncid, nvars, attlen, i
integer :: nx1, ny1
integer :: ncStatus
integer :: ndims, dimids(2)
character(len=128) :: attstr
integer :: start1(2), count1(2)
real(ESMF_KIND_R8), allocatable :: supercoord(:,:)
integer :: localrc
logical :: found_grid_tile_spec
logical :: found_geo_lon
logical :: found_geo_lat
if (present(rc)) rc=ESMF_SUCCESS
call ESMF_VMGetCurrent(vm, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! set up local pet info
call ESMF_VMGet(vm, localPet=PetNo, petCount=PetCnt, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
#ifdef ESMF_NETCDF
! Init variables to ensure that we find the correct things in file
found_grid_tile_spec = .false.
found_geo_lon=.false.
found_geo_lat=.false.
! Open file
ncStatus = nf90_open(path=filename, mode=nf90_nowrite, ncid=ncid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
filename, &
rc)) return
ncStatus = nf90_inquire(ncid, nVariables=nvars)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
filename, &
rc)) return
do i=1,nvars
! Check its standard_name attribute
ncStatus = nf90_inquire_attribute(ncid, i, 'standard_name', len=attlen)
if (ncStatus == nf90_noerr) then
ncStatus = nf90_get_att(ncid, i, 'standard_name', values=attstr)
if (attstr(1:attlen) .eq. 'grid_tile_spec') then
! skip checking the attributes -- not sure which one should be set to what
! but makesure this dummy variable exists
found_grid_tile_spec = .true.
#if 0
! check the projection attribute
ncStatus = nf90_inquire_attribute(ncid, i, 'projection', len=attlen)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"projection attribute does not exist", &
rc)) return
ncStatus = nf90_get_att(ncid, i, 'projection', values=attstr)
if (attstr(1:attlen) .ne. 'cube_gnomonic') then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- Only Cube Gnomonic projection is currently supported", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
#endif
elseif (attstr(1:attlen) .eq. 'geographic_longitude' .or. &
attstr(1:attlen) .eq. 'geographic_latitude') then
! read the longitude or latitude variable
! First find the dimension of this variable
ncStatus = nf90_inquire_variable(ncid, i, ndims=ndims, dimids=dimids)
if (ndims /= 2) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- The longitude variable should have dimension 2", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
! find out the dimenison size
ncStatus = nf90_inquire_dimension(ncid, dimids(1), len=nx1)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"contact dimension inquire", &
rc)) return
ncStatus = nf90_inquire_dimension(ncid, dimids(2), len=ny1)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"contact dimension inquire", &
rc)) return
if (nx1 /= (nx*2+1)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- The x dimension of the tile does not match with the supergrid dimension", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
if (ny1 /= (ny*2+1)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- The y dimension of the tile does not match with the supergrid dimension", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
if (present(start) .and. present(count)) then
! read a block instead of the entire array
count1 = count
if (staggerLoc == ESMF_STAGGERLOC_CENTER) then
start1=start*2
elseif (staggerLoc == ESMF_STAGGERLOC_CORNER) then
start1=start*2-1
elseif (staggerLoc == ESMF_STAGGERLOC_EDGE1) then
start1(1)=start(1)*2-1
start1(2)=start(2)*2
elseif (staggerLoc == ESMF_STAGGERLOC_EDGE2) then
start1(1)=start(1)*2
start1(2)=start(2)*2-1
endif
else
count1(1) = nx
count1(2) = ny
if (staggerLoc == ESMF_STAGGERLOC_CENTER) then
start1=2
elseif (staggerLoc == ESMF_STAGGERLOC_CORNER) then
start1=1
elseif (staggerLoc == ESMF_STAGGERLOC_EDGE1) then
start1(1)=1
start1(2)=2
elseif (staggerLoc == ESMF_STAGGERLOC_EDGE2) then
start1(1)=2
start1(2)=1
endif
endif
if (attstr(1:attlen) .eq. 'geographic_latitude') then
ncStatus = nf90_get_var(ncid, i, lat, start=start1, count=count1, stride=(/2,2/))
found_geo_lat=.true.
else
ncStatus = nf90_get_var(ncid, i, lon, start=start1, count=count1, stride=(/2,2/))
found_geo_lon=.true.
endif
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"error reading stagger coordinates", &
rc)) return
endif
endif
enddo
! Close file
ncStatus = nf90_close(ncid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"close tile file", &
rc)) return
! Error check that the correct variables were in file
if (.not. found_grid_tile_spec) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_NOT_FOUND, &
msg="No variable with a standard_name attribute set to grid_tile_spec found in file "//trim(filename), &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
if (.not. found_geo_lat) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_NOT_FOUND, &
msg="No variable with a standard_name attribute set to geographic_latitude found in file "//trim(filename), &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
if (.not. found_geo_lon) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_NOT_FOUND, &
msg="No variable with a standard_name attribute set to geographic_longitude found in file "//trim(filename), &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
! Leave
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 ESMF_GridspecReadStaggerR8