put_x_in_segments Subroutine

public subroutine put_x_in_segments(me, x_scaled)

Take the x vector from the solver, and use it to populate all the segments

Type Bound

mission_type

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x_scaled

opt var vector (scaled)


Calls

proc~~put_x_in_segments~~CallsGraph proc~put_x_in_segments mission_type%put_x_in_segments extract_vector extract_vector proc~put_x_in_segments->extract_vector proc~define_problem_size mission_type%define_problem_size proc~put_x_in_segments->proc~define_problem_size proc~set_segment_inputs segment%set_segment_inputs proc~put_x_in_segments->proc~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~~put_x_in_segments~~CalledByGraph proc~put_x_in_segments mission_type%put_x_in_segments 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

Source Code

    subroutine put_x_in_segments(me,x_scaled)

    implicit none

    class(mission_type),intent(inout) :: me
    real(wp),dimension(:),intent(in) :: x_scaled !! opt var vector (scaled)

    integer :: i,j,irev,iseg  !! counter
    real(wp),dimension(8) :: t0 !! [days]
    real(wp),dimension(8) :: tf !! [days]
    real(wp),dimension(6,8) :: x0_rotating  !! rotating frame
    real(wp),dimension(:),allocatable :: x !! opt var vector (unscaled)
    integer :: n_segs  !! number of segments

    call me%define_problem_size(n_segs=n_segs)

    ! unscale the x vector:
    allocate(x(size(x_scaled)))
    x = x_scaled * me%xscale

    ! first extract data from the opt var vector and put it into the segments:
    i = 0  ! initialize index, will be updated by extract_vector
    iseg = 0
    do irev = 1, me%n_revs

        ! extract the values:

        if (irev==1) then

            ! the first one has an extra opt point at the initial periapsis passage

            if (me%fix_initial_time) then
                t0(1) = me%segs(1)%data%t0    ! not in x, just keep the current value
            else
                call extract_vector(t0(1), x, i)
            end if
            if (me%fix_initial_r) then
                x0_rotating(1:3,1) = me%segs(1)%data%x0_rotating(1:3) ! not in x, just keep the current r value
                call extract_vector(x0_rotating(4:6,1), x, i) ! v only
            else
                call extract_vector(x0_rotating(:,1), x, i)
            end if
            tf(1)               = t0(1) + me%period8

            call extract_vector(t0(2)               , x, i)
            call extract_vector(x0_rotating(:,2)    , x, i)
            tf(2)               = tf(1)

            call extract_vector(t0(4)               , x, i)
            call extract_vector(x0_rotating(:,4)    , x, i)
            tf(4)               = t0(4) - me%period8

            t0(3)               = t0(2)
            x0_rotating(:,3)    = x0_rotating(:,2)
            tf(3)               = tf(4)

            t0(5)               = t0(4)
            x0_rotating(:,5)    = x0_rotating(:,4)
            tf(5)               = t0(5) + me%period8

            call extract_vector(t0(6)               , x, i)
            call extract_vector(x0_rotating(:,6)    , x, i)
            tf(6)               = tf(5)

            call extract_vector(t0(8)               , x, i)
            if (me%fix_ry_at_end_of_rev == 1) then
                x0_rotating(2,8) = me%segs(8)%data%x0_rotating(2) ! not in x, just keep the current ry value
                call extract_vector(x0_rotating(1,8) , x, i)
                call extract_vector(x0_rotating(3:6,8) , x, i)
            else
                call extract_vector(x0_rotating(:,8)    , x, i)
            end if
            tf(8)               = t0(8) - me%period8

            t0(7)               = t0(6)
            x0_rotating(:,7)    = x0_rotating(:,6)
            tf(7)               = tf(8)

            ! put the data for this rev into the segments:
            do j = 1, 8
                call me%segs(iseg+j)%set_input(t0(j),tf(j),x0_rotating(:,j))
            end do

            ! for next rev:
            iseg = iseg + 8

        else

            t0(1)               =  t0(8)                ! inherits from the previous rev
            x0_rotating(:,1)    =  x0_rotating(:,8)
            tf(1)               =  t0(1) + me%period8

            call extract_vector(t0(2)               , x, i)
            call extract_vector(x0_rotating(:,2)    , x, i)
            tf(2)               =  tf(1)

            call extract_vector(t0(4)               , x, i)
            call extract_vector(x0_rotating(:,4)    , x, i)
            tf(4)               =  t0(4) - me%period8

            t0(3)               =  t0(2)
            x0_rotating(:,3)    =  x0_rotating(:,2)
            tf(3)               =  tf(4)

            t0(5)               =  t0(4)
            x0_rotating(:,5)    =  x0_rotating(:,4)
            tf(5)               =  t0(5) + me%period8

            call extract_vector(t0(6)               , x, i)
            call extract_vector(x0_rotating(:,6)    , x, i)
            tf(6)               =  tf(5)

            call extract_vector(t0(8)               , x, i)

            if (me%fix_ry_at_end_of_rev == irev) then
                x0_rotating(2,8) = me%segs(irev*8)%data%x0_rotating(2) ! not in x, just keep the current ry value
                call extract_vector(x0_rotating(1,8)   , x, i)
                call extract_vector(x0_rotating(3:6,8) , x, i)
            else if (irev==me%n_revs .and. me%fix_final_ry_and_vx) then
                x0_rotating(2,8) = me%segs(n_segs)%data%x0_rotating(2) ! not in x, just keep the current ry value
                x0_rotating(4,8) = me%segs(n_segs)%data%x0_rotating(4) ! not in x, just keep the current vx value
                call extract_vector(x0_rotating(1,8)    , x, i)
                call extract_vector(x0_rotating(3,8)    , x, i)
                call extract_vector(x0_rotating(5:6,8)  , x, i)
            else
                call extract_vector(x0_rotating(:,8)    , x, i)
            end if
            tf(8)               =  t0(8) - me%period8

            t0(7)               =  t0(6)
            x0_rotating(:,7)    =  x0_rotating(:,6)
            tf(7)               =  tf(8)

            ! put the data for this rev into the segments:
            do j = 1, 8
                call me%segs(iseg+j)%set_input(t0(j),tf(j),x0_rotating(:,j))
            end do

            ! for next rev:
            iseg = iseg + 8

        end if

    end do

    end subroutine put_x_in_segments