get_c_cdot_two_body_rotating_pulsating Subroutine

private subroutine get_c_cdot_two_body_rotating_pulsating(me, eph, to_icrf, c, cdot, status_ok)

rotation matrix for ROTATING_PULSATING <-> ICRF

Type Bound

two_body_rotating_pulsating_frame

Arguments

Type IntentOptional Attributes Name
class(two_body_rotating_pulsating_frame), intent(inout) :: me
class(ephemeris_class), intent(inout) :: eph

for ephemeris computations (assumed to have already been initialized)

logical, intent(in) :: to_icrf
real(kind=wp), intent(out), dimension(3,3) :: c
real(kind=wp), intent(out), optional, dimension(3,3) :: cdot
logical, intent(out) :: status_ok

Calls

proc~~get_c_cdot_two_body_rotating_pulsating~~CallsGraph proc~get_c_cdot_two_body_rotating_pulsating transformation_module::two_body_rotating_pulsating_frame%get_c_cdot_two_body_rotating_pulsating proc~get_c_cdot_two_body_rotating transformation_module::two_body_rotating_frame%get_c_cdot_two_body_rotating proc~get_c_cdot_two_body_rotating_pulsating->proc~get_c_cdot_two_body_rotating get_rv get_rv proc~get_c_cdot_two_body_rotating->get_rv proc~cross vector_module::cross proc~get_c_cdot_two_body_rotating->proc~cross proc~cross_matrix vector_module::cross_matrix proc~get_c_cdot_two_body_rotating->proc~cross_matrix proc~unit vector_module::unit proc~get_c_cdot_two_body_rotating->proc~unit

Source Code

    subroutine get_c_cdot_two_body_rotating_pulsating(me,eph,to_icrf,c,cdot,status_ok)

    !! rotation matrix for ROTATING_PULSATING <-> ICRF

    implicit none

    class(two_body_rotating_pulsating_frame),intent(inout) :: me
    class(ephemeris_class),intent(inout)                   :: eph     !! for ephemeris computations (assumed to have already been initialized)
    logical,intent(in)                                     :: to_icrf
    real(wp),dimension(3,3),intent(out)                    :: c
    real(wp),dimension(3,3),intent(out),optional           :: cdot
    logical,intent(out)                                    :: status_ok

    ! local variables
    real(wp),dimension(3,3) :: cr, crdot
    real(wp),dimension(3) :: r12,v12
    real(wp) :: r12mag,factor
    logical :: need_cdot

    need_cdot = present(cdot)

    ! rotating frame transformation matrices:
    if (need_cdot) then
        call me%two_body_rotating_frame%get_c_cdot(eph,to_icrf,cr,crdot,status_ok=status_ok)
    else
        call me%two_body_rotating_frame%get_c_cdot(eph,to_icrf,cr,status_ok=status_ok)
    end if

    if (status_ok) then

        r12    = me%rv12(1:3) ! was computed in get_c_cdot_two_body_rotating
        v12    = me%rv12(4:6)
        r12mag = norm2(r12)

        if (to_icrf) then
            factor = r12mag/me%scale
            c = factor*cr
            if (need_cdot) cdot = factor*(dot_product(v12,r12)*cr/(r12mag**2)+crdot)
        else
            factor = me%scale/r12mag
            c = factor*cr
            if (need_cdot) cdot = factor*(-dot_product(v12,r12)*cr/(r12mag**2)+crdot)
        end if

    else  !error
        write(error_unit,'(A)') 'Error in get_c_cdot_two_body_rotating_pulsating: '//&
                                'Could not compute rotation matrix.'
        c = zero
        if (need_cdot) cdot = zero
    end if

    end subroutine get_c_cdot_two_body_rotating_pulsating