subroutine showit (what, u, file_name)

use cesrv_struct
use cesrv_interface
use nonlin_bpm_mod
use super_universe_com
use cesr_read_data_mod
use cesr_element_mod
use geometry_mod
use directory_mod
use bookkeeper_mod
use twiss_and_track_mod
use nr

implicit none

type (var_struct), pointer :: vp, v1, v2, v3
type (v1_var_struct), pointer :: var1
type (ele_struct), pointer :: ele, ip_design
type (ele_struct), save :: ave, runt, ele_p, ele_e
type (universe_struct), target :: u
type (but_value_struct) gain(0:120)
type (d2_data_struct), pointer :: d2_ptr
type (d1_data_struct), pointer :: d1_ptr
type (data_struct), pointer :: x_dat, y_dat, d_ptr
type (top10s) :: top_steer(150)
type (normal_modes_struct) global_model, global_e, global_p
type (butns_struct) butns
type (coord_struct) orb
type (coord_struct), pointer :: orbit
type (rad_int_all_ele_struct), target :: rad_int_by_ele_model
type (db_element_struct), pointer :: db_ele
type (tbt_all_data_struct), pointer :: tbt

real(rp) merit0, mr_model, mr_db, ff, f_phi, l_ring, a_max, ecc
real(rp) cu_model_min, cu_saved_min, d_target, s_pos, del_s
real(rp) sig_y_no, eps_x, eps_y, temp, theta_diff
real(rp) sig_x, sig_y, factor, sig_y_tot, model_val, design_val
real(rp) n_part, beta_a, beta_b, cbar12, dcbar22, f_xi_x, f_xi_y

real(rp), allocatable :: show_twiss_at(:)

real(rp) cbar_model(2,2), cbar_design(2,2), cbar_p(2,2), cbar_e(2,2), vec0(6)
real(rp) sig_a(6,6), sig_b(6,6), sig_c(6,6), mat6(6,6)

integer i, j, ix, ix1, ix2, ix_s2, show_what, nl, ix_99, ix_0, ix_ip
integer loc, ios, save_set, dcu, cu_model, n, ix_set
integer iu, dcu_model, this_species, opt_vars

character(*) what
character(*), optional :: file_name  

character(24) show_name, show2_name, mname
character(24) :: show_names(64) = [ &
      'QUADRUPOLE       ', 'SEXTUPOLE        ', 'DETECTORS        ', 'IP               ', 'WAVE             ', &
      'ALL              ', 'TOP10            ', 'STEERINGS        ', 'ORBIT            ', 'ASYMMETRY        ', &
      'SINGLE_CORRECTION', 'RING             ', 'VARIABLE         ', 'SEPARATORS       ', 'SKEW_QUAD        ', &
      'PHASE            ', 'ETA              ', 'ELEMENT          ', 'PARAMETERS       ', 'COMMENTS         ', &
      'TUNE             ', 'CALIBRATION      ', 'LUMINOSITY       ', 'COUPLING         ', 'CUSTOM           ', &
      'TWISS            ', 'DATA             ', 'GEOMETRY         ', 'WEIGHT           ', 'ENERGY_SHIFT     ', &
      'CUT_RING         ', 'OPT_VARIABLES    ', 'HORIZONTAL       ', 'VERTICAL         ', 'SETS             ', &
      'DMERIT           ', 'CHROMATICITY     ', 'LOGICALS         ', 'GROUP            ', 'ALL_TWISS        ', &
      'BPM_TILT         ', 'CBAR11           ', 'CBAR12           ', 'CBAR22           ', 'GLOBAL           ', &
      'OCTUPOLE         ', '_SKEW_SEX        ', 'LIGHTSOURCES     ', 'AC_ETA           ', 'COMMAND_FILES    ', &
      '2QX              ', '2QY              ', 'QX+QY            ', 'QX-QY            ', 'WOLSKI           ', &
      'MOD_ETA          ', 'MODE_ETA         ', 'RAD_INT          ', 'BPM_GAINS        ', 'BETA             ', &
      'XRAY             ', 'TRANSFER_MATRIX  ', 'FLOOR            ', 'TBT              ']

character(3) ew_name(0:120)
character(20) prefix
character(40) name, route_name, match, write_append_str
character(120) line, l0, l1, l2, l3
character(80) fmt, fmt2, fmt3, fmt4, fmt5, amt, imt, rmt, lmt
character(200), allocatable :: alloc_lines(:)
character(200) lines(2000), directory
character(16) ele_name
character(200) fname
character stub*16, date*24
character(140), save :: old_file_name = ' '

character data_date*20, lattice*40, comment*70, file_type*16

namelist / data_parameters / data_date, lattice, comment, &
        file_type, save_set, route_name

logical err_flag, name_found, init_needed / .true. /
logical ok, do_more, found, at_ends, match_route, match_lattice, show_all
logical :: take_pos = .false., take_ele = .false.

real(rp) epsilon_x, epsilon_y,  xi_y, i_tot, xi_x
integer n_bunch

namelist / lum_params / epsilon_x, epsilon_y, xi_x, xi_y, i_tot, n_bunch

! Init

rmt  = '(a, 9es16.8)'
imt  = '(a, 9i8)'
lmt  = '(a, 9(l3))'
amt  = '(9a)'

if (logic%lattice == 'INIT') then
  print *, 'Lattice not yet set...'
  return
endif  

! find what to show

if (init_needed) then
  do i = 0, 120
    if (i .le. 49) then
      write (ew_name(i), '(i2.2, a)') i, 'W'
    elseif (i .le. 99) then
      write (ew_name(i), '(i2.2, a)') 99-i, 'E'
    else
      write (ew_name(i), '(i3)') i
    endif
  enddo
  init_needed = .false.
endif

f_phi = logic%phase_conversion_factor

! add blank for parsing before a "/" for parsing.

line = what
ix1 = index(line, '/')
if (ix1 /= 0) line = line(:ix1-1) // ' ' // line(ix1:)

call string_trim (line, line, ix)

if (ix == 0) then
  show_name = 'DEFAULT'
else
  show_name = line(:ix)
  call string_trim (line(ix+1:), line, ix_s2)
  show2_name = line(1:ix_s2)

  if (index('CBAR', trim(show_name)) == 1) show_name = 'CBAR12'
  call match_word (show_name, show_names, show_what)

  if (show_what <= 0) then
    print *, 'ERROR: SHOW WHAT? I DO NOT UNDERSTAND: ', show_name
    return
  endif

  show_name = show_names(show_what)

endif

nl = 0

! select

select case (show_name)

!----------------------------------------------------------------------
! show eta

case ('AC_ETA')

  fmt = '(a, es12.2)'
  nl=nl+1; write (lines(nl), fmt) 'amp_fit   =', u%ac_eta%data_params%ac_z_amp_fit
  nl=nl+1; write (lines(nl), fmt) 'phase_fit =', u%ac_eta%data_params%ac_z_phase_fit
  nl=nl+1; write (lines(nl), fmt) 'chi^2     =', u%ac_eta%data_params%chisq
  nl=nl+1; lines(nl) = ''
  nl=nl+1; write (lines(nl), '(8x, a)') &
    '|          AC_Eta_X             |        AC_Eta_Y            |'
  nl=nl+1; write (lines(nl), '(8x, a)') &
    '|    Model    Design      Data  |  Model    Design      Data |'
  do i = lbound(u%phase%x%d, 1), ubound(u%phase%x%d, 1)
    if (u%ac_eta%x%d(i)%exists) then
      x_dat => u%ac_eta%x%d(i)
      y_dat => u%ac_eta%y%d(i)
      ele => u%ring%ele(x_dat%ix_ele)
      nl=nl+1
      write (lines(nl), '(1x, a7, 6f10.4, i5)') ele%name, &
          x_dat%model, x_dat%design, x_dat%meas, &
          y_dat%model, y_dat%design, y_dat%meas, i
    endif
  enddo
  nl=nl+1; write (lines(nl), '(8x, a)') &
    '|    Model    Design      Data  |  Model    Design      Data |'
  nl=nl+1; write (lines(nl), '(8x, a)') &
    '|          AC_Eta_X             |         AC_Eta_Y           |'

!----------------------------------------------------------------------
! ALL_TWISS

case ('ALL_TWISS')
  
  iu = lunget()
  open (iu, file = show2_name, status = 'old', action = 'read', iostat = ios)
  if (ios /= 0) then
    print *, "Could not open s positions file: ", show2_name
    return
  endif

  ! get s positions

  ! first find total number of positions
  i = 0
  do 
    read (iu, *, iostat = ios) temp
    if (ios > 0) then
      print *, "Read error in s positions file."
      return
    elseif (ios < 0) then
      exit
    else
      i = i + 1
    endif
  enddo

  allocate (show_twiss_at(i))

  ! header
  write (l1, '(21x, 2a)')  &
                '   |               X              (mm)  (mrad)', &
                   '|               Y              (mm)  (mrad)'
  write (l1, '(5x, a, 6x, 2a)') 'In Element', &
                'S  | Beta   Alpha   Eta    Phi     X      X_P ', &
       '| Beta   Alpha   Eta    Phi     Y      Y_P '
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  rewind(iu)
  do i = 1, size(show_twiss_at)
    read (iu, *) show_twiss_at(i)
    call twiss_and_track_at_s (u%ring, show_twiss_at(i), runt, u%orb, orb)
    nl=nl+1; write (lines(nl), '(1x, a15, f8.3, 2(2f7.2, f7.2, f8.3, 2f7.2))') &
      runt%name, show_twiss_at(i), &
      runt%a%beta, runt%a%alpha, runt%a%eta, f_phi*runt%a%phi, &
                          1000*orb%vec(1), 1000*orb%vec(2),&
      runt%b%beta, runt%b%alpha, runt%b%eta, f_phi*runt%b%phi, &
                          1000*orb%vec(3), 1000*orb%vec(4)
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1
  
  close(iu)

!----------------------------------------------------------------------
! show beta

case ('BETA')

   l1 = '|             Beta_x             |           Beta_y            |'
   d2_ptr => u%beta
   
   l2 = '|    Model    Design      Data  |  Model    Design      Data |'
   
   nl=nl+1; write (lines(nl), '(8x, a)') trim(l1)
   nl=nl+1; write (lines(nl), '(8x, a)') trim(l2)
   
   do i = lbound(u%phase%x%d, 1), ubound(u%phase%x%d, 1)
      if (d2_ptr%x%d(i)%exists) then
         x_dat => d2_ptr%x%d(i)
         y_dat => d2_ptr%y%d(i)
         ele => u%ring%ele(x_dat%ix_ele)
         nl=nl+1
         write (lines(nl), '(1x, a7, 6f10.4, i5)') ele%name, &
              x_dat%model, x_dat%design, x_dat%meas, &
              y_dat%model, y_dat%design, y_dat%meas, i
      endif
   enddo
   
   nl=nl+1; write (lines(nl), '(8x, a)') trim(l2)
   nl=nl+1; write (lines(nl), '(8x, a)') trim(l1)

!----------------------------------------------------------------------
! show calibration

case ('BPM_TILT')

  l1='        |      '
  l2='   Det  |   Tilt'  

  nl=nl+1; lines(nl) = ''
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(u%orbit%x%d, 1), ubound(u%orbit%x%d, 1)
    if (u%orbit%x%d(i)%exists) then
      ele => u%ring%ele(u%orbit%x%d(i)%ix_ele)

      nl=nl+1
      write (lines(nl), '(1x, a8, f10.4, f9.3, i5)') &
          ele%name, ele%value(tilt$), ele%s, i
    endif
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

!----------------------------------------------------------------------
! show calibration

case ('BPM_GAINS')

  call get_gain(gain)
  print *, 'Ix    Gain_1    Gain_2    Gain_3    Gain_4'
  do i = lbound(gain, 1), ubound(gain, 1)
    print '(i4, 4f10.6)', i, gain(i)%value
  enddo

!----------------------------------------------------------------------
! steerings

case ('CALIBRATION')


  L0 = '[dCU = Change in CU to make approx 1cm (@ beta = 20m) max amp ripple.]'
  L1 = "      | Current   dCU   Lower  Upper | z'(mrad)/              |"
  L2 = ' Name | Setting  (1cm)  Limit  Limit | 1000 CU   Beta    Phi  | Indx'
  L3 = '      |              [CU]            |                        |'

  nl=nl+1; lines(nl) = 'Horizontal:'
  nl=nl+1; lines(nl) = L0
  nl=nl+1; lines(nl) = L1
  nl=nl+1; lines(nl) = L2
  nl=nl+1; lines(nl) = L3

  if (f_phi > 10) then
    fmt = '(2x, a4, 3x, 4i7, f9.3, f8.2, f8.1, i6)'
  else
    fmt = '(2x, a4, 3x, 4i7, f9.3, f8.2, f8.3, i6)'
  endif

  factor = abs(sin(u%ring%a%tune/2))
  do i = lbound(u%hsteer_kick%v, 1), ubound(u%hsteer_kick%v, 1)
    vp => u%hsteer_kick%v(i)
    if (vp%good_var) then
      ix = vp%ix_ele
      call twiss_at_element (u%ring, ix, average = ave)
      dcu = 2 * 0.01 * factor / (sqrt(ave%a%beta * 20) * abs(vp%dvar_dcu))
      mr_model = 1000 * vp%model * sqrt(ave%a%beta)
      mr_db = 1000 * vp%cu_saved * vp%dvar_dcu * sqrt(ave%a%beta)
      nl = nl+1
      write (lines(nl), fmt) &
          vp%name(:4), vp%cu_saved, dcu, vp%cu_low_lim, &
          vp%cu_high_lim, 1.0e6*vp%dvar_dcu, ave%a%beta, f_phi*ave%a%phi, i
    endif
  enddo

  nl=nl+1; lines(nl) = L3
  nl=nl+1; lines(nl) = L2
  nl=nl+1; lines(nl) = L1
  nl=nl+1; lines(nl) = ' '
  nl=nl+1; lines(nl) = 'Vertical:'
  nl=nl+1; lines(nl) = L0
  nl=nl+1; lines(nl) = L1
  nl=nl+1; lines(nl) = L2
  nl=nl+1; lines(nl) = L3

  factor = abs(sin(u%ring%b%tune/2))
  do i = lbound(u%vsteer_kick%v, 1), ubound(u%vsteer_kick%v, 1)
    vp => u%vsteer_kick%v(i)
    if (vp%good_var) then
      ix = vp%ix_ele
      call twiss_at_element (u%ring, ix, average = ave)
      dcu = 2 * 0.01 * factor / (sqrt(ave%b%beta * 20) * abs(vp%dvar_dcu))
      mr_model = 1000 * vp%model * sqrt(ave%b%beta)
      mr_db = 1000 * vp%cu_saved * vp%dvar_dcu * sqrt(ave%b%beta)
      nl = nl+1
      write (lines(nl), fmt) &
          vp%name(:4), vp%cu_saved, dcu, vp%cu_low_lim, &
          vp%cu_high_lim, 1.0e6*vp%dvar_dcu, &
          ave%b%beta, f_phi*ave%b%phi, i
    endif
  enddo

  nl=nl+1; lines(nl) = L3
  nl=nl+1; lines(nl) = L2
  nl=nl+1; lines(nl) = L1

!----------------------------------------------------------------------
! show cbar

