export_trajectory_json_file Subroutine

public subroutine export_trajectory_json_file(me, filename, only_first_rev)

Export the trajectory JSON file.

Note

It is assumed that all the data is present in the segments needed to propagate.

Type Bound

mission_type

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
character(len=*), intent(in) :: filename

plot file name [without extension]

logical, intent(in), optional :: only_first_rev

to only do the first rev [Default is False]


Calls

proc~~export_trajectory_json_file~~CallsGraph proc~export_trajectory_json_file mission_type%export_trajectory_json_file add add proc~export_trajectory_json_file->add create_array create_array proc~export_trajectory_json_file->create_array create_object create_object proc~export_trajectory_json_file->create_object destroy destroy proc~export_trajectory_json_file->destroy initialize initialize proc~export_trajectory_json_file->initialize insert_after insert_after proc~export_trajectory_json_file->insert_after proc~compute_rdot_vecs compute_rdot_vecs proc~export_trajectory_json_file->proc~compute_rdot_vecs proc~destroy_trajectory trajectory%destroy_trajectory proc~export_trajectory_json_file->proc~destroy_trajectory proc~export_rp_ra_json_file mission_type%export_rp_ra_json_file proc~export_trajectory_json_file->proc~export_rp_ra_json_file proc~propagate_segment segment%propagate_segment proc~export_trajectory_json_file->proc~propagate_segment proc~export_rp_ra_json_file->add proc~export_rp_ra_json_file->destroy proc~export_rp_ra_json_file->initialize evaluate evaluate proc~export_rp_ra_json_file->evaluate solve solve proc~export_rp_ra_json_file->solve first_call first_call proc~propagate_segment->first_call icrf_frame icrf_frame proc~propagate_segment->icrf_frame integrate integrate proc~propagate_segment->integrate proc~set_segment_outputs segment%set_segment_outputs proc~propagate_segment->proc~set_segment_outputs transform transform proc~propagate_segment->transform two_body_rotating_frame two_body_rotating_frame proc~propagate_segment->two_body_rotating_frame

Called by

proc~~export_trajectory_json_file~~CalledByGraph proc~export_trajectory_json_file mission_type%export_trajectory_json_file proc~halo_solver_main halo_solver_main proc~halo_solver_main->proc~export_trajectory_json_file

Source Code

    subroutine export_trajectory_json_file(me,filename,only_first_rev)

    implicit none

    class(mission_type),intent(inout) :: me
    character(len=*),intent(in) :: filename !! plot file name [without extension]
    logical,intent(in),optional :: only_first_rev  !! to only do the first rev [Default is False]

    integer :: iseg  !! segment number counter
    integer :: nsegs_to_plot !! number of segments to plot
    type(json_core) :: json
    type(json_value),pointer :: p_root, p_segs, p_seg, p_current
    real(wp),dimension(:),allocatable :: rdot, rmag
    logical :: accumulate_rdot !! to generate the rdot/rp/ra file
    real(wp),dimension(:),allocatable :: cumulative_et,cumulative_rdot,cumulative_r
    integer :: n

    write(*,'(A)') ' * Export the trajectory JSON file.'

    ! optional arguments:
    nsegs_to_plot = size(me%segs) ! default export all the segments
    if (present(only_first_rev)) then
        if (only_first_rev) nsegs_to_plot = 8 ! only the first rev (8 segments)
    end if
    accumulate_rdot = me%generate_rp_ra_file

    ! the JSON file will contain an array of segments:
    call json%initialize(compress_vectors=.true.)
    call json%create_object(p_root, '')
    call json%create_array(p_segs, 'segs')
    call json%add(p_root, p_segs)
    p_current => null()

    do iseg = 1, nsegs_to_plot
        call destroy_traj(iseg)
    end do

    !====================================
    ! now propagate the segments:
!$OMP PARALLEL DO    !...FIRSTPRIVATE(me)
    do iseg = 1, nsegs_to_plot
        call me%segs(iseg)%propagate(mode=2)  ! [export points]
    end do
