ESMF_UtilCreateCSCoords Subroutine

public subroutine ESMF_UtilCreateCSCoords(npts, LonEdge, LatEdge, start, count, tile, LonCenter, LatCenter, schmidtTransform)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: npts
real(kind=ESMF_KIND_R8), intent(inout), optional :: LonEdge(:,:)
real(kind=ESMF_KIND_R8), intent(inout), optional :: LatEdge(:,:)
integer, intent(in), optional :: start(:)
integer, intent(in), optional :: count(:)
integer, intent(in), optional :: tile
real(kind=ESMF_KIND_R8), intent(inout), optional :: LonCenter(:,:)
real(kind=ESMF_KIND_R8), intent(inout), optional :: LatCenter(:,:)
type(ESMF_CubedSphereTransform_Args), intent(in), optional :: schmidtTransform

Source Code

subroutine ESMF_UtilCreateCSCoords(npts, LonEdge,LatEdge, start, count, tile, &
  LonCenter, LatCenter, schmidtTransform)

! !ARGUMENTS:
    integer,           intent(IN)     :: npts
!    integer,           intent(in)     :: petNo
    real(ESMF_KIND_R8), optional, intent(inout) :: LonEdge(:,:)
    real(ESMF_KIND_R8), optional, intent(inout) :: LatEdge(:,:)
    integer, optional, intent(in)     :: start(:)
    integer, optional, intent(in)     :: count(:)
    integer, optional, intent(in)     :: tile
    real(ESMF_KIND_R8), optional, intent(inout) :: LonCenter(:,:)
    real(ESMF_KIND_R8), optional, intent(inout) :: LatCenter(:,:)
    type(ESMF_CubedSphereTransform_Args), optional, intent(in)    :: schmidtTransform

! ErrLog variables
!-----------------

 integer                      :: STATUS

! Local variables
!-----------------
  integer, parameter            :: grid_type = 0
  integer                       :: ntiles=6
  integer                       :: ndims=2
  integer                       :: I, J, N
  integer                       :: IG, JG
  real(ESMF_KIND_R8)                          :: alocs(2)
  real(ESMF_KIND_R8), allocatable             :: tile1(:,:,:)
  integer                       :: rc
  real(ESMF_KIND_R8), allocatable             :: tile_local(:,:,:)
  real(ESMF_KIND_R8), allocatable, save       :: grid_global(:,:,:,:)

  if (allocated(grid_global)) then 
     if (size(grid_global,1) /= npts+1) then
        deallocate(grid_global)
     endif
  endif   
  if (.not. allocated(grid_global)) then
     allocate( grid_global(npts+1,npts+1,ndims,ntiles) )
     call gnomonic_grids(grid_type, npts, grid_global(:,:,1,1), grid_global(:,:,2,1))
     ! mirror_grid assumes that the tile=1 is centered on equator and greenwich meridian Lon[-pi,pi]
     call mirror_grid(grid_global, 0, npts+1, npts+1, 2, 6)
     do n=1,ntiles
        do j=1,npts+1
           do i=1,npts+1
!---------------------------------
! Shift the corner away from Japan
!---------------------------------
! This will result in the corner close to east coast of China
             if (.not.present(schmidtTransform)) grid_global(i,j,1,n) = grid_global(i,j,1,n) - pi/18.
             if ( grid_global(i,j,1,n) < 0. )              &
                grid_global(i,j,1,n) = grid_global(i,j,1,n) + 2.*pi
             if (ABS(grid_global(i,j,1,n)) < 1.e-10) grid_global(i,j,1,n) = 0.0
             if (ABS(grid_global(i,j,2,n)) < 1.e-10) grid_global(i,j,2,n) = 0.0
          enddo
       enddo
    enddo
    !---------------------------------
    ! Clean Up Corners
    !---------------------------------
    grid_global(  1,1:npts+1,:,2)=grid_global(npts+1,1:npts+1,:,1)
    grid_global(  1,1:npts+1,:,3)=grid_global(npts+1:1:-1,npts+1,:,1)
    grid_global(1:npts+1,npts+1,:,5)=grid_global(1,npts+1:1:-1,:,1)
    grid_global(1:npts+1,npts+1,:,6)=grid_global(1:npts+1,1,:,1)
    grid_global(1:npts+1,  1,:,3)=grid_global(1:npts+1,npts+1,:,2)
    grid_global(1:npts+1,  1,:,4)=grid_global(npts+1,npts+1:1:-1,:,2)
    grid_global(npts+1,1:npts+1,:,6)=grid_global(npts+1:1:-1,1,:,2)
    grid_global(  1,1:npts+1,:,4)=grid_global(npts+1,1:npts+1,:,3)
    grid_global(  1,1:npts+1,:,5)=grid_global(npts+1:1:-1,npts+1,:,3)
    !!!grid_global(npts+1,1:npts+1,:,3)=grid_global(1,1:npts+1,:,4)
    grid_global(1:npts+1,  1,:,5)=grid_global(1:npts+1,npts+1,:,4)
    grid_global(1:npts+1,  1,:,6)=grid_global(npts+1,npts+1:1:-1,:,4)
    grid_global(  1,1:npts+1,:,6)=grid_global(npts+1,1:npts+1,:,5)
  endif