case ('CBAR11', 'CBAR12', 'CBAR22')
  
  fmt = '(1x, a7, 3f9.4, es10.2, l11, i5)'

  if (show_name == 'CBAR11') then
    d1_ptr => u%cbar%m11
  elseif (show_name == 'CBAR12') then
    d1_ptr => u%cbar%m12
  elseif (show_name == 'CBAR22') then
    d1_ptr => u%cbar%m22
  endif

  write (l1, '(9x, a, 14x, a)') '|', d1_ptr%name
  write (l2, '(9x, a)') &
    '|  Model   Design     Data    Weight  Useit_Opt'

  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(d1_ptr%d, 1), ubound(d1_ptr%d, 1)
    if (.not. d1_ptr%d(i)%exists) cycle
    d_ptr => d1_ptr%d(i)
    nl=nl+1; write (lines(nl), fmt) d_ptr%ele_name, &
          d_ptr%model, d_ptr%design, d_ptr%meas, d_ptr%weight, &
          d_ptr%useit_opt, i
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

!----------------------------------------------------------------------
! show chromaticity

case ('CHROMATICITY')

  fmt2 = '(1x, a8, 4f11.3, 2x, a)'
  nl=nl+1; write (lines(nl), '(16x, a)') 'Meas        Ref      Model     Design '
  nl=nl+1; write (lines(nl), fmt2) 'Chrom_X:', u%chrom%x%d(1)%meas, u%chrom%x%d(1)%ref, &
                                     u%chrom%x%d(1)%model, u%chrom%x%d(1)%design, '! dQ/(dE/E)'
  nl=nl+1; write (lines(nl), fmt2) 'Chrom_Y:', u%chrom%y%d(1)%meas, u%chrom%y%d(1)%ref, &
                                     u%chrom%y%d(1)%model, u%chrom%y%d(1)%design, '! dQ/(dE/E)'


!----------------------------------------------------------------------
! show comment

case ('COMMAND_FILES')

  call fullfilename ('.', directory)
  ok = dir_open (directory)

  print *, 'In Current directory:'
  do
    ok = dir_read(fname)
    if (.not. ok) exit
    if (.not. match_wild(fname, '*.in')) cycle
    print *, '  ', trim(fname)
  enddo

  print *
  print *, 'In $CESR_ONLINE/acc_control/program_info/cesrv/command_files :'

  call fullfilename ('$CESR_ONLINE/acc_control/program_info/cesrv/command_files', directory)
  ok = dir_open (directory)
  do
    ok = dir_read(fname)
    if (.not. ok) exit
    if (.not. match_wild(fname, '*.in')) cycle
    print *, '  ', trim(fname)
  enddo

!----------------------------------------------------------------------
! show comment

case ('COMMENTS')

  match_lattice = .false.
  match_route = .false.

  if (line(1:1) == '/') then
    ix = index(line, '=')
    if (ix < 3) then
      print *, 'ERROR: CONFUSED /LATTICE OR /ROUTE SWITCH'
      return
    endif
    line(ix:ix) = ' '
    call string_trim(line, line, ix)
    ix1 = index('/LATTICE', line(1:ix))
    ix2 = index('/ROUTE', line(1:ix))
    if (ix1 == 1) then
      match_lattice = .true.
    elseif (ix2 == 1) then
      match_route = .true.
    else
      print *, 'ERROR: OPTION NOT "/LATTICE" OR "/ROUTE"'
      return
    endif
    call string_trim (line(ix+1:), line, ix)
    match = line(:ix)
    call string_trim (line(ix+1:), line, ix_s2)
    show2_name = line(:ix_s2)
  endif

  line = line(ix_s2+1:)

  ix = index(line, ':')
  if (ix /= 0) line(ix:ix) = ' '
  call string_trim (line, line, ix)

  if (ix == 0) then
    print *, 'ERROR IN SHOW COMMAND SYNTAX.'
    return
  endif

  read (line(:ix), *, iostat = ios) ix1
  if (ios /= 0) then
    print *, 'ERROR READING FIRST FILE NUMBER'
    return
  endif

  call string_trim (line(ix+1:), line, ix)
  if (ix == 0) then
    ix2 = 0
  else
    read (line, *, iostat = ios) ix2
    if (ios /= 0) then
      print *, 'ERROR READING SECOND FILE NUMBER'
      return
    endif
  endif

  do_more = .true.

  if (index('ORBIT', trim(show2_name)) == 1) then
    do i = ix1, ix2
      call number_to_file_name (i, 'orbit', fname, ix_set, err_flag)
      if (err_flag) return
      call read_butns_file (fname, butns, u%db, err_flag, .false., &
                              logic%nonlinear_calc, logic%gain_correction)
      if (err_flag) cycle
      nl = nl+1
      if (nl == size(lines)) then
        write (lines(nl), *) 'ARRAY OVERFLOW. TERMINATE SHOW.'
        exit
      endif
      write (lines(nl), '(1x, i6, 2(2x, a))') ix_set, butns%date(:15), &
                                                    trim(butns%comment(1))
    enddo
    do_more = .false.
  elseif (index('PHASE', trim(show2_name)) == 1) then
    call fullfilename('$CESR_ONLINE/machine_data/mach_meas/phase/phase.number',fname)
    stub = 'PHASE'
  elseif (index('MODEL', trim(show2_name)) == 1) then
    call fullfilename('$CESR_ONLINE/acc_control/program_info/cesrv/model/model.number',fname)
    stub = 'MODEL'
  elseif (index('ETA', trim(show2_name)) == 1) then
    call fullfilename('$CESR_ONLINE/machine_data/mach_meas/eta/eta.number',fname)
    stub = 'ETA'
  elseif (index('AC_ETA', trim(show2_name)) == 1) then
    call fullfilename('$CESR_ONLINE/machine_data/mach_meas/ac_eta/ac_eta.number',fname)
    stub = 'AC_ETA'
  else
    print *, 'ERROR: SHOW "PHASE", "MODEL", "ETA", OR "ORBIT" COMMENTS?'
    return
  endif

  if (do_more) then
    call calc_file_number (fname, ix1, ix1, err_flag)
    call calc_file_number (fname, ix2, ix2, err_flag)
    do i = ix1, ix2
      call form_file_name_with_number (trim(stub), i, fname, err_flag)
      if (err_flag) return
      open (1, file = fname, status = 'old', action = 'read', iostat = ios)
      if (ios /= 0) cycle
      read (1, nml = data_parameters, iostat = ios)
      close (1)
      if (match_lattice .and. lattice /= match) cycle
      if (match_route .and. route_name /= match) cycle
      nl = nl + 1
      if (nl == size(lines)) then
        write (lines(nl), *) 'ARRAY OVERFLOW. TERMINATE SHOW.'
        exit
      endif
      if (ios /= 0) then
        write (lines(nl), '(1x, i6, 2x, a)') i, '--- CANNOT READ HEADER ---'
      else
        write (lines(nl), '(1x, i6, 2(2x, a))') i, &
                                          data_date(3:17), trim(comment)
      endif
    enddo
  endif

!----------------------------------------------------------------------
! Coupling typeout

case ('COUPLING')

  l1 = '        |        Cbar_model     |        Cbar_data      | '
  l2 = '  Det   |    11      12      22 |    11      12      22 | Ix Used'

  nl=nl+1; lines(nl) = ' '
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(u%orbit%x%d, 1), ubound(u%orbit%x%d, 1)
    if (u%orbit%x%d(i)%exists) then
      nl=nl+1
      ix = u%orbit%x%d(i)%ix_ele
      write (lines(nl), '(1x, a7, 6f8.3, i5, l5)') u%ring%ele(ix)%name, &
          u%cbar%m11%d(i)%model, u%cbar%m12%d(i)%model, &
          u%cbar%m22%d(i)%model, &
          u%cbar%m11%d(i)%meas, u%cbar%m12%d(i)%meas, &
          u%cbar%m22%d(i)%meas, i, u%cbar%m12%d(i)%useit_plot
    endif
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

!----------------------------------------------------------------------
! show variables

case ('CUSTOM')

  a_max = 1.1
  do i = 1, n_custom_maxx
    if (.not. u%custom_var%v(i)%good_var) cycle
    a_max = max(a_max, abs(u%custom_var%v(i)%model))
    a_max = max(a_max, abs(u%custom_var%v(i)%saved))
  enddo
  n = max(0, 5 - int(log10(a_max)))

  nl=nl+1
  write (lines(nl), '(4x, a, 14x, a)') 'Ix  Name', &
                       'Attribute         Ix_ele     Model  Used'

  do i = 1, n_custom_maxx
    if (.not. u%custom_var%v(i)%good_var) cycle
    vp => u%custom_var%v(i)
    nl = nl+1
    write (lines(nl), '(i6, 2(2x, a16), i8, f10.<n>, l6)') i, &
          u%ring%ele(vp%ix_ele)%name, vp%attrib_name, vp%ix_ele, &
          vp%model, vp%useit
  enddo

  nl=nl+1; lines(nl) = ' '

!----------------------------------------------------------------------
! show init_orb

case ('CUT_RING')

  ix1 = logic%ix_cut_ring
  if (ix1 == -1) then
    write (lines(1), *) 'Ring has *NOT* been cut'
    nl = 1
  else
    write (lines(1), '(2a)') ' Ring cut at element: ', u%ring%ele(ix1)%name
    write (lines(2), *)
    do i = 1, 6
      write (lines(i+2), '(2x, a, i3, 3a, f10.6)') 'Init_Orb', i, &
              '  "',trim(u%init_orb%v(i)%name), '"', u%init_orb%v(i)%model
    enddo
    nl = 8
  endif

  nl=nl+1; lines(nl) = ''

  call type_twiss (u%ring%ele(ix1), cycles$, lines = lines(nl+1:), n_lines = n)
  nl = nl + n

!----------------------------------------------------------------------

case ('DATA')

  call match_word (show2_name, name$%data_type_name, ix1)

  if (show2_name == 'CBAR11') then
    d1_ptr => u%cbar%m11
    ix1 = cbar_data$
  elseif (show2_name == 'CBAR12') then
    d1_ptr => u%cbar%m12
    ix1 = cbar_data$
  elseif (show2_name == 'CBAR21') then
    d1_ptr => u%cbar%m21
    ix1 = cbar_data$
  elseif (show2_name == 'CBAR22') then
    d1_ptr => u%cbar%m22
    ix1 = cbar_data$
  elseif (show2_name == 'CMAT_A22') then
    d1_ptr => u%cmat_a%m22
    ix1 = cmat_a_data$
  elseif (show2_name == 'CMAT_A12') then
    d1_ptr => u%cmat_a%m12
    ix1 = cmat_a_data$
  elseif (show2_name == 'CMAT_B21') then
    d1_ptr => u%cmat_b%m21
    ix1 = cmat_b_data$
  elseif (show2_name == 'CMAT_B22') then
    d1_ptr => u%cmat_b%m22
    ix1 = cmat_b_data$
  elseif (ix1 == cbar_data$) then
    d1_ptr => u%cbar%m12
  elseif (ix1 == beta_data$) then
    d2_ptr => u%beta
  elseif (ix1 == orbit_data$) then
    d2_ptr => u%orbit
  elseif (ix1 == phase_data$) then
    d2_ptr => u%phase
  elseif (ix1 == eta_data$) then
    d2_ptr => u%eta
  elseif (ix1 == ac_eta_data$) then
    d2_ptr => u%ac_eta
  elseif (ix1 == mode_eta_data$) then
    d2_ptr => u%mode_eta
  elseif (ix1 == spline_data$) then
    d2_ptr => u%spline_beta
  elseif (ix1 == energy_data$) then
    d1_ptr => u%energy_data%d1
    d_ptr => d1_ptr%d(1)
  elseif (ix1 == tune_data$) then
    d2_ptr => u%tune
  elseif (ix1 == chrom_data$) then
    d2_ptr => u%chrom
  elseif (ix1 == q2x_data$) then
    d2_ptr => u%q2x
  elseif (ix1 == q2y_data$) then
    d2_ptr => u%q2y
  elseif (ix1 == qx_plus_qy_data$) then
    d2_ptr => u%qx_plus_qy
  elseif (ix1 == qx_minus_qy_data$) then
    d2_ptr => u%qx_minus_qy
  elseif (ix1 == e_xray_data$) then
    d2_ptr => u%e_xray
  elseif (ix1 == p_xray_data$) then
    d2_ptr => u%p_xray
  else
    print *, 'ERROR: WHAT TYPE OF DATA TO SHOW? (ORBIT, PHASE, ..._rp)'
    return
  endif

  call string_trim (line(ix_s2+1:), line, ix)

  if (ix1 /= cbar_data$ .and. ix1 /= energy_data$ .and. show2_name(1:4) /= 'CMAT') then
    call match_word (line, name$%plane, ix2)
    if (ix2 == x_plane$) then
      d1_ptr => d2_ptr%x
    elseif (ix2 == y_plane$) then
      d1_ptr => d2_ptr%y
    elseif (ix2 == z_plane$ .and. ix1 == tune_data$) then
      d1_ptr => d2_ptr%z
    elseif (ix2 == in_plane$) then
      d1_ptr => d2_ptr%a_in
    elseif (ix2 == out_plane$) then
      d1_ptr => d2_ptr%a_out
    else
      print *, 'ERROR: WHICH PLANE (X, Y, ...)?'
      return
    endif
    call string_trim (line(ix+1:), line, ix)
  endif

  if (ix1 == energy_data$) then
    i = 1
  else
    call cesr_locator (line, prefix, i, err_flag)
    if (prefix /= '' .or. err_flag .or. i < lbound(d1_ptr%d, 1) .or. i > ubound(d1_ptr%d, 1)) then
      print *, 'ERROR: CANNOT READ DATA NUMBER.'
      return
    endif
    d_ptr => d1_ptr%d(i)
  endif

  nl=nl+1; write (lines(nl), *) 'Index in %data array:', d1_ptr%ix_data + i
  nl=nl+1; write (lines(nl), *) 'Name:          ', d_ptr%name
  nl=nl+1; write (lines(nl), *) 'Ele_name:      ', d_ptr%ele_name
  nl=nl+1; write (lines(nl), *) 'Alias:         ', d_ptr%alias
  nl=nl+1; write (lines(nl), *) 'Ix_index:      ', d_ptr%ix_index
  nl=nl+1; write (lines(nl), '(1x, a, f15.2)') &
                                'ix_plot_index: ', d_ptr%ix_plot_index
  nl=nl+1; write (lines(nl), *) 'Ix_ele:        ', d_ptr%ix_ele
  nl=nl+1; write (lines(nl), *) 'Ix_dmeas:      ', d_ptr%ix_dmeas
  nl=nl+1; write (lines(nl), *) 'Meas:          ', d_ptr%meas
  nl=nl+1; write (lines(nl), *) 'Ref:           ', d_ptr%ref
  nl=nl+1; write (lines(nl), *) 'Model:         ', d_ptr%model
  nl=nl+1; write (lines(nl), *) 'Base_model:    ', d_ptr%base_model
  nl=nl+1; write (lines(nl), *) 'Delta:         ', d_ptr%delta
  nl=nl+1; write (lines(nl), *) 'Design:        ', d_ptr%design
  nl=nl+1; write (lines(nl), *) 'Design_nonlin: ', d_ptr%design_nonlin
  nl=nl+1; write (lines(nl), *) 'Old:           ', d_ptr%old
  nl=nl+1; write (lines(nl), *) 'Fit:           ', d_ptr%fit
  nl=nl+1; write (lines(nl), *) 'Merit:         ', d_ptr%merit
  nl=nl+1; write (lines(nl), *) 'Weight:        ', d_ptr%weight
  nl=nl+1; write (lines(nl), *) 'S:             ', d_ptr%s
  nl=nl+1; write (lines(nl), *) 'Exists:        ', d_ptr%exists
  nl=nl+1; write (lines(nl), *) 'Good_dat:      ', d_ptr%good_dat
  nl=nl+1; write (lines(nl), *) 'Good_ref:      ', d_ptr%good_ref
  nl=nl+1; write (lines(nl), *) 'Good_user:     ', d_ptr%good_user
  nl=nl+1; write (lines(nl), *) 'Good_opt:      ', d_ptr%good_opt
  nl=nl+1; write (lines(nl), *) 'Useit_plot:    ', d_ptr%useit_plot
  nl=nl+1; write (lines(nl), *) 'Useit_opt:     ', d_ptr%useit_opt

