ESMF_XGridValidate Subroutine

public subroutine ESMF_XGridValidate(xgrid, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_XGrid), intent(inout) :: xgrid
integer, intent(out), optional :: rc

Source Code

      subroutine ESMF_XGridValidate(xgrid, rc)
!
! !ARGUMENTS:
      type(ESMF_XGrid), intent(inout) :: xgrid 
      integer, intent(out), optional :: rc   
!
! !DESCRIPTION:
!      Validates that the {\tt xgrid} is internally consistent.
!      Currently this method determines if the {\tt xgrid} is uninitialized 
!      or already destroyed. 
!
!      The method returns an error code if problems are found.  
!
!     The arguments are:
!     \begin{description}
!     \item [xgrid]
!           {\tt ESMF\_XGrid} to validate.
!     \item [{[rc]}]
!           Return code; equals {\tt ESMF\_SUCCESS} if the {\tt xgrid} 
!           is valid.
!     \end{description}
!
!EOPI

      integer :: localrc, i

      type(ESMF_XGridType), pointer :: xgtypep
      type(ESMF_Status) :: xgridstatus
      integer :: ngridA, ngridB

      ! Initialize
      localrc = ESMF_RC_NOT_IMPL
      if (present(rc)) rc = ESMF_RC_NOT_IMPL

      ! Check variables
      ESMF_INIT_CHECK_DEEP(ESMF_XGridGetInit,xgrid,rc)

      if (.not.associated(xgrid%xgtypep)) then 
         call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
            msg="Uninitialized or already destroyed XGrid: xgtypep unassociated", &
             ESMF_CONTEXT, rcToReturn=rc)
         return
      endif 

      xgtypep => xgrid%xgtypep

      ! Make sure the xgrid is ready before trying to look at contents
      call ESMF_BaseGetStatus(xgtypep%base, xgridstatus, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
      if (xgridstatus .ne. ESMF_STATUS_READY) then
         call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
            msg="Uninitialized or already destroyed XGrid: xgridstatus not ready", &
             ESMF_CONTEXT, rcToReturn=rc)
         return
      endif 

      ! Verify all associated arrays are sized consistently
      if(associated(xgtypep%sideA)) ngridA = size(xgtypep%sideA, 1)
      if(associated(xgtypep%sideB)) ngridB = size(xgtypep%sideB, 1)

      if(associated(xgtypep%area) .and. associated(xgtypep%centroid)) then
        if(size(xgtypep%area, 1) /= size(xgtypep%centroid, 1)) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
             msg="number of area cells differs from number of centroid cells", &
             ESMF_CONTEXT, rcToReturn=rc)
           return
        endif
      endif

      if(associated(xgtypep%sparseMatA2X) .or. associated(xgtypep%sparseMatX2A)) then
        if(associated(xgtypep%distgridA)) then
          if(size(xgtypep%distgridA, 1) /= ngridA) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
               msg="number of distgrids on side A differs from number of Grids on side A", &
               ESMF_CONTEXT, rcToReturn=rc)
             return
          endif
        endif
      endif

      if(associated(xgtypep%sparseMatA2X)) then
        if(size(xgtypep%sparseMatA2X, 1) /= ngridA) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
             msg="number of sparseMat objects on side A differs from number of Grids on side A", &
             ESMF_CONTEXT, rcToReturn=rc)
           return
        endif
        do i = 1, ngridA
          if(size(xgtypep%sparseMatA2X(i)%factorIndexList, 2) /= &
            size(xgtypep%sparseMatA2X(i)%factorList, 1)) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
               msg="number of elements in sparseMatA2X is inconsistent", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
          endif
        enddo
      endif

      if(associated(xgtypep%sparseMatX2A)) then
        if(size(xgtypep%sparseMatX2A, 1) /= ngridA) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
             msg="number of sparseMat objects on side A differs from number of Grids on side A", &
             ESMF_CONTEXT, rcToReturn=rc)
           return
        endif
        do i = 1, ngridA
          if(size(xgtypep%sparseMatX2A(i)%factorIndexList, 2) /= &
            size(xgtypep%sparseMatX2A(i)%factorList, 1)) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
               msg="number of elements in sparseMatX2A is inconsistent", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
          endif
        enddo
      endif

      ! B side
      if(associated(xgtypep%sparseMatB2X) .or. associated(xgtypep%sparseMatX2B)) then
        if(associated(xgtypep%distgridB)) then
          if(size(xgtypep%distgridB, 1) /= ngridB) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
               msg="number of distgrids on side B differs from number of Grids on side B", &
               ESMF_CONTEXT, rcToReturn=rc)
             return
          endif
        endif
      endif

      if(associated(xgtypep%sparseMatB2X)) then
        if(size(xgtypep%sparseMatB2X, 1) /= ngridB) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
             msg="number of sparseMat objects on side B differs from number of Grids on side B", &
             ESMF_CONTEXT, rcToReturn=rc)
           return
        endif
        do i = 1, ngridB
          if(size(xgtypep%sparseMatB2X(i)%factorIndexList, 2) /= &
            size(xgtypep%sparseMatB2X(i)%factorList, 1)) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
               msg="number of elements in sparseMatB2X is inconsistent", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
          endif
        enddo
      endif

      if(associated(xgtypep%sparseMatX2B)) then
        if(size(xgtypep%sparseMatX2B, 1) /= ngridB) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
             msg="number of sparseMat objects on side B differs from number of Grids on side B", &
             ESMF_CONTEXT, rcToReturn=rc)
           return
        endif
        do i = 1, ngridB
          if(size(xgtypep%sparseMatX2B(i)%factorIndexList, 2) /= &
            size(xgtypep%sparseMatX2B(i)%factorList, 1)) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_OBJ_BAD, &
               msg="number of elements in sparseMatX2B is inconsistent", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
          endif
        enddo
      endif

      if (present(rc)) rc = ESMF_SUCCESS

      end subroutine ESMF_XGridValidate