next up previous contents
Next: Subroutine Showing a Beam-Beam Up: Small Example Programs Previous: Subroutines to test the   Contents


Subroutine displaying Radiation in PTC

!=================================================================
subroutine  Track_and_with_rad
use madx_ptc_module
use pointer_lattice
use gauss_dis
implicit none
type(layout), pointer :: ALS
real(dp) x(6),sig0(6),energy,deltap,sizes(6,6),bet(3),bet0,betf,t
type(probe) xs0,xs1
type(probe_8) xs
type(real_8) y(6)
type(damapspin) m,m1,m2
type(damap) id
type(internal_state) state,wrapped_state
type(fibre), pointer ::p
integer i,mf,j,n
type(normal_spin) nf
type(normalform) normal
real(dp) stoch(6,6),stochi(6,6),r(3),mat(6,6)

als=>M_U%start
state=default0+radiation0+envelope0

p=>als%start
do i=1,als%n
p%AG=.1d0    ! a_etienne
p=>p%next
enddo

sig0=1.d-6
sig0(6)=1.d-4

call MESS_UP_ALIGNMENT(als,SIG0,3.d0)




 CALL INIT(state,2,0)

 call alloc(m)
 call alloc(m1)
 call alloc(m2)
 call alloc(xs)
 call alloc(nf)
 call alloc(id)
 call alloc(y)
 call alloc(normal)
 nf%stochastic=my_true
 
!!!!  Tracking with a wrapped command     !!!!!!!!!!
 wrapped_state=default0+radiation0
 CALL FIND_ORBIT_x(als,X,wrapped_state,1.0e-5_dp,fibre1=1)
 
    WRITE(6,'(A)') " Closed orbit with Radiation "
    WRITE(6,'(6(1x,E15.8))') x
    call GET_loss(als,energy,deltap)
    write(6,'(a32,2(1x,E15.8))') "Energy loss: GEV and DeltaP/p0c ",energy,deltap

 id=1
 y=x+id
call TRACK_PROBE_X(ALS,y,wrapped_state, FIBRE1=1)
id=y

normal=id   ! Normal is a regular Normalform type
write(6,*) "Tunes "
write(6,*) normal%tune(1:3)
write(6,*) "Damping decrements"
write(6,*) normal%damping(1:3)
!!!!!!!!!!   end of Tracking with  a wrapped command    !!!!!!!!!!




!!!!!!!!!!   Tracking a probe_8 to get the map   !!!!!!!!!!
 state=default0+radiation0+envelope0
 xs0=x     ! closed orbit intializing a probe
 m=1       ! damapspin set to identity
 xs=xs0+m  ! Probe_8 = closed orbit probe + Identity

