propagate Subroutine

public pure subroutine propagate(mu, rv0, dt, rvf)

Basic two-body propagator using the Gooding universal element routines.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: mu

grav. parameter [km^3/s^2]

real(kind=wp), intent(in), dimension(6) :: rv0

initial state [km, km/s]

real(kind=wp), intent(in) :: dt

time step [sec]

real(kind=wp), intent(out), dimension(6) :: rvf

final state [km, km/s]


Calls

proc~~propagate~~CallsGraph proc~propagate gooding_module::propagate proc~els3pv gooding_module::els3pv proc~propagate->proc~els3pv proc~pv3els gooding_module::pv3els proc~propagate->proc~pv3els proc~els2pv gooding_module::els2pv proc~els3pv->proc~els2pv proc~pv2els gooding_module::pv2els proc~pv3els->proc~pv2els proc~dcbsol gooding_module::dcbsol proc~els2pv->proc~dcbsol proc~ekepl gooding_module::ekepl proc~els2pv->proc~ekepl proc~shkepl gooding_module::shkepl proc~els2pv->proc~shkepl proc~emkep gooding_module::emkep proc~pv2els->proc~emkep proc~shmkep gooding_module::shmkep proc~pv2els->proc~shmkep proc~dcubrt gooding_module::dcubrt proc~dcbsol->proc~dcubrt proc~ekepl->proc~dcbsol proc~ekepl->proc~emkep proc~shkepl->proc~shmkep proc~shkepl->proc~dcubrt

Source Code

    pure subroutine propagate(mu, rv0, dt, rvf)

    implicit none

    real(wp),intent(in)               :: mu    !! grav. parameter [km^3/s^2]
    real(wp),dimension(6),intent(in)  :: rv0   !! initial state [km, km/s]
    real(wp),intent(in)               :: dt    !! time step [sec]
    real(wp),dimension(6),intent(out) :: rvf   !! final state [km, km/s]

    real(wp),dimension(6) :: e

    !convert to elements, increment time,
    ! then convert back to cartesian:

    call pv3els(mu, rv0, e)
    e(6) = e(6) + dt
    call els3pv(mu, e, rvf)

    end subroutine propagate