test7d3_generic_fptr Subroutine

subroutine test7d3_generic_fptr(rc, minIndex, maxIndex, gridEdgeLWidth, gridEdgeUWidth, regDecomp, distgridToGridMap, datacopyflag, staggerloc, gridToFieldMap, totalLWidth, totalUWidth, fieldget)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: rc
integer, dimension(:) :: minIndex
integer, dimension(:) :: maxIndex
integer, optional, dimension(:) :: gridEdgeLWidth
integer, optional, dimension(:) :: gridEdgeUWidth
integer, optional, dimension(:) :: regDecomp
integer, optional, dimension(:) :: distgridToGridMap
type(ESMF_DataCopy_Flag), optional :: datacopyflag
type(ESMF_StaggerLoc), optional :: staggerloc
integer, optional, dimension(:) :: gridToFieldMap
integer, optional, dimension(:) :: totalLWidth
integer, optional, dimension(:) :: totalUWidth
logical, optional :: fieldget

Source Code

    subroutine test7d3_generic_fptr(rc, minindex, maxindex, &
        gridEdgeLWidth, gridEdgeUWidth, &
        regDecomp, &
        distgridToGridMap, &
        datacopyflag, &
        staggerloc, &
        gridToFieldMap, &
        totalLWidth, totalUWidth, &
        fieldget)
        integer, intent(out) :: rc
        integer, dimension(:)   :: minIndex
        integer, dimension(:)   :: maxIndex
        integer, dimension(:), optional   :: gridEdgeLWidth, gridEdgeUWidth
        integer, dimension(:), optional   :: regDecomp
        integer, dimension(:), optional   :: distgridToGridMap
        type(ESMF_DataCopy_Flag), optional     :: datacopyflag
        type(ESMF_StaggerLoc), optional   :: staggerloc
        integer, dimension(:), optional   :: gridToFieldMap
        integer, dimension(:), optional   :: totalLWidth, totalUWidth
        logical, optional                 :: fieldget

        real(ESMF_KIND_R8), dimension(:,:,:,:,:,:,:), pointer :: farray
        type(ESMF_Field)    :: field
        type(ESMF_Grid)     :: grid
        type(ESMF_DistGrid) :: distgrid
        integer             :: localrc
        integer             :: flb(7), fub(7)

        type(ESMF_Grid)         :: grid1
        type(ESMF_Array)        :: array
        type(ESMF_TypeKind_Flag)     :: typekind
        integer                 :: dimCount
        type(ESMF_StaggerLoc)   :: lstaggerloc
        integer, dimension(7) :: lgridToFieldMap
        integer, dimension(7,1) :: ltotalLWidth
        integer, dimension(7,1) :: ltotalUWidth

        integer                                     :: i, ii, ij, ik, il, im, io, ip
        integer, dimension(7)                       :: felb, feub, fclb, fcub, ftlb, ftub
        integer, dimension(7)                       :: fec, fcc, ftc
        integer, dimension(7)                       :: gclb, gcub, gcc, gelb, geub, gec
        real(ESMF_KIND_R8), dimension(:,:,:,:,:,:,:), pointer :: farray1
        real(ESMF_KIND_R8)                          :: n
        logical                                     :: t
        type(ESMF_STAGGERLOC)                       :: localStaggerLoc
        integer                 :: gminIndex(7), gmaxIndex(7), geleCount(7)
        integer                 :: lminIndex(7), lmaxIndex(7), leleCount(7)

        localrc = ESMF_SUCCESS
        rc = ESMF_SUCCESS
        nullify(farray1)

        if(present(staggerloc)) then
            localStaggerLoc=staggerloc
        else
            localStaggerLoc=ESMF_STAGGERLOC_CENTER
        endif

        distgrid = ESMF_DistGridCreate(minIndex=minIndex, maxIndex=maxIndex, &
            regDecomp=regDecomp, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        grid = ESMF_GridCreate(distgrid=distgrid, name="grid", &
            distgridToGridMap=distgridToGridMap, &
            gridEdgeLWidth=gridEdgeLWidth, gridEdgeUWidth=gridEdgeUWidth, &
            rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_GridGet(grid, localDe=0, staggerloc=localStaggerLoc, &
            exclusiveLBound=gelb, exclusiveUBound=geub, exclusiveCount=gec, &
            computationalLBound=gclb, computationalUBound=gcub, computationalCount=gcc, &
            rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_GridGetFieldBounds(grid, localDe=0, &
            staggerloc=localStaggerLoc, &
            totalLWidth=totalLWidth, totalUWidth=totalUWidth, &
            gridToFieldMap=gridToFieldMap, &
            totalLBound=flb, totalUBound=fub, &
            rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        allocate(farray(flb(1):fub(1), flb(2):fub(2), flb(3):fub(3), &
            flb(4):fub(4), flb(5):fub(5), flb(6):fub(6), flb(7):fub(7)) )

        if(present(fieldget)) then
          if(fieldget) then
            ! reverse looping order to make this a little faster by improving data locality
            do ip = flb(7), fub(7)
             do io = flb(6), fub(6)
              do im = flb(5), fub(5)
               do il = flb(4), fub(4)
                do ik = flb(3), fub(3)
                 do ij = flb(2), fub(2)
                  do ii = flb(1), fub(1)
                    farray(ii,ij,ik,il,im,io,ip) = ii+ij*2+ik+il*2+im+io*2+ip
                  enddo
                 enddo
                enddo
               enddo
              enddo
             enddo
            enddo
          endif
        endif

        field = ESMF_FieldEmptyCreate(name="field", rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldEmptyComplete(field, grid, farray, &
            totalLWidth=totalLWidth, totalUWidth=totalUWidth, &
            gridToFieldMap=gridToFieldMap, &
            datacopyflag=datacopyflag, &
            staggerloc=localStaggerLoc, &
            rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        if(present(fieldget)) then
          if(fieldget) then
            call ESMF_FieldGet(field, grid=grid1, array=array, typekind=typekind, &
                dimCount=dimCount, staggerloc=lstaggerloc, gridToFieldMap=lgridToFieldMap, &
                totalLWidth=ltotalLWidth, totalUWidth=ltotalUWidth, &
                rc=localrc)
            if (ESMF_LogFoundError(localrc, &
                ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
            call ESMF_FieldGet(field, localDe=0, farrayPtr=farray1, &
                exclusiveLBound=felb, exclusiveUBound=feub, exclusiveCount=fec, &
                computationalLBound=fclb, computationalUBound=fcub, computationalCount=fcc, &
                totalLBound=ftlb, totalUBound=ftub, totalCount=ftc, &
                rc=localrc)
            if (ESMF_LogFoundError(localrc, &
                ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

            ! test pointer equivalence
            t = associated(farray, farray1)
            do i = 1, 7
                t = t .and. (lbound(farray1, i) .eq. flb(i)) .and. (ubound(farray1, i) .eq. fub(i))
            enddo

            if(present(datacopyflag)) then
                if(datacopyflag==ESMF_DATACOPY_VALUE) t = .true.
            endif

            if(.not. t) then
              call ESMF_LogSetError(ESMF_RC_PTR_BAD, &
                msg="- pointer queried from object is not equivalent to the one passed in)", &
                ESMF_CONTEXT, rcToReturn=rc)
              return
            endif

            ! test field and grid bounds
            call ESMF_GridGet(grid, localDe=0, staggerloc=ESMF_STAGGERLOC_CENTER, &
                exclusiveLBound=gelb, exclusiveUBound=geub, &
                computationalLBound=gclb, computationalUBound=gcub, &
                rc=localrc) 
            if (ESMF_LogFoundError(localrc, &
                ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

            t = .true.
            do i = 1, 2
                t = t .and. (gelb(i) .eq. felb(i))
                t = t .and. (geub(i) .eq. feub(i))
                t = t .and. (gclb(i) .eq. fclb(i))
                t = t .and. (gcub(i) .eq. fcub(i))
            enddo
    
            if(.not. t) then
              call ESMF_LogSetError(ESMF_RC_PTR_BAD, &
                msg="- bounds queried from grid different from those queried from field)", &
                ESMF_CONTEXT, rcToReturn=rc)
              return
            endif

            ! reverse looping order to make this a little faster by improving data locality
            do ip = ftlb(7), ftub(7)
             do io = ftlb(6), ftub(6)
              do im = ftlb(5), ftub(5)
               do il = ftlb(4), ftub(4)
                do ik = ftlb(3), ftub(3)
                 do ij = ftlb(2), ftub(2)
                  do ii = ftlb(1), ftub(1)
                    n = ii+ij*2+ik+il*2+im+io*2+ip
                    if(farray1(ii,ij,ik,il,im,io,ip) .ne. n ) localrc = ESMF_FAILURE
                  enddo
                 enddo
                enddo
               enddo
              enddo
             enddo
            enddo
            if (ESMF_LogFoundError(localrc, &
                ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(field, minIndex = gminIndex, maxIndex = gmaxIndex, &
                                    elementCount = geleCount, &
                                    localMinIndex = lminIndex, &
                                    localMaxIndex = lmaxIndex, &
                                    localelementCount = leleCount, &
                                    rc=localrc)
            if (ESMF_LogFoundError(localrc, &
                ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
          endif ! fieldget = .true.
        endif ! present(fieldget) = .true.

        if(associated(farray1, farray)) then
            deallocate(farray1)
        else
            deallocate(farray)
        endif

        call ESMF_FieldDestroy(field, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_GridDestroy(grid, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_DistGridDestroy(distgrid, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    end subroutine test7d3_generic_fptr