!----------------------------------------------------------------------
case ('DEFAULT')

  opt_vars = logic%opt_vars

  if (opt_vars == opt_steering$ .or. opt_vars == opt_all_vars$ .or. &
             opt_vars == opt_skew_quad$) then
    call type_energy_shift (u, lines, nl)
  endif

  call merit_calc(merit0)

  fmt = '(1x, a, 2es11.2)'
  nl=nl+1; write (lines(nl), *)
  nl=nl+1; write (lines(nl), *) 'Optimizing Using: ', name$%opt_vars(opt_vars)
  nl=nl+1; write (lines(nl), *) 'Data to Optimize: ', name$%opt_data(logic%opt_data)
  nl=nl+1; write (lines(nl), *) 'Optimize Engine:  ', logic%engine
  nl=nl+1; write (lines(nl), *)
  nl=nl+1; write (lines(nl), *)   'Opt_loops     =', logic%opt_loops
  nl=nl+1; write (lines(nl), *)   'Opt_cycles    =', logic%opt_cycles
  nl=nl+1; write (lines(nl), fmt) 'Change_min    =', logic%change_min
  nl=nl+1; write (lines(nl), fmt) 'Opt_tolerance =', logic%opt_tolerance
  nl=nl+1; write (lines(nl), *)   'Merit:', merit0

! data wgt info

  nl=nl+1; write (lines(nl), *)   ''

  if (any(u%beta%x%d%good_opt) .or. any(u%phase%x%d%good_opt) .or. any(u%tune%x%d%good_opt) .or. &
      any(u%orbit%x%d%good_opt) .or. any(u%eta%x%d%good_opt) .or. any(u%ac_eta%x%d%good_opt) .or. &
      any(u%mode_eta%x%d%good_opt)) then
    nl=nl+1; write (lines(nl), *)   '                        X          Y'
  endif

  call this_weight_out (any(u%beta%x%d%good_opt),     'Beta_*_wgt     =', u%beta%x%d(1)%weight, u%beta%y%d(1)%weight)
  call this_weight_out (any(u%phase%x%d%good_opt),    'Phase_*_wgt    =', u%phase%x%d(1)%weight, u%phase%y%d(1)%weight)
  call this_weight_out (any(u%tune%x%d%good_opt),     'Tune_*_wgt     =', u%tune%x%d(1)%weight, u%tune%y%d(1)%weight)
  call this_weight_out (any(u%orbit%x%d%good_opt),    'Orbit_*_wgt    =', u%orbit%x%d(1)%weight, u%orbit%y%d(1)%weight)
  call this_weight_out (any(u%eta%x%d%good_opt),      'Eta_*_wgt      =', u%eta%x%d(1)%weight, u%eta%y%d(1)%weight)
  call this_weight_out (any(u%ac_eta%x%d%good_opt),   'AC_Eta_*_wgt   =', u%ac_eta%x%d(1)%weight, u%ac_eta%y%d(1)%weight)
  call this_weight_out (any(u%ac_eta%c12_a%d%good_opt),   'AC_Eta_c12_*_wgt   =', u%ac_eta%c12_a%d(1)%weight, u%ac_eta%c12_b%d(1)%weight)
  call this_weight_out (any(u%ac_eta%yxcos%d%good_opt),   'AC_Eta_xy_*_wgt   =', u%ac_eta%yxcos%d(1)%weight, u%ac_eta%yxsin%d(1)%weight)
  call this_weight_out (any(u%mode_eta%x%d%good_opt), 'Mode_Eta_*_wgt =', u%mode_eta%x%d(1)%weight, u%mode_eta%y%d(1)%weight)
  call this_weight_out (any(u%cbar%m12%d%good_opt),   'Cbar12_wgt     =    ', u%cbar%m12%d(1)%weight)
  call this_weight_out (any(u%cbar%m11%d%good_opt),   'Cbar11_wgt     =    ', u%cbar%m11%d(1)%weight)
  call this_weight_out (any(u%cbar%m22%d%good_opt),   'Cbar22_wgt     =    ', u%cbar%m22%d(1)%weight)
  call this_weight_out (any(u%cmat_a%m22%d%good_opt), 'Cmat_a22_wgt   =    ', u%cmat_a%m22%d(1)%weight)
  call this_weight_out (any(u%cmat_a%m12%d%good_opt), 'Cmat_a12_wgt   =    ', u%cmat_a%m12%d(1)%weight)
  call this_weight_out (any(u%cmat_b%m11%d%good_opt), 'Cmat_b11_wgt   =    ', u%cmat_b%m11%d(1)%weight)
  call this_weight_out (any(u%cmat_b%m12%d%good_opt), 'Cmat_b12_wgt   =    ', u%cmat_b%m12%d(1)%weight)
  call this_weight_out (any(u%chrom%x%d%good_opt),    'Chrom_wgt      =    ', u%chrom%x%d(1)%weight)
  call this_weight_out (any(u%energy_data%d1%d%good_opt), 'Energy_wgt     =    ', u%energy_data%d1%d(1)%weight)

  call this_weight_out (any(u%e_xray%y%d%good_opt),          'E_XRay_y_wgt   = ', u%e_xray%y%d(1)%weight)
  call this_weight_out (any(u%p_xray%y%d%good_opt),          'P_XRay_y_wgt   = ', u%p_xray%y%d(1)%weight)

  call this_weight_out (any(u%q2x%a_in%d%good_opt),          '2qx_in_wgt     = ', u%q2x%a_in%d(1)%weight)
  call this_weight_out (any(u%q2x%a_out%d%good_opt),         '2qx_out_wgt    = ', u%q2x%a_out%d(1)%weight)
  call this_weight_out (any(u%q2y%a_in%d%good_opt),          '2qy_in_wgt     = ', u%q2y%a_in%d(1)%weight)
  call this_weight_out (any(u%q2y%a_out%d%good_opt),         '2qy_out_wgt    = ', u%q2y%a_out%d(1)%weight)
  call this_weight_out (any(u%qx_plus_qy%a_in%d%good_opt),   'qx+qy_in_wgt   = ', u%qx_plus_qy%a_in%d(1)%weight)
  call this_weight_out (any(u%qx_plus_qy%a_out%d%good_opt),  'qx+qy_out_wgt  = ', u%qx_plus_qy%a_out%d(1)%weight)
  call this_weight_out (any(u%qx_minus_qy%a_in%d%good_opt),  'qx-qy_in_wgt   = ', u%qx_minus_qy%a_in%d(1)%weight)
  call this_weight_out (any(u%qx_minus_qy%a_out%d%good_opt), 'qx-qy_out_wgt  = ', u%qx_minus_qy%a_out%d(1)%weight)



  ! Var wgt info

  nl=nl+1; write (lines(nl), *) ' '

  call this_weight_out (any(u%quad_k1%v%good_opt),      'K1_wgt          =', u%quad_k1%v(1)%weight)
  call this_weight_out (any(u%skew_quad_k1%v%good_opt), 'Skew_wgt        =', u%skew_quad_k1%v(1)%weight)
  call this_weight_out (any(u%bpm_tilt%v%good_opt),     'bpm_tilt_wgt    =', u%bpm_tilt%v(1)%weight)

  call this_weight_out (any(u%hsteer_kick%v%good_opt), 'Horizontal_wgt  =', &
            sum(u%hsteer_kick%v%weight * u%hsteer_kick%v%dvar_dcu**2) / &
            count (u%hsteer_kick%v%dvar_dcu /= 0))
  call this_weight_out (any(u%vsteer_kick%v%good_opt), 'Vertical_wgt    =', &
            sum(u%vsteer_kick%v%weight * u%vsteer_kick%v%dvar_dcu**2) / &
            count (u%vsteer_kick%v%dvar_dcu /= 0))

  call this_weight_out (any(u%hsep_kick%v%good_opt),   'H_sep_wgt       =', u%hsep_kick%v(1)%weight)
  call this_weight_out (any(u%sex_k2%v%good_opt),      'K2_wgt          =', u%sex_k2%v(1)%weight)
  call this_weight_out (any(u%skew_sex_k2%v%good_opt), 'Skew_k2_wgt     =', u%skew_sex_k2%v(1)%weight)
  call this_weight_out (any(u%custom_var%v%good_opt),  'Var_wgt         =', u%custom_var%v(1)%weight)
  call this_weight_out (any(u%x_kick_quad%v%good_opt), 'x_kick_quad_wgt =', u%x_kick_quad%v(1)%weight)
  call this_weight_out (any(u%y_kick_quad%v%good_opt), 'y_kick_quad_wgt =', u%y_kick_quad%v(1)%weight)

! data used

  nl=nl+1; write (lines(nl), *)

  call data_type_useit (u%phase, lines, nl)
  call data_type_useit (u%orbit, lines, nl)
  call data_type_useit (u%e_xray, lines, nl)
  call data_type_useit (u%p_xray, lines, nl)
  call data_type_useit (u%cbar, lines, nl)
  call data_type_useit (u%beta, lines, nl)
  call data_type_useit (u%chrom, lines, nl)
  call data_type_useit (u%eta, lines, nl)
  call data_type_useit (u%mode_eta, lines, nl)
  call data_type_useit (u%q2x, lines, nl)
  call data_type_useit (u%q2y, lines, nl)
  call data_type_useit (u%qx_plus_qy, lines, nl)
  call data_type_useit (u%qx_minus_qy, lines, nl)

  ! optimization variables used

  call show_opt_vars

!----------------------------------------------------------------------
! Det typeout

case ('DETECTORS')

  l1='          |       Phi_model   |  Beta_model   | Orbit_model(mm)'
  l2='   Det    |      X         Y  |    X       Y  |    X       Y  |      Z  Det'
  nl=nl+1; write (lines(nl), *)
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l2

  do i = lbound(u%orbit%x%d, 1), ubound(u%orbit%x%d, 1)
    if (u%orbit%x%d(i)%exists) then
      ele => u%ring%ele(u%orbit%x%d(i)%ix_ele)

      nl=nl+1
      write (lines(nl), '(1x, a9, 2f10.4, 4f8.2, f9.3, i5)') &
          ele%name, f_phi*ele%a%phi, f_phi*ele%b%phi, ele%a%beta, ele%b%beta, &
          1000*u%orbit%x%d(i)%model, 1000*u%orbit%y%d(i)%model, ele%s, i
    endif
  enddo

  nl=nl+1; write (lines(nl), *) l2
  nl=nl+1; write (lines(nl), *) l1

!----------------------------------------------------------------------
! show dmerit

case ('DMERIT')

  if (line == '') then
    ix1 = 0
    ix2 = 0
  else
    read (line, *, iostat = ios) ix1, ix2
    if (ios /= 0) then
      print *, 'ERROR READING INDEXES.'
      return
    endif
  endif

  if (ix1 < 1 .and. ix2 < 1) then
    do i = 1, size(u%data)
      if (nl > 50) exit
      if (.not. u%data(i)%useit_opt) cycle
      ix1 = u%data(i)%ix_dmeas
      do j = 1, size(u%var)
        if (nl > 50) exit
        if (.not. u%var(j)%useit) cycle
        ix2 = u%var(j)%ix_dvar
        nl=nl+1; write (lines(nl), '(2(a, i3), a, es10.2, 2x, 2(a, 1x, i0, 5x))') &
          'dMerit(', ix1, ',', ix2, ') =', u%dm_dv(ix1, ix2), &
           trim(u%data(i)%name), u%data(i)%ix_index, trim(u%var(j)%name)
      enddo
    enddo

  elseif (ix1 < 1) then
    do i = 1, ubound(u%data, 1)
      if (.not. u%data(i)%useit_opt) cycle
      ix = u%data(i)%ix_dmeas
      nl=nl+1; write (lines(nl), '(2(a, i3), a, es10.2, 2x, a, 1x, i0)') &
          'dMerit(', ix, ',', ix2, ') =', u%dm_dv(ix, ix2), trim(u%data(i)%name), u%data(i)%ix_index
    enddo

  elseif (ix2 < 1) then
 
    do i = 1, size(u%var)
      if (.not. u%var(i)%useit) cycle
      ix = u%var(i)%ix_dvar
      nl=nl+1; write (lines(nl), '(2(a, i3), a, es10.2, 2x, a, 1x, i0)') &
          'dMerit(', ix1, ',', ix, ') =', u%dm_dv(ix1, ix), trim(u%var(i)%name)
    enddo

  else
    nl=nl+1; write (lines(nl), '(2(a, i3), a, es10.2, 2x, a, 1x, i0)') &
              'dMerit(', ix1, ',', ix2, ') =', &
              (super%u_(j)%dm_dv(ix1, ix2), j = 1, size(super%u_))
  endif

!----------------------------------------------------------------------
! show element

case ('ELEMENT')

  show_all = .false.
  if (line(:ix_s2) == '-ALL') then
    show_all = .true.
    call string_trim (line(ix_s2+1:), line, ix_s2)
    show2_name = line(1:ix_s2)
  endif

  name_found = .false.

  if (ix_s2 == 0) then
    print *, 'ERROR: PLEASE GIVE AN ELEMENT NAME OR NUMBER'
    return
  endif

  if (index(line(:ix_s2), '*') /= 0 .or. index(line(:ix_s2), '%') /= 0) then
    write (lines(1), *) 'Matches to name:'
    nl = 1
    do loc = 1, u%ring%n_ele_max
      if (match_wild(u%ring%ele(loc)%name, line(:ix_s2))) then
        nl = nl + 1
        write (lines(nl), '(i8, 2x, a)') loc, u%ring%ele(loc)%name
        name_found = .true.
      endif
    enddo
    if (.not. name_found) then
      nl = nl + 1
      write (lines(nl), *) '   *** No Matches to Name Found ***'
    endif
    goto 8000
  endif

  call locate_element (u%ring, line, loc, err_flag)
  if (err_flag) return

  nl=nl+1; write (lines(nl), *) 'Element #', loc

  if (show_all) then
    call type_ele (u%ring%ele(loc), .true., 6, .true., logic%phase_units, .true., lines = alloc_lines, n_lines = n)
  else
    call type_ele (u%ring%ele(loc), .true., 6, .false., logic%phase_units, &
                                          .true., .true., .true., .true., .true., lines = alloc_lines, n_lines = n)
  endif

  lines(nl+1:nl+n) = alloc_lines(1:n)
  nl = nl + n

  orb = u%orb(loc)
  fmt = '(2x, a, 3p2f15.8)'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *)   'Orbit:           [mm]         [mrad]'
  nl=nl+1; write (lines(nl), fmt) "X  X':", orb%vec(1), orb%vec(2)
  nl=nl+1; write (lines(nl), fmt) "Y  Y':", orb%vec(3), orb%vec(4)
  nl=nl+1; write (lines(nl), fmt) "Z  Z':", orb%vec(5), orb%vec(6)

  found = .false.
  do i = loc + 1, u%ring%n_ele_max
    if (u%ring%ele(i)%name == line(:ix_s2)) then
      if (found) then
        nl = nl + 1
        write (lines(nl), *)
        found = .true.
      endif
      nl = nl + 1
      write (lines(nl), *) 'Note: Found another element with same name at:', i
    endif
  enddo

