transform Subroutine

private subroutine transform(from, rv, to, et, eph, rv_out, status_ok)

Transform a Cartesian state from one reference frame to another at a specified epoch. The from and to reference_frames may each be defined at a different epoch. The et ephemeris time is the time the transformation is to be done, and accounts for the motion of the two frame centers from from%et and to%et to et.

Type Bound

reference_frame

Arguments

Type IntentOptional Attributes Name
class(reference_frame), intent(inout) :: from
real(kind=wp), intent(in), dimension(6) :: rv
class(reference_frame), intent(inout) :: to
real(kind=wp), intent(in) :: et

the time of the transformation [sec]

class(ephemeris_class), intent(inout) :: eph

for ephemeris computations (assumed to have already been initialized)

real(kind=wp), intent(out), dimension(6) :: rv_out
logical, intent(out) :: status_ok

Calls

proc~~transform~~CallsGraph proc~transform transformation_module::reference_frame%transform get_c_cdot get_c_cdot proc~transform->get_c_cdot proc~rvcto_rvcfrom_icrf transformation_module::rvcto_rvcfrom_icrf proc~transform->proc~rvcto_rvcfrom_icrf get_rv get_rv proc~rvcto_rvcfrom_icrf->get_rv proc~from_primary_to_center transformation_module::two_body_rotating_frame%from_primary_to_center proc~rvcto_rvcfrom_icrf->proc~from_primary_to_center proc~from_primary_to_center->get_rv

Called by

proc~~transform~~CalledByGraph proc~transform transformation_module::reference_frame%transform proc~transformation_module_test transformation_module::transformation_module_test proc~transformation_module_test->proc~transform

Source Code

    subroutine transform(from,rv,to,et,eph,rv_out,status_ok)

    implicit none

    class(reference_frame),intent(inout) :: from
    real(wp),dimension(6),intent(in)     :: rv
    class(reference_frame),intent(inout) :: to
    real(wp),intent(in)                  :: et      !! the time of the transformation [sec]
    class(ephemeris_class),intent(inout) :: eph     !! for ephemeris computations (assumed to have already been initialized)
    real(wp),dimension(6),intent(out)    :: rv_out
    logical,intent(out)                  :: status_ok

    ! local variables
    real(wp), dimension(3,3) :: c, cdot
    real(wp), dimension(3,3) :: rot1, rotd1
    real(wp), dimension(3,3) :: rot2, rotd2
    real(wp), dimension(3)   :: rc21_out, vc21_out
    real(wp), dimension(3)   :: rc21_icrf, vc21_icrf

    ! rotation matrix: input -> inertial
    call from%get_c_cdot(eph=eph,to_icrf=.true.,c=rot1,cdot=rotd1,status_ok=status_ok)
    if (.not. status_ok) then
        write(error_unit,'(A)') 'Error in transform: '//&
                        'Could not compute rotation matrix from FROM frame to inertial.'
        rv_out = zero
        return
    end if
    ! rotation matrix: inertial -> output
    call to%get_c_cdot(eph=eph,to_icrf=.false.,c=rot2,cdot=rotd2,status_ok=status_ok)
    if (.not. status_ok) then
        write(error_unit,'(A)') 'Error in transform: '//&
                        'Could not compute rotation matrix from inertial to TO frame.'
        rv_out = zero
        return
    end if

    ! rotation matrix: input -> output
    c    = matmul(rot2, rot1)
    cdot = matmul(rotd2, rot1) + matmul(rot2, rotd1)

    ! get the state of the `from` frame center w.r.t. the `to` frame center,
    ! at the transformation time:
    call rvcto_rvcfrom_icrf(from, to, eph, et, &
                rc21_icrf, vc21_icrf, status_ok)

    if (status_ok) then

        ! to->from frame center state: inertial -> output frame
        rc21_out = matmul(rot2, rc21_icrf)
        vc21_out = matmul(rotd2, rc21_icrf) + matmul(rot2, vc21_icrf)

        ! rotation + translation:
        rv_out(1:3) = matmul(c, rv(1:3)) + rc21_out
        rv_out(4:6) = matmul(c, rv(4:6)) + matmul(cdot,rv(1:3)) + vc21_out

    else !error
        rv_out = zero
    end if

    end subroutine transform