!+ ! subroutine MARK_LRBBI_ONLY(master_lat, master_lat_oppos, lat, crossings) ! ! Subroutine to insert markers at parasitic crossing point. ! ! Modules Needed: ! use bmad ! ! Input: ! master_lat -- Lat struct: Lat with markers at LRBBI locations. ! master_lat_oppos -- Lat struct: Lat for oppositely circulating ! particles, with markers at all LRBBI locations. ! lat -- Lat struct: Lattice ! crossings(:,:) -- Real(rp): First column is the position (as a ! fraction/percent of the total lat length) of a beam- ! beam interaction seen by a bunch. Second column is unused. ! (Third column should be an index ! indicating which crossing this is: 1rst, 2nd, etc, ! but it is not used or changed here.) The fourth and ! fifth columns should be empty and will hold, ! respectively, the index of the crossings in the lat, ! and the index in the master_lat. ! ! Output: ! lat -- Lat struct: lattice with markers placed where ! parasitic crossings are seen by that bunch. ! master_lat -- Lat struct: Master_lat with markers placed at every ! parasitic crossing (seen by any bunch). ! master_lat_oppos -- Lat struct: Master_lat_oppos with markers placed ! at every parasitic crossings (seen by any bunch). ! ! ! Note: The elements placed at the parasitic crossing sites are simply markers ! with unit 6x6 matrices. The fourth and fifth columns of crossings hold ! the indices of the location of all inserted markers (named ! 'LRBBI_MARKER'). lat_make_mat6 is called for all lattices. !- subroutine MARK_LRBBI_ONLY(master_lat, master_lat_oppos, lat, crossings) use bookkeeper_mod implicit none type (lat_struct) :: lat type (lat_struct) :: master_lat, master_lat_oppos type (ele_struct), allocatable, save :: insert_ele(:) real(rp), dimension(:,:) :: crossings real(rp) :: smallest, s_split integer :: i, j, k, m, ierr, loc_smallest, ix_split integer :: ix_ele, master_ix_split, master_ix_ele, total integer :: index1, index2, index3, index4 integer, dimension(1) :: minloc_array logical :: ok, split_done, split_done_1, split_done_2, split_done_3 ! total = size(crossings, 1) ierr = 0 if (.not. allocated(insert_ele)) then allocate(insert_ele(1:total), stat=ierr) elseif (size(insert_ele) /= total) then deallocate (insert_ele) allocate(insert_ele(1:total), stat=ierr) endif if (ierr .ne. 0) then print *, "INSERT_ELE: ALLOCATION REQUEST DENIED." call err_exit endif !--------------------------------------------------------------------------- do k = 1, total smallest = minval(crossings(k:total,1)) minloc_array = minloc(crossings(k:total,1)) loc_smallest = k - 1 + minloc_array(1) index1 = crossings(loc_smallest, 2) index2 = crossings(loc_smallest, 3) if (loc_smallest .ne. k) then crossings(loc_smallest,1) = crossings(k,1) crossings(loc_smallest, 2) = crossings(k,2) crossings(loc_smallest, 3) = crossings(k,3) crossings(k,1) = smallest crossings(k,2) = index1 crossings(k,3) = index2 endif enddo j_loop: do j = 1, total s_split = crossings(j,1) * master_lat%param%total_length if (s_split == 0) cycle if (j .gt. 1 ) then if (crossings(j,1) == crossings(j-1,1)) then call split_lat(lat, s_split, 0, ix_split, split_done) ! if (.not. split_done) cycle j_loop call init_ele(insert_ele(j)) insert_ele(j)%key = marker$ insert_ele(j)%name = "LRBBI_MARKER" call mat_make_unit(insert_ele(j)%mat6) insert_ele(j)%s = s_split ix_ele = ix_split + 1 crossings(j,4) = ix_ele crossings(j,5) = master_ix_ele call insert_element(lat, insert_ele(j), ix_ele) cycle j_loop endif endif ix_split = 0 call split_lat(master_lat_oppos, s_split, 0, master_ix_split,split_done_1) call split_lat(master_lat, s_split, 0, master_ix_split, split_done_2) call split_lat(lat, s_split, 0, ix_split, split_done_3) call init_ele(insert_ele(j)) insert_ele(j)%key = marker$ insert_ele(j)%name = "LRBBI_MARKER" call mat_make_unit(insert_ele(j)%mat6) insert_ele(j)%s = s_split ix_ele = ix_split + 1 crossings(j, 4) = ix_ele master_ix_ele = master_ix_split + 1 crossings(j,5) = master_ix_ele call insert_element(master_lat, insert_ele(j), master_ix_ele) call insert_element(master_lat_oppos, insert_ele(j), master_ix_ele) call insert_element(lat, insert_ele(j), ix_ele) ! print *,' MARK_LRBBI_ONLY: ',j,insert_ele(j)%name, s_split enddo j_loop ! call lattice_bookkeeper (master_lat) call lattice_bookkeeper (master_lat_oppos) call lattice_bookkeeper (lat) end subroutine MARK_LRBBI_ONLY