segment Derived Type

type, public, extends(ddeabm_with_event_class) :: segment

the main class for integrating a low-lunar orbit.


Inherits

type~~segment~~InheritsGraph type~segment segment jpl_ephemeris jpl_ephemeris type~segment->jpl_ephemeris eph geopotential_model_pines geopotential_model_pines type~segment->geopotential_model_pines grav ddeabm_with_event_class ddeabm_with_event_class type~segment->ddeabm_with_event_class

Contents

Source Code


Components

TypeVisibilityAttributesNameInitial
integer, public :: event =0

event function to use:

1 = integrate until minimum altitude 2 = integrate to apoapsis

real(kind=wp), public :: et_ref

reference ephemeris time (sec)

real(kind=wp), public :: nominal_altitude =100.0_wp

nominal altitude (km)

real(kind=wp), public :: deadband =10.0_wp

altitude below nominal to trigger maneuver (km)

real(kind=wp), public :: r_moon =1737.4_wp

radius of the moon (km)

logical, public :: include_third_bodies =.false.

to also include Earth and Sun in force model

type(geopotential_model_pines), public :: grav

central body geopotential model

type(jpl_ephemeris), public :: eph

the ephemeris

integer, public :: n_eoms =6

size of EOM derivative vector [x,y,z,vx,vy,vz]

real(kind=wp), public :: integrator_tol =1.0e-12_wp

integrator tols

integer, public :: maxsteps =1000000

integrator max steps

integer, public :: grav_n =8

max grav degree

integer, public :: grav_m =8

max grav order

real(kind=wp), public :: root_tol =1.0e-6_wp

event tolerance for deadband (km)


Type-Bound Procedures

procedure, public :: initialize_seg => initialize_segment

  • private subroutine initialize_segment(me, alt0, deadband_alt, grav_n, grav_m, ephemeris_file, gravfile)

    Initialize the segment for integration.

    Arguments

    TypeIntentOptionalAttributesName
    class(segment), intent(inout) :: me
    real(kind=wp), intent(in) :: alt0
    real(kind=wp), intent(in) :: deadband_alt
    integer, intent(in) :: grav_n
    integer, intent(in) :: grav_m
    character(len=*), intent(in) :: ephemeris_file
    character(len=*), intent(in) :: gravfile

procedure, public :: altitude_maintenance

  • private subroutine altitude_maintenance(seg, et0, inc0, ran0, tru0, dt_max, n_dvs, dv_total, xf)

    Altitude maintenance for a circular lunar orbit - periapsis only control.

    Arguments

    TypeIntentOptionalAttributesName
    class(segment), intent(inout) :: seg
    real(kind=wp), intent(in) :: et0

    initial ephemeris time (sec)

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

    initial inclination - IAU_MOON of date (deg)

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

    initial RAAN - IAU_MOON of date (deg)

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

    initial true anomaly - IAU_MOON of date (deg)

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

    how long to propagate (days)

    integer, intent(out) :: n_dvs

    number of DV maneuvers performed

    real(kind=wp), intent(out) :: dv_total

    total DV (km/s)

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

    final state - inertial frame (km, km/s)

Source Code

    type,extends(ddeabm_with_event_class),public :: segment

        !! the main class for integrating a low-lunar orbit.

        private

        integer :: event = 0 !! event function to use:
                             !!
                             !! 1 = integrate until minimum altitude
                             !! 2 = integrate to apoapsis

        real(wp) :: et_ref !! reference ephemeris time (sec)

        real(wp) :: nominal_altitude = 100.0_wp  !! nominal altitude (km)
        real(wp) :: deadband         = 10.0_wp   !! altitude below nominal to trigger maneuver (km)
        real(wp) :: r_moon           = 1737.4_wp !! radius of the moon (km)

        logical :: include_third_bodies = .false. !! to also include Earth and Sun in force model

        type(geopotential_model_pines) :: grav  !! central body geopotential model
        type(jpl_ephemeris)            :: eph   !! the ephemeris

        integer  :: n_eoms = 6                  !! size of EOM derivative vector [x,y,z,vx,vy,vz]
        real(wp) :: integrator_tol = 1.0e-12_wp !! integrator tols
        integer  :: maxsteps = 1000000          !! integrator max steps
        integer  :: grav_n = 8                  !! max grav degree
        integer  :: grav_m = 8                  !! max grav order
        real(wp) :: root_tol = 1.0e-6_wp        !! event tolerance for deadband (km)

        contains

        procedure,public :: initialize_seg => initialize_segment
        procedure,public :: altitude_maintenance

    end type segment