compare_redist_array Subroutine

private subroutine compare_redist_array(test_status, Array1, Array2, DistGrid1, DistGrid2, Memory, Grid, rc)

Arguments

Type IntentOptional Attributes Name
integer, intent(inout) :: test_status
type(ESMF_Array), intent(in) :: Array1
type(ESMF_Array), intent(in) :: Array2
type(ESMF_DistGrid), intent(in) :: DistGrid1
type(ESMF_DistGrid), intent(in) :: DistGrid2
type(memory_config), intent(in) :: Memory
type(grid_specification_record), intent(in) :: Grid
integer, intent(inout) :: rc

Source Code

  subroutine compare_redist_array(test_status, Array1, Array2,                 &
                                  DistGrid1, DistGrid2,                        &
                                  Memory, Grid, rc)
  !-----------------------------------------------------------------------------
  ! routine compares the contents of two arrays and returns if they agree
  !-----------------------------------------------------------------------------
  ! arguments
  integer, intent(inout) :: test_status
  type(ESMF_Array), intent(in   ) :: Array1, Array2
  type(ESMF_DistGrid), intent(in   ) :: DistGrid1, DistGrid2
  type(memory_config), intent(in   ) :: Memory
  type(grid_specification_record), intent(in   ) :: Grid
  integer, intent(inout) :: rc

  ! local ESMF types
  type(ESMF_LocalArray), allocatable :: larrayList1(:)
  type(ESMF_LocalArray), allocatable :: larrayList2(:)
  type(ESMF_Index_Flag) :: indexflag

  ! local integer variables
  integer :: de, i1, i2, i3, i4, i5, i6, i7, k
  integer :: localDeCount1, dimCount1, localDeCount2, dimCount2
  integer, allocatable ::  localDeToDeMap1(:), localDeToDeMap2(:)
  integer, allocatable :: LBnd(:,:), UBnd(:,:)
  integer, allocatable :: LBnd2(:,:), UBnd2(:,:)
  ! integer :: irank, tensorsize, fsize(7)
  ! integer, allocatable :: haloL(:), haloR(:)
  ! integer, allocatable :: top(:), bottom(:)
  integer :: allocRcToTest
  integer :: localrc ! local error status

  ! local logicals
  ! logical :: nohaloflag

  ! local real variables
  real(ESMF_KIND_R8), pointer :: farray1D(:), farray2D(:,:)
  real(ESMF_KIND_R8), pointer :: farray3D(:,:,:)
  real(ESMF_KIND_R8), pointer :: farray4D(:,:,:,:)
  real(ESMF_KIND_R8), pointer :: farray5D(:,:,:,:,:)
  real(ESMF_KIND_R8), pointer :: farray6D(:,:,:,:,:,:)
  real(ESMF_KIND_R8), pointer :: farray7D(:,:,:,:,:,:,:)

  real(ESMF_KIND_R8), pointer :: rarray1D(:), rarray2D(:,:)
  real(ESMF_KIND_R8), pointer :: rarray3D(:,:,:)
  real(ESMF_KIND_R8), pointer :: rarray4D(:,:,:,:)
  real(ESMF_KIND_R8), pointer :: rarray5D(:,:,:,:,:)
  real(ESMF_KIND_R8), pointer :: rarray6D(:,:,:,:,:,:)
  real(ESMF_KIND_R8), pointer :: rarray7D(:,:,:,:,:,:,:)


  ! initialize return flag
  localrc = ESMF_RC_NOT_IMPL
  rc = ESMF_RC_NOT_IMPL
  test_status = HarnessTest_SUCCESS

  !-----------------------------------------------------------------------------
  ! Sanity check - confirm that the two arrays have the same ranks and sizes
  !-----------------------------------------------------------------------------

  !-----------------------------------------------------------------------------
  ! get local array DE list from array1
  !-----------------------------------------------------------------------------
  call ESMF_ArrayGet(array1, localDeCount=localDeCount1, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local DE count from array", &
          rcToReturn=rc)) return

  allocate(localDeToDeMap1(localDeCount1), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="integer variable"//           &
     " localDeToDeMap1 in compare redist array", rcToReturn=rc)) then
  endif
  call ESMF_ArrayGet(array1, localDeToDeMap=localDeToDeMap1, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local DE list from array",  &
          rcToReturn=rc)) return

  allocate(larrayList1(localDeCount1), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="char variable"//              &
     " larrayList1 in compare redist array", rcToReturn=rc)) then
  endif
  call ESMF_ArrayGet(array1, localarrayList=larrayList1, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local array list",          &
          rcToReturn=rc)) return

  call ESMF_DistGridGet(DistGrid1, dimCount=dimCount1, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting dimCount from distGrid",    &
          rcToReturn=rc)) return

  !-----------------------------------------------------------------------------
  ! get local array DE list from array2
  !-----------------------------------------------------------------------------
  call ESMF_ArrayGet(array2, localDeCount=localDeCount2, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local DE count from array", &
          rcToReturn=rc)) return

  ! check localDeCount for agreement
  if( localDeCount1 /= localDeCount2 ) then
     localrc = ESMF_FAILURE
     call ESMF_LogSetError(ESMF_FAILURE, msg="local De Counts do not agree",     &
              rcToReturn=localrc)
     return
  endif

  allocate(localDeToDeMap2(localDeCount2), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="integer variable"//           &
     " localDeToDeMap2 in compare redist array", rcToReturn=rc)) then
  endif
  call ESMF_ArrayGet(array2, localDeToDeMap=localDeToDeMap2, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local DE list from array",  &
          rcToReturn=rc)) return

  ! check localDeToDeMap for agreement
  do de=1, localDeCount2
     if( localDeToDeMap2(de) /= localDeToDeMap1(de) ) then
        localrc = ESMF_FAILURE
        call ESMF_LogSetError(ESMF_FAILURE, msg="local De lists do not agree",   &
                 rcToReturn=localrc)
        return
     endif
  enddo

  ! get local De List
  allocate(larrayList2(localDeCount2), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="char variable"//              &
     " larrayList2 in compare redist array", rcToReturn=rc)) then
  endif
  call ESMF_ArrayGet(array2, localarrayList=larrayList2, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting local array list",          &
          rcToReturn=rc)) return

  !-----------------------------------------------------------------------------
  ! compare dimcounts for both arrays
  !-----------------------------------------------------------------------------
  call ESMF_DistGridGet(DistGrid2, dimCount=dimCount2, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting dimCount from distGrid",    &
          rcToReturn=rc)) return

  ! check localDeCount for agreement
  if( dimCount1 /= dimCount2 ) then
     localrc = ESMF_FAILURE
     call ESMF_LogSetError(ESMF_FAILURE, msg="array 1 and 2 dimCounts disagree", &
              rcToReturn=localrc)
     return
  endif

  !-----------------------------------------------------------------------------
  ! allocate bound arrays and extract exclusive bounds
  !-----------------------------------------------------------------------------
  allocate(UBnd(dimCount1, localDeCount1), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="integer variable"//           &
     " UBnd in compare redist array", rcToReturn=rc)) then
  endif
  allocate(LBnd(dimCount1, localDeCount1), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="integer variable"//           &
     " LBnd in compare redist array", rcToReturn=rc)) then
  endif

  call ESMF_ArrayGet(array=array1, indexflag=indexflag,                        &
           exclusiveLBound=LBnd, exclusiveUBound=UBnd, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting exclusive bound range",     &
          rcToReturn=rc)) return


  allocate(UBnd2(dimCount2, localDeCount2), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="integer variable"//           &
     " UBnd2 in compare redist array", rcToReturn=rc)) then
  endif
  allocate(LBnd2(dimCount2, localDeCount2), stat=allocRcToTest )
  if (ESMF_LogFoundAllocError(allocRcToTest, msg="integer variable"//           &
     " LBnd2 in compare redist array", rcToReturn=rc)) then
  endif

  call ESMF_ArrayGet(array=array2, indexflag=indexflag,                        &
           exclusiveLBound=LBnd2, exclusiveUBound=UBnd2, rc=localrc)
  if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error getting exclusive bound range",     &
          rcToReturn=rc)) return

  !-----------------------------------------------------------------------------
  ! check for exclusive bound agreement
  !-----------------------------------------------------------------------------
  do de=1, localDeCount1
     do k=1,dimCount1
        if( LBnd(k,de) /= LBnd2(k,de) ) then
           print*,'exclusive Lower bounds disagree',LBnd(k,de),LBnd2(k,de)
           localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE,msg="exclusive L bounds disagree",&
              rcToReturn=localrc)
           return
        endif
        if( UBnd(k,de) /= UBnd2(k,de) ) then
           print*,'exclusive Upper bounds disagree',UBnd(k,de),UBnd2(k,de)
           localrc = ESMF_FAILURE
           call ESMF_LogSetError(ESMF_FAILURE,msg="exclusive U bounds disagree",&
              rcToReturn=localrc)
           return
        endif
     enddo
  enddo
  !-----------------------------------------------------------------------------
  ! associate fortran pointers with the two local array objects and compare
  ! entries
  !-----------------------------------------------------------------------------

  !-----------------------------------------------------------------------------
  ! Memory Rank = Grid Rank, then there are no tensor dimensions
  !-----------------------------------------------------------------------------
  if( Memory%memRank ==  Memory%GridRank ) then

     select case(dimCount1)
     case(1)
     !--------------------------------------------------------------------------
     ! rank = 1
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray1D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray1D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
              if( farray1D(i1) /= rarray1D(i1) ) then
                 test_status = HarnessTest_FAILURE
                 print*,' arrays disagree ',i1,farray1D(i1),                   &
                        rarray1D(i1)
              endif
           enddo    !   i1
        enddo    ! de
     case(2)
     !--------------------------------------------------------------------------
     ! rank = 2
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray2D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray2D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
             do i2=LBnd(2,de), UBnd(2,de)
               if( farray2D(i1,i2) /= rarray2D(i1,i2) ) then
                 test_status = HarnessTest_FAILURE
                 print*,' arrays disagree ',i1,i2,farray2D(i1,i2),             &
                        rarray2D(i1,i2)
               endif
             enddo   !   i2
           enddo    !   i1
        enddo    ! de
     case(3)
     !--------------------------------------------------------------------------
     ! rank = 3
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray3D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray3D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
             do i2=LBnd(2,de), UBnd(2,de)
                do i3=LBnd(3,de), UBnd(3,de)
                   if( farray3D(i1,i2,i3) /= rarray3D(i1,i2,i3) ) then
                       test_status = HarnessTest_FAILURE
                       print*,' arrays disagree ',i1,i2,i3,farray3D(i1,i2,i3), &
                              rarray3D(i1,i2,i3)
                   endif
                enddo   !   i3
             enddo   !   i2
           enddo    !   i1
        enddo    ! de
     case(4)
     !--------------------------------------------------------------------------
     ! rank = 4
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray4D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray4D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
             do i2=LBnd(2,de), UBnd(2,de)
                do i3=LBnd(3,de), UBnd(3,de)
                   do i4=LBnd(4,de), UBnd(4,de)
                      if( farray4D(i1,i2,i3,i4) /= rarray4D(i1,i2,i3,i4) ) then
                       test_status = HarnessTest_FAILURE
                       print*,' arrays disagree ',i1,i2,i3,i4,                 &
                              farray4D(i1,i2,i3,i4), rarray4D(i1,i2,i3,i4)
                      endif
                   enddo   !   i4
                enddo   !   i3
             enddo   !   i2
           enddo    !   i1
        enddo    ! de
