Orbits_Swath Subroutine

public subroutine Orbits_Swath(slons, slats, Sat_name, nymd, nhms, deltat, SwathWidth, rc, wrap)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), pointer :: slons(:,:)
real(kind=dp), pointer :: slats(:,:)
character(len=*), intent(in) :: Sat_name
integer, intent(in) :: nymd(2)
integer, intent(in) :: nhms(2)
integer, intent(in) :: deltat
real(kind=dp), intent(in) :: SwathWidth(:)
integer, intent(out) :: rc
logical, intent(in), optional :: wrap

Calls

proc~~orbits_swath~~CallsGraph proc~orbits_swath Orbits_Swath proc~orbits_track Orbits_Track proc~orbits_swath->proc~orbits_track interface~orbits_track0 Orbits_Track0 proc~orbits_track->interface~orbits_track0

Source Code

    SUBROUTINE Orbits_Swath (slons, slats, Sat_name, nymd, nhms, deltat, &
                           SwathWidth, rc, wrap)
! SUBROUTINE Orbits_Swath (slons, slats, Sat_name, nymd, nhms, deltat, SwathWidth, wrapon, rc)

      IMPLICIT NONE
! !INPUT PARAMETERS:
      real(dp),    pointer         :: slons(:,:)      !Track longitude [degree]
      real(dp),    pointer         :: slats(:,:)      !Track latitude  [degree]
      character(len=*),intent(in)  :: Sat_name        !Satellite name Aqua etc ...
      integer,  intent(in)         :: nymd(2)         !Beginning/ending date: YYYYMMDD
      integer,  intent(in)         :: nhms(2)         !Beginning/ending time: HHMMSS
      REAL(dp), INTENT(IN)         :: SwathWidth(:)   !Right, left swath [km]
      integer, intent(in)          :: deltat          ! Time step [secs]
      LOGICAL, intent(in), optional:: wrap       ! if true then wrap to -180 to 180

      !LOGICAL                      :: wrapon       ! if true then wrap to -180 to 180


! !OUTPUT PARAMETERS:
      integer, intent(out)        :: rc

!     Other parameters!

      real(dp),    pointer  :: lons(:)      !Ground Track longitude [degree]
      real(dp),    pointer  :: lats(:)      !Ground Track latitude  [degree]
      REAL(dp),    DIMENSION(:,:,:), ALLOCATABLE  :: latlonshift
      real(dp)              :: temp1
      integer, parameter    :: iii = 3      !this is used for the number of swaths
                                            !here is 3 indicating mid left right
      integer               :: Sat          ! Satellite number
      integer               :: ierr1

      logical               :: wrapon

      if ( present(wrap) ) then
           wrapon = wrap
      else                 
           wrapon = .true.
      end if

      rc = 0 ! Initiate error code to 0
      CALL Orbits_Track(lons, lats, Sat_name, nymd, nhms, deltat, rc) !rc will return 3
      if ( rc /= 0 ) return
      CALL satname2int(Sat_name, Sat, rc)  ! name to sat number rc will return 99
      if ( rc /= 0 ) return

!     Check if the Swath width is given correctly
      if ( SwathWidth(1)<0.AND.SwathWidth(2) < 0 ) then
         rc = 90
         return
      else if (SwathWidth(1)<0.AND.SwathWidth(2) >= 0) then
         temp1 = SwathWidth(2) + SwathWidth(1)
         if (temp1<=0) then
            rc = 90
            return
         end if
      else if (SwathWidth(2)<0.AND.SwathWidth(1) >= 0) then
         temp1 = SwathWidth(1) + SwathWidth(2)
         if (temp1<=0) then
            rc = 90
            return
         end if
      end if

!     2- Find sweep points
      if (allocated(latlonshift)) deallocate(latlonshift) ! this is main
      ALLOCATE(latlonshift(iii, SIZE(lats) ,2), stat = ierr1 ) ! either is ok size of lat or lon
      if ( ierr1 /= 0 ) then
          rc = 3 !
          return
      endif
      latlonshift = 0
      !CALL find_sweeppoints(iii,lat_l(1:las),lon_l(1:los), real(swath(1),dp),  &
      !                  real(swath(2),dp), latlonshift)
      CALL find_sweeppoints(iii, lats, lons, -SwathWidth(1), SwathWidth(2), wrapon, latlonshift)
      ALLOCATE(slats(iii, SIZE(lats)), stat = ierr1 ) ! either is ok size of lat or lon
      ALLOCATE(slons(iii, SIZE(lons)), stat = ierr1 ) ! either is ok size of lat or lon

      slats(1:iii, 1:SIZE(lats)) = latlonshift(1:iii, 1:SIZE(lats) ,1)
      slons(1:iii, 1:SIZE(lons)) = latlonshift(1:iii, 1:SIZE(lons) ,2)
      slats(2,1:SIZE(lats))  = lats(1:SIZE(lats))  ! slats(1:3, :) 1- right 2- original ground track 3-left
      slons(2, 1:SIZE(lons)) = lons(1:SIZE(lons))

      DEALLOCATE(latlonshift, stat = ierr1)
      if (ierr1 /= 0) then
        rc = 3 !
        return
      endif       

  END SUBROUTINE Orbits_Swath