compute_halo_monodromy_matrix Subroutine

public subroutine compute_halo_monodromy_matrix(mu, rv, period, phi)

Uses

  • proc~~compute_halo_monodromy_matrix~~UsesGraph proc~compute_halo_monodromy_matrix halo_orbit_module::compute_halo_monodromy_matrix module~rk_module rk_module proc~compute_halo_monodromy_matrix->module~rk_module module~kind_module kind_module module~rk_module->module~kind_module module~numbers_module numbers_module module~rk_module->module~numbers_module iso_fortran_env iso_fortran_env module~kind_module->iso_fortran_env module~numbers_module->module~kind_module

Compute the halo orbit monodromy matrix (which is the state transition matrix propagated for one period) The input should be the result from the halo_to_rv_diffcorr routine.

Arguments

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

CRTBP parameter

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

halo orbit state vector (normalized)

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

halo orbit period (normalized)

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

monodromy matrix


Calls

proc~~compute_halo_monodromy_matrix~~CallsGraph proc~compute_halo_monodromy_matrix halo_orbit_module::compute_halo_monodromy_matrix proc~initialize~2 rk_module::rk_class%initialize proc~compute_halo_monodromy_matrix->proc~initialize~2 proc~integrate~2 rk_module::rk_class%integrate proc~compute_halo_monodromy_matrix->proc~integrate~2 proc~destroy~2 rk_module::rk_class%destroy proc~initialize~2->proc~destroy~2 step step proc~integrate~2->step

Source Code

    subroutine compute_halo_monodromy_matrix(mu,rv,period,phi)

    use rk_module

    implicit none

    real(wp),intent(in)                 :: mu     !! CRTBP parameter
    real(wp),dimension(6),intent(in)    :: rv     !! halo orbit state vector
                                                  !! (normalized)
    real(wp),intent(in)                 :: period !! halo orbit period
                                                  !! (normalized)
    real(wp),dimension(6,6),intent(out) :: phi    !! monodromy matrix

    real(wp),parameter :: t0 = zero  !! initial time (normalized)
                                     !! (epoch doesn't matter for cr3bp)
    integer,parameter  :: n_steps_per_rev = 100   !! number of integration steps
                                                  !! per orbit rev

    real(wp),dimension(42)  :: x0     !! initial normalized state and STM
    real(wp),dimension(42)  :: xf     !! final normalized state and STM
    real(wp),dimension(6,6) :: phi0   !! initial STM
    integer                 :: i      !! counter
    type(rk8_10_class)      :: prop   !! integrator
    real(wp)                :: dt     !! integration time step (normalized)

    ! initial state:
    x0(1:6) = rv

    ! initial stm is the identity matrix:
    phi0 = zero
    do i = 1, 6
       phi0(i,i) = one
    end do
    x0(7:42) = pack (phi0, mask=.true.)

    ! for now, use a fixed time step:
    ! (same as was used in [[halo_to_rv_diffcorr]])
    dt = period / real(n_steps_per_rev,wp)

    ! initialize the integrator:
    call prop%initialize(42,func)

    ! propagate for one period:
    call prop%integrate(t0,x0,dt,period,xf)

    ! extract the STM:
    phi = reshape(xf(7:42), shape=[6,6])

    contains
!*******************************************************************************

    !***************************************************************************
        subroutine func(me,t,x,xdot)
        !! CRTBP derivative function (with STM)
        implicit none
        class(rk_class),intent(inout)        :: me
        real(wp),intent(in)                  :: t
        real(wp),dimension(me%n),intent(in)  :: x
        real(wp),dimension(me%n),intent(out) :: xdot

        call crtbp_derivs_with_stm(mu,x,xdot)

        end subroutine func
    !**************************************************************************

    end subroutine compute_halo_monodromy_matrix