ESMF_NetCDFInquireDimension Function

private function ESMF_NetCDFInquireDimension(dimensionName, path, ncid, rc) result(n)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: dimensionName
character(len=*), intent(in), optional :: path
integer, intent(in), optional :: ncid
integer, intent(out), optional :: rc

Return Value integer


Source Code

  integer function ESMF_NetCDFInquireDimension(dimensionName, path, ncid, rc) result(n)

! ! ARGUMENTS:
    character(len=*), intent(in) :: dimensionName
    character(len=*), intent(in), optional :: path
    integer, intent(in), optional :: ncid
    integer, intent(out), optional :: rc

!-------------------------------------------------------------------------------
! !DESCRIPTION:
!
! Return the integer length of a dimension in a netCDF file using a file path or file
! identifier.
!
! The arguments are:
!
! \begin{description}
!
! \item [dimensionName]
!       The name of the target dimension.
!
! \item [{[path]}]
!       Path to the target file. Required if {\tt ncid} is not provided.
!
! \item [{[ncid]}]
!       NetCDF file identifier of the target file. Required if {\tt path} is not provided.
!
! \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!
! \end{description}
!
!EOPI
!-------------------------------------------------------------------------------

    ! LOCAL VARIABLES:
    integer :: dimid, ncidV, ncStatus

    ! --------------------------------------------------------------------------

#ifdef ESMF_NETCDF
    if (present(rc)) rc = ESMF_RC_NOT_IMPL

    if (present(ncid)) ncidV = ncid

    ! Open the file and update the file identifier variable.
    if (.not. present(ncid)) then
      ncStatus = nf90_open(path, NF90_NOWRITE, ncidV)
      if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, path, __LINE__, rc)) return
    endif

    ! Get the dimension identifier.
    ncStatus = nf90_inq_dimid(ncidV, dimensionName, dimid)
    if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, path, __LINE__, rc)) return
    ! Get the length/size of the dimension.
    ncStatus = nf90_inquire_dimension(ncidV, dimid, len=n)
    if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, path, __LINE__, rc)) return

    ! Close the file.
    if (.not. present(ncid)) then
      ncStatus = nf90_close(ncidV)
      if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, path, __LINE__, rc)) return
    endif

    if (present(rc)) rc = ESMF_SUCCESS
#else
    n=-1  ! indicate invalid
    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 function ESMF_NetCDFInquireDimension