subroutine ESMF_FactorRead(filename, factorList, factorIndexList, rc)
! ! ARGUMENTS:
character(len=*), intent(in) :: filename
real(ESMF_KIND_R8), dimension(:), allocatable, intent(out) :: factorList
integer, dimension(:, :), allocatable, intent(out) :: factorIndexList
integer, intent(out), optional :: rc
!-------------------------------------------------------------------------------
! !DESCRIPTION:
!
! The arguments are:
!
! \begin{description}
!
! \item [filename]
! Path to the file containing weights for creating an {\tt ESMF\_RouteHandle}.
! See ~(\ref{sec:weightfileformat}) for a description of the SCRIP weight
! file format. Only "row", "col", and "S" variables are required. They
! must be one-dimensional with dimension "n\_s".
!
! \item [factorList]
! The weight factors / interpolation weights to be read from file.
!
! \item [factorIndexList]
! The indices into the source and destination arrays to be read from file. The
! first dimension are the source indices. The second dimension are the
! destination indices.
!
! \item [{[rc]}]
! Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!
! \end{description}
!
!EOPI
!-------------------------------------------------------------------------------
! LOCAL VARIABLES ----------------------------------------------------------
! Variable name for index into the destination array.
character(len=3), parameter :: vfactorindexlist_dst = "row"
! Variable name for index into the source array.
character(len=3), parameter :: vfactorindexlist_src = "col"
! Variable name for the factor value.
character(len=1), parameter :: vfactorlist = "S"
! Dimensions for the factor variables.
character(len=3), parameter :: dfactor = "n_s"
integer :: ncid, varid, dimid, localPet, petCount, nElements, esplit, lb, ub, remainder
integer, dimension(1) :: nElementsArray, startArray
integer :: ncStatus, theSize
type(ESMF_VM) :: vm
! --------------------------------------------------------------------------
#ifdef ESMF_NETCDF
! Initialize return code; assume routine not implemented
if (present(rc)) rc = ESMF_RC_NOT_IMPL
! Get all VM information
call ESMF_VMGetCurrent(vm, rc=rc)
if (ESMF_LogFoundError(rc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Set up local pet info
call ESMF_VMGet(vm, localPet=localPet, petCount=petCount, rc=rc)
if (ESMF_LogFoundError(rc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! --------------------------------------------------------------------------
! Open the netCDF file.
ncStatus = nf90_open(filename, NF90_NOWRITE, ncid)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) return
! Get the number of factors from the target netCDF file.
nElements = ESMF_NetCDFInquireDimension(dfactor, ncid=ncid, rc=rc)
if (ESMF_LogFoundError(rc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Compute the lower and upper bounds for the current PET -------------------
esplit = floor(real(nElements) / real(petCount))
remainder = nElements - esplit*petCount
if (localPet<remainder) then
lb = localPet * (esplit+1) + 1
theSize = esplit+1
else
lb = localPet * esplit + remainder + 1
theSize = esplit
endif
ub = lb + theSize - 1
! --------------------------------------------------------------------------
! Convert to arrays as expected by the netCDF Fortran interface. Set the size
! of the factor arrays to zero if the localPet will not read any values.
startArray = (/lb/)
nElementsArray = (/theSize/)
! Allocate our factor arrays now that we know the size of the dimension.
allocate(factorList(theSize))
allocate(factorIndexList(2, theSize))
!--- Fill variables from the netCDF file -----------------------------------
! Only read from file if the lower bound is reasonable - falls inside the size
! of the target variable.
if (lb <= nElements) then
ncStatus = nf90_inq_varid(ncid, vfactorindexlist_dst, varid)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) &
return
ncStatus = nf90_get_var(ncid, varid, factorIndexList(2, :), &
start=startArray, count=nElementsArray)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) &
return
ncStatus = nf90_inq_varid(ncid, vfactorindexlist_src, varid)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) &
return
ncStatus = nf90_get_var(ncid, varid, factorIndexList(1, :), &
start=startArray, count=nElementsArray)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) &
return
ncStatus = nf90_inq_varid(ncid, vfactorlist, varid)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) &
return
ncStatus = nf90_get_var(ncid, varid, factorList, start=startArray, &
count=nElementsArray)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) &
return
endif
!---------------------------------------------------------------------------
! Close the file.
ncStatus = nf90_close(ncid)
if (ESMF_NetCDFCheckError(ncStatus, ESMF_METHOD, filename, __LINE__, rc)) return
if (present(rc)) 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)
return
#endif
end subroutine ESMF_FactorRead