Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Added calc_twiss argument to twiss_and_track_from_s_to_s and twiss_and_track_intra_ele. #1358

Merged
merged 1 commit into from
Jan 10, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
68 changes: 37 additions & 31 deletions bmad/code/create_element_slice.f90
Original file line number Diff line number Diff line change
Expand Up @@ -164,37 +164,43 @@ recursive subroutine create_element_slice (sliced_ele, ele_in, l_slice, offset,

! Use a speedier tracking method.

select case (sliced_ele%tracking_method)
case (taylor$, symp_lie_ptc$)
if (sliced_ele%field_calc == fieldmap$) then
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%tracking_method = symp_lie_bmad$
case default; sliced_ele%tracking_method = symp_lie_ptc$
end select
else
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%tracking_method = symp_lie_bmad$
case (em_field$); sliced_ele%tracking_method = symp_lie_ptc$
case default; sliced_ele%tracking_method = bmad_standard$
end select
endif
end select

select case (sliced_ele%mat6_calc_method)
case (taylor$, symp_lie_ptc$)
if (sliced_ele%field_calc == fieldmap$) then
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%mat6_calc_method = symp_lie_bmad$
case default; sliced_ele%mat6_calc_method = symp_lie_ptc$
end select
else
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%mat6_calc_method = symp_lie_bmad$
case (em_field$); sliced_ele%mat6_calc_method = symp_lie_ptc$
case default; sliced_ele%mat6_calc_method = bmad_standard$
end select
endif
end select
if (sliced_ele%tracking_method == taylor$ .and. ptc_private%taylor_order_ptc == 1) then
sliced_ele%tracking_method = linear$
sliced_ele%tracking_method = auto$

else
select case (sliced_ele%tracking_method)
case (taylor$, symp_lie_ptc$)
if (sliced_ele%field_calc == fieldmap$) then
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%tracking_method = symp_lie_bmad$
case default; sliced_ele%tracking_method = symp_lie_ptc$
end select
else
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%tracking_method = symp_lie_bmad$
case (em_field$); sliced_ele%tracking_method = symp_lie_ptc$
case default; sliced_ele%tracking_method = bmad_standard$
end select
endif
end select

select case (sliced_ele%mat6_calc_method)
case (taylor$, symp_lie_ptc$)
if (sliced_ele%field_calc == fieldmap$) then
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%mat6_calc_method = symp_lie_bmad$
case default; sliced_ele%mat6_calc_method = symp_lie_ptc$
end select
else
select case (sliced_ele%key)
case (wiggler$, undulator$); sliced_ele%mat6_calc_method = symp_lie_bmad$
case (em_field$); sliced_ele%mat6_calc_method = symp_lie_ptc$
case default; sliced_ele%mat6_calc_method = bmad_standard$
end select
endif
end select
endif

sliced_ele%field_calc = refer_to_lords$

Expand Down
59 changes: 31 additions & 28 deletions bmad/code/twiss_and_track_from_s_to_s.f90
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
!+
! Subroutine twiss_and_track_from_s_to_s (branch, orbit_start, s_end, orbit_end, ele_start, ele_end,
! err, compute_floor_coords)
! err, compute_floor_coords, compute_twiss)
!
! Routine to track a particle from one location to another
!
Expand All @@ -18,16 +18,17 @@
! Note: Use element_at_s to properly initialize orbit_start.
!
! Input:
! branch -- branch_struct: Lattice branch to track through.
! orbit_start -- Coord_struct: Starting phase space coordinates at s_start.
! %s -- Starting position.
! %ix_ele -- Starting element.
! %location -- Location relative element.
! s_end -- real(rp): Ending position.
! ele_start -- Ele_struct, optional: Holds the starting parameters at s_start.
! compute_floor_coords
! -- logical, optional: If present and True then the global "floor" coordinates will be
! calculated and put in ele_end%floor.
! branch -- branch_struct: Lattice branch to track through.
! orbit_start -- Coord_struct: Starting phase space coordinates at s_start.
! %s -- Starting position.
! %ix_ele -- Starting element.
! %location -- Location relative element.
! s_end -- real(rp): Ending position.
! ele_start -- Ele_struct, optional: Holds the starting parameters at s_start.
! compute_floor_coords -- logical, optional: If present and True then the global "floor" coordinates will be
! calculated and put in ele_end%floor.
! compute_twiss -- logical, optional: Default True. If False, to save a little time, do not
! compute Twiss parameters.
!
! Output:
! orbit_end -- Coord_struct: End phase space coordinates.
Expand All @@ -38,7 +39,7 @@
!-