call track_probe(als,xs,state,fibre1=1)
m=xs       ! damapspin = Probe_8 
nf=m       ! normal_spin = damapspin (Normalization including spin (if present) or radiation
           ! envelope if present. (Spin without radiation)
write(6,*) nf%n%tune(1:3)
write(6,*) nf%emittance
write(6,*) nf%s_ij0(1,1),nf%s_ij0(5,5)
!!!!!!!!!!  end Tracking a probe_8  to get the map     !!!!!!!!!!
mat=m%m
!!!!!!!!!!   Tracking a probe_8 to beam sizes around the ring  !!!!!!!!!!
 call kanalnummer(mf,"beam_xx_inf.dat")
 
 state=default0+radiation0+envelope0
 xs0=x     ! closed orbit intializing a probe
 m=1       ! damapspin set to identity
 xs=xs0+m  ! Probe_8 = closed orbit probe + Identity
 
  p=>als%start 

  call extract_beam_sizes(xs,nf%s_ij0,sizes)
  write(mf,*) p%mag%name, sizes(1,1)

do i=1,als%n
 call track_probe(als,xs,state,fibre1=i,fibre2=i+1) 
 call extract_beam_sizes(xs,nf%s_ij0,sizes)
 p=>p%next
 write(mf,*) p%mag%name, sizes(1,1)
enddo
close(mf)

!!!!!!!!!!   Tracking a probe_8 to beam sizes around the ring using emittances  !!!!!!!!!!

 call kanalnummer(mf,"beam_xx_inf_using_emittances.dat")
 
 state=default0+radiation0+envelope0
 xs0=x     ! closed orbit intializing a probe
 m=nf%a_t  ! damapspin set to A
 xs=xs0+m  ! Probe_8 = closed orbit probe + Identity
 
  p=>als%start 

 bet(1)=(m%m%v(1).sub.'1')**2+(m%m%v(1).sub.'01')**2
 bet(2)=(m%m%v(1).sub.'001')**2+(m%m%v(1).sub.'0001')**2
 bet(3)=(m%m%v(1).sub.'00001')**2+(m%m%v(1).sub.'000001')**2
 sizes(1,1)=(bet(1)*nf%emittance(1)+bet(2)*nf%emittance(2)+bet(3)*nf%emittance(3)) 
 bet0=sizes(1,1)
 write(mf,*) p%mag%name, sizes(1,1)

do i=1,als%n
 call track_probe(als,xs,state,fibre1=i,fibre2=i+1) 
 m=xs
 bet(1)=(m%m%v(1).sub.'1')**2+(m%m%v(1).sub.'01')**2
 bet(2)=(m%m%v(1).sub.'001')**2+(m%m%v(1).sub.'0001')**2
 bet(3)=(m%m%v(1).sub.'00001')**2+(m%m%v(1).sub.'000001')**2
 sizes(1,1)=(bet(1)*nf%emittance(1)+bet(2)*nf%emittance(2)+bet(3)*nf%emittance(3)) 
 p=>p%next
 write(mf,*) p%mag%name, sizes(1,1)
enddo
 betf=sizes(1,1)
 write(mf,*) " log of change after one turn and damping "
 write(mf,*)  log(betf/bet0)/2.d0,nf%n%damping(1)
close(mf)
!!!!!!!!!!   Tracking a probe_8 to beam sizes around the ring using emittances  !!!!!!!!!!

 call kanalnummer(mf,"beam_xx_inf_using_emittances_ripken.dat")
 
 state=default0+radiation0+envelope0
 xs0=x     ! closed orbit intializing a probe
 m=nf%a_t  ! damapspin set to A
 xs=xs0+m  ! Probe_8 = closed orbit probe + Identity
 
  p=>als%start 
 
 id=(m%m.sub.1)**(-1)
 bet(1)=(m%m%v(1).sub.'1')*(id%v(2).sub.'01')-(m%m%v(1).sub.'01')*(id%v(1).sub.'01')
 bet(2)=(m%m%v(1).sub.'001')*(id%v(4).sub.'01')-(m%m%v(1).sub.'0001')*(id%v(3).sub.'01')
 bet(3)=(m%m%v(1).sub.'00001')*(id%v(6).sub.'01')-(m%m%v(1).sub.'000001')*(id%v(5).sub.'01')
 sizes(1,1)=(bet(1)*nf%emittance(1)+bet(2)*nf%emittance(2)+bet(3)*nf%emittance(3)) 
 write(mf,*) p%mag%name, sizes(1,1)

do i=1,als%n
 call track_probe(als,xs,state,fibre1=i,fibre2=i+1) 
 m=xs
 id=(m%m.sub.1)**(-1)
 bet(1)=(m%m%v(1).sub.'1')*(id%v(2).sub.'01')-(m%m%v(1).sub.'01')*(id%v(1).sub.'01')
 bet(2)=(m%m%v(1).sub.'001')*(id%v(4).sub.'01')-(m%m%v(1).sub.'0001')*(id%v(3).sub.'01')
 bet(3)=(m%m%v(1).sub.'00001')*(id%v(6).sub.'01')-(m%m%v(1).sub.'000001')*(id%v(5).sub.'01')
 sizes(1,1)=(bet(1)*nf%emittance(1)+bet(2)*nf%emittance(2)+bet(3)*nf%emittance(3)) 
 p=>p%next
 write(mf,*) p%mag%name, sizes(1,1)
enddo
close(mf)

!!!!!!!!!!   End of Tracking a probe_8 to beam sizes around the ring using emittances and super-lattices functions  !!!!!!!!!!
stoch=nf%stoch
id=stoch
id=id**(-1)
stochi=id

x=0.d0
sizes=0.d0
write(6,*) " Give the number of turns > 1 million"
read(5,*) n
do i=1,n
 x=matmul(mat,x)
 x=matmul(stoch,x)
  do j=1,3
       t=RANF()
       if(t>half) then
          t=one
       else
          t=-one
       endif
   x(2*j-1)= x(2*j-1)+ t*nf%kick(j)
   x(2*j )  = x(2*j)+  t*nf%kick(j)
 enddo
 x=matmul(stochi,x)
 sizes(1,1)=sizes(1,1)+x(1)**2
 sizes(2,2)=sizes(2,2)+x(2)**2
 sizes(1,3)=sizes(1,3)+x(1)*x(3)

enddo

sizes=sizes/n
write(6,*)  sizes(1,1), sizes(2,2), sizes(1,3)
write(6,*)  nf%s_ij0(1,1), nf%s_ij0(2,2), nf%s_ij0(1,3)
 call kill(m)
 call kill(m1)
 call kill(m2)
 call kill(xs)
 call kill(nf)
 call kill(id)
 call kill(y)
 call kill(normal)

end subroutine  Track_and_with_rad



Frank Schmidt 2010-10-15