!----------------------------------------------------------------------

case ('ENERGY_SHIFT')

  call type_energy_shift (u, lines, nl)

!----------------------------------------------------------------------
! show eta

case ('ETA', 'MODE_ETA')

  if (show_name == 'ETA') then
    l1 = '|             Eta_X             |           Eta_Y            |'
    d2_ptr => u%eta
  else
    l1 = '|        Mode_Eta_X             |      Mode_Eta_Y            |'
    d2_ptr => u%mode_eta
  endif

  l2 = '|    Model    Design      Data  |  Model    Design      Data |'

  nl=nl+1; write (lines(nl), '(8x, a)') trim(l1)
  nl=nl+1; write (lines(nl), '(8x, a)') trim(l2)

  do i = lbound(u%phase%x%d, 1), ubound(u%phase%x%d, 1)
    if (d2_ptr%x%d(i)%exists) then
      x_dat => d2_ptr%x%d(i)
      y_dat => d2_ptr%y%d(i)
      ele => u%ring%ele(x_dat%ix_ele)
      nl=nl+1
      write (lines(nl), '(1x, a7, 6f10.4, i5)') ele%name, &
          x_dat%model, x_dat%design, x_dat%meas, &
          y_dat%model, y_dat%design, y_dat%meas, i
    endif
  enddo

  nl=nl+1; write (lines(nl), '(8x, a)') trim(l2)
  nl=nl+1; write (lines(nl), '(8x, a)') trim(l1)

!----------------------------------------------------------------------

case ('GEOMETRY', 'FLOOR')

  call lat_geometry (u%ring)
  call find_range_endpoints (line, loc, ix1, ix2, err_flag)
  if (err_flag) return

  nl=nl+1; write (lines(nl), '(36x, a)') 'Geometry at End of Element'
  nl=nl+1; write (lines(nl), '(5x, a, 18x, a)') 'Element',  'S         Len        Z          X                Theta'
  nl=nl+1; write (lines(nl), '(30x, a)')                    '                                           (rad)       (deg)'

  do ix = ix1, ix2
    ele => u%ring%ele(ix)
    nl = nl + 1
    write (lines(nl), '(i4, 1x, a20, 2f10.5, 3f12.7, f12.5)') ix, &
            ele%name, ele%s, ele%value(l$), &
            ele%floor%r(3), ele%floor%r(1), ele%floor%theta, ele%floor%theta * 180 / pi
  enddo

!----------------------------------------------------------------------
! Global

case ('GLOBAL')

  call type_energy_shift (u, lines, nl)

  l_ring = u%ring%param%total_length
  nl=nl+1; write (lines(nl), *) ''
  nl=nl+1; write (lines(nl), *) 'Lattice:          ', u%ring%lattice
  nl=nl+1; write (lines(nl), *) 'Lattice File:     ', trim(u%ring%input_file_name)
  nl=nl+1; write (lines(nl), *) 'Circumference:    ', l_ring
  nl=nl+1; write (lines(nl), *) 'Energy (GeV):     ', u%ring%ele(0)%value(E_TOT$)/1e9
  nl=nl+1; write (lines(nl), *) 'Species:          ', particle_name(u%ring%param%particle)
  nl=nl+1; write (lines(nl), *) 'Reverse_tracking: ', logic%reverse_tracking
  nl=nl+1; write (lines(nl), *) 'Radiation_on:     ', logic%radiation_on
  nl=nl+1; write (lines(nl), *) 'RF_on:            ', on_off(logic%rf_on)
  nl=nl+1; write (lines(nl), *) '# Particles in a bunch:   ', u%ring%param%n_part
  nl=nl+1; write (lines(nl), *) '# Tracking ring elements: ', u%ring%n_ele_track
  nl=nl+1; write (lines(nl), *) 'Total # of elements:      ', u%ring%n_ele_max

  call set_on_off (rfcavity$, u%ring, on$, u%orb)
  call radiation_integrals (u%ring, u%orb, global_model)
  call set_on_off (rfcavity$, u%ring, on_off_int(logic%rf_on), u%orb)
  
  call chrom_calc (u%ring, 1.0e-4_rp, global_model%a%chrom, global_model%b%chrom)

  nl=nl+1; write (lines(nl), *)
  nl=nl+1; write (lines(nl), '(17x, a)') '       X          |            Y'
  nl=nl+1; write (lines(nl), '(17x, a)') 'Model     Design  |     Model     Design'
  fmt = '(1x, a10, 1p 2e11.3, 2x, 2e11.3, 2x, a)'
  fmt2 = '(1x, a10, 2f11.3, 2x, 2f11.3, 2x, a)'
  fmt3 = '(1x, a10, 2f11.4, 2x, 2f11.4, 2x, a)'
  f_phi = 1 / twopi
  nl=nl+1; write (lines(nl), fmt2) 'Q', f_phi*u%tune%x%d(1)%model, &
          f_phi*u%tune%x%d(1)%design, f_phi*u%tune%y%d(1)%model, &
          f_phi*u%tune%y%d(1)%design, '! Tune'
  nl=nl+1; write (lines(nl), fmt2) 'Chrom', global_model%a%chrom, &
        u%global_design%a%chrom, global_model%b%chrom, &
        u%global_design%b%chrom, '! dQ/(dE/E)'
  nl=nl+1; write (lines(nl), fmt2) 'J_damp', global_model%a%j_damp, &
        u%global_design%a%j_damp, global_model%b%j_damp, &
        u%global_design%b%j_damp, '! Damping Partition #'
  nl=nl+1; write (lines(nl), fmt) 'Emittance', global_model%a%emittance, &
        u%global_design%a%emittance, global_model%b%emittance, &
        u%global_design%b%emittance, '! Meters'
  nl=nl+1; write (lines(nl), fmt) 'Alpha_damp', global_model%a%alpha_damp, &
        u%global_design%a%alpha_damp, global_model%b%alpha_damp, &
        u%global_design%b%alpha_damp, '! Damping per turn'
  nl=nl+1; write (lines(nl), fmt) 'I4', global_model%a%synch_int(4), &
        u%global_design%a%synch_int(4), global_model%b%synch_int(4), &
        u%global_design%b%synch_int(4), '! Radiation Integral'
  nl=nl+1; write (lines(nl), fmt) 'I5', global_model%a%synch_int(5), &
        u%global_design%a%synch_int(5), global_model%b%synch_int(5), &
        u%global_design%b%synch_int(5), '! Radiation Integral'

  nl=nl+1; write (lines(nl), *)
  nl=nl+1; write (lines(nl), '(19x, a)') 'Model     Design'
  fmt = '(1x, a12, 2es11.3, 3x, a)'
  nl=nl+1; write (lines(nl), fmt) 'Sig_E/E:', &
              global_model%sigE_E, u%global_design%sigE_E, '! Energy spread'
  nl=nl+1; write (lines(nl), fmt) 'Sig_z:', &
              global_model%sig_z, u%global_design%sig_z, '! Bunch length'
  nl=nl+1; write (lines(nl), fmt) 'Voltage:', &
              global_model%rf_voltage, u%global_design%rf_voltage, '! Total RF voltage (Volt)'
  nl=nl+1; write (lines(nl), fmt) 'Energy Loss:', global_model%e_loss, &
              u%global_design%e_loss, '! eV / Turn'
  nl=nl+1; write (lines(nl), fmt) 'J_damp:', global_model%z%j_damp, &
              u%global_design%z%j_damp, '! Longitudinal Damping Partition #'
  nl=nl+1; write (lines(nl), fmt) 'Alpha_damp:', global_model%z%alpha_damp, &
              u%global_design%z%alpha_damp, '! Longitudinal Damping per turn'
  nl=nl+1; write (lines(nl), fmt) 'Alpha_p:', global_model%synch_int(1)/l_ring, &
              u%global_design%synch_int(1)/l_ring, '! Momentum Compaction'
  nl=nl+1; write (lines(nl), fmt) 'I1:', global_model%synch_int(1), &
              u%global_design%synch_int(1), '! Radiation Integral'
  nl=nl+1; write (lines(nl), fmt) 'I2:', global_model%synch_int(2), &
              u%global_design%synch_int(2), '! Radiation Integral'
  nl=nl+1; write (lines(nl), fmt) 'I3:', global_model%synch_int(3), &
              u%global_design%synch_int(3), '! Radiation Integral'
  nl=nl+1; write (lines(nl), *) 'Note: Radiation integrals are calculated with RF on'

!----------------------------------------------------------------------