subroutine twiss_and_track_from_s_to_s (branch, orbit_start, s_end, orbit_end, ele_start, ele_end, &
err, compute_floor_coords)
err, compute_floor_coords, compute_twiss)

use bookkeeper_mod, dummy => twiss_and_track_from_s_to_s

Expand All @@ -58,10 +59,10 @@ subroutine twiss_and_track_from_s_to_s (branch, orbit_start, s_end, orbit_end, e
integer ix_ele

logical, optional, intent(inout) :: err
logical, optional :: compute_floor_coords
logical, optional :: compute_floor_coords, compute_twiss
logical track_upstream_end, err_flag

character(40), parameter :: r_name = 'twiss_and_track_from_s_to_s'
character(*), parameter :: r_name = 'twiss_and_track_from_s_to_s'

! Easy case & error check

Expand Down Expand Up @@ -114,17 +115,17 @@ subroutine twiss_and_track_from_s_to_s (branch, orbit_start, s_end, orbit_end, e
! Track within a single element case

if (s_true_end > s_start .and. ix_start == ix_end) then
call twiss_and_track_intra_ele (ele0, branch%param, s_start-s0, s_true_end-s0, track_upstream_end, &
.true., orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords)
call twiss_and_track_intra_ele (ele0, branch%param, s_start-s0, s_true_end-s0, track_upstream_end, .true., &
orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, compute_twiss)
return
endif

!--------------------------
! Track through multiple elements...
! First track to end of current element

call twiss_and_track_intra_ele (ele0, branch%param, s_start-s0, ele0%value(l$), track_upstream_end, &
.true., orbit_start, orbit_end, ele_start, ele_end, err_flag, compute_floor_coords)
call twiss_and_track_intra_ele (ele0, branch%param, s_start-s0, ele0%value(l$), track_upstream_end, .true., &
orbit_start, orbit_end, ele_start, ele_end, err_flag, compute_floor_coords, compute_twiss)

if (present(err)) err = err_flag
if (err_flag) return
Expand All @@ -150,17 +151,19 @@ subroutine twiss_and_track_from_s_to_s (branch, orbit_start, s_end, orbit_end, e

if (orbit_start%species == photon$) cycle

call transfer_twiss (ele_end, ele_here)
ele_here%mat6 = ele_end%mat6
ele_here%vec0 = ele_end%vec0
ele_end%key = ele_track%key
ele_end%value = ele_track%value
ele_end%mat6 = ele_track%mat6
ele_end%map_ref_orb_in = ele_track%map_ref_orb_in ! Needed for dispersion calc.
ele_end%map_ref_orb_out = ele_track%map_ref_orb_out ! Needed for dispersion calc.
call twiss_propagate1 (ele_here, ele_end, err_flag)
if (present(err)) err = err_flag
if (err_flag) return
if (logic_option(.true., compute_twiss)) then
call transfer_twiss (ele_end, ele_here)
ele_here%mat6 = ele_end%mat6
ele_here%vec0 = ele_end%vec0
ele_end%key = ele_track%key
ele_end%value = ele_track%value
ele_end%mat6 = ele_track%mat6
ele_end%map_ref_orb_in = ele_track%map_ref_orb_in ! Needed for dispersion calc.
call twiss_propagate1 (ele_here, ele_end, err_flag)
if (present(err)) err = err_flag
if (err_flag) return
endif
ele_end%vec0 = matmul(ele_end%mat6, ele_here%vec0) + ele_track%vec0
ele_end%mat6 = matmul(ele_end%mat6, ele_here%mat6)
endif
Expand Down
16 changes: 9 additions & 7 deletions bmad/code/twiss_and_track_intra_ele.f90
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
!+
! Subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, track_upstream_end, track_downstream_end,
! orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, reuse_ele_end)
! orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, compute_twiss, reuse_ele_end)
!
! Routine to track a particle within an element.
!
Expand All @@ -26,6 +26,8 @@
! instead of recomputing ele_end from scratch. This can save time.
! compute_floor_coords -- logical, optional: If present and True then the global "floor" coordinates
! (without misalignments) will be calculated and put in ele_end%floor.
! compute_twiss -- logical, optional: Default True. If False, to save a little time, do not
! compute Twiss parameters.
! reuse_ele_end -- logical, optional: If present and True, and if ele_end has the correct
! lonigitudianal length and key type, reuse ele_end from trancking instead of
! recomputing ele_end from scratch. This can save time.
Expand All @@ -42,8 +44,8 @@
! the particle gets lost in tracking
!-

recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, track_upstream_end, &
track_downstream_end, orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, reuse_ele_end)
recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, track_upstream_end, track_downstream_end, &
orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, compute_twiss, reuse_ele_end)

