PV2S Subroutine

public subroutine PV2S(pv, theta, phi, r, td, pd, rd)

Convert position/velocity from Cartesian to spherical coordinates.

Status: vector/matrix support routine.

Notes

  1. If the position part of PV is null, THETA, PHI, TD and PD are indeterminate. This is handled by extrapolating the position through unit time by using the velocity part of PV. This moves the origin without changing the direction of the velocity component. If the position and velocity components of PV are both null, zeroes are returned for all six results.

  2. If the position is a pole, THETA, TD and PD are indeterminate. In such cases zeroes are returned for all three.

History

  • IAU SOFA revision: 2008 May 10

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in), dimension(3,2):: pv

pv-vector

real(kind=wp), intent(out) :: theta

longitude angle (radians)

real(kind=wp), intent(out) :: phi

latitude angle (radians)

real(kind=wp), intent(out) :: r

radial distance

real(kind=wp), intent(out) :: td

rate of change of THETA

real(kind=wp), intent(out) :: pd

rate of change of PHI

real(kind=wp), intent(out) :: rd

rate of change of R


Called by

proc~~pv2s~~CalledByGraph proc~pv2s PV2S proc~fk425 FK425 proc~fk425->proc~pv2s proc~hfk5z HFK5Z proc~hfk5z->proc~pv2s proc~pvstar PVSTAR proc~pvstar->proc~pv2s proc~fk524 FK524 proc~fk524->proc~pv2s proc~fk52h FK52H proc~fk52h->proc~pvstar proc~starpm STARPM proc~starpm->proc~pvstar proc~h2fk5 H2FK5 proc~h2fk5->proc~pvstar proc~fk54z FK54Z proc~fk54z->proc~fk524 proc~pmsafe PMSAFE proc~pmsafe->proc~starpm

Contents

Source Code


Source Code

    subroutine PV2S ( pv, theta, phi, r, td, pd, rd )

    implicit none

    real(wp),dimension(3,2),intent(in) :: pv !! pv-vector
    real(wp),intent(out) :: theta !! longitude angle (radians)
    real(wp),intent(out) :: phi !! latitude angle (radians)
    real(wp),intent(out) :: r !! radial distance
    real(wp),intent(out) :: td !! rate of change of THETA
    real(wp),intent(out) :: pd !! rate of change of PHI
    real(wp),intent(out) :: rd !! rate of change of R

    real(wp) :: x, y, z, xd, yd, zd, rxy2, rxy, r2, &
                rtrue, rw, xyp

    !  Components of position/velocity vector.
    x =  pv(1,1)
    y =  pv(2,1)
    z =  pv(3,1)
    xd = pv(1,2)
    yd = pv(2,2)
    zd = pv(3,2)

    !  Component of R in XY plane squared.
    rxy2 = x*x + y*y

    !  Modulus squared.
    r2 = rxy2 + z*z

    !  Modulus.
    rtrue = sqrt(r2)

    !  If null vector, move the origin along the direction of movement.
    rw = rtrue
    if ( rtrue == 0.0_wp ) then
       x = xd
       y = yd
       z = zd
       rxy2 = x*x + y*y
       r2 = rxy2 + z*z
       rw = sqrt(r2)
    end if

    !  Position and velocity in spherical coordinates.
    rxy = sqrt(rxy2)
    xyp = x*xd + y*yd
    if ( rxy2 /= 0.0_wp ) then
       theta = atan2(y,x)
       phi = atan2(z,rxy)
       td = ( x*yd - y*xd ) / rxy2
       pd = ( zd*rxy2 - z*xyp ) / ( r2*rxy )
    else
       theta = 0.0_wp
       if ( z/=0.0_wp ) then
          phi = atan2(z,rxy)
       else
          phi = 0.0_wp
       end if
       td = 0.0_wp
       pd = 0.0_wp
    end if
    r = rtrue
    if ( rw/=0.0_wp ) then
       rd = ( xyp + z*zd ) / rw
    else
       rd = 0.0_wp
    end if

    end subroutine PV2S