!........................................................................ ! ! Subroutine : beta_v_amplitude(LAT, nonlin) ! ! Description: Compute beta at each element for an orbit that starts ! with zero angle at IP (cos) and for another that starts ! with zero displacement at IP (sin). This gives amplitude ! dependence of beta ! ! Arguments : INPUT: ! LAT -- lat_struct : Lat ! OUTPUT: ! NONLIN -- Nonlin_ele_struct : ! Nonlin.non_ele(i).x.dbeta_dcos, etc ! ! Mod/Commons: ! ! Calls : ! ! Author : ! ! Modified : ! !........................................................................ ! ! ! $Log$ ! Revision 1.11 2007/01/30 16:15:13 dcs ! merged with branch_bmad_1. ! ! Revision 1.8 2006/02/08 19:09:46 mjf7 ! Replaced obsolete bmad subroutine calls with new versions, and found variables which needed a save attribute. ! ! Revision 1.7 2004/11/08 19:16:26 dlr ! replace n_ele_maxx with n_ele_max ! ! Revision 1.6 2003/07/22 12:12:17 mjf7 ! Synchronized some orbit allocation sizes. - mjf ! ! Revision 1.5 2003/07/17 20:46:50 mjf7 ! Fixed a bug which occurred at the end of an optimization loop. ! Changed many allocatable arrays to have the save parameter for faster ! performance. - mjf ! ! Revision 1.4 2003/07/16 19:04:27 dcs ! memory leak fix. ! ! Revision 1.3 2003/07/08 19:27:48 mjf7 ! ! ! Modified all subroutines to correctly use allocatable lat elements. n_ele_maxx is no longer global, but a member variable of the lat struct. -mjf ! ! Revision 1.2 2003/04/30 17:14:47 cesrulib ! dlr's changes since last import ! ! Revision 1.1.1.1 2002/12/13 19:23:27 cesrulib ! import bmadz ! ! !........................................................................ ! subroutine beta_v_amplitude(lat, nonlin) use bmad use nonlin_mod implicit none type (lat_struct) lat type (lat_struct), save :: lat_amp type (nonlin_ele_struct) nonlin type (coord_struct), allocatable, save :: co_a(:) type(coord_struct) co real(rp) :: amp = 0.002 integer i call reallocate_coord( co_a, lat%n_ele_max ) ! reference co_a(0)%vec(:) = 0. call closed_orbit_calc (lat, co_a, 4) co%vec(:) = co_a(0)%vec(:) call track_all (lat, co_a) call lat_make_mat6(lat,-1,co_a) call twiss_at_start (lat) call twiss_propagate_all (lat) lat_amp = lat ! horizontal and vertical cos co_a(0)%vec(:) = co%vec(:) co_a(0)%vec(1) = co_a(0)%vec(1) + amp / sqrt(lat%ele(0)%a%beta) call track_all (lat, co_a) call lat_make_mat6(lat_amp,-1,co_a) call twiss_at_start (lat_amp) call twiss_propagate_all (lat_amp) do i = 1, lat%n_ele_track nonlin%non_ele(i)%a%dbeta_dcos = & (lat_amp%ele(i)%a%beta - lat%ele(i)%a%beta)/lat%ele(i)%a%beta/amp nonlin%non_ele(i)%b%dbeta_dcos = & (lat_amp%ele(i)%b%beta - lat%ele(i)%b%beta)/lat%ele(i)%b%beta/amp end do ! horizontal and vertical sin co_a(0)%vec(:) = co%vec(:) co_a(0)%vec(2) = co%vec(2) + amp * sqrt(lat%ele(0)%a%beta) call track_all (lat, co_a) call lat_make_mat6(lat_amp,-1,co_a) call twiss_at_start (lat_amp) call twiss_propagate_all (lat_amp) do i = 1, lat%n_ele_track nonlin%non_ele(i)%a%dbeta_dsin = & (lat_amp%ele(i)%a%beta - lat%ele(i)%a%beta)/lat%ele(i)%a%beta/amp nonlin%non_ele(i)%b%dbeta_dsin = & (lat_amp%ele(i)%b%beta - lat%ele(i)%b%beta)/lat%ele(i)%b%beta/amp end do return end subroutine beta_v_amplitude