case ('GROUP')

  call match_var_type ('GROUP ' // line, var1, u, err_flag)
  if (err_flag) return

  print *, '            Name     Model     Saved'
  do i = lbound(var1%v, 1), ubound(var1%v, 1)
    if (.not. var1%v(i)%exists) cycle
    print '(1x, a, f10.1, i10)', var1%v(i)%name, &
                                     var1%v(i)%model, var1%v(i)%cu_saved
  enddo
  print *, '            Name     Model     Saved'

!----------------------------------------------------------------------
! Horizontal

case ('HORIZONTAL')
  call show_quad_etc (u%hsteer_kick, 'Horizontal', 'K.')

!----------------------------------------------------------------------
!  LIGHTSOURCES

case ('LIGHTSOURCES')

  call show_lightsources( u%ring, lines, nl )
    
!----------------------------------------------------------------------
! Logicals

case ('LOGICALS')

  nl=nl+1; write (lines(nl), amt) 'logic%btns_gain_correction_file  = ', trim(logic%btns_gain_correction_file)
  nl=nl+1; write (lines(nl), amt) 'logic%engine                     = ', trim(logic%engine)
  nl=nl+1; write (lines(nl), amt) 'logic%lattice                    = ', trim(logic%lattice)
  nl=nl+1; write (lines(nl), amt) 'logic%plot_what                  =' , trim(logic%plot_what)
  nl=nl+1; write (lines(nl), amt) 'logic%queue                      = ', trim(logic%queue)
  nl=nl+1; write (lines(nl), amt) 'logic%tilt_correction_file       = ', trim(logic%tilt_correction_file)

  nl=nl+1; write (lines(nl), imt) 'logic%ac_eta_freq_range          =', logic%ac_eta_freq_range
  nl=nl+1; write (lines(nl), imt) 'logic%beam_bunch                 =', logic%beam_bunch
  nl=nl+1; write (lines(nl), imt) 'logic%beam_species               =', logic%beam_species
  nl=nl+1; write (lines(nl), imt) 'logic%biggrp_set                 =', logic%biggrp_set
  nl=nl+1; write (lines(nl), imt) 'logic%bpm_calib_max_cycles       =', logic%bpm_calib_max_cycles
  nl=nl+1; write (lines(nl), imt) 'logic%csr_set                    =', logic%csr_set
  nl=nl+1; write (lines(nl), imt) 'logic%iu_remember                =', logic%iu_remember
  nl=nl+1; write (lines(nl), imt) 'logic%ix_cut_ring                =', logic%ix_cut_ring
  nl=nl+1; write (lines(nl), imt) 'logic%ix_bpm_fft                 =', logic%ix_bpm_fft
  nl=nl+1; write (lines(nl), imt) 'logic%last_read                  =', logic%last_read
  nl=nl+1; write (lines(nl), imt) 'logic%n_1dim_calls               =', logic%n_1dim_calls
  nl=nl+1; write (lines(nl), imt) 'logic%n_db_group_max             =', logic%n_db_group_max
  nl=nl+1; write (lines(nl), imt) 'logic%opt_base                   =', logic%opt_base
  nl=nl+1; write (lines(nl), imt) 'logic%opt_cycles                 =', logic%opt_cycles
  nl=nl+1; write (lines(nl), imt) 'logic%opt_data                   =', logic%opt_data
  nl=nl+1; write (lines(nl), imt) 'logic%opt_loops                  =', logic%opt_loops
  nl=nl+1; write (lines(nl), imt) 'logic%opt_vars                   =', logic%opt_vars
  nl=nl+1; write (lines(nl), imt) 'logic%phase_units                =', logic%phase_units
  nl=nl+1; write (lines(nl), amt) 'logic%pretzel                    =', pretzel_names(logic%pretzel)
  nl=nl+1; write (lines(nl), imt) 'logic%ref_particle               =', logic%ref_particle
  nl=nl+1; write (lines(nl), imt) 'logic%sex_calib_delta_cu         =', logic%sex_calib_delta_cu
  nl=nl+1; write (lines(nl), imt) 'logic%sex_calib_n_bump_set       =', logic%sex_calib_n_bump_set
  nl=nl+1; write (lines(nl), imt) 'logic%tune_units                 =', logic%tune_units
  nl=nl+1; write (lines(nl), imt) 'logic%u_num                      =', logic%u_num
  nl=nl+1; write (lines(nl), imt) 'logic%u_view                     =', logic%u_view

  nl=nl+1; write (lines(nl), rmt) 'logic%bpm_calib_max_dorbit       =', logic%bpm_calib_max_dorbit
  nl=nl+1; write (lines(nl), rmt) 'logic%change_min                 =', logic%change_min
  nl=nl+1; write (lines(nl), rmt) 'logic%dQ_max_quad_calib          =', logic%dQ_max_quad_calib
  nl=nl+1; write (lines(nl), rmt) 'logic%dQ_max_sex_calib           =', logic%dQ_max_sex_calib
  nl=nl+1; write (lines(nl), rmt) 'logic%h_scale                    =', logic%h_scale
  nl=nl+1; write (lines(nl), rmt) 'logic%l_times_alpha              =', logic%l_times_alpha
  nl=nl+1; write (lines(nl), rmt) 'logic%opt_tolerance              =', logic%opt_tolerance
  nl=nl+1; write (lines(nl), rmt) 'logic%phase_conversion_factor    =', logic%phase_conversion_factor
  nl=nl+1; write (lines(nl), rmt) 'logic%tune_conversion_factor     =', logic%tune_conversion_factor

  nl=nl+1; write (lines(nl), lmt) 'logic%auto_measurement           =', logic%auto_measurement
  nl=nl+1; write (lines(nl), lmt) 'logic%bpm_calib_load_golden      =', logic%bpm_calib_load_golden
  nl=nl+1; write (lines(nl), lmt) 'logic%calc_twiss_with_cut_ring   =', logic%calc_twiss_with_cut_ring
  nl=nl+1; write (lines(nl), lmt) 'logic%command_file_open          =', logic%command_file_open
  nl=nl+1; write (lines(nl), lmt) 'logic%debug                      =', logic%debug
  nl=nl+1; write (lines(nl), lmt) 'logic%digital_shaker_on          =', logic%digital_shaker_on
  nl=nl+1; write (lines(nl), lmt) 'logic%dmerit_calc_on             =', logic%dmerit_calc_on
  nl=nl+1; write (lines(nl), lmt) 'logic%gain_correction            =', logic%gain_correction
  nl=nl+1; write (lines(nl), lmt) 'logic%gui                        =', logic%gui
  nl=nl+1; write (lines(nl), lmt) 'logic%limit_on                   =', logic%limit_on
  nl=nl+1; write (lines(nl), lmt) 'logic%lrbbi_inserted             =', logic%lrbbi_inserted
  nl=nl+1; write (lines(nl), lmt) 'logic%new_cbar_tilt_calc         =', logic%new_cbar_tilt_calc
  nl=nl+1; write (lines(nl), lmt) 'logic%nonlinear_calc             =', logic%nonlinear_calc
  nl=nl+1; write (lines(nl), lmt) 'logic%offset_correction          =', logic%offset_correction
  nl=nl+1; write (lines(nl), lmt) 'logic%ok_to_read_dmerit_file     =', logic%ok_to_read_dmerit_file
  nl=nl+1; write (lines(nl), lmt) 'logic%opt_base_locked            =', logic%opt_base_locked
  nl=nl+1; write (lines(nl), lmt) 'logic%opt_cbar11                 =', logic%opt_cbar11
  nl=nl+1; write (lines(nl), lmt) 'logic%opt_cbar22                 =', logic%opt_cbar22
  nl=nl+1; write (lines(nl), lmt) 'logic%opt_running                =', logic%opt_running
  nl=nl+1; write (lines(nl), lmt) 'logic%opt_vars_locked            =', logic%opt_vars_locked
  nl=nl+1; write (lines(nl), lmt) 'logic%plot_ip_at_center          =', logic%plot_ip_at_center
  nl=nl+1; write (lines(nl), lmt) 'logic%plot_locked                =', logic%plot_locked
  nl=nl+1; write (lines(nl), lmt) 'logic%plot_single_bottom         =', logic%plot_single_bottom
  nl=nl+1; write (lines(nl), lmt) 'logic%plot_single_top            =', logic%plot_single_top
  nl=nl+1; write (lines(nl), lmt) 'logic%plot_special               =', logic%plot_special
  nl=nl+1; write (lines(nl), lmt) 'logic%plotit                     =', logic%plotit
  nl=nl+1; write (lines(nl), lmt) 'logic%radiation_on               =', logic%radiation_on
  nl=nl+1; write (lines(nl), lmt) 'logic%read_twiss_with_tbt        =', logic%read_twiss_with_tbt
  nl=nl+1; write (lines(nl), lmt) 'logic%remembering                =', logic%remembering
  nl=nl+1; write (lines(nl), lmt) 'logic%restricted_cbar_plot       =', logic%restricted_cbar_plot
  nl=nl+1; write (lines(nl), lmt) 'logic%rf_on                      =', logic%rf_on
  nl=nl+1; write (lines(nl), lmt) 'logic%ring_initialized           =', logic%ring_initialized
  nl=nl+1; write (lines(nl), lmt) 'logic%tilt_correction            =', logic%tilt_correction
  nl=nl+1; write (lines(nl), lmt) 'logic%use_old_in_loading         =', logic%use_old_in_loading
  nl=nl+1; write (lines(nl), lmt) 'logic%use_bpm_quad_offsets       =', logic%use_bpm_quad_offsets
  nl=nl+1; write (lines(nl), lmt) 'logic%wolski_normal_mode_calc_on =', logic%wolski_normal_mode_calc_on

!----------------------------------------------------------------------
! Luminosity

case ('LUMINOSITY')

  fname = 'cesrv_luminosity.in'
  open (1, file = fname, status = 'old', action = 'read', iostat = ios)
  if (ios /= 0) then
    call fullfilename(&
         'CESR_ONLINE/acc_control/program_info/cesrv/init/cesrv_luminosity.in',&
         fname)
    open (1, file = fname, status = 'old', action = 'read', iostat = ios)
    if (ios /= 0) then
      print *, 'ERROR OPENING "CESRV_LUMINOSITY.IN" IN CURRENT DIRECTORY'
      print *, '      OR IN $CESR_ONLINE/acc_control/program_info/cesrv/init'
      return
    endif
  endif

  read (1, nml = lum_params, iostat = ios)
  close (1)
  if (ios /= 0) then
    print *, 'ERROR READING "LUM_PARAMS" NAMELIST FROM: ', trim(fname)
    return
  endif

! positron calc

  this_species = u%ring%param%particle

  u%ring%param%particle = +1
  call ring_calc (u)
  call radiation_integrals (u%ring, u%orb, global_p)

  ele_p = u%ring%ele(0)
  call c_to_cbar (u%ring%ele(0), cbar_p)

! electron calc

  u%ring%param%particle = -1
  call ring_calc (u)
  call radiation_integrals (u%ring, u%orb, global_e)

  ele_e = u%ring%ele(0)
  call c_to_cbar (u%ring%ele(0), cbar_e)

  u%ring%param%particle = this_species

! average

  if (epsilon_x == 0) then
    eps_x = (global_p%a%emittance + global_e%a%emittance) / 2
  else
    eps_x = epsilon_x
  endif

  n_part = i_tot * u%ring%param%total_length / (n_bunch * e_charge * c_light)
  beta_a = (ele_p%a%beta + ele_e%a%beta) / 2
  beta_b = (ele_p%b%beta + ele_e%b%beta) / 2
  cbar12 = (cbar_p(1,2) + cbar_e(1,2)) / 2
  dcbar22 = abs(cbar_p(2,2) - cbar_e(2,2)) / 2
  sig_x = sqrt(beta_a * eps_x)

  f_xi_x = beta_a * r_e * n_part * e_mass / (2 * pi * &
                                        ele_e%value(E_TOT$) / 1e9)
  f_xi_y = beta_b * r_e * n_part * e_mass / (2 * pi * &
                                        ele_e%value(E_TOT$) / 1e9)

  if (epsilon_y == 0) then
    sig_y_no = f_xi_y / (sig_x * xi_y)
    sig_y_no = sig_y_no * (1 - sig_y_no/sig_x)
    eps_y = (f_xi_y / (xi_y * (sig_x + sig_y_no)))**2 / beta_b
  else
    eps_y = epsilon_y
    sig_y_no = sqrt(eps_y * beta_b)
  endif

  sig_y = sqrt(sig_y_no**2 + cbar12**2 * eps_x * beta_b)

  factor = (1e-4 * i_tot**2 * u%ring%param%total_length) / &
                      (4 * pi * n_bunch * sig_x * c_light * e_charge**2)

  theta_diff = sqrt(beta_b/beta_a) * dcbar22
  sig_y_tot = sqrt(sig_y**2 + (sig_x*theta_diff)**2)

!

  write (lines( 1), *) 'Read luminosity parameters from: ', trim(fname)
  write (lines( 2), *) 'N_bunch =', n_bunch
  write (lines( 3), *) 'I_tot   =', i_tot
  write (lines( 4), *)
  write (lines( 5), '(13x, a)') 'Lum_Calc    E+ Lat    E- Lat'
  fmt = '(a, 1p3e10.2)'
  write (lines( 6), fmt) 'Epsilon_x  ', eps_x, global_p%a%emittance, &
                                                  global_p%a%emittance
  write (lines( 7), fmt) 'Epsilon_y  ', eps_y, global_p%b%emittance, &
                                                  global_p%b%emittance
  fmt = '(a, 3f10.3)'
  write (lines( 8), fmt) 'Beta_a     ', beta_a, ele_p%a%beta, ele_e%a%beta
  write (lines( 9), fmt) 'Beta_b     ', beta_b, ele_p%a%beta, ele_e%a%beta
  write (lines(10), fmt) 'Cbar12     ', cbar12, cbar_p(1,2), cbar_e(1,2)
  fmt = '(a, 10x, 2f10.3)'
  write (lines(11), fmt) 'Cbar22     ', cbar_p(2,2), cbar_e(2,2)
  fmt = '(a, 3f10.3)'
  write (lines(12), fmt) 'd_Cbar22   ', dcbar22
  write (lines(13), *)
  write (lines(14), '(21x, a)') 'Lum     Sig_x     Sig_y      Xi_x      Xi_y'
  fmt = '(a, es10.2, 6p2f10.2, 0p2f10.3)'
  ff = 1 / (1 + sig_y_no / sig_x)
  write (lines(15), fmt) 'No   Coupling:', factor / sig_y_no, &
      sig_x, sig_y_no, ff*f_xi_x/(sig_x**2), ff*f_xi_y/(sig_x*sig_y_no)
  ff = 1 / (1 + sig_y_tot / sig_x)
  write (lines(16), fmt) 'With Coupling:', factor / sig_y_tot, &
      sig_x, sig_y_tot, ff*f_xi_x/(sig_x**2), ff*f_xi_y/(sig_x*sig_y_tot)
  fmt = '(a, f10.1, 2x, a)'
  write (lines(17), fmt) '%Reduction', 100 * (1 - sig_y_no / sig_y), &
                     '! due to Cbar12 (sig_y blowup)'
  write (lines(18), fmt) '%Reduction', 100 * (1 - sig_y / sig_y_tot), &
                     '! due to Cbar22 (differential angle)'
  write (lines(19), fmt) '%Reduction', 100 * (1 - sig_y_no / sig_y_tot), &
                     '! Total due to Coupling'
  nl = 19


!----------------------------------------------------------------------
! Octupole

case ('OCTUPOLE')
  call show_quad_etc (u%oct_k3, 'Octupole', 'K3')

!----------------------------------------------------------------------
! orbit

case ('ORBIT')

  l1 = '     |      Data       |       Ref       |      Model      |'
  l2 = '  ix | X (mm)   Y (mm) | X (mm)   Y (mm) | X (mm)   Y (mm) | Used      ix_plot'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l2
  do i = lbound(u%orbit%x%d, 1), ubound(u%orbit%x%d, 1)
    if (u%orbit%x%d(i)%exists) then
      ele_name = u%ring%ele(u%orbit%x%d(i)%ix_ele)%name(5:)
      nl=nl+1
      write (lines(nl), '(i5, 3p6f9.3, 0p, l7, 2x, a, f8.2)') i, &
           u%orbit%x%d(i)%meas, u%orbit%y%d(i)%meas, &
           u%orbit%x%d(i)%ref, u%orbit%y%d(i)%ref, &
           u%orbit%x%d(i)%model, u%orbit%y%d(i)%model, &
           u%orbit%x%d(i)%useit_plot, trim(ele_name), u%orbit%x%d(i)%ix_plot_index
    endif
  enddo
  nl=nl+1; write (lines(nl), *) l2
  nl=nl+1; write (lines(nl), *) l1

!----------------------------------------------------------------------
! opt_variables

case ('OPT_VARIABLES')

  call show_opt_vars

!----------------------------------------------------------------------
! show parameters

case ('PARAMETERS')

  write (lines(1), *) 'Use "SHOW GLOBAL" instead'
  nl = 1

!----------------------------------------------------------------------
! show phase

case ('PHASE')

  if (f_phi > 10) then
    fmt = '(1x, a7, 6f10.1, i5)'
  else
    fmt = '(1x, a7, 6f10.4, i5)'
  endif

  write (l1, '(9x, a)') &
    '|              X               |             Y              |'
  write (l2, '(9x, a)') &
    '|   Model    Design      Data  |  Model    Design      Data |'

  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(u%phase%x%d, 1), ubound(u%phase%x%d, 1)
    if (u%phase%x%d(i)%exists) then
      x_dat => u%phase%x%d(i)
      y_dat => u%phase%y%d(i)
      ele => u%ring%ele(x_dat%ix_ele)
      nl=nl+1
      write (lines(nl), fmt) ele%name, &
          f_phi*x_dat%model, f_phi*x_dat%design, f_phi*x_dat%meas, &
          f_phi*y_dat%model, f_phi*y_dat%design, f_phi*y_dat%meas, i
    endif
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

!----------------------------------------------------------------------
! Sextupole resonances

case ('2QX', '2QY', 'QX+QY', 'QX-QY')

  select case (show_name)
  case ('2QX');  d2_ptr => u%q2x
  case('2QY');   d2_ptr => u%q2y
  case('QX+QY'); d2_ptr => u%qx_plus_qy
  case('QX-QY'); d2_ptr => u%qx_minus_qy
  end select

  fmt = '(1x, a7, 6f10.2, i5)'

  write (l1, '(9x, a)') '|            In_Phase          |        Out_Phase           |'
  write (l2, '(9x, a)') '|   Model    Design      Data  |  Model    Design      Data |'

  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(d2_ptr%a_in%d, 1), ubound(d2_ptr%a_in%d, 1)
    if (.not. d2_ptr%a_in%d(i)%exists) cycle
    x_dat => d2_ptr%a_in%d(i)
    y_dat => d2_ptr%a_out%d(i)
    ele => u%ring%ele(x_dat%ix_ele)
    nl=nl+1; write (lines(nl), fmt) ele%name, &
        x_dat%model, x_dat%design, x_dat%meas, y_dat%model, y_dat%design, y_dat%meas, i
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

  nl=nl+1; lines(nl) = ''
  nl=nl+1; write (lines(nl), '(a, f10.6)') 'Qx Amplitude:    ', u%tune%data_params%x_amp
  nl=nl+1; write (lines(nl), '(a, f10.6)') 'Qy Amplitude:    ', u%tune%data_params%y_amp
  nl=nl+1; write (lines(nl), '(a, f10.6)') '2Qx Amplitude:   ', u%q2x%data_params%x_amp
  nl=nl+1; write (lines(nl), '(a, f10.6)') '2Qy Amplitude:   ', u%q2y%data_params%x_amp
  nl=nl+1; write (lines(nl), '(a, f10.6)') 'Qx+Qy Amplitude: ', u%qx_plus_qy%data_params%x_amp
  nl=nl+1; write (lines(nl), '(a, f10.6)') 'Qx-Qy Amplitude: ', u%qx_minus_qy%data_params%x_amp

!----------------------------------------------------------------------
! Quad

case ('QUADRUPOLE')
  call show_quad_etc (u%quad_k1, 'Quadrupole', 'K1')

!----------------------------------------------------------------------
! show rad_int

case ('RAD_INT')

  call set_on_off (rfcavity$, u%ring, on$, u%orb)
  call radiation_integrals (u%ring, u%orb, global_model, rad_int_by_ele = rad_int_by_ele_model)
  call set_on_off (rfcavity$, u%ring, on_off_int(logic%rf_on), u%orb)
  
  call match_word (show2_name, &
                    ['I0 ', 'I1 ', 'I2 ', 'I3 ', 'I4A', 'I4B', 'I5A', 'I5B', 'I6 '], &
                    ix1, matched_name = mname)

  l1 = '           Model           Design              Diff'
  nl=nl+1; write (lines(nl), '(21x, a)') trim(l1)

  do i = 1, ubound(rad_int_by_ele_model%ele, 1)
    ele => u%ring%ele(i)

    select case (mname)
    case ('I0')
      model_val  = rad_int_by_ele_model%ele(i)%i0
      design_val = u%rad_int_by_ele_design%ele(i)%i0
    case ('I1')
      model_val  = rad_int_by_ele_model%ele(i)%i1
      design_val = u%rad_int_by_ele_design%ele(i)%i1
    case ('I2')
      model_val  = rad_int_by_ele_model%ele(i)%i2
      design_val = u%rad_int_by_ele_design%ele(i)%i2
    case ('I3')
      model_val  = rad_int_by_ele_model%ele(i)%i3
      design_val = u%rad_int_by_ele_design%ele(i)%i3
    case ('I4A')
      model_val  = rad_int_by_ele_model%ele(i)%i4a
      design_val = u%rad_int_by_ele_design%ele(i)%i4a
    case ('I4B')
      model_val  = rad_int_by_ele_model%ele(i)%i4b
      design_val = u%rad_int_by_ele_design%ele(i)%i4b
    case ('I5A')
      model_val  = rad_int_by_ele_model%ele(i)%i5a
      design_val = u%rad_int_by_ele_design%ele(i)%i5a
    case ('I5B')
      model_val  = rad_int_by_ele_model%ele(i)%i5b
      design_val = u%rad_int_by_ele_design%ele(i)%i5b
    case ('I6B')
      model_val  = rad_int_by_ele_model%ele(i)%i6b
      design_val = u%rad_int_by_ele_design%ele(i)%i6b
    case default
      print *, 'BAD RAD_INT'
      return
    end select

    if (model_val == 0 .and. design_val == 0) cycle

    nl=nl+1; write (lines(nl), '(i4, 1x, a16, 3es18.6)') i, &
      ele%name, model_val, design_val, model_val - design_val
  enddo

  nl=nl+1; write (lines(nl), '(21x, a)') trim(l1)

!----------------------------------------------------------------------
! show u%ring info

case ('RING')

  ix = index('ENDS', line(:ix_s2))
  if (ix == 1) then
    at_ends = .true.
    call string_trim (line(ix_s2+1:), line, ix_s2)
    nl=nl+1; write (lines(nl), '(36x, a)') 'Model values at End of Element:'
  else
    at_ends = .false.
    nl=nl+1; write (lines(nl), '(36x, a)') 'Model values at Center of Element:'
  endif

  call find_range_endpoints (line, loc, ix1, ix2, err_flag)
  if (err_flag) return

  nl=nl+1; write (lines(nl), '(40x, a)')    '|              X         (mm) |             Y         (mm)'
  nl=nl+1; write (lines(nl), '(5x, a, 20x, a)') &
                         'Name', 'S       L  |  Beta     Phi   Eta     Orb | Beta     Phi   Eta     Orb'

  do ix = ix1, ix2
    ele => u%ring%ele(ix)
    if (ix == 0 .or. at_ends) then
      runt = ele
      orb = u%orb(ix)
    else
      call twiss_and_track_intra_ele (ele, u%ring%param, 0.0_rp, ele%value(l$)/2, .true., .false., &
                                                             u%orb(ix-1), orb, u%ring%ele(ix-1), runt)
      runt%s = ele%s - ele%value(l$) / 2
    endif
    nl=nl+1
    write (lines(nl), '(i4, 1x, a19, 2f8.3, 2(f7.2, f8.3, f6.2, f8.3))') &
        ix, ele%name, runt%s, ele%value(l$), &
        runt%a%beta, f_phi*runt%a%phi, runt%a%eta, 1000*orb%vec(1), &
        runt%b%beta, f_phi*runt%b%phi, runt%b%eta, 1000*orb%vec(3)
  enddo

!----------------------------------------------------------------------
! show separators

case ('SEPARATORS')

  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *) &
     '                              Kick (mr)        |           CU'
  nl=nl+1; write (lines(nl), *) &
     'Index     Name       Model    Design     Saved | Model  Design   Saved'
  nl=nl+1; write (lines(nl), *) ' '

  do i = 1, ubound(u%hsep_kick%v, 1)
    vp => u%hsep_kick%v(i)
    if (vp%exists) then
        nl=nl+1
        write (lines(nl), '(i5, 2x, a10, 3f10.3, 3i8)') i, vp%name, &
            vp%model*1000, vp%design*1000,  vp%saved*1000, &
            nint(vp%model/vp%dvar_dcu), vp%cu_design, vp%cu_saved
    endif
  enddo

  nl=nl+1
  write (lines(nl), *) ' '

  nl=nl+1
  write (lines(nl), *) ' '