#ifndef ESMF_NO_GREATER_THAN_4D
     case(5)
     !--------------------------------------------------------------------------
     ! rank = 5
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray5D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray5D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
             do i2=LBnd(2,de), UBnd(2,de)
                do i3=LBnd(3,de), UBnd(3,de)
                   do i4=LBnd(4,de), UBnd(4,de)
                      do i5=LBnd(5,de), UBnd(5,de)
                         if( farray5D(i1,i2,i3,i4,i5) /=                       &
                             rarray5D(i1,i2,i3,i4,i5) ) then
                               test_status = HarnessTest_FAILURE
                               print*,' arrays disagree ',i1,i2,i3,i4,i5,      &
                              farray5D(i1,i2,i3,i4,i5),                        &
                              rarray5D(i1,i2,i3,i4,i5)
                         endif
                      enddo   !   i5
                   enddo   !   i4
                enddo   !   i3
             enddo   !   i2
           enddo    !   i1
        enddo    ! de
     case(6)
     !--------------------------------------------------------------------------
     ! rank = 6
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray6D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray6D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
             do i2=LBnd(2,de), UBnd(2,de)
                do i3=LBnd(3,de), UBnd(3,de)
                   do i4=LBnd(4,de), UBnd(4,de)
                      do i5=LBnd(5,de), UBnd(5,de)
                      do i6=LBnd(6,de), UBnd(6,de)
                         if( farray6D(i1,i2,i3,i4,i5,i6) /=                    &
                             rarray6D(i1,i2,i3,i4,i5,i6) ) then
                               test_status = HarnessTest_FAILURE
                               print*,' arrays disagree ',i1,i2,i3,i4,i5,i6,   &
                              farray6D(i1,i2,i3,i4,i5,i6),                     &
                              rarray6D(i1,i2,i3,i4,i5,i6)
                         endif
                      enddo   !   i6
                      enddo   !   i5
                   enddo   !   i4
                enddo   !   i3
             enddo   !   i2
           enddo    !   i1
        enddo    ! de
     case(7)
     !--------------------------------------------------------------------------
     ! rank = 7
     !--------------------------------------------------------------------------
        do de=1, localDeCount1
           call ESMF_LocalArrayGet(larrayList1(de), farrayPtr=farray7D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                  "array1 list", rcToReturn=rc)) return

           call ESMF_LocalArrayGet(larrayList2(de), farrayPtr=rarray7D,             &
                    datacopyflag=ESMF_DATACOPY_REFERENCE, rc=localrc)

           if (CheckError(checkpoint, __LINE__, __FILE__, localrc,"error connecting pointer to " // &
                   "array2 list", rcToReturn=rc)) return

           do i1=LBnd(1,de), UBnd(1,de)
             do i2=LBnd(2,de), UBnd(2,de)
                do i3=LBnd(3,de), UBnd(3,de)
                   do i4=LBnd(4,de), UBnd(4,de)
                      do i5=LBnd(5,de), UBnd(5,de)
                      do i6=LBnd(6,de), UBnd(6,de)
                      do i7=LBnd(7,de), UBnd(7,de)
                         if( farray7D(i1,i2,i3,i4,i5,i6,i7) /=                 &
                             rarray7D(i1,i2,i3,i4,i5,i6,i7) ) then
                              test_status = HarnessTest_FAILURE
                              print*,' arrays disagree ',i1,i2,i3,i4,i5,i6,i7, &
                              farray7D(i1,i2,i3,i4,i5,i6,i7),                  &
                              rarray7D(i1,i2,i3,i4,i5,i6,i7)
                         endif
                      enddo   !   i7
                      enddo   !   i6
                      enddo   !   i5
                   enddo   !   i4
                enddo   !   i3
             enddo   !   i2
           enddo    !   i1
        enddo    ! de
