user_run Subroutine

public subroutine user_run(comp, importState, exportState, clock, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_GridComp) :: comp
type(ESMF_State) :: importState
type(ESMF_State) :: exportState
type(ESMF_Clock) :: clock
integer, intent(out) :: rc

Source Code

    subroutine user_run(comp, importState, exportState, clock, rc)
      type(ESMF_GridComp) :: comp
      type(ESMF_State) :: importState, exportState
      type(ESMF_Clock) :: clock
      integer, intent(out) :: rc

!   ! Local variables
      type(ESMF_Field) :: dstField
      type(ESMF_Grid) :: dstGrid
      integer :: lDE, localDECount, localrc,i1,i2
      real(ESMF_KIND_R8), pointer :: farrayPtrXC(:,:), farrayPtr1D(:)
      real(ESMF_KIND_R8), pointer :: farrayPtrYC(:,:)
      real(ESMF_KIND_R8), pointer :: farrayPtr(:,:),farrayPtr2(:,:)
      integer :: clbnd(2),cubnd(2)

      rc = ESMF_SUCCESS

      ! Get information from the component.
      call ESMF_StateGet(importState, "dst", dstField, rc=rc)
      if(rc/=ESMF_SUCCESS) return


      ! Get Grid from field
      call ESMF_FieldGet(dstField, grid=dstGrid, rc=localrc)
      if (localrc /=ESMF_SUCCESS) then
          rc=ESMF_FAILURE
          return
      endif

      ! Get number of local DEs
      call ESMF_GridGet(dstGrid, localDECount=localDECount, rc=localrc)
      if (localrc /=ESMF_SUCCESS) then
          rc=ESMF_FAILURE
          return
      endif

      ! Check error
      do lDE=0,localDECount-1

         !! get coords
         call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1, &
                                farrayPtr=farrayPtrXC, rc=localrc)
         if (localrc /=ESMF_SUCCESS) then
            rc=ESMF_FAILURE
            return
         endif

         call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=2, &
                               farrayPtr=farrayPtrYC, rc=localrc)
         if (localrc /=ESMF_SUCCESS) then
             rc=ESMF_FAILURE
             return
         endif

         call ESMF_FieldGet(dstField, lDE, farrayPtr, computationalLBound=clbnd, &
                            computationalUBound=cubnd,  rc=localrc)
         if (localrc /=ESMF_SUCCESS) then
             rc=ESMF_FAILURE
            return
         endif

         !! check error
         do i1=clbnd(1),cubnd(1)
          do i2=clbnd(2),cubnd(2)

             !! if error is too big report an error
             if (abs(farrayPtr(i1,i2)-(20.0+farrayPtrXC(i1,i2)+farrayPtrYC(i1,i2))) > 0.0001) then
!                 write(*,*) farrayPtr(i1,i2),".ne.",(20.0+farrayPtrXC(i1,i2)+farrayPtrYC(i1,i2))
                 rc=ESMF_FAILURE
                 return
             endif      
         enddo
       enddo

       ! RESET DESTINATION BACK TO 0
       farrayPtr=0.0

      enddo    ! lDE


    end subroutine user_run