!----------------------------------------------------------------------
! Sets

case ('SETS')

  fmt = '(1x, a, i7, 4x, a)'

  call get_save_set_date ('CSR', logic%csr_set, date)
  nl=nl+1; write (lines(nl), fmt) 'CSR    Save Set:', logic%csr_set, date
  call get_save_set_date ('BIGGRP', logic%biggrp_set, date)
  nl=nl+1; write (lines(nl), fmt) 'BIGGRP Save Set:', logic%biggrp_set, date
  nl=nl+1; write (lines(nl), *)
  nl=nl+1; write (lines(nl), fmt) 'Orbit  Data  #', u%orbit%ix_meas, u%orbit%date
  nl=nl+1; write (lines(nl), fmt) 'Orbit  Ref   #', u%orbit%ix_ref, u%orbit%ref_date
  nl=nl+1; write (lines(nl), fmt) 'Phase  Data  #', u%phase%ix_meas, u%phase%date
  nl=nl+1; write (lines(nl), fmt) 'Phase  Ref   #', u%phase%ix_ref, u%phase%ref_date
  nl=nl+1; write (lines(nl), fmt) 'Eta    Data  #', u%eta%ix_meas, u%eta%date
  nl=nl+1; write (lines(nl), fmt) 'Eta    Ref   #', u%eta%ix_ref, u%eta%ref_date
  nl=nl+1; write (lines(nl), fmt) 'AC_Eta Data  #', u%ac_eta%ix_meas, u%ac_eta%date
  nl=nl+1; write (lines(nl), fmt) 'AC_Eta Ref   #', u%ac_eta%ix_ref, u%ac_eta%ref_date

!----------------------------------------------------------------------
! Sextupole

case ('SEXTUPOLE')
  call show_quad_etc (u%sex_k2, 'Sextupole', 'K2')

!----------------------------------------------------------------------
! single corrector info
! merit function is of the form: Merit = a*x^2 + b*x + c

case ('SINGLE_CORRECTION')

  if (logic%opt_vars == opt_steering$) then

    call to_top10(top_steer, 0.0_rp, 'INIT_TOP10', 0, 'max')

    call single_corrector_calc (u%hsteer_kick, top_steer, 'Horiz', u)
    call single_corrector_calc (u%vsteer_kick, top_steer, 'Vert',  u)

    nl=nl+1; write (lines(nl), *) '                         |    @minimum'
    nl=nl+1; write (lines(nl), *) 'Steering  ix  Merit_min  |Model  Saved-Model'
    do i = 1, 150
      if (top_steer(i)%valid) then
        ix = top_steer(i)%index
        if (top_steer(i)%name(1:1) == 'H') then
          cu_model_min = u%hsteer_kick%v(ix)%old / u%hsteer_kick%v(ix)%dvar_dcu
          cu_saved_min = u%hsteer_kick%v(ix)%cu_saved - cu_model_min
        else
          cu_model_min = u%vsteer_kick%v(ix)%old / u%vsteer_kick%v(ix)%dvar_dcu
          cu_saved_min = u%vsteer_kick%v(ix)%cu_saved - cu_model_min
        endif
        nl=nl+1
        write (lines(nl), '(2x, a7, i4, f11.3, i8, i11)') top_steer(i)%name, &
               top_steer(i)%index, top_steer(i)%value, &
               nint(cu_model_min), nint(cu_saved_min)

      endif
    enddo
    nl=nl+1; write (lines(nl), *) 'Steering  ix  Merit_min  |Model  Saved-Model'
    nl=nl+1; write (lines(nl), *) '                         |    @minimum'
    nl=nl+1; write (lines(nl), *) '[Note: Merit_min calculated assuming steering_wgt = 0]'
    nl=nl+1; write (lines(nl), *) '[Note: -Model = Change in variable to make a correcten.]'
    nl=nl+1; write (lines(nl), *) '[Note: Saved-Model = What you want to set to make a correction]'
  endif

!----------------------------------------------------------------------
! Skew_quad

case ('SKEW_QUAD')
  call show_quad_etc (u%skew_quad_k1, 'Skew_Quad', 'K1')

!----------------------------------------------------------------------
! Skew_sex

case ('_SKEW_SEX')
  call show_quad_etc (u%skew_sex_k2, '_Skew_sex', 'K2')

!----------------------------------------------------------------------
! Steerings

case ('STEERINGS')
  call show_quad_etc (u%hsteer_kick, 'Horizontal', 'K.')
  nl=nl+1; write (lines(nl), *) ' '
  call show_quad_etc (u%vsteer_kick, 'Vertical', 'K.')

!----------------------------------------------------------------------
! TBT

case ('TBT')

  ix = index('REFERENCE', show2_name(:ix_s2))
  if (ix_s2 == 0) then
    tbt => u%tbt_meas
  elseif (ix == 1) then
    tbt => u%tbt_ref
  else
    nl=nl+1; write (lines(nl), amt) 'I DO NOT UNDERSTAND: "', show2_name(:ix_s2), '"'
    goto 8000
  endif

  if (.not. allocated(tbt%bunch)) then
    nl=nl+1; write (lines(nl), amt) 'No data read in.'
    goto 8000
  endif

  nl=nl+1; write (lines(nl), amt)          'Time:               ', tbt%time_stamp
  nl=nl+1; write (lines(nl), amt)          'Data file:          ', trim(tbt%data_file)
  nl=nl+1; write (lines(nl), '(3(a, i0))') 'Data file num:      ', tbt%file_number, ' (and ', tbt%file_number-1, ')'
  nl=nl+1; write (lines(nl), imt)          'CONDX set num:      ', tbt%condx_number
  nl=nl+1; write (lines(nl), imt)          'Number of turns:    ', size(tbt%shaker)
  nl=nl+1; write (lines(nl), lmt)          'Shake_at_reflection_tune_a: ', tbt%shake_at_reflection_tune_a
  nl=nl+1; write (lines(nl), lmt)          'Shake_at_reflection_tune_b: ', tbt%shake_at_reflection_tune_b
  nl=nl+1; write (lines(nl), imt)          'logic%ix_bunch_fft:         ', logic%ix_bunch_fft
  nl=nl+1; write (lines(nl), lmt)          'Bunches:'
  nl=nl+1; write (lines(nl), lmt)          '        Ix Bunch_Index   RF_Bucket'
  do i = 1, size(tbt%bunch)
    if (i == logic%ix_bunch_fft) then
      nl=nl+1; write (lines(nl), '(a, i6, 2i12)') '  **', i, tbt%bunch(i)%ix_bunch, tbt%bunch(i)%ix_rf_bucket
    else
      nl=nl+1; write (lines(nl), '(i10, 2i12)') i, tbt%bunch(i)%ix_bunch, tbt%bunch(i)%ix_rf_bucket
    endif
  enddo

!----------------------------------------------------------------------
! top10

case ('TOP10')

  call top10_type (u, lines, nl)

!----------------------------------------------------------------------
! TRANSFER_MATRIX

case ('TRANSFER_MATRIX')

  loc = 0
  if (show2_name /= '') then
    call locate_element (u%ring, show2_name, loc, err_flag)
    if (err_flag) return
  endif

  call transfer_matrix_calc (u%ring, .true., mat6, vec0, loc, loc, one_turn = .true.)
  nl=nl+1; write (lines(nl), '(a, es12.3, a)') 'Transfer Matrix : Kick  [Matrix symplectic error:', &
                                             mat_symp_error(mat6), ']'
  if (any(abs(mat6) >= 1000)) then
    do i = 1, 6
      nl=nl+1; write (lines(nl), '(6es11.3, a, es11.3)') (mat6(i, j), j = 1, 6), '   : ', vec0(i)
    enddo
  else
    do i = 1, 6
      nl=nl+1; write (lines(nl), '(6f10.5, a, es11.3)') (mat6(i, j), j = 1, 6), '   : ', vec0(i)
    enddo
  endif

!----------------------------------------------------------------------
! Tune

case ('TUNE')

  call show_tune (u, lines, nl)

!----------------------------------------------------------------------

case ('VARIABLE')

  call match_var_name (line, vp, ix, u, err_flag)
  if (err_flag) return

  nl=nl+1; write (lines(nl), *) 'Index in %var array: ', ix
  nl=nl+1; write (lines(nl), *) 'Name:           ', vp%name
  nl=nl+1; write (lines(nl), *) 'Ele_name:       ', vp%ele_name
  nl=nl+1; write (lines(nl), *) 'Alias:          ', vp%alias
  nl=nl+1; write (lines(nl), *) 'Db_node_name:   ', vp%db_node_name
  nl=nl+1; write (lines(nl), *) 'Attrib_name:    ', vp%attrib_name
  nl=nl+1; write (lines(nl), *) 'Ix_index:       ', vp%ix_index
  nl=nl+1; write (lines(nl), *) 'Ix_db:          ', vp%ix_db
  nl=nl+1; write (lines(nl), *) 'Ix_dvar:        ', vp%ix_dvar
  nl=nl+1; write (lines(nl), *) 'Ix_ele:         ', vp%ix_ele
  nl=nl+1; write (lines(nl), *) 'Ix_attrib:      ', vp%ix_attrib
  nl=nl+1; write (lines(nl), *) 'Model:          ', vp%model
  nl=nl+1; write (lines(nl), *) 'Base_model:     ', vp%base_model
  nl=nl+1; write (lines(nl), *) 'Design:         ', vp%design
  nl=nl+1; write (lines(nl), *) 'Old:            ', vp%old
  nl=nl+1; write (lines(nl), *) 'Base_cu0:       ', vp%base_cu0
  nl=nl+1; write (lines(nl), *) 'Saved:          ', vp%saved
  nl=nl+1; write (lines(nl), *) 'Saved_ref:      ', vp%saved_ref
  nl=nl+1; write (lines(nl), *) 'Target:         ', vp%target
  nl=nl+1; write (lines(nl), *) 'High_target_lim:', vp%high_target_lim
  nl=nl+1; write (lines(nl), *) 'Low_target_lim: ', vp%low_target_lim
  nl=nl+1; write (lines(nl), *) 'Delta:          ', vp%delta
  nl=nl+1; write (lines(nl), *) 'Merit:          ', vp%merit
  nl=nl+1; write (lines(nl), *) 'Dmerit:         ', vp%dmerit
  nl=nl+1; write (lines(nl), *) 'Weight:         ', vp%weight
  nl=nl+1; write (lines(nl), *) 'Step:           ', vp%step
  nl=nl+1; write (lines(nl), *) 'Dvar_dcu:       ', vp%dvar_dcu
  nl=nl+1; write (lines(nl), *) 'Cu_saved:       ', vp%cu_saved
  nl=nl+1; write (lines(nl), *) 'Cu_saved_ref:   ', vp%cu_saved_ref
  nl=nl+1; write (lines(nl), *) 'Cu_target:      ', vp%cu_target
  nl=nl+1; write (lines(nl), *) 'Cu_design:      ', vp%cu_design
  nl=nl+1; write (lines(nl), *) 'Cu_high_lim:    ', vp%cu_high_lim
  nl=nl+1; write (lines(nl), *) 'Cu_low_lim:     ', vp%cu_low_lim
  nl=nl+1; write (lines(nl), *) 'Cu_zero_lim:    ', vp%cu_zero_lim
  nl=nl+1; write (lines(nl), *) 'S:              ', vp%s
  nl=nl+1; write (lines(nl), *) 'Exists:         ', vp%exists
  nl=nl+1; write (lines(nl), *) 'Good_var:       ', vp%good_var
  nl=nl+1; write (lines(nl), *) 'Good_user:      ', vp%good_user
  nl=nl+1; write (lines(nl), *) 'Good_opt:       ', vp%good_opt
  nl=nl+1; write (lines(nl), *) 'Useit:          ', vp%useit
  nl=nl+1; write (lines(nl), *) 'Makeup_method:  ', name$%makeup_method(vp%makeup_method)
  nl=nl+1; write (lines(nl), *) 'Do_limit_calc:  ', vp%do_limit_calc
  nl=nl+1; write (lines(nl), *) 'Dm_dv_computed: ', vp%dm_dv_computed
  nl=nl+1; write (lines(nl), *) 'Units:          ', var_units_name(vp%v1%units)

