test_locstreambkgsph Subroutine

subroutine test_locstreambkgsph(rc)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: rc

Source Code

      subroutine test_locstreambkgsph(rc)
        integer, intent(out)  :: rc
  logical :: correct
  integer :: localrc
  type(ESMF_Grid) :: gridA
  type(ESMF_VM) :: vm
  real(ESMF_KIND_R8), pointer :: farrayPtrLonC(:,:)
  real(ESMF_KIND_R8), pointer :: farrayPtrLatC(:,:)
  integer :: clbnd(2),cubnd(2)
  integer :: i1,i2,i3, index(2)
  real(ESMF_KIND_R8) :: coord(2)
  integer A_nlon, A_nlat
  integer :: localPet, petCount
  real(ESMF_KIND_R8) :: de_minlon, de_maxlon
  real(ESMF_KIND_R8) :: de_minlat, de_maxlat
  integer :: pntCount
  real(ESMF_KIND_R8), allocatable :: Lon(:),Lat(:)
  real(ESMF_KIND_R8), pointer :: tstLon(:),tstLat(:)
  type(ESMF_LocStream) :: locstream,  newlocstream
  real(ESMF_KIND_R8) :: tmpLonC, tmpLatC
  real(ESMF_KIND_R8) :: A_dlon, A_dlat

  ! result code
  integer :: finalrc
  
  ! init success flag
  correct=.true.
  rc=ESMF_SUCCESS

  ! get pet info
  call ESMF_VMGetGlobal(vm, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
    rc=ESMF_FAILURE
    return
  endif
 
  call ESMF_VMGet(vm, petCount=petCount, localPet=localpet, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
    rc=ESMF_FAILURE
    return
  endif

  ! Establish the resolution of the grids
  A_nlon = 16
  A_nlat = 16

  
  ! setup source grid
  gridA=ESMF_GridCreateNoPeriDim(minIndex=(/1,1/),maxIndex=(/A_nlon,A_nlat/),regDecomp=(/petCount,1/), &
                              coordSys=ESMF_COORDSYS_SPH_DEG, &
                              indexflag=ESMF_INDEX_GLOBAL, &
                              rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
    rc=ESMF_FAILURE
    return
  endif


  ! Allocate coordinates
  call ESMF_GridAddCoord(gridA, staggerloc=ESMF_STAGGERLOC_CORNER, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
    rc=ESMF_FAILURE
    return
  endif


  ! Construct Grid A
  ! (Get memory and set coords for src)
  !! get coord 1
  call ESMF_GridGetCoord(gridA, localDE=0, staggerLoc=ESMF_STAGGERLOC_CORNER, coordDim=1, &
                            computationalLBound=clbnd, computationalUBound=cubnd, farrayPtr=farrayPtrLonC, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
     rc=ESMF_FAILURE
     return
  endif

  call ESMF_GridGetCoord(gridA, localDE=0, staggerLoc=ESMF_STAGGERLOC_CORNER, coordDim=2, &
                            computationalLBound=clbnd, computationalUBound=cubnd, farrayPtr=farrayPtrLatC, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
        rc=ESMF_FAILURE
        return
  endif

  A_dlon=360.0/A_nlon
  A_dlat=180.0/A_nlat

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

    ! Set source coordinates as 0 to 360
    farrayPtrLonC(i1,i2) = REAL(i1-1)*A_dlon

    ! Set source coordinates as -90 to 90
    farrayPtrLatC(i1,i2) = -90. + (REAL(i2-1)*A_dlat + 0.5*A_dlat)

  enddo
  enddo


  !!!!!!!! Create LocStream !!!!!!!!!!!
  ! Set number of points
  if (localPet .eq. 0) then  
      pntCount=2
  else if (localPet .eq. 1) then  
      pntCount=3
  else if (localPet .eq. 2) then  
      pntCount=1
  else if (localPet .eq. 3) then  
      pntCount=1
  endif


  ! Create LocStream
  locstream=ESMF_LocStreamCreate(localCount=pntCount, rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE
     

  ! Allocate Lon array
  allocate(Lon(pntCount))

  ! allocate Lat array
  allocate(Lat(pntCount))


 ! Fill in points
  if (localPet .eq. 0) then  
      Lon(1)=190.0
      Lat(1)=10.0
      Lon(2)=280.0
      Lat(2)=-15.0
  else if (localPet .eq. 1) then  
      Lon(1)=10.0
      Lat(1)=50.0
      Lon(2)=20.0
      Lat(2)=40.0
      Lon(3)=100.0
      Lat(3)=0.0
  else if (localPet .eq. 2) then  
      Lon(1)=355.0
      Lat(1)=-60.0
  else if (localPet .eq. 3) then  
      Lon(1)=170.0
      Lat(1)=-80.0
  endif


  ! Add key Lon
  call ESMF_LocStreamAddKey(locstream, keyName="ESMF:Lon", farray=Lon,  rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE

  ! Add key Lat
  call ESMF_LocStreamAddKey(locstream, keyName="ESMF:Lat", farray=Lat,  rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE

  ! Do locStream create from background mesh
  newLocstream=ESMF_LocStreamCreate(locstream, &
                 background=gridA, rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) then
     rc=ESMF_FAILURE
     return
  endif

  
  call ESMF_LocStreamDestroy(locstream,rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE

  ! deallocate array
  deallocate(Lon)
  deallocate(Lat)



  ! Test Newly Created LocStream
  ! Since the grid is setup to be rectilinear checking the min/max should be fine
 
  !! get coord 1
  call ESMF_GridGetCoord(gridA, localDE=0, staggerLoc=ESMF_STAGGERLOC_CORNER, coordDim=1, &
                            farrayPtr=farrayPtrLonC, rc=localrc)
     if (localrc /=ESMF_SUCCESS) then
        rc=ESMF_FAILURE
        return
     endif

  call ESMF_GridGetCoord(gridA, localDE=0, staggerLoc=ESMF_STAGGERLOC_CORNER, coordDim=2, &
                            computationalLBound=clbnd, computationalUBound=cubnd, farrayPtr=farrayPtrLatC, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
        rc=ESMF_FAILURE
        return
  endif

  !! Init min/max of DE
  de_minlon=360.0
  de_minlat=90.0
  de_maxlon=0.0
  de_maxlat=-90.0

  !! Adjust loop to cover min-max of cells not just nodes
  !! i.e. extend by 1 to cover other point of cell.
  if (localPet .lt. PetCount-1) cubnd(1) =cubnd(1)+1

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

     ! Set source coordinates
    ! Set source coordinates as 0 to 360
    tmpLonC = REAL(i1-1)*A_dlon

    ! Set source coordinates as -90 to 90
    tmpLatC = -90. + (REAL(i2-1)*A_dlat + 0.5*A_dlat)
 
     ! Min/max off coordinates
     if (tmpLonC < de_minlon) de_minlon=tmpLonC 
     if (tmpLonC > de_maxlon) de_maxlon=tmpLonC

     if (tmpLatC < de_minlat) de_minlat=tmpLatC
     if (tmpLatC > de_maxlat) de_maxlat=tmpLatC

  enddo
  enddo

 !  write(*,*) localPet," [",de_minlon,de_maxlon,"]"

  !!!!!!!!! Check locstream points vs Grid min max !!!!!!!!!!!!!!!!!
  call ESMF_LocStreamGetKey(newlocStream,keyName="ESMF:Lon", &
         farray=tstLon, &
         exclusiveLBound=el, exclusiveUBound=eu, rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE

  call ESMF_LocStreamGetKey(newLocStream,keyName="ESMF:Lat", &
         farray=tstLat, rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE


  ! Test points
  do i=el,eu
     if ((tstLon(i) <  de_minlon) .or. (tstLon(i) >  de_maxlon)) then
        write(*,*) tstLon(i),"not in [",de_minlon,de_maxlon,"]"
        correct=.false.
     endif

     if ((tstLat(i) <  de_minlat) .or. (tstLat(i) >  de_maxlat)) then
        write(*,*) tstLat(i),"not in [",de_minlat,de_maxlat,"]"
        correct=.false.
     endif
   enddo


  ! Get rid of the new locstream
  call ESMF_LocStreamDestroy(newLocstream,rc=localrc)
  if (localrc .ne. ESMF_SUCCESS) rc=ESMF_FAILURE
  
  ! Free the grids
  call ESMF_GridDestroy(gridA, rc=localrc)
  if (localrc /=ESMF_SUCCESS) then
      rc=ESMF_FAILURE
      return
   endif


  ! return answer based on correct flag
  if (correct) then
    rc=ESMF_SUCCESS
  else
    rc=ESMF_FAILURE
  endif

 end subroutine test_locstreambkgsph