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

Calls

proc~~esmf_xgridvalidate~~CallsGraph proc~esmf_xgridvalidate ESMF_XGridValidate proc~esmf_basegetstatus ESMF_BaseGetStatus proc~esmf_xgridvalidate->proc~esmf_basegetstatus proc~esmf_imerr ESMF_IMErr proc~esmf_xgridvalidate->proc~esmf_imerr proc~esmf_logfounderror ESMF_LogFoundError proc~esmf_xgridvalidate->proc~esmf_logfounderror proc~esmf_logseterror ESMF_LogSetError proc~esmf_xgridvalidate->proc~esmf_logseterror proc~esmf_xgridgetinit ESMF_XGridGetInit proc~esmf_xgridvalidate->proc~esmf_xgridgetinit proc~esmf_basegetstatus->proc~esmf_logfounderror c_esmc_basegetstatus c_esmc_basegetstatus proc~esmf_basegetstatus->c_esmc_basegetstatus proc~esmf_imerr->proc~esmf_logfounderror proc~esmf_initcheckdeep ESMF_InitCheckDeep proc~esmf_imerr->proc~esmf_initcheckdeep esmf_breakpoint esmf_breakpoint proc~esmf_logfounderror->esmf_breakpoint proc~esmf_logrc2msg ESMF_LogRc2Msg proc~esmf_logfounderror->proc~esmf_logrc2msg proc~esmf_logwrite ESMF_LogWrite proc~esmf_logfounderror->proc~esmf_logwrite proc~esmf_logseterror->esmf_breakpoint proc~esmf_logseterror->proc~esmf_logrc2msg proc~esmf_logseterror->proc~esmf_logwrite c_esmc_loggeterrormsg c_esmc_loggeterrormsg proc~esmf_logrc2msg->c_esmc_loggeterrormsg c_esmc_vmwtime c_esmc_vmwtime proc~esmf_logwrite->c_esmc_vmwtime proc~esmf_logclose ESMF_LogClose proc~esmf_logwrite->proc~esmf_logclose proc~esmf_logflush ESMF_LogFlush proc~esmf_logwrite->proc~esmf_logflush proc~esmf_logopenfile ESMF_LogOpenFile proc~esmf_logwrite->proc~esmf_logopenfile proc~esmf_utiliounitflush ESMF_UtilIOUnitFlush proc~esmf_logwrite->proc~esmf_utiliounitflush proc~esmf_utilstring2array ESMF_UtilString2Array proc~esmf_logwrite->proc~esmf_utilstring2array proc~esmf_logclose->proc~esmf_logflush proc~esmf_logflush->proc~esmf_utiliounitflush proc~esmf_utilarray2string ESMF_UtilArray2String proc~esmf_logflush->proc~esmf_utilarray2string proc~esmf_logopenfile->proc~esmf_utiliounitflush proc~esmf_utiliounitget ESMF_UtilIOUnitGet proc~esmf_logopenfile->proc~esmf_utiliounitget

Called by

proc~~esmf_xgridvalidate~~CalledByGraph proc~esmf_xgridvalidate ESMF_XGridValidate proc~checkproxy checkProxy proc~checkproxy->proc~esmf_xgridvalidate proc~esmf_geomvalidate ESMF_GeomValidate proc~esmf_geomvalidate->proc~esmf_xgridvalidate proc~esmf_xgridcreatefromsparsemat ESMF_XGridCreateFromSparseMat proc~esmf_xgridcreatefromsparsemat->proc~esmf_xgridvalidate proc~esmf_fieldvalidate ESMF_FieldValidate proc~esmf_fieldvalidate->proc~esmf_geomvalidate proc~test3 test3 proc~test3->proc~esmf_xgridcreatefromsparsemat proc~test_regrid2xg test_regrid2xg proc~test_regrid2xg->proc~esmf_xgridcreatefromsparsemat proc~test_regrid2xg_clip test_regrid2xg_clip proc~test_regrid2xg_clip->proc~esmf_xgridcreatefromsparsemat proc~test_regrid2xg_contain test_regrid2xg_contain proc~test_regrid2xg_contain->proc~esmf_xgridcreatefromsparsemat proc~test_regrid2xg_half test_regrid2xg_half proc~test_regrid2xg_half->proc~esmf_xgridcreatefromsparsemat proc~test_regrid2xg_online test_regrid2xg_online proc~test_regrid2xg_online->proc~checkproxy proc~test_regrid2xgsph test_regrid2xgSph proc~test_regrid2xgsph->proc~esmf_xgridcreatefromsparsemat proc~test_regrid2xgsph~2 test_regrid2xgSph proc~test_regrid2xgsph~2->proc~checkproxy proc~test_regridxg test_regridxg proc~test_regridxg->proc~esmf_xgridcreatefromsparsemat proc~test_regridxg_const test_regridxg_const proc~test_regridxg_const->proc~esmf_xgridcreatefromsparsemat program~esmf_xgridsparsematex ESMF_XGridSparseMatEx program~esmf_xgridsparsematex->proc~esmf_xgridcreatefromsparsemat proc~esmf_statevalidate ESMF_StateValidate proc~esmf_statevalidate->proc~esmf_fieldvalidate proc~test7d4_generic test7d4_generic proc~test7d4_generic->proc~esmf_fieldvalidate proc~test_mesh_create_gt_1localde test_mesh_create_gt_1localde proc~test_mesh_create_gt_1localde->proc~esmf_fieldvalidate program~esmf_meshutest ESMF_MeshUTest program~esmf_meshutest->proc~esmf_fieldvalidate program~esmf_meshutest->proc~test_mesh_create_gt_1localde program~esmf_xgridutest ESMF_XGridUTest program~esmf_xgridutest->proc~test3 program~esmf_statereconcileutest ESMF_StateReconcileUTest program~esmf_statereconcileutest->proc~esmf_statevalidate

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