!----------------------------------------------------------------------
! Vertical

case ('VERTICAL')
  call show_quad_etc (u%vsteer_kick, 'Vertical', 'K.')

!----------------------------------------------------------------------
! Weight

case ('WEIGHT')

  if (index('TWISS', show2_name(:ix_s2)) == 1) then
    L1 = '   |         Phase       |  Cbar12    Cbar11    Cbar22   |' // &
                                                    '         Beta      |'
    L2 = 'ix |      X         Y    |                               |' // &
                                                    '      X         Y  |'

    nl=nl+1; write (lines(nl), *)
    nl=nl+1; write (lines(nl), *) L1
    nl=nl+1; write (lines(nl), *) L2

    do i = 0, 99
      nl = nl + 1
      write (lines(nl), '(i4)') i
      if (u%phase%x%d(i)%exists) write (lines(nl)(5:), '(2es10.2)') &
                                u%phase%x%d(i)%weight, u%phase%y%d(i)%weight
      if (u%cbar%m12%d(i)%exists) write (lines(nl)(27:), '(1p3e10.2)') &
                  u%cbar%m12%d(i)%weight, u%cbar%m11%d(i)%weight, &
                  u%cbar%m22%d(i)%weight
      if (u%beta%x%d(i)%exists) write (lines(nl)(59:), '(2es10.2)') &
                                u%beta%x%d(i)%weight, u%beta%y%d(i)%weight
    enddo

    nl=nl+1; write (lines(nl), *) L2
    nl=nl+1; write (lines(nl), *) L1

  elseif (index('VARIABLES', show2_name(:ix_s2)) == 1) then
    L1 = ' ix | Horizontal    Vertical  Quadrupole   Sextupole |'

    nl=nl+1; write (lines(nl), *)
    nl=nl+1; write (lines(nl), *) L1

    do i = 0, 120
      if (.not. (u%quad_k1%v(i)%exists .or. u%sex_k2%v(i)%exists .or. &
            u%hsteer_kick%v(i)%exists .or. u%vsteer_kick%v(i)%exists)) cycle
      nl = nl + 1
      write (lines(nl), '(i4)') i
      if (u%hsteer_kick%v(i)%exists) write (lines(nl)(6:), '(es12.2)') &
               u%hsteer_kick%v(i)%weight * u%hsteer_kick%v(i)%dvar_dcu**2
      if (u%vsteer_kick%v(i)%exists) write (lines(nl)(18:), '(es12.2)') &
               u%vsteer_kick%v(i)%weight * u%vsteer_kick%v(i)%dvar_dcu**2
      if (u%quad_k1%v(i)%exists) write (lines(nl)(30:), '(es12.2)') &
                                                    u%quad_k1%v(i)%weight
      if (u%sex_k2%v(i)%exists) write (lines(nl)(42:), '(es12.2)') &
                                                    u%sex_k2%v(i)%weight
    enddo

    nl=nl+1; write (lines(nl), *) L1

  elseif (index('ORBIT_AND_DISPERSION', show2_name(:ix_s2)) == 1) then
    L1 = '    |         Orbit       |          Eta       |'
    L2 = ' ix |      X         Y    |       X         Y  |'

    nl=nl+1; write (lines(nl), *)
    nl=nl+1; write (lines(nl), *) L1
    nl=nl+1; write (lines(nl), *) L2

    do i = 0, 99
      nl = nl + 1
      write (lines(nl), '(i4)') i
      if (u%orbit%x%d(i)%exists) write (lines(nl)(6:), '(2es10.2)') &
                                u%orbit%x%d(i)%weight, u%orbit%y%d(i)%weight
      if (u%eta%x%d(i)%exists) write (lines(nl)(29:), '(2es10.2)') &
                                u%eta%x%d(i)%weight, u%eta%y%d(i)%weight
    enddo

    nl=nl+1; write (lines(nl), *) L2
    nl=nl+1; write (lines(nl), *) L1

  else
    print *, 'Show weight of WHAT? (Twiss, Variables, or Orbit_and_Dispersion)'
    return
  endif

!----------------------------------------------------------------------
! show wolski

case ('WOLSKI')

  read (line, *, iostat = ios) i
  if (ios /= 0) then
    print *, 'ERROR: CANNOT READ BUTTON NUMBER!'
    return
  endif

  nl=nl+1; write(lines(nl), '(5x, 4es12.4)') u%but_to_mode(i)%mat(1,:)
  nl=nl+1; write(lines(nl), '(5x, 4es12.4)') u%but_to_mode(i)%mat(2,:)

!----------------------------------------------------------------------
! show xray

case ('XRAY')

  fmt = '(i3, 3p5f10.4, 0p2f10.3, l8, f8.2, l7, i7, 2x a)'

  write (l1, '(a)') '    |            Vertical         |      At Source Pt |         source_to                 |'
  write (l2, '(a)') ' ix |    Meas       Ref     Model |       y        py |      S    monitor   Mirror    Mag | Used  ix_db  Name'

  nl=nl+1; lines(nl) = 'Electron:'
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(u%e_xray%y%d, 1), ubound(u%e_xray%y%d, 1)
    if (.not. u%e_xray%y%d(i)%exists) cycle
    y_dat => u%e_xray%y%d(i)
    ele => u%ring%ele(y_dat%ix_ele)
    db_ele => u%db%e_chess_monitor_source(i)
    orbit => u%orb(y_dat%ix_ele)
    nl=nl+1
    write (lines(nl), fmt) i, y_dat%meas, y_dat%ref, y_dat%model, orbit%vec(1), orbit%vec(2), &
            ele%s, db_ele%val(dist_to_detec$), db_ele%logic(has_mirror$), db_ele%val(magnification$), &
            y_dat%useit_opt, y_dat%ix_db, ele%name
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

  nl=nl+1; lines(nl) = ''
  nl=nl+1; lines(nl) = 'Positron:'
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(u%p_xray%y%d, 1), ubound(u%p_xray%y%d, 1)
    if (.not. u%p_xray%y%d(i)%exists) cycle
    y_dat => u%p_xray%y%d(i)
    ele => u%ring%ele(y_dat%ix_ele)
    db_ele => u%db%p_chess_monitor_source(i)
    orbit => u%orb(y_dat%ix_ele)
    nl=nl+1
    write (lines(nl), fmt) i, y_dat%meas, y_dat%ref, y_dat%model, orbit%vec(1), orbit%vec(2), &
            ele%s, db_ele%val(dist_to_detec$), db_ele%logic(has_mirror$), db_ele%val(magnification$), &
            y_dat%useit_opt, y_dat%ix_db, ele%name
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

!----------------------------------------------------------------------
! show twiss

case ('TWISS')

  read (line, *, iostat = ios) s_pos
  if (ios /= 0) then
    print *, 'ERROR: CANNOT READ S OFFSET'
    return
  endif

  s_pos = modulo (s_pos, u%ring%param%total_length)

  do i = 1, u%ring%n_ele_track
    if (s_pos < u%ring%ele(i)%s) then
      del_s = s_pos - u%ring%ele(i-1)%s
      call twiss_and_track_intra_ele (u%ring%ele(i), u%ring%param, 0.0_rp, del_s, &
                               .true., .false., u%orb(i-1), orb, u%ring%ele(i-1), runt)
      exit
    endif
  enddo

  call radiation_integrals (u%ring, u%orb, global_model)
  call beam_sigma_calc (runt, global_model, sig_a, sig_b, sig_c)

  write (lines(1), '(a, f10.5)') 'At S =', s_pos
  write (lines(2), '(2a)')            'In Element: ', u%ring%ele(i)%name

  call type_twiss (runt, logic%phase_units, lines = lines(3:), n_lines = nl)
  nl = nl + 2

  call c_to_cbar (runt, cbar_model)
  fmt = '(2f11.5, 4x, 2f11.5)'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), '(1x, a, 23x, a)') 'Cbar:', 'C:'
  nl=nl+1; write (lines(nl), fmt) &
          cbar_model(1,1), cbar_model(1,2), runt%c_mat(1,1), runt%c_mat(1,2)
  nl=nl+1; write (lines(nl), fmt) &
          cbar_model(2,1), cbar_model(2,2), runt%c_mat(2,1), runt%c_mat(2,2)

  fmt = '(2x, a, 3p2f11.4)'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *)   'Orbit: [mm, mrad]'
  nl=nl+1; write (lines(nl), fmt) "X  X':", orb%vec(1), orb%vec(2)
  nl=nl+1; write (lines(nl), fmt) "Y  Y':", orb%vec(3), orb%vec(4)
  nl=nl+1; write (lines(nl), fmt) "Z  Z':", orb%vec(5), orb%vec(6)

  nl=nl+1; lines(nl)= ''
  nl=nl+1; lines(nl) = 'Beam sigmas from radiation damping and excitation:'
  nl=nl+1; lines(nl) = '            A-mode      B-mode      Z-mode       Total'
  nl=nl+1; write (lines(nl), '(a, 4es12.3)') 'sig_x ', sqrt(sig_a(1,1)), &
             sqrt(sig_b(1,1)), sqrt(sig_c(1,1)), sqrt(sig_a(1,1)+sig_b(1,1)+sig_c(1,1))
  nl=nl+1; write (lines(nl), '(a, 4es12.3)') 'sig_px', sqrt(sig_a(2,2)), &
             sqrt(sig_b(2,2)), sqrt(sig_c(2,2)), sqrt(sig_a(2,2)+sig_b(2,2)+sig_c(2,2))
  nl=nl+1; write (lines(nl), '(a, 4es12.3)') 'sig_y ', sqrt(sig_a(3,3)), &
             sqrt(sig_b(3,3)), sqrt(sig_c(3,3)), sqrt(sig_a(3,3)+sig_b(3,3)+sig_c(3,3))
  nl=nl+1; write (lines(nl), '(a, 4es12.3)') 'sig_py', sqrt(sig_a(4,4)), &
             sqrt(sig_b(4,4)), sqrt(sig_c(4,4)), sqrt(sig_a(4,4)+sig_b(4,4)+sig_c(4,4))
  nl=nl+1; write (lines(nl), '(a, 4es12.3)') 'sig_z ', sqrt(sig_a(5,5)), &
             sqrt(sig_b(5,5)), sqrt(sig_c(5,5)), sqrt(sig_a(5,5)+sig_b(5,5)+sig_c(5,5))
  nl=nl+1; write (lines(nl), '(a, 4es12.3)') 'sig_pz', sqrt(sig_a(6,6)), &
             sqrt(sig_b(6,6)), sqrt(sig_c(6,6)), sqrt(sig_a(6,6)+sig_b(6,6)+sig_c(6,6))

!----------------------------------------------------------------------
! show wave

case ('WAVE')

  if (u%wave%wave_what /= '') then
    nl=nl+1; write (lines(nl), *) 'Wave data: ', trim(u%wave%wave_what), ' ', trim(name$%plane(u%wave%plane))
  endif
  nl=nl+1; write (lines(nl), *) 'IX_A1 =', u%wave%ix_a1, '    ! Left  edge of A region'
  nl=nl+1; write (lines(nl), *) 'IX_A2 =', u%wave%ix_a2, '    ! Right edge of A region'
  nl=nl+1; write (lines(nl), *) 'IX_B1 =', u%wave%ix_b1, '    ! Left  edge of B region'
  nl=nl+1; write (lines(nl), *) 'IX_B2 =', u%wave%ix_b2, '    ! Right edge of B region'

!----------------------------------------------------------------------

case default

  print *, 'Eh? What to show?'

end select

!----------------------------------------------------------------------
! end stuff

8000 continue

do i=1,nl
 if (logic%gui) then
  call cesrv_show(lines(i))
 else
   print '(1x, a)', trim(lines(i))
 endif
enddo

if (present(file_name)) then

  if (file_name == '*REM') then
    iu = logic%iu_remember
    if (iu == 0) then
      print *, 'ERROR: NO "REMEMBER" FILE. USE THE "REMEMBER" COMMAND FIRST.'
      return
    endif
    old_file_name = logic%remember_file_name
    ios = 0
    write_append_str = 'Appended to: '
  elseif (file_name == ' ' .and. old_file_name /= ' ') then
    iu = lunget()
    open (iu, file = old_file_name, status = 'old', position = 'append', iostat = ios)
    write_append_str = 'Appended to: '
  else
    if (old_file_name == ' ') old_file_name = 'CESRV.DAT'
    if (file_name /= ' ') old_file_name = file_name
    iu = lunget()
    open (iu, file = old_file_name, iostat = ios)
    write_append_str = 'Written:'
  endif

  if (ios /= 0) then
    print *, 'ERROR OPENING: ', trim(file_name)
    return
  endif

  do i = 1, nl
    write (iu, *) trim(lines(i)(2:))
  enddo

  close(iu)

  print *, trim(write_append_str)," ", trim(old_file_name)

endif

!------------------------------------------------------------
contains

subroutine find_range_endpoints (line, loc, ix1, ix2, err_flag)

implicit none

character(*) line
integer ix1, ix2, loc, ix, i1, i2
logical err_flag

! if in the form "nnn:mmm" then just set ix1 = nnn, ix2 = mmm

ix = index(line, ':')
if (ix /= 0) then
  err_flag = .true.
  read (line(:ix-1), *, iostat = ios) i1
  if (ios /= 0) then
    print *, 'ERROR READING 1ST LOCATION'
    return
  endif
  read (line(ix+1:), *, iostat = ios) i2
  if (ios /= 0) then
    print *, 'ERROR READING 2nd LOCATION'
    return
  endif
  ix1 = max(min(i1, i2, u%ring%n_ele_max), 0)
  ix2 = min(max(i1, i2, 0), u%ring%n_ele_max)
  err_flag = .false.
  return
endif

! otherwise input is, for example, something like: "45E"
! The idea now is to choose ix1 and ix2 to bracket the location with
! ix2-ix1 = 30 (approximately)

call cesr_locator(line, prefix, loc, err_flag)

if (err_flag .or. loc < 0 .or. loc > 99) then
  print *, 'ERROR: BAD RING LOCATION'
  err_flag = .true.
  return
endif

if (loc == 0) then
  ix1 = 1
else
  ix1 = u%quad_k1%v(loc-1)%ix_ele
  if (ix1 > u%ring%n_ele_track) then
    ix = u%ring%ele(ix1)%ix1_slave
    ix1 = u%ring%control(ix)%ix_slave
  endif
endif

if (loc == 99) then
  ix2 = u%ring%n_ele_track
else
  ix2 = u%quad_k1%v(loc+1)%ix_ele
  if (ix2 > u%ring%n_ele_track) then
    ix = u%ring%ele(ix2)%ix2_slave
    ix2 = u%ring%control(ix)%ix_slave
  endif
endif

if (ix1 == 0 .and. ix2 == 0) then
  ix1 = u%ring%n_ele_track / loc - 15
  ix2 = u%ring%n_ele_track / loc + 15
