LDN Subroutine

public subroutine LDN(n, b, ob, sc, sn)

For a star, apply light deflection by multiple solar-system bodies, as part of transforming coordinate direction into natural direction.

Status: support routine.

  1. The array B contains N entries, one for each body to be considered. If N = 0, no gravitational light deflection will be applied, not even for the Sun.

  2. The array B should include an entry for the Sun as well as for any planet or other body to be taken into account. The entries should be in the order in which the light passes the body.

  3. In the entry in the B array for body I, the mass parameter B(1,I) can, as required, be adjusted in order to allow for such effects as quadrupole field.

  4. The deflection limiter parameter B(2,I) is phi^2/2, where phi is the angular separation (in radians) between star and body at which limiting is applied. As phi shrinks below the chosen threshold, the deflection is artificially reduced, reaching zero for phi = 0. Example values suitable for a terrestrial observer, together with masses, are as follows:

    body I     B(1,I)         B(2,I)
    
    Sun        1D0            6D-6
    Jupiter    0.00095435D0   3D-9
    Saturn     0.00028574D0   3D-10
    
  5. For cases where the starlight passes the body before reaching the observer, the body is placed back along its barycentric track by the light time from that point to the observer. For cases where the body is "behind" the observer no such shift is applied. If a different treatment is preferred, the user has the option of instead using the LD routine. Similarly, LD can be used for cases where the source is nearby, not a star.

  6. The returned vector SN is not normalized, but the consequential departure from unit magnitude is always negligible.

  7. For efficiency, validation is omitted. The supplied masses must be greater than zero, the position and velocity vectors must be right, and the deflection limiter greater than zero.

Reference

  • Urban, S. & Seidelmann, P. K. (eds), Explanatory Supplement to the Astronomical Almanac, 3rd ed., University Science Books (2013), Section 7.2.4.

History

  • IAU SOFA revision: 2017 March 16

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

number of bodies (Note 1)

real(kind=wp), intent(in), dimension(8,n):: b

data for each of the N bodies (Notes 1,2): (1,I) mass of the body (solar masses, Note 3) (2,I) deflection limiter (Note 4) (3-5,I) barycentric position of the body (au) (6-8,I) barycentric velocity of the body (au/day)

real(kind=wp), intent(in), dimension(3):: ob

barycentric position of the observer (au)

real(kind=wp), intent(in), dimension(3):: sc

observer to star coordinate direction (unit vector)

real(kind=wp), intent(out), dimension(3):: sn

observer to deflected star (unit vector)


Calls

proc~~ldn~~CallsGraph proc~ldn LDN proc~pmp PMP proc~ldn->proc~pmp proc~pdp PDP proc~ldn->proc~pdp proc~ld LD proc~ldn->proc~ld proc~ppsp PPSP proc~ldn->proc~ppsp proc~ld->proc~pdp proc~pxp PXP proc~ld->proc~pxp

Called by

proc~~ldn~~CalledByGraph proc~ldn LDN proc~atciqn ATCIQN proc~atciqn->proc~ldn proc~aticqn ATICQN proc~aticqn->proc~ldn

Contents

Source Code

LDN

Source Code

    subroutine LDN ( n, b, ob, sc, sn )

    implicit none

    integer,intent(in) :: n !! number of bodies (Note 1)
    real(wp),dimension(8,n),intent(in) :: b !! data for each of the N bodies (Notes 1,2):
                                            !! (1,I)     mass of the body (solar masses, Note 3)
                                            !! (2,I)     deflection limiter (Note 4)
                                            !! (3-5,I)   barycentric position of the body (au)
                                            !! (6-8,I)   barycentric velocity of the body (au/day)
    real(wp),dimension(3),intent(in) :: ob !! barycentric position of the observer (au)
    real(wp),dimension(3),intent(in) :: sc !! observer to star coordinate direction (unit vector)
    real(wp),dimension(3),intent(out) :: sn !! observer to deflected star (unit vector)

    !  Astronomical unit (m, IAU 2012)
    real(wp),parameter :: aum = 149597870.7d3

    !  Light time for 1 au (day)
    real(wp),parameter :: cr = aum/cmps/d2s

    integer :: i
    real(wp) :: s(3), v(3), d , dt, ev(3), em, e(3)

    !  Star direction prior to deflection.
    call CP ( sc, s )

    !  Body by body.
    do i=1,n

       !  Body to observer vector at epoch of observation (au).
       call PMP ( ob, b(3,i), v )

       !  Minus the time since the light passed the body (days).
       call PDP ( s, v, d )
       dt = d * cr

       !  Neutralize if the star is "behind" the observer.
       dt = min ( dt, 0.0_wp )

       !  Backtrack the body to the time the light was passing the body.
       call PPSP ( v, -dt, b(6,i), ev )

       !  Separate the body to observer vector into magnitude and direction.
       call PN ( ev, em, e )

       !  Apply light deflection for this body.
       call LD ( b(1,i), s, s, e, em, b(2,i), v )

       !  Update the star direction.
       call CP ( v, s )

       !  Next body.
    end do

    !  Return the deflected star direction.
    call CP ( s, sn )

    end subroutine LDN