user_init Subroutine

public subroutine user_init(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_init(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_VM) :: vm
      type(ESMF_grid) :: dstGrid
      type(ESMF_ArraySpec) :: arrayspec
      integer :: localPET, petCount
      integer dst_nx, dst_ny, i1,i2
      real(ESMF_KIND_R8) :: dst_minx,dst_miny
      real(ESMF_KIND_R8) :: dst_maxx,dst_maxy
      integer :: lDE, localDECount, localrc
      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
      ! Initially import state contains a field with a grid but no data.

      ! Query component for VM and create a layout with the right breakdown
      call ESMF_GridCompGet(comp, vm=vm, rc=rc)
      if(rc/=ESMF_SUCCESS) return
      call ESMF_VMGet(vm, localPet=localPet, petCount=petCount, rc=rc)
      if(rc/=ESMF_SUCCESS) return


      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      ! Create Destination grid
      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

      ! Establish the resolution of the grids
      dst_nx = 10
      dst_ny = 10

      ! Establish the coordinates of the grids
      dst_minx = 0.1
      dst_miny = 0.1

      dst_maxx = 1.9
      dst_maxy = 1.9

      ! Create Grid
      dstGrid=ESMF_GridCreateNoPeriDim(minIndex=(/1,1/),maxIndex=(/dst_nx,dst_ny/), &
                coordSys=ESMF_COORDSYS_CART, &
                regDecomp=(/2,2/), indexflag=ESMF_INDEX_GLOBAL, rc=localrc)
      if (localrc /=ESMF_SUCCESS) then
         rc=ESMF_FAILURE
         return
      endif

      ! Create source/destination fields
      call ESMF_ArraySpecSet(arrayspec, 2, ESMF_TYPEKIND_R8, rc=rc)

      dstField = ESMF_FieldCreate(dstGrid, arrayspec, &
                         staggerloc=ESMF_STAGGERLOC_CENTER, name="dst", rc=localrc)
      if (localrc /=ESMF_SUCCESS) then
          rc=ESMF_FAILURE
         return
      endif

      call ESMF_GridAddCoord(dstGrid, staggerloc=ESMF_STAGGERLOC_CENTER, 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

      ! Get memory and set coords for dst
      do lDE=0,localDECount-1

         !! get coords
         call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1, &
                           computationalLBound=clbnd, computationalUBound=cubnd, 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, &
                            computationalLBound=clbnd, computationalUBound=cubnd, farrayPtr=farrayPtrYC, rc=localrc)
         if (localrc /=ESMF_SUCCESS) then
            rc=ESMF_FAILURE
            return
         endif


         call ESMF_FieldGet(dstField, lDE, farrayPtr, rc=localrc)
         if (localrc /=ESMF_SUCCESS) then
             rc=ESMF_FAILURE
             return
         endif

        !! set coords
        do i1=clbnd(1),cubnd(1)
        do i2=clbnd(2),cubnd(2)

           ! Set source coordinates
           farrayPtrXC(i1,i2) = ((dst_maxx-dst_minx)*REAL(i1-1)/REAL(dst_nx-1))+dst_minx
           farrayPtrYC(i1,i2) = ((dst_maxy-dst_miny)*REAL(i2-1)/REAL(dst_ny-1))+dst_miny

           ! initialize destination field
           farrayPtr(i1,i2)=0.0

        enddo
        enddo
      enddo    ! lDE

      ! Set Field Into State
      call ESMF_StateAdd(importState, (/dstField/), rc=rc)
      if(rc/=ESMF_SUCCESS) return

    end subroutine user_init