subroutine ESMF_UtilCreateCSCoordsPar(npts, LonEdge,LatEdge, start, count, tile, &
     LonCenter, LatCenter, schmidtTransform, local_algorithm)
! !ARGUMENTS:
    integer,           intent(IN)     :: npts
    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
    logical, optional, intent(in) :: local_algorithm
 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       :: global_tile1(:,:,:)
  integer                       :: shapLon(2), shapLat(2)
  logical :: local_algorithm_
  local_algorithm_ = .false.
  if (present(local_algorithm)) local_algorithm_ = local_algorithm
  allocate(tile_local(count(1)+1,count(2)+1,ndims) )
  if (local_algorithm_) then
     call get_gnomonic_local_coords(grid_type, npts, start, tile_local(:,:,1), tile_local(:,:,2))
     call mirror_grid_local_new(tile_local, tile)
  else ! legacy algorithm using global sized tile
     allocate(global_tile1(npts+1,npts+1,ndims))
     call gnomonic_grids(grid_type, npts, global_tile1(:,:,1), global_tile1(:,:,2))
     call mirror_grid_local(tile_local, global_tile1, start, count, 2, tile)
     deallocate(global_tile1)
  end if
    ! mirror_grid assumes that the tile=1 is centered on equator and greenwich meridian Lon[-pi,pi]
!---------------------------------
! Shift the corner away from Japan for global tile #1
!---------------------------------
! This will result in the corner close to east coast of China
    ! fix the values in the local tile
    do j=1,count(2)+1
       do i=1,count(1)+1
           if (.not.present(schmidtTransform)) tile_local(i,j,1) = tile_local(i,j,1) - JAPAN_SHIFT
           if ( tile_local(i,j,1) < 0. )              &
                tile_local(i,j,1) = tile_local(i,j,1) + 2.*pi
           if (ABS(tile_local(i,j,1)) < 1.e-10) tile_local(i,j,1) = 0.0
           if (ABS(tile_local(i,j,2)) < 1.e-10) tile_local(i,j,2) = 0.0
       enddo
    enddo
    if (present(schmidtTransform)) then
     call direct_transform(schmidtTransform%stretch_factor,schmidtTransform%target_lon,&
          schmidtTransform%target_lat,tile_local(:,:,1),tile_local(:,:,2))
    end if
    
    if (present(LonEdge) .and. present(LatEdge)) then
       shapLon=shape(LonEdge)
       shapLat=shape(LatEdge)
       LonEdge=tile_local(1:shapLon(1),1:shapLon(2),1)
       LatEdge=tile_local(1:shapLat(1),1:shapLat(2),2)
    endif
    if (present(LonCenter) .and. present(LatCenter)) then
        do j=1, count(2)
           do i=1,count(1)
              call cell_center2(tile_local(i,j,  1), tile_local(i,j,    2),   &
                                tile_local(i+1,j,1), tile_local(i+1,j,  2),   &
                                tile_local(i,j+1,1), tile_local(i,j+1,  2),   &
                                tile_local(i+1,j+1,1), tile_local(i+1,j+1,2),   &
                                alocs)
              LonCenter(i,j) = alocs(1)
              LatCenter(i,j) = alocs(2)
           enddo
        enddo
     end if
     deallocate(tile_local)
     return
  end subroutine ESMF_UtilCreateCSCoordsPar