!$OMP END PARALLEL DO
    !====================================

    if (accumulate_rdot) then
        allocate(cumulative_et(0))
        allocate(cumulative_rdot(0))
        allocate(cumulative_r(0))
    end if

    do iseg = 1, nsegs_to_plot

        associate (seg => me%segs(iseg))
            ! generate the trajectory for this segment:
            !call destroy_traj(iseg)
            !call seg%propagate(mode=2)  ! [export points]

            ! create the segment object for exporting the trajectory:
            call json%create_object(p_seg, '')
            if (associated(p_current)) then
                call json%insert_after(p_current, p_seg) ! next one
            else
                call json%add(p_segs, p_seg) ! first one
            end if
            p_current => p_seg ! update for next seg

            call json%add(p_seg, 'iseg', iseg)
            call json%add(p_seg, 'et', seg%traj_inertial%et)

            call json%add(p_seg, 'x_inertial',  seg%traj_inertial%x)
            call json%add(p_seg, 'y_inertial',  seg%traj_inertial%y)
            call json%add(p_seg, 'z_inertial',  seg%traj_inertial%z)
            call json%add(p_seg, 'vx_inertial', seg%traj_inertial%vx)
            call json%add(p_seg, 'vy_inertial', seg%traj_inertial%vy)
            call json%add(p_seg, 'vz_inertial', seg%traj_inertial%vz)

            if (accumulate_rdot) then
                call compute_rdot_vecs(seg%traj_inertial%x,seg%traj_inertial%y,seg%traj_inertial%z,&
                                       seg%traj_inertial%vx,seg%traj_inertial%vy,seg%traj_inertial%vz,&
                                       rmag, rdot)
                call append_traj_to_arrays(seg)
                !call json%add(p_seg, 'rdot_inertial', rdot)
                ! don't include the last time step since that will overlap the next segment
                ! n = size(seg%traj_inertial%x)
                ! cumulative_et   = [cumulative_et,   seg%traj_inertial%et]
                ! cumulative_rdot = [cumulative_rdot, rdot]
                ! cumulative_r    = [cumulative_r,    rmag]
            end if

            call json%add(p_seg, 'x_rotating',  seg%traj_rotating%x)
            call json%add(p_seg, 'y_rotating',  seg%traj_rotating%y)
            call json%add(p_seg, 'z_rotating',  seg%traj_rotating%z)
            ! call json%add(p_seg, 'vx_rotating', seg%traj_rotating%vx)
            ! call json%add(p_seg, 'vy_rotating', seg%traj_rotating%vy)
            ! call json%add(p_seg, 'vz_rotating', seg%traj_rotating%vz)
            ! call compute_rdot_vecs(seg%traj_rotating%x,seg%traj_rotating%y,seg%traj_rotating%z,&
            !                        seg%traj_rotating%vx,seg%traj_rotating%vy,seg%traj_rotating%vz,&
            !                        rdot)
            ! call json%add(p_seg, 'rdot_rotating', rdot)

            call json%add(p_seg, 'x_se_rotating',  seg%traj_se_rotating%x)
            call json%add(p_seg, 'y_se_rotating',  seg%traj_se_rotating%y)
            call json%add(p_seg, 'z_se_rotating',  seg%traj_se_rotating%z)
            ! call json%add(p_seg, 'vx_se_rotating', seg%traj_se_rotating%vx) ! don't need these
            ! call json%add(p_seg, 'vy_se_rotating', seg%traj_se_rotating%vy)
            ! call json%add(p_seg, 'vz_se_rotating', seg%traj_se_rotating%vz)

            !TODO:
            !  - maybe also earth & sun ephemeris for plotting
            !  - solar fraction to color the trajectory [but really that should go in the eclipse file?]

            !call destroy_traj(iseg) ! keep them for the eclipse file generation ...

        end associate
    end do

    if (accumulate_rdot) then
        call me%export_rp_ra_json_file(cumulative_et, cumulative_r, cumulative_rdot, filename)
    end if

    do iseg = 1, nsegs_to_plot
        call destroy_traj(iseg)
    end do

    call json%print(p_root, trim(filename)//'.json')
    call json%destroy(p_root)

    contains

        subroutine destroy_traj(iseg)
            integer,intent(in) :: iseg !! segment number
            call me%segs(iseg)%traj_inertial%destroy()
            call me%segs(iseg)%traj_rotating%destroy()
            call me%segs(iseg)%traj_se_rotating%destroy()
        end subroutine destroy_traj

        subroutine append_traj_to_arrays(seg)
            type(segment),intent(in) :: seg

            integer :: n, i

            n = size(seg%traj_inertial%et)
            if (seg%traj_inertial%et(1)<=seg%traj_inertial%et(2)) then
                ! forward propagated segment
                cumulative_et   = [cumulative_et,   seg%traj_inertial%et(1:n-1)]
                cumulative_rdot = [cumulative_rdot, rdot(1:n-1)]
                cumulative_r    = [cumulative_r,    rmag(1:n-1)]
            else
                ! some of the segments are propagated backwards, so we have to reverse them
                ! [note this is very inefficient... need to fix this]  -TODO
                cumulative_et   = [cumulative_et,   seg%traj_inertial%et(n:2:-1)]
                cumulative_rdot = [cumulative_rdot, rdot(n:2:-1)]
                cumulative_r    = [cumulative_r,    rmag(n:2:-1)]
            end if

        end subroutine append_traj_to_arrays

    end subroutine export_trajectory_json_file