use bmad_interface, dummy => twiss_and_track_intra_ele

Expand All @@ -62,7 +64,7 @@ recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, trac

logical track_upstream_end, track_downstream_end, do_upstream, do_downstream, err_flag
logical track_up, track_down, length_ok
logical, optional :: err, compute_floor_coords, reuse_ele_end
logical, optional :: err, compute_floor_coords, reuse_ele_end, compute_twiss

character(*), parameter :: r_name = 'twiss_and_track_intra_ele'

Expand Down Expand Up @@ -91,8 +93,8 @@ recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, trac
l1 = min(slave%value(l$), s_end - slave%s_start)
track_up = (track_upstream_end .or. s_start < slave%s_start)
track_down = (track_downstream_end .or. s_end > slave%s)
call twiss_and_track_intra_ele (slave, param, l0, l1, track_up, &
track_down, orbit_end, orbit_end, ele_end, ele_end, err, compute_floor_coords)
call twiss_and_track_intra_ele (slave, param, l0, l1, track_up, track_down, orbit_end, &
orbit_end, ele_end, ele_end, err, compute_floor_coords, compute_twiss = compute_twiss)
if (s_end <= slave%s + bmad_com%significant_length) exit
enddo
return
Expand Down Expand Up @@ -190,7 +192,7 @@ recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, trac
call make_mat6 (ele_p, param, err_flag = err_flag)
endif
if (err_flag) return
call twiss_propagate1 (ele_start, ele_p, err_flag)
if (logic_option(.true., compute_twiss)) call twiss_propagate1 (ele_start, ele_p, err_flag)
if (logic_option(.false., compute_floor_coords)) call ele_geometry (ele_start%floor, ele_p, ele_p%floor)
if (.not. associated(ele_p, ele_end)) call transfer_ele(ele_p, ele_end, .true.)
if (err_flag) return
Expand Down
10 changes: 5 additions & 5 deletions bmad/modules/bmad_routine_interface.f90
Original file line number Diff line number Diff line change
Expand Up @@ -3101,19 +3101,19 @@ subroutine transfer_wall3d (wall3d_in, wall3d_out)
end subroutine

subroutine twiss_and_track_from_s_to_s (branch, orbit_start, s_end, orbit_end, &
ele_start, ele_end, err, compute_floor_coords)
ele_start, ele_end, err, compute_floor_coords, compute_twiss)
import
implicit none
type (coord_struct) :: orbit_start, orbit_end
type (ele_struct), optional :: ele_start, ele_end
type (branch_struct), target :: branch
real(rp) s_end
logical, optional, intent(inout) :: err
logical, optional :: compute_floor_coords
logical, optional :: compute_floor_coords, compute_twiss
end subroutine

recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, track_upstream_end, &
track_downstream_end, orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, reuse_ele_end)
recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, track_upstream_end, track_downstream_end, &
orbit_start, orbit_end, ele_start, ele_end, err, compute_floor_coords, compute_twiss, reuse_ele_end)
import
implicit none
type (coord_struct), optional :: orbit_start, orbit_end
Expand All @@ -3122,7 +3122,7 @@ recursive subroutine twiss_and_track_intra_ele (ele, param, l_start, l_end, trac
type (lat_param_struct) param
real(rp) l_start, l_end
logical track_upstream_end, track_downstream_end
logical, optional :: err, compute_floor_coords, reuse_ele_end
logical, optional :: err, compute_floor_coords, reuse_ele_end, compute_twiss
end subroutine

recursive subroutine twiss_at_element (ele, start_ele, end_ele, average)
Expand Down
2 changes: 1 addition & 1 deletion tao/code/tao_graph_setup_mod.f90
Original file line number Diff line number Diff line change
Expand Up @@ -2447,7 +2447,7 @@ subroutine tao_calc_data_at_s_pts (tao_lat, curve, comp_sign, good)

else
call twiss_and_track_from_s_to_s (branch, orbit, s_now, orbit_end, ele, ele, err_flag, compute_floor_coords = .true.)
mat6 = matmul(ele%mat6, mat6)
mat6 = matmul(ele%mat6, mat6) ! Matrix from beginning of branch.
vec0 = matmul(ele%mat6, vec0) + ele%vec0
orbit = orbit_end
endif
Expand Down
Loading
Loading