set_segment_inputs Subroutine

public subroutine set_segment_inputs(me, t0, tf, x0_rotating)

Sets all the info in a segment for it to be propagated.

Type Bound

segment

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me
real(kind=wp), intent(in) :: t0

[days]

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

[days]

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

state in rotating frame


Calls

proc~~set_segment_inputs~~CallsGraph proc~set_segment_inputs segment%set_segment_inputs icrf_frame icrf_frame proc~set_segment_inputs->icrf_frame transform transform proc~set_segment_inputs->transform two_body_rotating_frame two_body_rotating_frame proc~set_segment_inputs->two_body_rotating_frame

Called by

proc~~set_segment_inputs~~CalledByGraph proc~set_segment_inputs segment%set_segment_inputs proc~initialize_the_mission mission_type%initialize_the_mission proc~initialize_the_mission->proc~set_segment_inputs proc~put_x_in_segments mission_type%put_x_in_segments proc~put_x_in_segments->proc~set_segment_inputs proc~constraint_violations mission_type%constraint_violations proc~constraint_violations->proc~put_x_in_segments proc~halo_solver_main halo_solver_main proc~halo_solver_main->proc~put_x_in_segments proc~halo_solver_main->proc~constraint_violations proc~initialize_the_solver my_solver_type%initialize_the_solver proc~halo_solver_main->proc~initialize_the_solver proc~initialize_the_solver->proc~initialize_the_mission

Source Code

    subroutine set_segment_inputs(me,t0,tf,x0_rotating)

    implicit none

    class(segment),intent(inout) :: me
    real(wp),intent(in) :: t0 !! [days]
    real(wp),intent(in) :: tf !! [days]
    real(wp),dimension(6),intent(in) :: x0_rotating  !! state in rotating frame

    real(wp),dimension(6) :: x0_inertial !! inertial frame
    type(icrf_frame) :: inertial
    type(two_body_rotating_frame) :: rotating
    real(wp) :: et0
    logical :: status_ok

    ! note that the position vectors are in the rotating frame,
    ! have to transform to inertial to integrate.
    ! create the frames for the transformation:
    et0 = me%data%et_ref + t0 * day2sec  ! convert to ephemeris time [sec]
    rotating = two_body_rotating_frame(primary_body=body_earth,&
                                        secondary_body=body_moon,&
                                        center=center_at_secondary_body,&
                                        et=et0)
    inertial = icrf_frame(b=body_moon)

    ! from rotating to inertial:
    call rotating%transform(x0_rotating,inertial,et0,me%eph,x0_inertial,status_ok)
    if (.not. status_ok) then
        write(*,*) ''
        write(*,'(A/,*(E30.16/))') 't0:         ', t0
        write(*,'(A/,*(E30.16/))') 'et0:        ', et0
        write(*,'(A/,*(E30.16/))') 'x0_rotating:', x0_rotating
        write(*,*) ''
        error stop 'transformation error in set_segment_inputs'
    end if

    ! set the inputs (needed by the propagator)
    me%data%t0          = t0
    me%data%x0_rotating = x0_rotating
    me%data%x0          = x0_inertial
    me%data%tf          = tf

    end subroutine set_segment_inputs