test_gather_1d Subroutine

subroutine test_gather_1d(totalLWidth, totalUWidth, dgCase, rc)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: totalLWidth(:)
integer, intent(in) :: totalUWidth(:)
character(len=*), intent(in) :: dgCase
integer, intent(out) :: rc

Calls

proc~~test_gather_1d~~CallsGraph proc~test_gather_1d test_gather_1d esmf_arraycreate esmf_arraycreate proc~test_gather_1d->esmf_arraycreate esmf_arraydestroy esmf_arraydestroy proc~test_gather_1d->esmf_arraydestroy esmf_arraygather esmf_arraygather proc~test_gather_1d->esmf_arraygather esmf_arrayget esmf_arrayget proc~test_gather_1d->esmf_arrayget interface~esmf_distgridcreate ESMF_DistGridCreate proc~test_gather_1d->interface~esmf_distgridcreate interface~esmf_vmget ESMF_VMGet proc~test_gather_1d->interface~esmf_vmget proc~esmf_arrayspecset ESMF_ArraySpecSet proc~test_gather_1d->proc~esmf_arrayspecset proc~esmf_distgriddestroy ESMF_DistGridDestroy proc~test_gather_1d->proc~esmf_distgriddestroy proc~esmf_logfounderror ESMF_LogFoundError proc~test_gather_1d->proc~esmf_logfounderror proc~esmf_logseterror ESMF_LogSetError proc~test_gather_1d->proc~esmf_logseterror proc~esmf_vmgetcurrent ESMF_VMGetCurrent proc~test_gather_1d->proc~esmf_vmgetcurrent

Called by

proc~~test_gather_1d~~CalledByGraph proc~test_gather_1d test_gather_1d program~esmf_arraygatherutest ESMF_ArrayGatherUTest program~esmf_arraygatherutest->proc~test_gather_1d

Source Code

    subroutine test_gather_1d(totalLWidth, totalUWidth, dgCase, rc)
        integer, intent(in)       :: totalLWidth(:), totalUWidth(:)
        character(*), intent(in)  :: dgCase
        integer, intent(out)      :: rc

        ! local arguments used to create field etc
        type(ESMF_DistGrid)                         :: distgrid
        type(ESMF_VM)                               :: vm
        type(ESMF_Array)                            :: array
        type(ESMF_ArraySpec)                        :: arrayspec
        integer                                     :: localrc, localPet, i, j
        integer, allocatable                        :: deBlockList(:,:,:)

        integer, pointer                            :: farray(:)
        integer, pointer                            :: farrayDst(:)

        rc = ESMF_SUCCESS
        localrc = ESMF_SUCCESS

        call ESMF_VMGetCurrent(vm, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_VMGet(vm, localPet=localPet, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
          
        if (trim(dgCase)=="regDecomp") then
          ! default DistGrid with regDecomp
          distgrid = ESMF_DistGridCreate(minIndex =(/1/), maxIndex=(/16/), &
            rc=localrc)
          if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        else if (trim(dgCase)=="deBlockList") then
          ! DistGrid with deBlockList
          allocate(deBlockList(1, 2, 4))  ! dimCount, 2, deCount
          deBlockList(1, 1, 1) = 5        ! 1st DE minIndex
          deBlockList(1, 2, 1) = 8        ! 1st DE maxIndex
          deBlockList(1, 1, 2) = 9        ! 2nd DE minIndex
          deBlockList(1, 2, 2) = 12       ! 2nd DE maxIndex
          deBlockList(1, 1, 3) = 13       ! 3rd DE minIndex
          deBlockList(1, 2, 3) = 16       ! 3rd DE maxIndex
          deBlockList(1, 1, 4) = 1        ! 4th DE minIndex
          deBlockList(1, 2, 4) = 4        ! 4th DE maxIndex
          distgrid = ESMF_DistGridCreate(minIndex =(/1/), maxIndex=(/16/), &
            deBlockList=deBlockList, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
          deallocate(deBlockList)
        else
          call ESMF_LogSetError(ESMF_RC_ARG_BAD, &
            msg="An invalid 'option' argument was provided.", &
            ESMF_CONTEXT, rcToReturn=rc)
          return
        endif

        call ESMF_ArraySpecSet(arrayspec, typekind=ESMF_TYPEKIND_I4, rank=1, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

        array = ESMF_ArrayCreate(distgrid, arrayspec, &
          totalLWidth=totalLWidth, totalUWidth=totalUWidth, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_ArrayGet(array, farrayPtr=farray, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

        farray = 0  ! initialize the entire local array
        do i=1, 4
          if (trim(dgCase)=="regDecomp") then
            farray(lbound(farray,1)+totalLWidth(1)-1+i) = localPet * 10 + i
          else if (trim(dgCase)=="deBlockList") then
            farray(lbound(farray,1)+totalLWidth(1)-1+i) = mod(localPet+1,4) * 10 + i
          endif
        enddo

        if(localPet .eq. 0) then
          allocate(farrayDst(16))  ! rootPet
        else
          allocate(farrayDst(1))
        end if
        call ESMF_ArrayGather(array, farrayDst, rootPet=0, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
            
        ! check that the values gathered on rootPet are correct
        if(localPet .eq. 0) then
          do j = 1, 4
            do i = 1, 4
              if(farrayDst((j-1)*4+i) .ne. (j-1)*10+i) then
                localrc=ESMF_FAILURE
              endif
            enddo
          enddo
          if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        endif

        call ESMF_ArrayDestroy(array, 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
        
        deallocate(farrayDst)

        rc = ESMF_SUCCESS
    end subroutine test_gather_1d