subroutine f_esmf_xgridgetsparsematx2a(xgrid, &
sideAIndex_base0, &
factorListCount, &
factorList, &
factorIndexList, &
rc)
use ESMF_XGridMod
use ESMF_XGridGetMod
use ESMF_XGridCreateMod
use ESMF_XGridGeomBaseMod
use ESMF_UtilTypesMod
use ESMF_BaseMod
use ESMF_LogErrMod
use iso_c_binding
implicit none
type(ESMF_XGrid) :: xgrid
integer, intent(in) :: sideAIndex_base0
integer, intent(out) :: factorListCount
type(C_PTR) :: factorList
type(C_PTR) :: factorIndexList
integer, intent(out) :: rc
! Local Variables
real(ESMF_KIND_R8), pointer :: factorListFPtr(:)
integer(ESMF_KIND_I4), pointer :: factorIndexListFPtr(:,:)
integer :: sideAGridCount, sideAMeshCount
type(ESMF_XGridSpec), allocatable :: sparseMatX2A(:)
integer :: sideAIndex
! Init rc
rc = ESMF_RC_NOT_IMPL
! Make sideA index base 1
sideAIndex=sideAIndex_base0+1
! Get number of side A Grids
call ESMF_XGridGet(xgrid, sideAGridCount=sideAGridCount, rc=rc)
if (ESMF_LogFoundError(rc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Get number of side A Meshes
call ESMF_XGridGet(xgrid, sideAMeshCount=sideAMeshCount, rc=rc)
if (ESMF_LogFoundError(rc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Make sure the index is not too big
if (sideAIndex > sideAGridCount+sideAMeshCount) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
msg="sideAIndex bigger than the number of Grids and Meshes on sideA", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
! Make sure the index is not too small
if (sideAIndex < 1) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
msg="sideAIndex below 1 in Fortran or 0 in C.", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
! Allocate XGridSpec array
allocate(sparseMatX2A(sideAGridCount+sideAMeshCount))
! Get info
call ESMF_XGridGet(xgrid, sparseMatX2A=sparseMatX2A, rc=rc)
if (ESMF_LogFoundError(rc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Get sparse matrix information
factorListCount = size(sparseMatX2A(sideAIndex)%factorList, 1)
! Associate the Fortran pointers with C pointers. Only do this if
! factors were created during the regrid store call.
if (factorListCount > 0) then
factorList = C_LOC(sparseMatX2A(sideAIndex)%factorList(1))
factorIndexList = C_LOC(sparseMatX2A(sideAIndex)%factorIndexList(1,1))
endif
! set rc to success
rc = ESMF_SUCCESS
end subroutine f_esmf_xgridgetsparsematx2a