#if 0
  ! check if they are identical
  if (PetNo == 0) then
  do n=1,npts+1
    do j=1,2
    if (grid_global(1,n,j,2) /= grid_global(npts+1,n,j,1)) then
       print *, "mismatch tile 1&2", n, grid_global(1,n,j,2),grid_global(npts+1,n,j,1)
    endif
    if (grid_global(1,n,j,3) /= grid_global(npts+2-n,npts+1,j,1)) then
       print *, "mismatch tile 1&3", n, grid_global(1,n,j,3),grid_global(npts+2-n,npts+1,j,1)
    endif
    if (grid_global(n,npts+1,j,5) /= grid_global(1,npts+2-n,j,1)) then
       print *, "mismatch tile 1&5", n, grid_global(n,npts+1,j,5),grid_global(1,npts+2-n,j,1)
    endif
    if (grid_global(n,npts+1,j,6) /= grid_global(n,1,j,1)) then
       print *, "mismatch tile 1&6", n, grid_global(n,npts+1,j,6),grid_global(n,1,j,1)
    endif
    if (grid_global(n,1,j,3) /= grid_global(n,npts+1,j,2)) then
       print *, "mismatch tile 2&3", n, grid_global(n,1,j,3),grid_global(n,npts+1,j,2)
    endif
    if (grid_global(n,1,j,4) /= grid_global(npts+1,npts+2-n,j,2)) then
       print *, "mismatch tile 2&4", n, grid_global(n,1,j,4),grid_global(npts+1,npts+2-n,j,2)
    endif
    if (grid_global(npts+1,n,j,6) /= grid_global(npts+2-n,1,j,2)) then
       print *, "mismatch tile 2&6", n, grid_global(npts+1,n,j,6),grid_global(npts+2-n,1,j,2)
    endif
    if (grid_global(1,n,j,4) /= grid_global(npts+1,n,j,3)) then
       print *, "mismatch tile 3&4", n, grid_global(i,n,j,4),grid_global(npts+1,n,j,3)
    endif
    if (grid_global(1,n,j,5) /= grid_global(npts+2-n,npts+1,j,3)) then
       print *, "mismatch tile 5&3", n, grid_global(1,n,j,5),grid_global(npts+2-n,npts+1,j,3)
    endif
    if (grid_global(n,1,j,5) /= grid_global(n,npts+1,j,4)) then
       print *, "mismatch tile 5&4", n, grid_global(n,1,j,5),grid_global(n,npts+1,j,4)
    endif
    if (grid_global(n,1,j,6) /= grid_global(npts+1,npts+2-n,j,4)) then
       print *, "mismatch tile 6&4", n, grid_global(n,1,j,6),grid_global(npts+1,npts+2-n,j,4)
    endif
    if (grid_global(1,n,j,6) /= grid_global(npts+1,n,j,5)) then
       print *, "mismatch tile 1&2", n, grid_global(1,n,j,6),grid_global(npts+1,n,j,5)
    endif
  enddo
  enddo
  endif
#endif

  if (present(schmidtTransform)) then
     do n=1,ntiles
       call direct_transform(schmidtTransform%stretch_factor,schmidtTransform%target_lon,&
          schmidtTransform%target_lat,grid_global(:,:,1,n),grid_global(:,:,2,n))
     enddo
  end if

  if (present(LonEdge) .and. present(LatEdge)) then
    do n=1,ntiles
       do j=1,npts+1
          do i=1,npts+1
            jg=(n-1)*(npts+1)+j
            LonEdge(i,jg) = grid_global(i,j,1,n)
            LatEdge(i,jg) = grid_global(i,j,2,n)
          enddo
       enddo
    enddo
  endif
 
  if (present(LonCenter) .and. present(LatCenter)) then
    if (present(start) .and. present(count) .and. present(tile)) then
        n=tile
        do j=start(2), start(2)+count(2)-1
           do i=start(1),start(1)+count(1)-1
              call cell_center2(grid_global(i,j,  1,n), grid_global(i,j,  2,n),   &
                                grid_global(i+1,j,1,n), grid_global(i+1,j,2,n),   &
                                grid_global(i,j+1,1,n), grid_global(i,j+1,2,n),   &
                                grid_global(i+1,j+1,1,n), grid_global(i+1,j+1,2,n),   &
                                alocs)
              LonCenter(i-start(1)+1,j-start(2)+1) = alocs(1)
              LatCenter(i-start(1)+1,j-start(2)+1) = alocs(2)
           enddo
        enddo
     else
        do n=1,ntiles
           do j=1,npts
              do i=1,npts
              call cell_center2(grid_global(i,j,  1,n), grid_global(i,j,  2,n),   &
                                grid_global(i+1,j,1,n), grid_global(i+1,j,2,n),   &
                                grid_global(i,j+1,1,n), grid_global(i,j+1,2,n),   &
                                grid_global(i+1,j+1,1,n), grid_global(i+1,j+1,2,n),   &
                                alocs)
                  jg = (n-1)*npts + j
                  LonCenter(i,jg) = alocs(1)
                  LatCenter(i,jg) = alocs(2)
              enddo
           enddo
        enddo
     end if
  end if
  deallocate( grid_global )

  return

end subroutine ESMF_UtilCreateCSCoords