elseif (ix1 == 0) then
  ix1 = ix2 - 30
elseif (ix2 == 0) then
  ix2 = ix1 + 30
endif

ix1 = max(min(ix1, u%ring%n_ele_track), 1)
ix2 = max(min(ix2, u%ring%n_ele_track), 1)

if (ix2-ix1 < 28) then
  ix1 = max (1, (ix1+ix2)/2 - 15)
  ix2 = min (u%ring%n_ele_track, (ix1+ix2)/2 + 15)
endif

end subroutine find_range_endpoints

!-------------------------------------------------
! contains

subroutine get_save_set_date (kind, set_num, date)

character(*) kind, date
integer set_num, iu, ix, ios
character(80) line
logical err

!

call form_file_name_with_number (kind, set_num, fname, err_flag)

if (set_num <= 0 .or. err) then
  date = '!! Does not Exist'
  return
endif

iu = lunget()
open (iu, file = fname, status = 'old', iostat = ios)

if (ios == 0) then
  read (iu, '(a)') line
  call string_trim (line(32:), date, ix)
  close (iu)
else
  date = '!! Does Not Exist'
endif

end subroutine get_save_set_date

!-------------------------------------------------
! contains

subroutine show_opt_vars

integer i, j

nl=nl+1; write (lines(nl), *)

opt_vars = logic%opt_vars

if (opt_vars == opt_custom$ .or. opt_vars == opt_all_vars$) then
  call var_type_useit (u%custom_var, .false., lines, nl)
endif

if (opt_vars == opt_quad$ .or. opt_vars == opt_all_vars$ .or. &
  opt_vars == opt_skew_quad$) then
  call var_type_useit (u%quad_k1, .true., lines, nl)
  call var_type_useit (u%skew_quad_k1, .true., lines, nl)
endif

if (opt_vars == opt_sex$ .or. opt_vars == opt_all_vars$) then
  call var_type_useit (u%sex_k2, .true., lines, nl)
endif
if (opt_vars == opt_sex$ .or. opt_vars == opt_all_vars$) then
  call var_type_useit (u%skew_sex_k2, .false., lines, nl)
endif

if (opt_vars == opt_steering$ .or. opt_vars == opt_all_vars$) then
  call var_type_useit (u%hsteer_kick, .true., lines, nl)
  call var_type_useit (u%vsteer_kick, .true., lines, nl)
  call var_type_useit (u%hsep_kick, .false., lines, nl)
  if (u%init_orb%v(1)%exists) call var_type_useit (u%init_orb, .false., lines, nl)
endif

if (opt_vars == opt_all_vars$) then
  call var_type_useit (u%bpm_tilt, .true., lines, nl)
  call var_type_useit (u%x_kick_quad, .true., lines, nl)
  call var_type_useit (u%y_kick_quad, .true., lines, nl)
endif

do i = 1, logic%n_db_group_max
  do j = 1, size(u%db_group(i)%v1%v)
    if (u%db_group(i)%v1%v(j)%useit) then
      call var_type_useit (u%db_group(i)%v1, .false., lines, nl)
      exit
    endif
  enddo
enddo

end subroutine show_opt_vars

!-------------------------------------------------
! contains

subroutine show_quad_etc (v1_mag, this_name, k_name)

type (v1_var_struct) v1_mag
type (var_struct), pointer :: vp
type (ele_struct) ele_twiss

real(rp) target, saved, l_ele, saved_ref
real(rp) beta, phi, ff
integer i, j
integer, allocatable :: this_indexx(:)
logical is_steering, is_quad
character(*) this_name, k_name

! init

call cu_target_calc (v1_mag%v, u, 1.0_rp)

is_quad = .false.  
if (v1_mag%type == quad_k1$) is_quad = .true.

if (v1_mag%type == hsteer_kick$ .or. v1_mag%type == vsteer_kick$) then
  if (show2_name == 'KL') then
    nl=nl+1; write (lines(nl), *) 'WARNING "KL" DOES *NOT* MAKE SENSE FOR A STEERING!'
    return
  endif
  is_steering = .true.
  ff = 1000
else
  is_steering = .false.
  ff = 1
endif

! quad blank

if (len_trim(show2_name) == 0 .and. is_steering) then

  if (f_phi > 10) then
    fmt = '(1x, a9, 3f8.3, 3i8, f8.2, f7.1, 3x, i3)'
  else
    fmt = '(1x, a9, 3f8.3, 3i8, f8.2, f7.3, 3x, i3)'
  endif

  L1 = '   Name   | Model   dBase dB-Modl | Model   dBase dB-Modl |  Beta    Phi |  Ix'
  L2 = '          |  [mrad*sqrt(u%beta)]  |          [CU]         |              |'
  nl=nl+1; write (lines(nl), *) this_name, ':'
  nl=nl+1; write (lines(nl), *) L1
  nl=nl+1; write (lines(nl), *) L2

  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (vp%dvar_dcu == 0 .or. .not. vp%exists) cycle
    ix = vp%ix_ele
    call twiss_at_element (u%ring, ix, average = ave)
    cu_model = nint(vp%model / vp%dvar_dcu)
    if (v1_mag%type == hsteer_kick$) then
      beta = ave%a%beta
      phi = f_phi*ave%a%phi
    elseif (v1_mag%type == vsteer_kick$) then
      beta = ave%b%beta
      phi = f_phi*ave%b%phi
    else
      print *, 'INTERNAL ERROR IN: SHOW_STEERINGS_ETC'
      call err_exit
    endif
    mr_model = 1000 * vp%model * sqrt(beta)
    mr_db = 1000 * vp%cu_saved * vp%dvar_dcu * sqrt(beta)
    nl=nl+1; write (lines(nl), fmt) &
            trim(vp%name), mr_model, mr_db, mr_db-mr_model, &
            cu_model, vp%cu_saved, vp%cu_saved-cu_model,  &
            beta, phi, i
  enddo

  nl=nl+1; write (lines(nl), *) L2
  nl=nl+1; write (lines(nl), *) L1

elseif (len_trim(show2_name) == 0 .and. .not. is_steering) then

  L0='         |                                 |       Model      |          |'
  L1='         |    ............' // k_name // '.............. |    Beta      Eta |          |'
  L2='Name     |    Model     Design      Saved  |  X      Y     X  |  d' // k_name // '_dCU | Ix'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *) l0
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l2

  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (.not. vp%exists) cycle
    ix = vp%ix_ele
    nl=nl+1
    if (is_quad) then
      write (lines(nl), '(a11, f9.6, 2f11.6, 2f7.2, f6.2, 1p, e11.2, 0p, i5)') &
        vp%name, vp%model, vp%design, vp%saved, &
        u%beta%x%d(i)%model, u%beta%y%d(i)%model, u%eta%x%d(i)%model, &
        vp%dvar_dcu, i
    else
      call twiss_at_element (u%ring, ix, average = ele_twiss)
      write (lines(nl), '(a11, f9.6, 2f11.6, 2f7.2, f6.2, 1p, e11.2, 0p, i5)') &
        vp%name, vp%model, vp%design, vp%saved, &
        ele_twiss%a%beta, ele_twiss%b%beta, ele_twiss%a%eta, &
        vp%dvar_dcu, i
    endif
  enddo

  nl=nl+1; write (lines(nl), *) l2
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l0

! LIMIT

elseif (show2_name(1:2) == 'LI') then

  l1='           ..............................CU................................         |'
  l2=' Name      Lim_diff  Target  Hi_lim  Lo_lim  Change   Saved   Model  Design Useit | Indx'
  nl=nl+1; lines(nl) = ' '
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  v1_mag%v%cu_lim_diff = 0
  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (.not. vp%good_var) cycle
    if (vp%dvar_dcu == 0) cycle
    vp%cu_lim_diff = min(vp%cu_high_lim-vp%cu_target, vp%cu_target-vp%cu_low_lim)
  enddo

  call re_allocate(this_indexx, size(v1_mag%v))
  call indexx(v1_mag%v%cu_lim_diff, this_indexx)

  do i = size(v1_mag%v), 1, -1
    j = this_indexx(i) - lbound(v1_mag%v, 1) - 1
    vp => v1_mag%v(j)
    if (.not. vp%good_var) cycle
    if (vp%dvar_dcu == 0) cycle
    cu_model = nint((vp%model-vp%design)/vp%dvar_dcu + vp%cu_design)
    nl=nl+1; write (lines(nl), '(1x, a10, i8, i8, 6i8, l6, i6)') &
        vp%name, vp%cu_lim_diff, vp%cu_target, vp%cu_high_lim, vp%cu_low_lim, &
        vp%cu_target-vp%cu_saved, vp%cu_saved, &
        cu_model, vp%cu_design, vp%useit, j
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

! CU

elseif (trim(show2_name) == 'CU') then

  l1='             ..............................CU.............................       |'
  l2=' Name        Target Change   Saved Ref_Sav   Model  Design  Hi_lim  Lo_lim Useit | Indx'
  nl=nl+1; lines(nl) = ' '
  nl=nl+1; lines(nl) = l1
  nl=nl+1; lines(nl) = l2

  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (.not. vp%good_var) cycle
    if (vp%dvar_dcu == 0) then
      cu_model = 0
    else
      cu_model = nint((vp%model-vp%design)/vp%dvar_dcu + vp%cu_design)
    endif
    nl=nl+1; write (lines(nl), '(1x, a10, i8, i7, 6i8, l6, i6)') &
        vp%name, vp%cu_target, vp%cu_target-vp%cu_saved, vp%cu_saved, vp%cu_saved_ref, &
        cu_model, vp%cu_design, vp%cu_high_lim, vp%cu_low_lim, vp%useit, i
  enddo

  nl=nl+1; lines(nl) = l2
  nl=nl+1; lines(nl) = l1

! KL

elseif (trim(show2_name) == 'KL') then

  l1='             ....................' // k_name // '*L...................|'
  l2=' Name        Target   Change    Saved    Model   Design | d' // k_name // 'L_dCU | Indx'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l2

  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (.not. vp%exists) cycle
    target = vp%cu_target * vp%dvar_dcu + vp%base_cu0
    saved = vp%cu_saved * vp%dvar_dcu + vp%base_cu0
    l_ele = u%ring%ele(vp%ix_ele)%value(l$)
    nl=nl+1
    write (lines(nl), '(1x, a10, 5f9.4, es11.2, i7)') &
        vp%name, target*l_ele, target-saved*l_ele, saved*l_ele, &
        vp%model*l_ele, vp%design*l_ele, vp%dvar_dcu*l_ele, i
  enddo

  nl=nl+1; write (lines(nl), *) l2
  nl=nl+1; write (lines(nl), *) l1

! K

elseif (trim(show2_name) == 'K') then

  if (is_steering) then
    l1='             ...................Kick (mrad).............|'
  else
    l1='             ......................K....................|'
  endif
  l2=' Name        Target    Saved   Change    Model   Design |   dK_dCU | Indx'
  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l2

  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (.not. vp%exists) cycle
    target = (vp%cu_target * vp%dvar_dcu + vp%base_cu0) * ff
    saved = (vp%cu_saved * vp%dvar_dcu + vp%base_cu0) * ff
    nl=nl+1
    write (lines(nl), '(1x, a10, 5f9.4, es11.2, i7)') &
        vp%name, target, saved, target-saved, ff*vp%model, &
        ff*vp%design, ff*vp%dvar_dcu, i
  enddo

  nl=nl+1; write (lines(nl), *) l2
  nl=nl+1; write (lines(nl), *) l1

! difference

elseif (index('REF_DIFFERENCE', trim(show2_name)) == 1) then

  if (is_steering) then
    l1='              ............Kick (mrad).......|...........CU..........|'
  else
    l1='              ................K.............|...........CU..........|'
  endif
  l2=' Name            Meas        Ref       Diff |  Meas     Ref    Diff | Indx'

  nl=nl+1; write (lines(nl), *) ' '
  nl=nl+1; write (lines(nl), *) l1
  nl=nl+1; write (lines(nl), *) l2
  nl=nl+3

  do i = lbound(v1_mag%v, 1), ubound(v1_mag%v, 1)
    vp => v1_mag%v(i)
    if (.not. vp%exists) cycle
    saved = ff*(vp%cu_saved * vp%dvar_dcu + vp%base_cu0)
    saved_ref = ff*(vp%cu_saved_ref * vp%dvar_dcu + vp%base_cu0)
    nl=nl+1
    write (lines(nl), '(1x, a10, 3f11.6, 3i8, i7)') &
        vp%name, saved, saved_ref, saved-saved_ref, vp%cu_saved, &
        vp%cu_saved_ref, vp%cu_saved-vp%cu_saved_ref, i
  enddo
  nl=nl+1; write (lines(nl), *) l2
  nl=nl+1; write (lines(nl), *) l1

else
  print *, 'Show ', this_name, ' WHAT?'
endif

end subroutine

!-------------------------------------------------
! contains

subroutine this_weight_out (good_opt, who, weight1, weight2)

real(rp) weight1
real(rp), optional :: weight2
character(*) who
logical good_opt
 
!

if (.not. good_opt) return

nl = nl + 1
if (present(weight2)) then
  write (lines(nl), '(1x, a, 2es11.2)') who, weight1, weight2
else
  write (lines(nl), '(1x, a, 2es11.2)') who, weight1
endif  

end subroutine

end subroutine

!-------------------------------------------------
!+
! Subroutine beam_sigma_calc (ele, modes, sig_a, sig_b, sig_c)
!
! Routine to calculate the contribution to the beam sigma matrix at a point
! due to each mode.
! This calculation is appropriate for circular lattices.
!
! Input:
!   ele   -- Ele_struct: Element defining the point in the lattice.
!   modes -- Normal_modes_struct: Structure holding the normal mode info.
!
! Output:
!-

subroutine beam_sigma_calc (ele, modes, sig_a, sig_b, sig_c)

use bmad_struct
use bmad_interface

implicit none

type (ele_struct) ele
type (normal_modes_struct) modes

real(rp) sig_a(6,6), sig_b(6,6), sig_c(6,6)
real(rp) v_mat(6,6), v_mat_trans(6,6)

!

call mat_make_unit (v_mat)
call make_v_mats (ele, v_mat(1:4,1:4), v_mat_trans(1:4,1:4))
v_mat_trans = transpose(v_mat)

sig_a = 0
sig_a(1, 1:2) = (/  ele%a%beta, -ele%a%alpha /) * modes%a%emittance
sig_a(2, 1:2) = (/ -ele%a%alpha, ele%a%gamma /) * modes%a%emittance
sig_a = matmul(matmul(v_mat, sig_a), v_mat_trans)

sig_b = 0
sig_b(3, 3:4) = (/  ele%b%beta, -ele%b%alpha /) * modes%b%emittance
sig_b(4, 3:4) = (/ -ele%b%alpha, ele%b%gamma /) * modes%b%emittance
sig_b = matmul(matmul(v_mat, sig_b), v_mat_trans)

sig_c = 0
sig_c(1,1) = (ele%x%eta  * modes%sigE_E) ** 2
sig_c(2,2) = (ele%x%etap * modes%sigE_E) ** 2
sig_c(3,3) = (ele%y%eta  * modes%sigE_E) ** 2
sig_c(4,4) = (ele%y%etap * modes%sigE_E) ** 2
sig_c(5,5) = modes%sig_z**2
sig_c(6,6) = modes%sigE_E**2

end subroutine