#endif
     case default
        ! error
        localrc = ESMF_FAILURE
        call ESMF_LogSetError(ESMF_FAILURE, msg="DimCount out of range",  &
                 rcToReturn=localrc)
        return
     end select

  !-----------------------------------------------------------------------------
  ! Memory Rank > Grid Rank, then there are MemRank-GridRank tensor dimensions
  !-----------------------------------------------------------------------------
  elseif (Memory%memRank >  Memory%GridRank) then
    call ESMF_LogSetError (ESMF_FAILURE, msg="tensor dimensions not supported", &
         rcToReturn=rc)
    return

  !-----------------------------------------------------------------------------
  ! Memory Rank < Grid Rank, then parameters don't make sense
  !-----------------------------------------------------------------------------
  else
    call ESMF_LogSetError (ESMF_FAILURE, msg="memory rank < grid rank", &
         rcToReturn=rc)
    return
  endif

  !-----------------------------------------------------------------------------
  ! clean up allocated arrays
  !-----------------------------------------------------------------------------
  deallocate(larrayList1, larrayList2)
  deallocate(localDeToDeMap1, localDeToDeMap2)
  deallocate(LBnd, UBnd)
  deallocate( LBnd2, UBnd2 )

  !-----------------------------------------------------------------------------
  rc = ESMF_SUCCESS
  !-----------------------------------------------------------------------------

  !-----------------------------------------------------------------------------
  end subroutine compare_redist_array