subroutine ACCEL(fpar,xyz,tt,e,b) c c Compute fields in RF cavity c implicit double precision (a-h,o-z) include 'icommon.inc' c save jsold save first,jt,zadi save tstart,gap,nbox,v0,v1,v2,v3,v4,v5,v6,v7,v8 ! real*8 xyz(3),e(3),b(3) real*8 k,kc,fpar(15) logical first character*10 fname ! data jsold/-1/,tstart/0./,zadi/0./ data first/.true./ c do i=1,3 e(i) = 0. b(i) = 0. end do ! model = fpar(1) ! --------------------------------------------------------- ! ! Initialization of cavity models ! ! --------------------------------------------------------- ! Only do this at start of region if( jsrg.ne.jsold .and. jpar.gt.-1 ) then ! if( phasemodel.eq.1 .or. phasemodel.eq.5 ) then call ACCEL_PHASE(fpar) end if ! ! Following is initialized at start of every region ! if( model .eq. 4 ) then call CAV4_SET(fpar) ! ! ......................................................... else if( model .eq. 5 ) then call SFISH_SET(fpar) ! Superfish RF field ! ! ........................................................ else if( model .eq. 6 ) then ! Induction linac (poly fit) dt = fpar(2) ! time offset in [s] gap = fpar(3) v0 = fpar(5) v1 = fpar(6) v2 = fpar(7) v3 = fpar(8) v4 = fpar(9) v5 = fpar(10) v6 = fpar(11) v7 = fpar(12) v8 = fpar(13) ! if( fpar(4) .eq. 1. ) then ! tpref is the time the ref particle exits the region ! We want tstart to refer to the time it enters tstart = tpref -2.*(tpref-trefmean) + dt end if ! ......................................................... else if( model .eq. 7 ) then ! Induction linac model 2 call WAVEFORM(fpar) ! ......................................................... else if( model .eq. 8 ) then ! Induction linac (table) ! if( fpar(4) .eq. 1. ) then dt = fpar(2) ! time offset in [s] gap = fpar(3) mfile = fpar(5) ! ! tpref is the time the ref particle exits the region ! We want tstart to refer to the time it enters tstart = tpref -2.*(tpref-trefmean) + dt ! if( mfile .gt. 19 ) then ! input waveform file write(fname,80) mfile 80 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open induction linac wave table' iflag = -1 go to 900 end if read(mfile,*) nbox if( nbox.le.0 .or. nbox.gt.mxnbx ) then write(2,*)'Error reading induction linac wave table' iflag = -1 go to 900 end if do i=1,nbox read(mfile,*)tinduc(i),vinduc(i) end do close(mfile) end if ! end of test on mfile ! kfile = fpar(7) ! output diagnostics file if( kfile.gt.19 .and. first ) then write(fname,80) kfile open(unit=kfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open induction linac diagnostic file' iflag = -1 go to 900 end if write(kfile,90) 90 format(t3,'reg',t8,'par',t19,'z',t32,'t',t45,'Ez') dt_il = fpar(8) t0_il = tstart - dt_il first = .false. end if end if ! end of test on fpar(4) ! ! ......................................................... else if( model .eq. 10 ) then ! variable frequency pillbox if( fpar(6) .eq. 1. ) then ! reset parameter ! The first call to this new region comes from a reference particle ! However, the ref particle is not actually tracked. It's z value has ! already been computed as the end of the region. So take region length away ! before computing start position for real particles. zadi = xp(3) - slen end if ! ......................................................... else if( model .eq. 11 ) then ! straight pillbox plus dipole call RFDIP_INIT(fpar) ! ! ......................................................... ! Rectangular cross-section sector RF cavity else if( model .eq. 12 ) then if( htrk .ne. 0.d0 ) then rcurv = 1.d0 / htrk else write(2,*) 'Called ACCEL model 12 with 0 curvature' stop end if ! call INITSC(fpar,phase,rcurv) ! end if ! end of test on model ! ......................................................... jsold = jsrg ! flag initialization done ! end if ! end of test on new jsrg ! --------------------------------------------------------- ! ! Normal tracking starts here ! ! --------------------------------------------------------- rp = SQRT( xyz(1)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then cphi = xyz(1) / rp sphi = xyz(2) / rp else cphi = 1.d0 sphi = 0.d0 end if ! ! ........................................................ if ( model .eq. 1 ) then ! transverse-independent gap mode = fpar(8) ! Set the 5th parameter to be the radius of curvature ! to approximate a rectangular cavity in a circular ICOOL geometry rc = fpar(5) if( mode .eq. 0 ) then ! time-independent e(3) = emagrf else ! time-dependent targ = omrf * tt + phase e(3) = emagrf * COS(targ) end if if( rc .ne. 0.d0 ) then fac = 1.d0 - xyz(1)/rc e(3) = e(3) * fac end if ! ! ........................................................ else if( model .eq. 2 ) then ! TM01p cylindrical resonator mode = fpar(8) ! Set the 5th parameter to be the radius of curvature ! to approximate a rectangular cavity in a circular ICOOL geometry rc = fpar(5) xo = fpar(6) yo = fpar(7) rp = SQRT( (xyz(1)-xo)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then cphi = (xyz(1)-xo) / rp sphi = (xyz(2)-yo) / rp else cphi = 1.d0 sphi = 0.d0 end if ! if( mode .eq. 0 ) then ! TM010 mode call PILLBOX(rp,sphi,cphi,tt,e,b) ! if( rc .ne. 0.d0 ) then fac = 1.d0 - xyz(1)/rc e(3) = e(3) * fac b(1) = b(1) * fac b(2) = b(2) * fac end if ! else ! TM011 mode arg = (omrf/cc)**2 - (pi/slen)**2 ! if( arg .lt. 0. ) then ! non-physical write(2,*)'Non-physical parameters for TM011 cavity' iflag = -1 go to 900 ! else if( arg .eq. 0. ) then ! TM011 degenerate mode rarg = rp / 2. zarg = pi * xyz(3) / slen targ = omrf * tt + phase sz = SIN(zarg) ct = COS(targ) ! e(3) = emagrf * sz * ct er = -pi*emagrf/slen * rarg * COS(zarg) * ct e(1) = er * cphi e(2) = er * sphi bphi = -emagrf*omrf/(cc*cc) * rarg * sz * SIN(targ) b(1) = -bphi * sphi b(2) = bphi * cphi ! else ! TM011 non-degenerate mode kc = SQRT( arg ) rcav = j01 / kc if( rp .gt. rcav ) then iflag = -39 go to 900 end if rarg = kc * rp zarg = pi * xyz(3) / slen targ = omrf * tt + phase r1 = BESSJ1(rarg) sz = SIN(zarg) ct = COS(targ) ! e(3) = emagrf * BESSJ0(rarg) * sz * ct er = -pi*emagrf/kc/slen * r1 * COS(zarg) * ct e(1) = er * cphi e(2) = er * sphi bphi = -emagrf*omrf/(cc*cc*kc) * r1 * sz * SIN(targ) b(1) = -bphi * sphi b(2) = bphi * cphi end if ! end of test on arg end if ! end of test on mode ! ! ........................................................ else if( model .eq. 3 ) then ! travelling wave cavity ! Ref: Handbuch der Physik, V. 44, 1959, p. 345 betw = fpar(8) if( betw.le.0. .or. betw.ge.1. ) then write(2,*)'Phase velocity not in range 0= R2 upot3=1.-(R1/R2)**3 upot4=1.-(R1/R2)**4 upot5=1.-(R1/R2)**5 upot6=1.-(R1/R2)**6 upot7=1.-(R1/R2)**7 upot8=1.-(R1/R2)**8 upot9=1.-(R1/R2)**9 upot10=1.-(R1/R2)**10 upot11=1.-(R1/R2)**11 c tec1=BL_Q1(r,s,XL)-BL_Q1(r,s,-XL) uec1=BL_Z1(r,s,XL)-BL_Z1(r,s,-XL) c tec2=BL_Q2(r,s,XL)-BL_Q2(r,s,-XL) uec2=BL_Z2(r,s,XL)-BL_Z2(r,s,-XL) c tec3=BL_Q3(r,s,XL)-BL_Q3(r,s,-XL) uec3=BL_Z3(r,s,XL)-BL_Z3(r,s,-XL) c tec4=BL_Q4(r,s,XL)-BL_Q4(r,s,-XL) uec4=BL_Z4(r,s,XL)-BL_Z4(r,s,-XL) c tec5=BL_Q5(r,s,XL)-BL_Q5(r,s,-XL) uec5=BL_Z5(r,s,XL)-BL_Z5(r,s,-XL) c tec6=BL_Q6(r,s,XL)-BL_Q6(r,s,-XL) uec6=BL_Z6(r,s,XL)-BL_Z6(r,s,-XL) c tec7=BL_Q7(r,s,XL)-BL_Q7(r,s,-XL) uec7=BL_Z7(r,s,XL)-BL_Z7(r,s,-XL) c tec8=BL_Q8(r,s,XL)-BL_Q8(r,s,-XL) uec8=BL_Z8(r,s,XL)-BL_Z8(r,s,-XL) c tec9=BL_Q9(r,s,XL)-BL_Q9(r,s,-XL) uec9=BL_Z9(r,s,XL)-BL_Z9(r,s,-XL) c VECTOR POTENTIAL ! vector(k,i)=b00/r*R2**3*(tec1*upot3/6. ! + +R2**2*(-tec3*upot5/80. ! + +R2**2*(tec5*upot7/2688. ! + +R2**2*(-tec7*upot9/165888. ! + +R2**2*tec9*upot11/8110080.)))) c there is factor of 1/r in the functions tec* which is explicitly taken out c RADIAL MAGNETIC FIELD br8 = b00/r*R2**3*(tec2*upot3/6. + +R2**2*(-tec4*upot5/80. + +R2**2*(tec6*upot7/2688. + +R2**2*(-tec8*upot9/165888.)))) c LONG MAGNETIC FIELD bz8 = -b00*R2**3*(upot3*(uec1)/6. + +R2**2*(-upot5*(uec3)/80. + +R2**2*(upot7*(uec5)/2688. + +R2**2*(-upot9*(uec7)/165888. + +R2**2*upot11*(uec9)/8110080.)))) c endif ! end testing radial position ! ! br4 = br8 ! sign flip needed for compatibility with SHEET ! bz4 = -bz8 br4 = -qr ! sign flip needed for compatibility with SHEET bz4 = -qz ! end ! ********************************************************* function BBLOCU(a) ! ! Bu current sheet function for radial current block integration ! implicit double precision (a-h,o-z) common/BBLOCDAT/bbu,bbv,bbclen ! call BSHEET(bbu,bbclen,bbv,a,bu,bv) ! bblocu = bu ! end ! ********************************************************* function BBLOCV(a) ! ! Bv current sheet function for radial current block integration ! implicit double precision (a-h,o-z) common/BBLOCDAT/bbu,bbv,bbclen ! call BSHEET(bbu,bbclen,bbv,a,bu,bv) ! bblocv = bv ! end c ********************************************************* subroutine BENTSOL(fpar,xx,e,b) c c Compute field of bent solenoid ! (or region with superimposed solenoid and dipole fields) c implicit double precision (a-h,o-z) include 'icommon.inc' c save jsold,jcold,npts,jt,mxord ! save s0,bs0,by0,by1,by2,by3,by4,by5,bx0,h0 save bsc,bsd,b0c,b0d,a0c,a0d,b1c,b1d,a1c,a1d save b2c,b2d,a2c,a2d,b3c,b3d,a3c,a3d save b4c,b4d,a4c,a4d,b5c,b5d,a5c,a5d save pers,fns,fna0,fnb0,fna1,fnb1,fna2,fnb2, & fna3,fnb3,fna4,fnb4,fna5,fnb5 save fnss,fna0s,fnb0s,fna1s,fnb1s,fna2s,fnb2s, & fna3s,fnb3s,fna4s,fnb4s,fna5s,fnb5s save cs,sn,sold,pref,iord,kf save h,hp,hpp,h3p,h4p,g,gp,gpp,g3p ! real*8 xx(3),e(3),b(3) real*8 lams,lamd,fpar(15),c(5) real*8 ds0(5) ! real*8 s0(2000),bs0(2000),ds0(5),h0(2000) ! dimension by0(2000),by1(2000),by2(2000),by3(2000),by4(2000) ! dimension by5(2000),bx0(2000) ! real*8 cs(0:199),sn(0:199) real*8 bsc(0:199),bsd(0:199) real*8 b0c(0:199),b0d(0:199),a0c(0:199),a0d(0:199) real*8 b1c(0:199),b1d(0:199),a1c(0:199),a1d(0:199) real*8 b2c(0:199),b2d(0:199),a2c(0:199),a2d(0:199) real*8 b3c(0:199),b3d(0:199),a3c(0:199),a3d(0:199) real*8 b4c(0:199),b4d(0:199),a4c(0:199),a4d(0:199) real*8 b5c(0:199),b5d(0:199),a5c(0:199),a5d(0:199) ! character*80 title character*10 fname data jcold/-1/,jsold/-1/,sold/-1./ ! do i=1,3 e(i) = 0. b(i) = 0. end do ! model = fpar(1) ! ! --------------------------------------------------------- if( model .eq. 1 ) then ! first order Bs and Bd bsp0 = fpar(2) ! solenoid strength b0p0 = fpar(3) ! dipole strength kf = fpar(4) ! curvature factor pref = fpar(5) b1p0 = fpar(6) ! quad strength in dipole field ! bsp1 = 0d0 a0p0 = 0d0 a0p1 = 0d0 a1p0 = 0d0 b0p1 = 0d0 if( kf .eq. 2 ) then h = pcharge * b0p0 / pref else h = fpar(7) end if htrk = h hptrk = 0d0 gtrk = 0d0 gptrk = 0d0 ! call BS1ORD(xx(1),xx(2),h,gtrk,b) ! ! --------------------------------------------------------- else if( model .eq. 2 ) then ! Bs and Bd components with hyperbolic s-dependence bsol = fpar(2) ! solenoid strength dip = fpar(3) ! dipole strength quad = fpar(4) ! quad strength in dipole field pref = fpar(5) clens = fpar(6) ! length of solenoid elens = fpar(7) lams = fpar(8) clend = fpar(9) ! length of dipole and quad elend = fpar(10) lamd = fpar(11) bcen = fpar(12) kf = fpar(13) ! curvature factor iord = fpar(15) ! order of expansion ! ! if( iord .ge. 1 ) then ! first order coefficients z = xx(3) - elens bsp0 = bcen + bsol * FTANH(z,clens,lams,0) bsp1 = bsol * FTANH(z,clens,lams,1) z = xx(3) - elend b0p0 = dip * FTANH(z,clend,lamd,0) b0p1 = dip * FTANH(z,clend,lamd,1) b1p0 = quad * FTANH(z,clend,lamd,0) a0p0 = 0d0 a0p1 = 0d0 a1p0 = 0d0 end if ! if( iord .ge. 2 ) then ! second order coefficients z = xx(3) - elens bsp2 = bsol * FTANH(z,clens,lams,2) z = xx(3) - elend b0p2 = dip * FTANH(z,clend,lamd,2) b1p1 = quad * FTANH(z,clend,lamd,1) b2p0 = 0d0 a0p2 = 0d0 a1p1 = 0d0 a2p0 = 0d0 end if ! if( iord .ge. 3 ) then ! third order coefficients z = xx(3) - elens bsp3 = bsol * FTANH(z,clens,lams,3) z = xx(3) - elend b0p3 = dip * FTANH(z,clend,lamd,3) b1p2 = quad * FTANH(z,clend,lamd,2) b2p1 = 0d0 b3p0 = 0d0 a0p3 = 0d0 a1p2 = 0d0 a2p1 = 0d0 a3p0 = 0d0 end if ! if( kf .eq. 2 ) then h = pcharge * b0p0 / pref else h = fpar(14) end if ! htrk = h hptrk = 0d0 gtrk = 0d0 gptrk = 0d0 ! x = xx(1) y = xx(2) ! if( iord .eq. 1 ) then call BS1ORD(x,y,h,gtrk,b) else if( iord .eq. 2 ) then call BS2ORDH(x,y,h,b) else call BS3ORDH(x,y,h,b) end if ! ! --------------------------------------------------------- else if( model .eq. 3 ) then ! Compute field components from tables of on-axis fields mfile = fpar(2) pref = fpar(3) kf = fpar(5) iord = fpar(15) ! if( (kbcell.eq.0.and.jsrg.ne.jsold) .or. & (kbcell.eq.1.and.jcel.ne.jcold) ) then sclbs = fpar(6) sclb0 = fpar(7) sclb1 = fpar(8) sclb2 = fpar(9) sclb3 = fpar(10) sclb4 = fpar(11) sclb5 = fpar(12) scla0 = fpar(13) sclh = fpar(14) write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open bent solensoid field input file' iflag = -1 go to 900 end if write(2,*) ' Opened file ',mfile,' for bent solenoid' ! read(mfile,'(a80)') title read(mfile,*) npts ! do i=1,npts read(mfile,*) soa(i),bsoa(i),b0oa(i),b1oa(i),b2oa(i), & b3oa(i),b4oa(i),b5oa(i),a0oa(i),hoa(i) bsoa(i) = bsoa(i) * sclbs b0oa(i) = b0oa(i) * sclb0 b1oa(i) = b1oa(i) * sclb1 b2oa(i) = b2oa(i) * sclb2 b3oa(i) = b3oa(i) * sclb3 b4oa(i) = b4oa(i) * sclb4 b5oa(i) = b5oa(i) * sclb5 a0oa(i) = a0oa(i) * scla0 hoa(i) = hoa(i) * sclh end do ! a1p0 = 0.d0 ! unused skew multipoles a1p1 = 0.d0 a1p2 = 0.d0 a1p3 = 0d0 a1p4 = 0d0 a2p0 = 0.d0 a2p1 = 0.d0 a2p2 = 0d0 a2p3 = 0d0 a3p0 = 0.d0 a3p1 = 0d0 a3p2 = 0d0 a4p0 = 0d0 a4p1 = 0d0 a5p0 = 0d0 ! close( mfile ) jsold = jsrg jcold = jcel end if ! end of test on kbcell ! ......................................................... s = xx(3) call HUNT(soa,npts,s,jt) if( jt.lt.0 .or. jt.gt.npts ) then ! outside range go to 900 else if( jt .eq. 0 ) then jt = 1 else if( jt .eq. npts ) then jt = npts - 1 end if ! nord = 5 kk = MIN( MAX(jt-nord/2,1),npts-nord ) do i=1,5 ds0(i) = soa(kk+i-1) - soa(jt) end do ds = s - soa(jt) ! if( sclbs .ne. 0d0 ) then call POLCOF(ds0,bsoa(kk),nord,c) bsp0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 bsp1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 bsp2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 bsp3 = 6.*c(4) + 24.*c(5)*ds bsp4 = 24.*c(5) bsp5 = 0d0 else bsp0 = 0.d0 bsp1 = 0.d0 bsp2 = 0.d0 bsp3 = 0.d0 bsp4 = 0d0 bsp5 = 0d0 end if ! if( sclb0 .ne. 0d0 ) then call POLCOF(ds0,b0oa(kk),nord,c) b0p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b0p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 b0p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 b0p3 = 6.*c(4) + 24.*c(5)*ds b0p4 = 24.*c(5) b0p5 = 0d0 else b0p0 = 0.d0 b0p1 = 0.d0 b0p2 = 0.d0 b0p3 = 0.d0 b0p4 = 0d0 b0p5 = 0d0 end if ! if( sclb1 .ne. 0d0 ) then call POLCOF(ds0,b1oa(kk),nord,c) b1p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b1p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 b1p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 b1p3 = 6.*c(4) + 24.*c(5)*ds b1p4 = 24.*c(5) else b1p0 = 0.d0 b1p1 = 0.d0 b1p2 = 0.d0 b1p3 = 0d0 b1p4 = 0d0 end if ! if( sclb2 .ne. 0d0 ) then call POLCOF(ds0,b2oa(kk),nord,c) b2p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b2p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 b2p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 b2p3 = 6.*c(4) + 24.*c(5)*ds else b2p0 = 0.d0 b2p1 = 0.d0 b2p2 = 0d0 b2p3 = 0d0 end if ! if( sclb3 .ne. 0d0 ) then call POLCOF(ds0,b3oa(kk),nord,c) b3p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b3p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 b3p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 else b3p0 = 0.d0 b3p1 = 0d0 b3p2 = 0d0 end if ! if( sclb4 .ne. 0d0 ) then call POLCOF(ds0,b4oa(kk),nord,c) b4p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b4p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 else b4p0 = 0d0 b4p1 = 0d0 end if ! if( sclb5 .ne. 0d0 ) then call POLCOF(ds0,b5oa(kk),nord,c) b5p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 else b5p0 = 0d0 end if ! if( scla0 .ne. 0d0 ) then call POLCOF(ds0,a0oa(kk),nord,c) a0p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 a0p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 a0p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 a0p3 = 6.*c(4) + 24.*c(5)*ds a0p4 = 24.*c(5) a0p5 = 0d0 g = -pcharge * a0p0 / pref gp = 0d0 gpp = 0.d0 g3p = 0d0 else a0p0 = 0.d0 a0p1 = 0.d0 a0p2 = 0.d0 a0p3 = 0.d0 a0p4 = 0d0 a0p5 = 0d0 g = 0d0 gp = 0d0 gpp = 0.d0 g3p = 0d0 end if ! ......................................................... if( kf .eq. 1 ) then h = fpar(4) hp = 0d0 hpp = 0.d0 h3p = 0d0 h4p = 0d0 else if( kf .eq. 2 ) then h = pcharge * b0p0 / pref hp = pcharge * b0p1 / pref hpp = pcharge * b0p2 / pref h3p = pcharge * b0p3 / pref h4p = pcharge * b0p4 / pref else if( kf .eq. 3 ) then call POLCOF(ds0,hoa(kk),nord,c) h = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 hp = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 hpp = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 h3p = 6.*c(4) + 24.*c(5)*ds h4p = 24.*c(5) end if ! htrk = h hptrk = hp gtrk = g gptrk = gp ! ......................................................... x = xx(1) y = xx(2) ! if( iord .eq. 1 ) then ! call BS1ORD(x,y,h,g,b) ! ......................................................... else if( iord .eq. 2 ) then ! if( kf .eq. 1 ) then call BS2ORDH(x,y,h,b) else call BS2ORD(x,y,h,hp,g,gp,b) end if ! ......................................................... ! else if( iord .eq. 3 ) then if( kf .eq. 1 ) then call BS3ORDH(x,y,h,b) else call BS3ORD(x,y,h,hp,hpp,g,gp,gpp,b) end if ! ......................................................... else if( iord .eq. 4 ) then ! if( kf .eq. 1 ) then call BS4ORDH(x,y,h,b) else call BS4ORD(x,y,h,hp,hpp,h3p,g,gp,gpp,g3p,b) end if ! ......................................................... else if( iord .eq. 5 ) then ! call BS5ORD(x,y,h,hp,hpp,h3p,h4p,b) end if ! ! --------------------------------------------------------- else if( model .eq. 4 ) then ! Compute field components from tables of on-axis Fourier coefficients pref = fpar(3) iord = fpar(4) kf = fpar(5) ! if( (kbcell.eq.0.and.jsrg.ne.jsold.and.jsec.eq.1) .or. & (kbcell.eq.1.and.jcel.ne.jcold.and.jsec.eq.1) ) then ! OPEN(UNIT=66,FILE='ICMULT.DAT',STATUS='UNKNOWN') mfile = fpar(2) sclbs = fpar(6) sclb0 = fpar(7) scla0 = fpar(8) sclb1 = fpar(9) scla1 = fpar(10) sclb2 = fpar(11) scla2 = fpar(12) sclb3 = fpar(13) scla3 = fpar(14) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open BSOL multipole input file' iflag = -1 go to 900 end if write(2,*) ' Opened file ',mfile,' for bent solenoid' ! read(mfile,'(a80)') title read(mfile,*) per,fns,fnb0,fna0,fnb1,fna1,fnb2,fna2, & fnb3,fna3,fnb4,fna4,fnb5,fna5 read(mfile,*) mxord read(mfile,'(a80)') title ! do i=0,mxord read(mfile,*) idum,bsc(i),bsd(i), & b0c(i),b0d(i),a0c(i),a0d(i), & b1c(i),b1d(i),a1c(i),a1d(i), & b2c(i),b2d(i),a2c(i),a2d(i), & b3c(i),b3d(i),a3c(i),a3d(i), & b4c(i),b4d(i),a4c(i),a4d(i), & b5c(i),b5d(i),a5c(i),a5d(i) end do ! fns = fns * sclbs fna0 = fna0 * scla0 fnb0 = fnb0 * sclb0 fna1 = fna1 * scla1 fnb1 = fnb1 * sclb1 fna2 = fna2 * scla2 fnb2 = fnb2 * sclb2 fna3 = fna3 * scla3 fnb3 = fnb3 * sclb3 ! close( mfile ) ! jsold = jsrg jcold = jcel end if ! end of test on kbcell ! ! --------------------------------------------------------- else if( model.eq.5 .and. kbcell.eq.1 ) then pref = fpar(3) iord = fpar(4) kf = fpar(5) ! if( jcel.ne.jcold .and. jsec.eq.1 ) then ! at start of cell structure ! Read in data on tapered scaling function mfile = fpar(6) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open TAPER file' iflag = -1 go to 900 end if ! write(2,*)' TAPER for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) ntaprec if( ioc .ne. 0 ) then write(2,*)'Error reading TAPER file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NCELLS = ',ntaprec if( ntaprec .ne. ncelltap ) then write(2,*)'Taper data mismatch' iflag = -1 go to 900 end if ! do i=1,ntaprec read(mfile,*,iostat=ioc)idm,ftap(i) if( ioc .ne. 0 ) then write(2,*)'Error reading TAPER file' iflag = -1 go to 900 end if write(2,'(a,i4,2e12.4)') ' #, ftap:',idm,ftap(i) end do ! close(mfile) ! ! Read in tables of on-axis Fourier coefficients mfile = fpar(2) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open BSOL multipole input file' iflag = -1 go to 900 end if write(2,*) ' Opened file ',mfile,' for BSOL' ! read(mfile,'(a80)') title read(mfile,*) pers,fnss,fnb0s,fna0s,fnb1s,fna1s,fnb2s, & fna2s,fnb3s,fna3s,fnb4s,fna4s,fnb5s,fna5s read(mfile,*) mxord read(mfile,'(a80)') title ! do i=0,mxord read(mfile,*) idum,bsc(i),bsd(i), & b0c(i),b0d(i),a0c(i),a0d(i), & b1c(i),b1d(i),a1c(i),a1d(i), & b2c(i),b2d(i),a2c(i),a2d(i), & b3c(i),b3d(i),a3c(i),a3d(i), & b4c(i),b4d(i),a4c(i),a4d(i), & b5c(i),b5d(i),a5c(i),a5d(i) end do ! close( mfile ) ! ! Open taper diagnostic file mfile = fpar(7) kdiagtap = mfile if( mfile.ge.20 .and. mfile.le.99 ) then write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open TAPER diagnostic file' iflag = -1 go to 900 end if write(mfile,*)' TAPER diagnostic file ',mfile write(mfile,150) 'reg','is','ic','cfld','rfld', & 'fsc','slen','stot','rap','Bz','per','freq', & 'uw','zw','ww','hw' 150 format(t3,a,t8,a,t12,a,t17,a,t24,a,t32,a,t41,a, & t53,a,t64,a,t75,a,t85,a,t95,a,t105,a,t115,a, & t125,a,t135,a) end if ! jctap = 0 jcold = jcel acctap = 0d0 end if ! end of test on jcel ! ! Scale values of all multipole normalizations if( kcel .ne. kcold ) then ! for each cell in structure jctap = jctap + 1 per = pers / ftap(jctap) pertap = per fns = fnss * ftap(jctap) fna0 = fna0s * ftap(jctap) fnb0 = fnb0s * ftap(jctap) fna1 = fna1s * ftap(jctap) fnb1 = fnb1s * ftap(jctap) fna2 = fna2s * ftap(jctap) fnb2 = fnb2s * ftap(jctap) fna3 = fna3s * ftap(jctap) fnb3 = fnb3s * ftap(jctap) fna4 = fna4s * ftap(jctap) fnb4 = fnb4s * ftap(jctap) fna5 = fna5s * ftap(jctap) fnb5 = fnb5s * ftap(jctap) ! kcold = kcel end if ! end of test on kcel ! end if ! end of test on model ! ! --------------------------------------------------------- if( model.eq.4 .or. model.eq.5 ) then ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! Enter here for tracking ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ x = xx(1) y = xx(2) s = xx(3) ! ! OPEN(UNIT=100,FILE='BSOL2D.MAP',STATUS='UNKNOWN') ! X = 0d0 ! DO 1000 IS=1,51 ! S = (IS-1)*0.010006554d0 ! DO 1000 IY=1,18 ! Y = (IY-1)*0.01d0 ! ! Compute multipole moments for this s location if( s .ne. sold ) then do n=0,mxord arg = twopi * n * s / per cs(n) = COS(arg) sn(n) = SIN(arg) end do ! ! Define the longitudinal functions and their derivatives if( fns .ne. 0d0 ) then bsp0 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,0) bsp1 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,1) bsp2 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,2) bsp3 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,3) else bsp0 = 0d0 bsp1 = 0d0 bsp2 = 0d0 bsp3 = 0d0 end if ! if( fnb0 .ne. 0d0 ) then b0p0 = FOURIER(cs,sn,b0c,b0d,mxord,fnb0,per,0) b0p1 = FOURIER(cs,sn,b0c,b0d,mxord,fnb0,per,1) b0p2 = FOURIER(cs,sn,b0c,b0d,mxord,fnb0,per,2) b0p3 = FOURIER(cs,sn,b0c,b0d,mxord,fnb0,per,3) else b0p0 = 0d0 b0p1 = 0d0 b0p2 = 0d0 b0p3 = 0d0 end if ! if( fna0 .ne. 0d0 ) then a0p0 = FOURIER(cs,sn,a0c,a0d,mxord,fna0,per,0) a0p1 = FOURIER(cs,sn,a0c,a0d,mxord,fna0,per,1) a0p2 = FOURIER(cs,sn,a0c,a0d,mxord,fna0,per,2) a0p3 = FOURIER(cs,sn,a0c,a0d,mxord,fna0,per,3) else a0p0 = 0d0 a0p1 = 0d0 a0p2 = 0d0 a0p3 = 0d0 end if ! if( fnb1 .ne. 0.d0 ) then ! add gradient terms b1p0 = FOURIER(cs,sn,b1c,b1d,mxord,fnb1,per,0) b1p1 = FOURIER(cs,sn,b1c,b1d,mxord,fnb1,per,1) b1p2 = FOURIER(cs,sn,b1c,b1d,mxord,fnb1,per,2) else b1p0 = 0.d0 b1p1 = 0.d0 b1p2 = 0.d0 end if ! if( fna1 .ne. 0.d0 ) then ! add gradient terms a1p0 = FOURIER(cs,sn,a1c,a1d,mxord,fna1,per,0) a1p1 = FOURIER(cs,sn,a1c,a1d,mxord,fna1,per,1) a1p2 = FOURIER(cs,sn,a1c,a1d,mxord,fna1,per,2) else a1p0 = 0.d0 a1p1 = 0.d0 a1p2 = 0.d0 end if ! if( fnb2 .ne. 0.d0 ) then ! add sextupole terms b2p0 = FOURIER(cs,sn,b2c,b2d,mxord,fnb2,per,0) b2p1 = FOURIER(cs,sn,b2c,b2d,mxord,fnb2,per,1) else b2p0 = 0.d0 b2p1 = 0.d0 end if ! if( fna2 .ne. 0.d0 ) then ! add sextupole terms a2p0 = FOURIER(cs,sn,a2c,a2d,mxord,fna2,per,0) a2p1 = FOURIER(cs,sn,a2c,a2d,mxord,fna2,per,1) else a2p0 = 0.d0 a2p1 = 0.d0 end if ! if( fnb3 .ne. 0.d0 ) then ! add octupole terms b3p0 = FOURIER(cs,sn,b3c,b3d,mxord,fnb3,per,0) else b3p0 = 0.d0 end if ! if( fna3 .ne. 0.d0 ) then ! add octupole terms a3p0 = FOURIER(cs,sn,a3c,a3d,mxord,fna3,per,0) else a3p0 = 0.d0 end if ! ......................................................... if( iord .ge. 4 ) then ! 4th order terms if( fns .ne. 0.d0 ) then bsp4 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,4) else bsp4 = 0d0 end if ! if( fnb0 .ne. 0.d0 ) then b0p4 = FOURIER(cs,sn,b0c,b0d,mxord,fnb0,per,4) else b0p4 = 0d0 end if ! if( fna0 .ne. 0.d0 ) then a0p4 = FOURIER(cs,sn,a0c,a0d,mxord,fna0,per,4) else a0p4 = 0d0 end if if( fnb1 .ne. 0.d0 ) then b1p3 = FOURIER(cs,sn,b1c,b1d,mxord,fnb1,per,3) else b1p3 = 0.d0 end if ! if( fna1 .ne. 0.d0 ) then a1p3 = FOURIER(cs,sn,a1c,a1d,mxord,fna1,per,3) else a1p3 = 0.d0 end if ! if( fnb2 .ne. 0.d0 ) then b2p2 = FOURIER(cs,sn,b2c,b2d,mxord,fnb2,per,2) b2p3 = FOURIER(cs,sn,b2c,b2d,mxord,fnb2,per,3) else b2p2 = 0.d0 b2p3 = 0.d0 end if ! if( fna2 .ne. 0.d0 ) then a2p2 = FOURIER(cs,sn,a2c,a2d,mxord,fna2,per,2) a2p3 = FOURIER(cs,sn,a2c,a2d,mxord,fna2,per,3) else a2p2 = 0.d0 a2p3 = 0.d0 end if ! if( fnb3 .ne. 0.d0 ) then ! add octupole terms b3p1 = FOURIER(cs,sn,b3c,b3d,mxord,fnb3,per,1) else b3p1 = 0.d0 end if ! if( fna3 .ne. 0.d0 ) then ! add octupole terms a3p1 = FOURIER(cs,sn,a3c,a3d,mxord,fna3,per,1) else a3p1 = 0.d0 end if ! if( fnb4 .ne. 0.d0 ) then ! add n=4 terms b4p0 = FOURIER(cs,sn,b4c,b4d,mxord,fnb4,per,0) b4p1 = FOURIER(cs,sn,b4c,b4d,mxord,fnb4,per,1) else b4p0 = 0.d0 b4p1 = 0.d0 end if ! if( fna4 .ne. 0.d0 ) then ! add n=4 terms a4p0 = FOURIER(cs,sn,a4c,a4d,mxord,fna4,per,0) a4p1 = FOURIER(cs,sn,a4c,a4d,mxord,fna4,per,1) else a4p0 = 0.d0 a4p1 = 0.d0 end if end if ! end of test on iord = 4 ! ......................................................... if( iord .ge. 5 ) then ! 5th order bent fields, h(s) only if( fns .ne. 0.d0 ) then bsp5 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,5) else bsp5 = 0d0 end if ! if( fnb0 .ne. 0.d0 ) then b0p5 = FOURIER(cs,sn,b0c,b0d,mxord,fnb0,per,5) else b0p5 = 0d0 end if ! if( fna0 .ne. 0.d0 ) then a0p5 = FOURIER(cs,sn,a0c,a0d,mxord,fna0,per,5) else a0p5 = 0d0 end if ! if( fnb1 .ne. 0.d0 ) then b1p4 = FOURIER(cs,sn,b1c,b1d,mxord,fnb1,per,4) else b1p4 = 0.d0 end if ! if( fna1 .ne. 0.d0 ) then a1p4 = FOURIER(cs,sn,a1c,a1d,mxord,fna1,per,4) else a1p4 = 0.d0 end if ! if( fnb3 .ne. 0.d0 ) then ! add octupole terms b3p2 = FOURIER(cs,sn,b3c,b3d,mxord,fnb3,per,2) else b3p2 = 0.d0 end if ! if( fna3 .ne. 0.d0 ) then ! add octupole terms a3p2 = FOURIER(cs,sn,a3c,a3d,mxord,fna3,per,2) else a3p2 = 0.d0 end if ! if( fnb5 .ne. 0.d0 ) then ! add n=5 terms b5p0 = FOURIER(cs,sn,b5c,b5d,mxord,fnb5,per,0) else b5p0 = 0.d0 end if ! if( fna5 .ne. 0.d0 ) then ! add n=5 terms a5p0 = FOURIER(cs,sn,a5c,a5d,mxord,fna5,per,0) else a5p0 = 0.d0 end if ! end if ! end of test on iord = 5 ! ......................................................... if( kf .eq. 2 ) then h = pcharge * b0p0 / pref hp = pcharge * b0p1 / pref hpp = pcharge * b0p2 / pref h3p = pcharge * b0p3 / pref h4p = pcharge * b0p4 / pref g = -pcharge * a0p0 / pref gp = -pcharge * a0p1 / pref gpp = -pcharge * a0p2 / pref g3p = -pcharge * a0p3 / pref else h = fpar(15) hp = 0.d0 hpp = 0.d0 h3p = 0.d0 h4p = 0.d0 g = 0.d0 gp = 0.d0 gpp = 0.d0 g3p = 0.d0 end if ! htrk = h ! htrk must have same sign as b0 hptrk = hp gtrk = g gptrk = gp ! sold = s end if ! end of test on S .NE. SOLD ! ......................................................... if( iord .eq. 1 ) then ! call BS1ORD(x,y,h,g,b) ! ! ......................................................... else if( iord .eq. 2 ) then ! if( kf .eq. 2 ) then call BS2ORD(x,y,h,hp,g,gp,b) else call BS2ORDH(x,y,h,b) end if ! ! ......................................................... else if( iord .eq. 3 ) then ! 3rd order bent fields ! if( kf .eq. 2 ) then call BS3ORD(x,y,h,hp,hpp,g,gp,gpp,b) else call BS3ORDH(x,y,h,b) end if ! ! ......................................................... else if( iord .eq. 4 ) then ! 4th order bent fields ! if( kf .eq. 2 ) then call BS4ORD(x,y,h,hp,hpp,h3p,g,gp,gpp,g3p,b) else call BS4ORDH(x,y,h,b) end if ! ! ......................................................... else if( iord .eq. 5 ) then ! 5th order bent fields (x-s only) ! call BS5ORD(x,y,h,hp,hpp,h3p,h4p,b) ! end if ! end of test on iord ! ! WRITE(66,'(14f14.6)') S,BSP0,B0P0,B1P0,B2P0,B3P0,B4P0,B5P0, ! & A0P0,A1P0,A2P0,A3P0,A4P0,A5P0 ! WRITE(100,'(2f7.2,f16.8,3E16.8)') X,Y,S,B !1000 CONTINUE ! end if ! end of test on model ! --------------------------------------------------------- 900 continue end c ********************************************************* subroutine BFIELDINT(jz,jr,kk,z,r,jj,fz,fr) ! ! Explicit 3rd order polynomial interpolation for bz and br ! on 2-dimensional z-r grid, modified from BIPOLYINT ! implicit double precision (a-h,o-z) include 'icommon.inc' ! dimension frtmp(4),fztmp(4),frhold(4),fzhold(4) ! kz = MIN( MAX(jz-kk/2,1), nzgr(jj)-kk ) kr = MIN( MAX(jr-kk/2,1), nrgr(jj)-kk ) ! delta = delzgr(jj) if (kk .eq. 2) then fact = 1d0/(2d0*(delta**2)) else fact = 1d0/(6d0*(delta**3)) end if xtmp = z - zgr(kz+kk/2,jj) do j=1,kk+1 ! interp grid to current z for each r gridpoint do k=1,kk+1 ! get the z grid points for fixed r for 1D interp fzhold(k) = bzgr(kz-1+k,kr-1+j,jj) frhold(k) = brgr(kz-1+k,kr-1+j,jj) end do if (kk .eq. 2) then fztmp(j) = fact*(-xtmp*(delta-xtmp)*fzhold(1) & +2d0*(delta+xtmp)*(delta-xtmp)*fzhold(2) & +xtmp*(delta+xtmp)*fzhold(3)) frtmp(j) = fact*(-xtmp*(delta-xtmp)*frhold(1) & +2d0*(delta+xtmp)*(delta-xtmp)*frhold(2) & +xtmp*(delta+xtmp)*frhold(3)) else fztmp(j) = fact*(-xtmp*(delta-xtmp) & *((2d0*delta-xtmp)*fzhold(1)+(delta+xtmp)*fzhold(4)) & + 3d0*(delta+xtmp)*(2d0*delta-xtmp) & *((delta-xtmp)*fzhold(2)+xtmp*fzhold(3))) frtmp(j) = fact*(-xtmp*(delta-xtmp) & *((2d0*delta-xtmp)*frhold(1)+(delta+xtmp)*frhold(4)) & + 3d0*(delta+xtmp)*(2d0*delta-xtmp) & *((delta-xtmp)*frhold(2)+xtmp*frhold(3))) end if end do ! delta = delrgr(jj) if (kk .eq. 2) then fact = 1d0/(2d0*(delta**2)) else fact = 1d0/(6d0*(delta**3)) end if xtmp = r - rgr(kr+kk/2,jj) if (kk .eq. 2) then fz = fact*(-xtmp*(delta-xtmp)*fztmp(1) & +2d0*(delta+xtmp)*(delta-xtmp)*fztmp(2) & +xtmp*(delta+xtmp)*fztmp(3)) fr = fact*(-xtmp*(delta-xtmp)*frtmp(1) & +2d0*(delta+xtmp)*(delta-xtmp)*frtmp(2) & +xtmp*(delta+xtmp)*frtmp(3)) else fz = fact*(-xtmp*(delta-xtmp) & *((2d0*delta-xtmp)*fztmp(1)+(delta+xtmp)*fztmp(4)) & + 3d0*(delta+xtmp)*(2d0*delta-xtmp) & *((delta-xtmp)*fztmp(2)+xtmp*fztmp(3))) fr = fact*(-xtmp*(delta-xtmp) & *((2d0*delta-xtmp)*frtmp(1)+(delta+xtmp)*frtmp(4)) & + 3d0*(delta+xtmp)*(2d0*delta-xtmp) & *((delta-xtmp)*frtmp(2)+xtmp*frtmp(3))) end if ! return end c ********************************************************* subroutine BILINEAR(jzlo,jrlo,z,rp,kk,bz,br) ! ! Make bilinear interpolation on r-z magnetic field grid ! implicit double precision (a-h,o-z) include 'icommon.inc' ! u = (z - zgr(jzlo,kk)) / delzgr(kk) v = (rp - rgr(jrlo,kk)) / delrgr(kk) br = (1d0-u)*(1d0-v)*brgr(jzlo,jrlo,kk) & + u*(1d0-v)*brgr(jzlo+1,jrlo,kk) & + u*v*brgr(jzlo+1,jrlo+1,kk) & + (1d0-u)*v*brgr(jzlo,jrlo+1,kk) bz = (1d0-u)*(1d0-v)*bzgr(jzlo,jrlo,kk) & + u*(1d0-v)*bzgr(jzlo+1,jrlo,kk) & + u*v*bzgr(jzlo+1,jrlo+1,kk) & + (1d0-u)*v*bzgr(jzlo,jrlo+1,kk) end c ********************************************************* subroutine BLOCK(fpar,xyz,e,b) ! ! Get field from thick solenoidal current block ! implicit double precision (a-h,o-z) include 'icommon.inc' c save jsold,jcold,kcold,jzlo,jrlo,fjf,fjc real*8 xyz(3),e(3),b(3),fpar(15) real*8 fjf(200),fjc(200) character*80 title character*10 fname data bcnorm/4d-7/ ! mu0 / pi data jcold/-1/,jsold/-1/ c model = fpar(1) ! do i=1,3 e(i) = 0. b(i) = 0. end do ! rp = SQRT( xyz(1)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) else phi = 0. end if sphi = sin(phi) cphi = cos(phi) ! ! --------------------------------------------------------- if( model .eq. 1 ) then ! exact field from single block a1 = fpar(3) ! inner radius of block a2 = fpar(4) ! outer radius of block clen = fpar(5) ! length of block cur = fpar(6) ! current density z = xyz(3) - fpar(2) ! if( fpar(7) .eq. 1.d0 ) then call BBLOCS(z,rp,clen,a1,a2,cur,bz,br) else call BBLOCI(z,rp,clen,a1,a2,bu,bv) br = bv * bcnorm * cur * 1.d6 bz = bu * bcnorm * cur * 1.d6 end if if( iflag .ne. 0 ) go to 900 ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! go to 900 end if ! -------------------------------------------------------- ! Initialize for model 2 ! if( model.eq.2 ) then if((kbcell.eq.0 .and. jsrg.ne.jsold ) .or. & (kbcell.eq.1 .and. jcel.ne.jcold ) ) then ! need data file ! mfile = fpar(2) cscale = fpar(11) if( mfile .lt. 0 ) then ! allow short cut for polarity flip mfile = -mfile fcsign = -1. else fcsign = 1. end if write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open BLOCK data file' iflag = -1 go to 900 end if ! write(2,*)' BLOCK SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NBLOCKS = ',nelgr ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),zlelgr(i), & relgr(i),r2elgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK data file' iflag = -1 go to 900 end if curelgr(i) = curelgr(i) * fcsign * cscale write(2,50) idm,dzelgr(i),zlelgr(i), & relgr(i),r2elgr(i),curelgr(i) 50 format(' #,DZ,L,R1,R2,I:',i4,5e12.4) end do ! close(mfile) ! jcold = jcel jsold = jsrg end if ! end of test on kbcell ... end if ! end of test on model = 2 ! End of initialization ! ! -------------------------------------------------------- if( model .eq. 2 ) then ! Exact field from sum of blocks in data file ! do i=1,nelgr clen = zlelgr(i) a1 = relgr(i) a2 = r2elgr(i) cur = curelgr(i) z = xyz(3) - dzelgr(i) - clen/2. ! if( fpar(3) .eq. 1.d0 ) then call BBLOCS(z,rp,clen,a1,a2,cur,bz,br) else call BBLOCI(z,rp,clen,a1,a2,bu,bv) br = bv * bcnorm * cur * 1.d6 bz = bu * bcnorm * cur * 1.d6 end if if( iflag .ne. 0 ) go to 900 ! b(1) = b(1) + br * cphi b(2) = b(2) + br * sphi b(3) = b(3) + bz end do ! go to 900 end if ! ! -------------------------------------------------------- if( model .eq. 3 ) then ! Interpolated field from predefined grid kk = fpar(2) ! grid number interp = fpar(3) ! interpolation model z = xyz(3) call HUNT(zgr(1,kk),nzgr(kk),z,jzlo) ! Watch out for initial point exactly on left edge of grid if( jzlo .eq. 0 ) jzlo = 1 call HUNT(rgr(1,kk),nrgr(kk),rp,jrlo) if( jrlo .eq. 0 ) jrlo = 1 if( z .ge. zgr(nzgr(kk),kk) ) go to 900 ! B=0 beyond grid if( rp .ge. rgr(nrgr(kk),kk) ) go to 900 ! ! if( interp .le. 1 ) then ! bi-linear interpolation call BILINEAR(jzlo,jrlo,z,rp,kk,bz,br) ! else if( interp.eq.2 .or. interp.eq.3 ) then ! Polynomial interpolation ! ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,bzgr,mxzgr,mxrgr,interp, ! & z,rp,bz,df) ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,brgr,mxzgr,mxrgr,interp, ! & z,rp,br,df) ! call BFIELDINT(jzlo,jrlo,interp,z,rp,kk,bz,br) ! end if ! end of test on interpolation ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! go to 900 end if ! ! -------------------------------------------------------- if( model.eq.4 .and. ncelltap.gt.0 ) then if( kbcell.eq.1 .and. jcel.ne.jcold ) then ! Read in data on current blocks mfile = fpar(2) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open BLOCK data file' iflag = -1 go to 900 end if ! write(2,*)' BLOCK SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NBLOCKS = ',nelgr if( nelgr .ne. 15 ) then write(2,*)'Must use 15 blocks for taper calculation' iflag = -1 go to 900 end if ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),zlelgr(i), & relgr(i),r2elgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK data file' iflag = -1 go to 900 end if write(2,150) idm,dzelgr(i),zlelgr(i), & relgr(i),r2elgr(i),curelgr(i) 150 format(' #,DZ,L,R1,R2,I:',i4,5e12.4) end do ! close(mfile) ! ! Read in data on tapered current scaling functions mfile = fpar(9) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open BLOCK taper file' iflag = -1 go to 900 end if ! write(2,*)' BLOCK TAPER for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) ntaprec if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK taper file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NCELLS = ',ntaprec if( ntaprec .ne. ncelltap ) then write(2,*)'Taper data mismatch' iflag = -1 go to 900 end if ! do i=1,ntaprec read(mfile,*,iostat=ioc)idm,fjf(i),fjc(i) if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK taper file' iflag = -1 go to 900 end if write(2,160) idm,fjf(i),fjc(i) 160 format(' #, Jf, Jc:',i4,2e12.4) end do ! close(mfile) ntap = 0 jcold = jcel ! end if ! end of test on jcel ! ! ........................................................ if( kbcell.eq.1 .and. kcel.ne.kcold ) then ! set up scaled grid ntap = ntap + 1 ! current cell location ! kk = fpar(10) ! grid number delzgr(kk) = fpar(3) delrgr(kk) = fpar(4) zglen = fpar(5) rglen = fpar(6) nzgr(kk) = ANINT( zglen / delzgr(kk) + 1. ) nrgr(kk) = ANINT( rglen / delrgr(kk) + 1. ) ! do iz=1,nzgr(kk) if( termout .and. MOD(iz,20).eq.1 ) then write(*,180)' Computing grid:',iz,' /',nzgr(kk),' done' 180 format(a,i5,a,i5,a) end if zgr(iz,kk) = (iz-1) * delzgr(kk) ! do ir=1,nrgr(kk) rgr(ir,kk) = (ir-1) * delrgr(kk) bzgr(iz,ir,kk) = 0. brgr(iz,ir,kk) = 0. ! do 200 ic=1,nelgr clen = zlelgr(ic) a1 = relgr(ic) a2 = r2elgr(ic) cur = curelgr(ic) ! jm2 = ntap - 2 if( jm2 .lt. 1 ) jm2 = 1 jm1 = ntap -1 if( jm1 .lt. 1 ) jm1 = 1 jp1 = ntap + 1 if( jp1 .gt. ntaprec ) jp1 = ntaprec jp2 = ntap + 2 if( jp2 .gt. ntaprec ) jp2 = ntaprec ! if( ic.eq.1 .or. ic.eq.3 ) cur = cur * fjf(jm2) if( ic .eq. 2 ) cur = cur * fjc(jm2) if( ic.eq.4 .or. ic.eq.6 ) cur = cur * fjf(jm1) if( ic .eq. 5 ) cur = cur * fjc(jm1) if( ic.eq.7 .or. ic.eq.9 ) cur = cur * fjf(ntap) if( ic .eq. 8 ) cur = cur * fjc(ntap) if( ic.eq.10 .or. ic.eq.12 ) cur = cur * fjf(jp1) if( ic .eq. 11 ) cur = cur * fjc(jp1) if( ic.eq.13 .or. ic.eq.15 ) cur = cur * fjf(jp2) if( ic .eq. 14 ) cur = cur * fjc(jp2) ! z = zgr(iz,kk) - dzelgr(ic) - clen/2. ! if( fpar(8) .eq. 1.d0 ) then call BBLOCS(z,rgr(ir,kk),clen,a1,a2,cur,bz,br) else call BBLOCI(z,rgr(ir,kk),clen,a1,a2,bu,bv) br = bv * bcnorm * cur * 1.d6 bz = bu * bcnorm * cur * 1.d6 end if if( iflag .ne. 0 ) go to 900 ! bzgr(iz,ir,kk) = bzgr(iz,ir,kk) + bz brgr(iz,ir,kk) = brgr(iz,ir,kk) + br 200 continue ! end do end do ! ncelltap = ncelltap - 1 ! remaining cells kcold = kcel end if ! end of test on kcel end if ! end of test on model=4 ! --------------------------------------------------------- if( model .eq. 4 ) then ! Interpolated field from predefined grid kk = fpar(10) ! grid number interp = fpar(7) ! interpolation model z = xyz(3) call HUNT(zgr(1,kk),nzgr(kk),z,jzlo) ! Watch out for initial point exactly on left edge of grid if( jzlo .eq. 0 ) jzlo = 1 call HUNT(rgr(1,kk),nrgr(kk),rp,jrlo) if( jrlo .eq. 0 ) jrlo = 1 if( z .ge. zgr(nzgr(kk),kk) ) go to 900 ! B=0 beyond grid if( rp .ge. rgr(nrgr(kk),kk) ) go to 900 ! if( interp .le. 1 ) then ! bi-linear interpolation call BILINEAR(jzlo,jrlo,z,rp,kk,bz,br) ! else if( interp.eq.2 .or. interp.eq.3 ) then call BFIELDINT(jzlo,jrlo,interp,z,rp,kk,bz,br) ! end if ! end of test on interpolation ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! go to 900 end if ! 900 continue end c ********************************************************* subroutine BLOOP(a,z,rho,bz,br) ! ! Get field from a single circular current loop (unnormalized) ! c cf. W.R. Smythe, Static and Dynamic Electricity, p.266-7 ! ! a radius of coil ! z distance of observation point to center plane of coil ! rho radius of observation point ! bz longitudinal field component ! br radial field component ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 k ! data p2/1.570796327d0/ ! ! Dont allow observation points on the coil if( z.eq.0.d0 .and. rho.eq.a ) then if(termout) then write(*,*) 'Asking for B on a coil, terminate' write(*,'(a,2f10.4)')'Zo,Rsh: ',z,a write(2,*) 'Asking for B on a coil, terminate' write(2,'(a,2f10.4)')'Zo,Rsh: ',z,a end if stop ! bz = 0.d0 ! br = 0.d0 ! iflag = -59 ! go to 900 end if ! f1 = (a + rho)**2 + z*z k = 4.d0*a*rho / f1 k = SQRT(k) f2 = (a*a+rho**2+z*z) / ((a-rho)**2+z*z) f3 = (a*a-rho**2-z*z) / ((a-rho)**2+z*z) ! call ELLTHREE(p2,c2,k,ele,elk,elp) ele = ELLE(piov2,k) elk = ELLF(piov2,k) if( rho .ne. 0.d0 ) then br = z / (rho*SQRT(f1)) * ( -elk + f2*ele ) else br = 0.d0 end if bz = 1.d0 / SQRT(f1) * ( elk + f3*ele ) ! 900 continue return end ! ********************************************************* double precision function BL_C0(ss,R,S) ! function ! ! This and following BL_ subroutines are used by BBLOCS ! The BL_Ci routines were updated 31 July 2000 ! implicit real*8 (a-h,o-z) c yh1=R**2 + (-S + ss)**2 aux0=(ss - S)*Log(R + Sqrt(yh1)) bl_c0=aux0 end ! ********************************************************* real*8 function BL_C00(ss,R,T,S) ! function ! ! This and following BL_ subroutines are used by BBLOCS ! The BL_Ci routines were updated 31 July 2000 ! implicit real*8 (a-h,o-z) c aux0=(ss - S)*Log((R + Sqrt(R**2 + (ss - S )**2))/ & (T + Sqrt(T**2 + (ss - S )**2))) bl_c00=aux0 end ! ********************************************************* real*8 function BL_C1(ss,R,S) ! 1st derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + ss)**2 aux0=((ss - S)**2/ - ((R + Sqrt(yh1))* - Sqrt(yh1)))+ - Log(R + Sqrt(yh1)) bl_c1=aux0 end ! ********************************************************* real*8 function BL_C1x(ss,R,T,S) ! 1st derv implicit real*8 (a-h,o-z) ! aux0=((ss - S)**2/ - ((R + Sqrt(R**2 + (ss - S)**2))* - Sqrt(R**2 + (ss - S)**2))) + - Log((R + Sqrt(R**2 + (ss - S)**2))/ & (T + Sqrt(T**2 + (ss - S)**2))) - & ((ss - S)**2/ - ((T + Sqrt(T**2 + (ss - S)**2))* - Sqrt(T**2 + (ss - S)**2))) bl_c1x=aux0 end ! ********************************************************* real*8 function BL_C2(ss,R,S) ! 2nd derv implicit real*8 (a-h,o-z) c yh1=R**2 + (-S + ss)**2 aux0= -((-S + ss)**3/ - ((yh1)*(R + Sqrt(yh1))**2))- - (-S + ss)**3/ - ((yh1)**1.5*(R + Sqrt(yh1)))+ - (3.d0*(-S + ss))/ - (Sqrt(yh1)*(R + Sqrt(yh1))) bl_c2=aux0 end ! ********************************************************* real*8 function BL_C3(ss,R,s) ! 3th derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + ss)**2 aux0= (2.d0*(-S + ss)**4)/ - ((yh1)**1.5*(R + Sqrt(yh1))**3)+ - (3.*(-S + ss)**4)/ - ((yh1)**2*(R + Sqrt(yh1))**2) - - (6.*(-S + ss)**2)/ - ((yh1)*(R + Sqrt(yh1))**2) + - (3.*(-S + ss)**4)/ - ((yh1)**2.5*(R + Sqrt(yh1))) - - (6.*(-S + ss)**2)/ - ((yh1)**1.5*(R + Sqrt(yh1))) + - 3./(Sqrt(yh1)*(R + Sqrt(yh1))) bl_c3=aux0 end ! ********************************************************* real*8 function BL_C4(ss,R,s) ! 4th derv implicit real*8 (a-h,o-z) c yh1=R**2 + (-S + ss)**2 aux0= (-6.*(-S + ss)**5)/ - ((yh1)**2*(R + Sqrt(yh1))**4) - - (12.*(-S + ss)**5)/ - ((yh1)**2.5*(R + Sqrt(yh1))**3)+ - (20.*(-S + ss)**3)/ - ((yh1)**1.5*(R + Sqrt(yh1))**3)- - (15.*(-S + ss)**5)/ - ((yh1)**3*(R + Sqrt(yh1))**2)+ - (30.*(-S + ss)**3)/ - ((yh1)**2*(R + Sqrt(yh1))**2) - - (15.*(-S + ss))/ - ((yh1)*(R + Sqrt(yh1))**2) - - (15.*(-S + ss)**5)/ - ((yh1)**3.5*(R + Sqrt(yh1))) + - (30.*(-S + ss)**3)/ - ((yh1)**2.5*(R + Sqrt(yh1))) - - (15.*(-S + ss))/ - ((yh1)**1.5*(R + Sqrt(yh1))) bl_c4=aux0 end ! ********************************************************* real*8 function BL_C5(ss,R,S) ! 5th derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + ss)**2 aux0= (24.*(-S + ss)**6)/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)+ - (60.*(-S + ss)**6)/ - ((yh1)**3*(R + Sqrt(yh1))**4) - - (90.*(-S + ss)**4)/ - ((yh1)**2*(R + Sqrt(yh1))**4) + - (90.*(-S + ss)**6)/ - ((yh1)**3.5*(R + Sqrt(yh1))**3)- - (180.*(-S + ss)**4)/ - ((yh1)**2.5*(R + Sqrt(yh1))**3)+ - (90.*(-S + ss)**2)/ - ((yh1)**1.5*(R + Sqrt(yh1))**3)+ - (105.*(-S + ss)**6)/ - ((yh1)**4*(R + Sqrt(yh1))**2) - - (225.*(-S + ss)**4)/ - ((yh1)**3*(R + Sqrt(yh1))**2) + - (135.*(-S + ss)**2)/ - ((yh1)**2*(R + Sqrt(yh1))**2) - - 15./((yh1)*(R + Sqrt(yh1))**2) + - (105.*(-S + ss)**6)/ - ((yh1)**4.5*(R + Sqrt(yh1))) - - (225.*(-S + ss)**4)/ - ((yh1)**3.5*(R + Sqrt(yh1))) + - (135.*(-S + ss)**2)/ - ((yh1)**2.5*(R + Sqrt(yh1))) - - 15./((yh1)**1.5*(R + Sqrt(yh1))) bl_c5=aux0 end ! ********************************************************* real*8 function BL_C6(s,R,XL) ! 6th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 aux0= (-120.*(-XL + s)**7)/ - ((yh1)**3*(R + Sqrt(yh1))**6) - - (360.*(-XL + s)**7)/ - ((yh1)**3.5*(R + Sqrt(yh1))**5)+ - (504.*(-XL + s)**5)/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)- - (630.*(-XL + s)**7)/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - (1260.*(-XL + s)**5)/ - ((yh1)**3*(R + Sqrt(yh1))**4) - - (630.*(-XL + s)**3)/ - ((yh1)**2*(R + Sqrt(yh1))**4) - - (840.*(-XL + s)**7)/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)+ - (1890.*(-XL + s)**5)/ - ((yh1)**3.5*(R + Sqrt(yh1))**3)- - (1260.*(-XL + s)**3)/ - ((yh1)**2.5*(R + Sqrt(yh1))**3)+ - (210.*(-XL + s))/ - ((yh1)**1.5*(R + Sqrt(yh1))**3)- - (945.*(-XL + s)**7)/ - ((yh1)**5*(R + Sqrt(yh1))**2) + - (2205.*(-XL + s)**5)/ - ((yh1)**4*(R + Sqrt(yh1))**2) - - (1575.*(-XL + s)**3)/ - ((yh1)**3*(R + Sqrt(yh1))**2) + - (315.*(-XL + s))/ - ((yh1)**2*(R + Sqrt(yh1))**2) - - (945.*(-XL + s)**7)/ - ((yh1)**5.5*(R + Sqrt(yh1))) + - (2205.*(-XL + s)**5)/ - ((yh1)**4.5*(R + Sqrt(yh1))) - - (1575.*(-XL + s)**3)/ - ((yh1)**3.5*(R + Sqrt(yh1))) + - (315.*(-XL + s))/ - ((yh1)**2.5*(R + Sqrt(yh1))) bl_c6=aux0 end ! ********************************************************* real*8 function BL_C7(s,R,XL) ! 7th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 aux0= (720.*(-XL + s)**8)/ - ((yh1)**3.5*(R + Sqrt(yh1))**7)+ - (2520.*(-XL + s)**8)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - (3360.*(-XL + s)**6)/ - ((yh1)**3*(R + Sqrt(yh1))**6) + - (5040.*(-XL + s)**8)/ - ((yh1)**4.5*(R + Sqrt(yh1))**5)- - (10080.*(-XL + s)**6)/ - ((yh1)**3.5*(R + Sqrt(yh1))**5)+ - (5040.*(-XL + s)**4)/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)+ - (7560.*(-XL + s)**8)/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - (17640.*(-XL + s)**6)/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - (12600.*(-XL + s)**4)/ - ((yh1)**3*(R + Sqrt(yh1))**4) - - (2520.*(-XL + s)**2)/ - ((yh1)**2*(R + Sqrt(yh1))**4) + - (9450.*(-XL + s)**8)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3)- - (23520.*(-XL + s)**6)/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)+ - (18900.*(-XL + s)**4)/ - ((yh1)**3.5*(R + Sqrt(yh1))**3)- - (5040.*(-XL + s)**2)/ - ((yh1)**2.5*(R + Sqrt(yh1))**3)+ - 210./((yh1)**1.5* - (R + Sqrt(yh1))**3) + - (10395.*(-XL + s)**8)/ - ((yh1)**6*(R + Sqrt(yh1))**2) - - (26460.*(-XL + s)**6)/ - ((yh1)**5*(R + Sqrt(yh1))**2) + - (22050.*(-XL + s)**4)/ - ((yh1)**4*(R + Sqrt(yh1))**2) - - (6300.*(-XL + s)**2)/ - ((yh1)**3*(R + Sqrt(yh1))**2) + - 315./((yh1)**2* - (R + Sqrt(yh1))**2) + - (10395.*(-XL + s)**8)/ - ((yh1)**6.5*(R + Sqrt(yh1))) - - (26460.*(-XL + s)**6)/ - ((yh1)**5.5*(R + Sqrt(yh1))) + - (22050.*(-XL + s)**4)/ - ((yh1)**4.5*(R + Sqrt(yh1))) - - (6300.*(-XL + s)**2)/ - ((yh1)**3.5*(R + Sqrt(yh1))) + - 315./((yh1)**2.5*(R + Sqrt(yh1))) bl_c7=aux0 end ! ********************************************************* real*8 function BL_C8(s,R,XL) ! 8th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 aux0= (-5040.*(-XL + s)**9)/ - ((yh1)**4*(R + Sqrt(yh1))**8) - - (20160.*(-XL + s)**9)/ - ((yh1)**4.5*(R + Sqrt(yh1))**7)+ - (25920.*(-XL + s)**7)/ - ((yh1)**3.5*(R + Sqrt(yh1))**7)- - (45360.*(-XL + s)**9)/ - ((yh1)**5*(R + Sqrt(yh1))**6) + - (90720.*(-XL + s)**7)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - (45360.*(-XL + s)**5)/ - ((yh1)**3*(R + Sqrt(yh1))**6) - - (75600.*(-XL + s)**9)/ - ((yh1)**5.5*(R + Sqrt(yh1))**5)+ - (181440.*(-XL + s)**7)/ - ((yh1)**4.5*(R + Sqrt(yh1))**5)- - (136080.*(-XL + s)**5)/ - ((yh1)**3.5*(R + Sqrt(yh1))**5)+ - (30240.*(-XL + s)**3)/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)- - (103950.*(-XL + s)**9)/ - ((yh1)**6*(R + Sqrt(yh1))**4) + - (272160.*(-XL + s)**7)/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - (238140.*(-XL + s)**5)/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - (75600.*(-XL + s)**3)/ - ((yh1)**3*(R + Sqrt(yh1))**4) - - (5670.*(-XL + s))/ - ((yh1)**2*(R + Sqrt(yh1))**4) - - (124740.*(-XL + s)**9)/ - ((yh1)**6.5*(R + Sqrt(yh1))**3)+ - (340200.*(-XL + s)**7)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3)- - (317520.*(-XL + s)**5)/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)+ - (113400.*(-XL + s)**3)/ - ((yh1)**3.5*(R + Sqrt(yh1))**3)- - (11340.*(-XL + s))/ - ((yh1)**2.5*(R + Sqrt(yh1))**3)- - (135135.*(-XL + s)**9)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + - (374220.*(-XL + s)**7)/ - ((yh1)**6*(R + Sqrt(yh1))**2) - - (357210.*(-XL + s)**5)/ - ((yh1)**5*(R + Sqrt(yh1))**2) + - (132300.*(-XL + s)**3)/ - ((yh1)**4*(R + Sqrt(yh1))**2) - - (14175.*(-XL + s))/ - ((yh1)**3*(R + Sqrt(yh1))**2) - - (135135.*(-XL + s)**9)/ - ((yh1)**7.5*(R + Sqrt(yh1))) + - (374220.*(-XL + s)**7)/ - ((yh1)**6.5*(R + Sqrt(yh1))) - - (357210.*(-XL + s)**5)/ - ((yh1)**5.5*(R + Sqrt(yh1))) + - (132300.*(-XL + s)**3)/ - ((yh1)**4.5*(R + Sqrt(yh1))) - - (14175.*(-XL + s))/ - ((yh1)**3.5*(R + Sqrt(yh1))) bl_c8=aux0 end ! ********************************************************* real*8 function BL_C9(s,R,XL) ! 9th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 fun9= (40320.*(-XL + s)**10)/ - ((yh1)**4.5*(R + Sqrt(yh1))**9)+ - (181440.*(-XL + s)**10)/ - ((yh1)**5*(R + Sqrt(yh1))**8) - - (226800.*(-XL + s)**8)/ - ((yh1)**4*(R + Sqrt(yh1))**8) + - (453600.*(-XL + s)**10)/ - ((yh1)**5.5*(R + Sqrt(yh1))**7)- - (907200.*(-XL + s)**8)/ - ((yh1)**4.5*(R + Sqrt(yh1))**7)+ - (453600.*(-XL + s)**6)/ - ((yh1)**3.5*(R + Sqrt(yh1))**7)+ - (831600.*(-XL + s)**10)/ - ((yh1)**6*(R + Sqrt(yh1))**6) - - (2041200.*(-XL + s)**8)/ - ((yh1)**5*(R + Sqrt(yh1))**6) + - (1587600.*(-XL + s)**6)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - (378000.*(-XL + s)**4)/ - ((yh1)**3*(R + Sqrt(yh1))**6) + - (1247400.*(-XL + s)**10)/ - ((yh1)**6.5*(R + Sqrt(yh1))**5)- - (3402000.*(-XL + s)**8)/ - ((yh1)**5.5*(R + Sqrt(yh1))**5)+ - (3175200.*(-XL + s)**6)/ - ((yh1)**4.5*(R + Sqrt(yh1))**5)- - (1134000.*(-XL + s)**4)/ - ((yh1)**3.5*(R + Sqrt(yh1))**5)+ - (113400.*(-XL + s)**2)/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)+ - (1621620.*(-XL + s)**10)/ - ((yh1)**7*(R + Sqrt(yh1))**4) - - (4677750.*(-XL + s)**8)/ - ((yh1)**6*(R + Sqrt(yh1))**4) + - (4762800.*(-XL + s)**6)/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - (1984500.*(-XL + s)**4)/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - (283500.*(-XL + s)**2)/ - ((yh1)**3*(R + Sqrt(yh1))**4) - - 5670/((yh1)**2* - (R + Sqrt(yh1))**4) + - (1891890.*(-XL + s)**10)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3)- - (5613300.*(-XL + s)**8)/ - ((yh1)**6.5*(R + Sqrt(yh1))**3)+ - (5953500.*(-XL + s)**6)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3)- - (2646000.*(-XL + s)**4)/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)+ - (425250.*(-XL + s)**2)/ - ((yh1)**3.5*(R + Sqrt(yh1))**3)- - 11340./((yh1)**2.5* - (R + Sqrt(yh1))**3) + - (2027025.*(-XL + s)**10)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - - (6081075.*(-XL + s)**8)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + - (6548850.*(-XL + s)**6)/ - ((yh1)**6*(R + Sqrt(yh1))**2)- - (2976750.*(-XL + s)**4)/ - ((yh1)**5*(R + Sqrt(yh1))**2)+ - (496125.*(-XL + s)**2)/ - ((yh1)**4*(R + Sqrt(yh1))**2)- - 14175./((yh1)**3* - (R + Sqrt(yh1))**2) + - (2027025.*(-XL + s)**10)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - - (6081075.*(-XL + s)**8)/ - ((yh1)**7.5*(R + Sqrt(yh1))) + - (6548850.*(-XL + s)**6)/ - ((yh1)**6.5*(R + Sqrt(yh1))) - - (2976750.*(-XL + s)**4)/ - ((yh1)**5.5*(R + Sqrt(yh1))) + - (496125.*(-XL + s)**2)/ - ((yh1)**4.5*(R + Sqrt(yh1))) - - 14175./((yh1)**3.5* - (R + Sqrt(yh1))) bl_c9=fun9 end ! ********************************************************* real*8 function BL_C10(s,R,XL) ! 10th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 fun10= (-362880.*(-XL + s)**11)/ - ((yh1)**5*(R + Sqrt(yh1))**10)- - (1814400.*(-XL + s)**11)/ - ((yh1)**5.5*(R + Sqrt(yh1))**9)+ - (2217600.*(-XL + s)**9)/ - ((yh1)**4.5*(R + Sqrt(yh1))**9)- - (4989600.*(-XL + s)**11)/ - ((yh1)**6*(R + Sqrt(yh1))**8) + - (9979200.*(-XL + s)**9)/ - ((yh1)**5*(R + Sqrt(yh1))**8) - - (4989600.*(-XL + s)**7)/ - ((yh1)**4*(R + Sqrt(yh1))**8) - - (9979200.*(-XL + s)**11)/ - ((yh1)**6.5*(R + Sqrt(yh1))**7)+ - (24948000.*(-XL + s)**9)/ - ((yh1)**5.5*(R + Sqrt(yh1))**7)- - (19958400.*(-XL + s)**7)/ - ((yh1)**4.5*(R + Sqrt(yh1))**7)+ - (4989600.*(-XL + s)**5)/ - ((yh1)**3.5*(R + Sqrt(yh1))**7)- - (16216200.*(-XL + s)**11)/ - ((yh1)**7*(R + Sqrt(yh1))**6) + - (45738000.*(-XL + s)**9)/ - ((yh1)**6*(R + Sqrt(yh1))**6) - - (44906400.*(-XL + s)**7)/ - ((yh1)**5*(R + Sqrt(yh1))**6) + - (17463600.*(-XL + s)**5)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - (2079000.*(-XL + s)**3)/ - ((yh1)**3*(R + Sqrt(yh1))**6) - - (22702680.*(-XL + s)**11)/ - ((yh1)**7.5*(R + Sqrt(yh1))**5)+ - (68607000.*(-XL + s)**9)/ - ((yh1)**6.5*(R + Sqrt(yh1))**5)- - (74844000.*(-XL + s)**7)/ - ((yh1)**5.5*(R + Sqrt(yh1))**5)+ - (34927200.*(-XL + s)**5)/ - ((yh1)**4.5*(R + Sqrt(yh1))**5)- - (6237000.*(-XL + s)**3)/ - ((yh1)**3.5*(R + Sqrt(yh1))**5)+ - (249480.*(-XL + s))/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)- - (28378350.*(-XL + s)**11)/ - ((yh1)**8*(R + Sqrt(yh1))**4) + - (89189100.*(-XL + s)**9)/ - ((yh1)**7*(R + Sqrt(yh1))**4) - - (102910500.*(-XL + s)**7)/ - ((yh1)**6*(R + Sqrt(yh1))**4) + - (52390800.*(-XL + s)**5)/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - (10914750.*(-XL + s)**3)/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - (623700.*(-XL + s))/ - ((yh1)**3*(R + Sqrt(yh1))**4) - - (32432400.*(-XL + s)**11)/ - ((yh1)**8.5*(R + Sqrt(yh1))**3)+ - (104053950.*(-XL + s)**9)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3)- - (123492600.*(-XL + s)**7)/ - ((yh1)**6.5*(R + Sqrt(yh1))**3)+ - (65488500.*(-XL + s)**5)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3)- - (14553000.*(-XL + s)**3)/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)+ - (935550.*(-XL + s))/ - ((yh1)**3.5*(R + Sqrt(yh1))**3)- - (34459425.*(-XL + s)**11)/ - ((yh1)**9*(R + Sqrt(yh1))**2) + - (111486375.*(-XL + s)**9)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - - (133783650.*(-XL + s)**7)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + - (72037350.*(-XL + s)**5)/ - ((yh1)**6*(R + Sqrt(yh1))**2)- - (16372125.*(-XL + s)**3)/ - ((yh1)**5*(R + Sqrt(yh1))**2) + - (1091475.*(-XL + s))/ - ((yh1)**4*(R + Sqrt(yh1))**2) - - (34459425.*(-XL + s)**11)/ - ((yh1)**9.5*(R + Sqrt(yh1))) + - (111486375.*(-XL + s)**9)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - - (133783650.*(-XL + s)**7)/ - ((yh1)**7.5*(R + Sqrt(yh1))) + - (72037350.*(-XL + s)**5)/ - ((yh1)**6.5*(R + Sqrt(yh1))) - - (16372125.*(-XL + s)**3)/ - ((yh1)**5.5*(R + Sqrt(yh1))) + - (1091475.*(-XL + s))/ - ((yh1)**4.5*(R + Sqrt(yh1))) bl_c10=fun10 end ! ! ********************************************************* real*8 function BL_C11(s,R,XL) ! 11th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 fun11= (3628800*(-XL + s)**12)/ -((yh1)**5.5*(R + Sqrt(yh1))**11)+ - (19958400.*(-XL + s)**12)/ - ((yh1)**6*(R + Sqrt(yh1))**10)- - (23950080*(-XL + s)**10)/ - ((yh1)**5*(R + Sqrt(yh1))**10)+ - (59875200.*(-XL + s)**12)/ - ((yh1)**6.5*(R + Sqrt(yh1))**9)- - (119750400.*(-XL + s)**10)/ - ((yh1)**5.5*(R + Sqrt(yh1))**9)+ - (59875200.*(-XL + s)**8)/ - ((yh1)**4.5*(R + Sqrt(yh1))**9)+ - (129729600.*(-XL + s)**12)/ - ((yh1)**7*(R + Sqrt(yh1))**8) - - (329313600.*(-XL + s)**10)/ - ((yh1)**6*(R + Sqrt(yh1))**8) + - (269438400.*(-XL + s)**8)/ - ((yh1)**5*(R + Sqrt(yh1))**8) - - (69854400.*(-XL + s)**6)/ - ((yh1)**4*(R + Sqrt(yh1))**8) + - (227026800.*(-XL + s)**12)/ - ((yh1)**7.5*(R + Sqrt(yh1))**7)- - (658627200.*(-XL + s)**10)/ - ((yh1)**6.5*(R + Sqrt(yh1))**7)+ - (673596000.*(-XL + s)**8)/ - ((yh1)**5.5*(R + Sqrt(yh1))**7)- - (279417600.*(-XL + s)**6)/ - ((yh1)**4.5*(R + Sqrt(yh1))**7)+ - (37422000.*(-XL + s)**4)/ - ((yh1)**3.5*(R + Sqrt(yh1))**7)+ - (340540200.*(-XL + s)**12)/ - ((yh1)**8*(R + Sqrt(yh1))**6) - - (1070269200.*(-XL + s)**10)/ - ((yh1)**7*(R + Sqrt(yh1))**6) + - (1234926000.*(-XL + s)**8)/ - ((yh1)**6*(R + Sqrt(yh1))**6) - - (628689600.*(-XL + s)**6)/ - ((yh1)**5*(R + Sqrt(yh1))**6) + - (130977000.*(-XL + s)**4)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - (7484400.*(-XL + s)**2)/ - ((yh1)**3*(R + Sqrt(yh1))**6) + - (454053600.*(-XL + s)**12)/ - ((yh1)**8.5*(R + Sqrt(yh1))**5)- - (1498376880*(-XL + s)**10)/ - ((yh1)**7.5*(R + Sqrt(yh1))**5)+ - (1852389000.*(-XL + s)**8)/ - ((yh1)**6.5*(R + Sqrt(yh1))**5)- - (1047816000.*(-XL + s)**6)/ - ((yh1)**5.5*(R + Sqrt(yh1))**5)+ - (261954000.*(-XL + s)**4)/ - ((yh1)**4.5*(R + Sqrt(yh1))**5)- - (22453200.*(-XL + s)**2)/ - ((yh1)**3.5*(R + Sqrt(yh1))**5)+ - 249480/ - ((yh1)**2.5*(R + Sqrt(yh1))**5)+ - (551350800.*(-XL + s)**12)/ - ((yh1)**9*(R + Sqrt(yh1))**4) - - (1872971100.*(-XL + s)**10)/ - ((yh1)**8*(R + Sqrt(yh1))**4) + - (2408105700.*(-XL + s)**8)/ - ((yh1)**7*(R + Sqrt(yh1))**4) - - (1440747000.*(-XL + s)**6)/ - ((yh1)**6*(R + Sqrt(yh1))**4) + - (392931000.*(-XL + s)**4)/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - (39293100.*(-XL + s)**2)/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - 623700/((yh1)**3* _ (R + Sqrt(yh1))**4) + - (620269650.*(-XL + s)**12)/ - ((yh1)**9.5*(R + Sqrt(yh1))**3)- - (2140538400.*(-XL + s)**10)/ - ((yh1)**8.5*(R + Sqrt(yh1))**3)+ - (2809456650.*(-XL + s)**8)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3)- - (1728896400.*(-XL + s)**6)/ - ((yh1)**6.5*(R + Sqrt(yh1))**3)+ - (491163750.*(-XL + s)**4)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3)- - (52390800.*(-XL + s)**2)/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)+ - 935550./ - ((yh1)**3.5*(R + Sqrt(yh1))**3)+ - (654729075.*(-XL + s)**12)/ - ((yh1)**10*(R + Sqrt(yh1))**2)- - (2274322050.*(-XL + s)**10)/ - ((yh1)**9*(R + Sqrt(yh1))**2) + - (3010132125.*(-XL + s)**8)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - - (1872971100.*(-XL + s)**6)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + - (540280125.*(-XL + s)**4)/ - ((yh1)**6*(R + Sqrt(yh1))**2) - - (58939650.*(-XL + s)**2)/ - ((yh1)**5*(R + Sqrt(yh1))**2) + - 1091475/ - ((yh1)**4*(R + Sqrt(yh1))**2) + - (654729075.*(-XL + s)**12)/ - ((yh1)**10.5*(R + Sqrt(yh1))) - - (2274322050.*(-XL + s)**10)/ - ((yh1)**9.5*(R + Sqrt(yh1))) + - (3010132125.*(-XL + s)**8)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - - (1872971100.*(-XL + s)**6)/ - ((yh1)**7.5*(R + Sqrt(yh1))) + - (540280125.*(-XL + s)**4)/ - ((yh1)**6.5*(R + Sqrt(yh1))) - - (58939650.*(-XL + s)**2)/ - ((yh1)**5.5*(R + Sqrt(yh1))) + - 1091475./((yh1)**4.5* - (R + Sqrt(yh1))) bl_c11=fun11 end ! ********************************************************* real*8 function BL_C12(s,R,XL) ! 12th derv implicit real*8 (a-h,o-z) real*8 XL ! yh1=R**2 + (-XL + s)**2 fun12= (-39916800.*(-XL + s)**13)/ - ((yh1)**6*(R + Sqrt(yh1))**12) - - (239500800.*(-XL + s)**13)/ - ((yh1)**6.5*(R + Sqrt(yh1))**11)+ - (283046400.*(-XL + s)**11)/ -((yh1)**5.5*(R + Sqrt(yh1))**11)- - (778377600.*(-XL + s)**13)/ - ((yh1)**7*(R + Sqrt(yh1))**10)+ -(1556755200.*(-XL + s)**11)/ - ((yh1)**6*(R + Sqrt(yh1))**10)- -(778377600.*(-XL + s)**9)/ - ((yh1)**5*(R + Sqrt(yh1))**10)- -(1816214400.*(-XL + s)**13)/ -((yh1)**7.5*(R + Sqrt(yh1))**9)+ -(4670265600.*(-XL + s)**11)/ - ((yh1)**6.5*(R + Sqrt(yh1))**9)- -(3891888000.*(-XL + s)**9)/ - ((yh1)**5.5*(R + Sqrt(yh1))**9)+ - (1037836800.*(-XL + s)**7)/ - ((yh1)**4.5*(R + Sqrt(yh1))**9)- -(3405402000.*(-XL + s)**13)/ -((yh1)**8*(R + Sqrt(yh1))**8) + - (10118908800.*(-XL + s)**11)/ - ((yh1)**7*(R + Sqrt(yh1))**8) - - (10702692000.*(-XL + s)**9)/ - ((yh1)**6*(R + Sqrt(yh1))**8) + -(4670265600.*(-XL + s)**7)/ -((yh1)**5*(R + Sqrt(yh1))**8) - - (681080400.*(-XL + s)**5)/ - ((yh1)**4*(R + Sqrt(yh1))**8) - -(5448643200.*(-XL + s)**13)/ -((yh1)**8.5*(R + Sqrt(yh1))**7)+ -(17708090400.*(-XL + s)**11)/ -((yh1)**7.5*(R + Sqrt(yh1))**7)- -(21405384000.*(-XL + s)**9)/ -((yh1)**6.5*(R + Sqrt(yh1))**7)+ -(11675664000.*(-XL + s)**7)/ -((yh1)**5.5*(R + Sqrt(yh1))**7)- -(2724321600.*(-XL + s)**5)/ -((yh1)**4.5*(R + Sqrt(yh1))**7)+ -(194594400.*(-XL + s)**3)/ -((yh1)**3.5*(R + Sqrt(yh1))**7)- -(7718911200.*(-XL + s)**13)/ -((yh1)**9*(R + Sqrt(yh1))**6) + -(26562135600.*(-XL + s)**11)/ -((yh1)**8*(R + Sqrt(yh1))**6) - -(34783749000.*(-XL + s)**9)/ -((yh1)**7*(R + Sqrt(yh1))**6) + -(21405384000.*(-XL + s)**7)/ -((yh1)**6*(R + Sqrt(yh1))**6) - -(6129723600.*(-XL + s)**5)/ -((yh1)**5*(R + Sqrt(yh1))**6) + -(681080400.*(-XL + s)**3)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - -(16216200.*(-XL + s))/ - ((yh1)**3*(R + Sqrt(yh1))**6) - -(9924314400.*(-XL + s)**13)/ -((yh1)**9.5*(R + Sqrt(yh1))**5)+ -(35416180800.*(-XL + s)**11)/ -((yh1)**8.5*(R + Sqrt(yh1))**5)- -(48697248600.*(-XL + s)**9)/ -((yh1)**7.5*(R + Sqrt(yh1))**5)+ -(32108076000.*(-XL + s)**7)/ -((yh1)**6.5*(R + Sqrt(yh1))**5)- -(10216206000.*(-XL + s)**5)/ -((yh1)**5.5*(R + Sqrt(yh1))**5)+ -(1362160800.*(-XL + s)**3)/ -((yh1)**4.5*(R + Sqrt(yh1))**5)- -(48648600.*(-XL + s))/ -((yh1)**3.5*(R + Sqrt(yh1))**5)- -(11785123350.*(-XL + s)**13)/ -((yh1)**10*(R + Sqrt(yh1))**4)+ -(43005362400.*(-XL + s)**11)/ -((yh1)**9*(R + Sqrt(yh1))**4) - -(60871560750.*(-XL + s)**9)/ -((yh1)**8*(R + Sqrt(yh1))**4) + -(41740498800.*(-XL + s)**7)/ -((yh1)**7*(R + Sqrt(yh1))**4) - -(14047283250.*(-XL + s)**5)/ -((yh1)**6*(R + Sqrt(yh1))**4) + -(2043241200.*(-XL + s)**3)/ -((yh1)**5*(R + Sqrt(yh1))**4) - -(85135050.*(-XL + s))/ -((yh1)**4*(R + Sqrt(yh1))**4) - -(13094581500.*(-XL + s)**13)/ -((yh1)**10.5*(R + Sqrt(yh1))**3)+ -(48381032700.*(-XL + s)**11)/ -((yh1)**9.5*(R + Sqrt(yh1))**3)- -(69567498000.*(-XL + s)**9)/ -((yh1)**8.5*(R + Sqrt(yh1))**3)+ - (48697248600.*(-XL + s)**7)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3)- -(16856739900.*(-XL + s)**5)/ -((yh1)**6.5*(R + Sqrt(yh1))**3)+ - (2554051500.*(-XL + s)**3)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3)- -(113513400.*(-XL + s))/ - ((yh1)**4.5*(R + Sqrt(yh1))**3)- - (13749310575.*(-XL + s)**13)/ - ((yh1)**11*(R + Sqrt(yh1))**2)+ - (51068867850.*(-XL + s)**11)/ - ((yh1)**10*(R + Sqrt(yh1))**2)- - (73915466625.*(-XL + s)**9)/ - ((yh1)**9*(R + Sqrt(yh1))**2) + - (52175623500.*(-XL + s)**7)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - -(18261468225.*(-XL + s)**5)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + -(2809456650.*(-XL + s)**3)/ - ((yh1)**6*(R + Sqrt(yh1))**2) - - (127702575.*(-XL + s))/ - ((yh1)**5*(R + Sqrt(yh1))**2) - -(13749310575.*(-XL + s)**13)/ - ((yh1)**11.5*(R + Sqrt(yh1))) + - (51068867850.*(-XL + s)**11)/ -((yh1)**10.5*(R + Sqrt(yh1))) - -(73915466625.*(-XL + s)**9)/ - ((yh1)**9.5*(R + Sqrt(yh1))) + -(52175623500.*(-XL + s)**7)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - -(18261468225.*(-XL + s)**5)/ -((yh1)**7.5*(R + Sqrt(yh1)))+ - (2809456650.*(-XL + s)**3)/ - ((yh1)**6.5*(R + Sqrt(yh1)))- - (127702575.*(-XL + s))/ - ((yh1)**5.5*(R + Sqrt(yh1))) bl_c12=fun12 end ! ********************************************************* real*8 function BL_C13(x,R,S) ! 13th derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + x)**2 fun13= (479001600.*(-S + x)**14)/ - ((yh1)**6.5*(R + Sqrt(yh1))**13) + - (3113510400.*(-S + x)**14)/ - ((yh1)**7*(R + Sqrt(yh1))**12) - - (3632428800.*(-S + x)**12)/ - ((yh1)**6*(R + Sqrt(yh1))**12) + - (10897286400.*(-S + x)**14)/ - ((yh1)**7.5*(R + Sqrt(yh1))**11) - - (21794572800.*(-S + x)**12)/ - ((yh1)**6.5*(R + Sqrt(yh1))**11) + - (10897286400.*(-S + x)**10)/ - ((yh1)**5.5*(R + Sqrt(yh1))**11) + - (27243216000.*(-S + x)**14)/ - ((yh1)**8*(R + Sqrt(yh1))**10) - - (70832361600.*(-S + x)**12)/ - ((yh1)**7*(R + Sqrt(yh1))**10) + - (59935075200.*(-S + x)**10)/ - ((yh1)**6*(R + Sqrt(yh1))**10) - - (16345929600.*(-S + x)**8)/ - ((yh1)**5*(R + Sqrt(yh1))**10) + - (54486432000.*(-S + x)**14)/ - ((yh1)**8.5*(R + Sqrt(yh1))**9) - - (165275510400.*(-S + x)**12)/ - ((yh1)**7.5*(R + Sqrt(yh1))**9) + - (179805225600.*(-S + x)**10)/ - ((yh1)**6.5*(R + Sqrt(yh1))**9) - - (81729648000.*(-S + x)**8)/ - ((yh1)**5.5*(R + Sqrt(yh1))**9) + - (12713500800.*(-S + x)**6)/ - ((yh1)**4.5*(R + Sqrt(yh1))**9) + - (92626934400.*(-S + x)**14)/ - ((yh1)**9*(R + Sqrt(yh1))**8) - - (309891582000.*(-S + x)**12)/ - ((yh1)**8*(R + Sqrt(yh1))**8) + - (389577988800.*(-S + x)**10)/ - ((yh1)**7*(R + Sqrt(yh1))**8) - - (224756532000.*(-S + x)**8)/ - ((yh1)**6*(R + Sqrt(yh1))**8) + - (57210753600.*(-S + x)**6)/ - ((yh1)**5*(R + Sqrt(yh1))**8) - - (4767562800.*(-S + x)**4)/ - ((yh1)**4*(R + Sqrt(yh1))**8) + - (138940401600.*(-S + x)**14)/ - ((yh1)**9.5*(R + Sqrt(yh1))**7) - - (495826531200.*(-S + x)**12)/ - ((yh1)**8.5*(R + Sqrt(yh1))**7) + - (681761480400.*(-S + x)**10)/ - ((yh1)**7.5*(R + Sqrt(yh1))**7) - - (449513064000.*(-S + x)**8)/ - ((yh1)**6.5*(R + Sqrt(yh1))**7) + - (143026884000.*(-S + x)**6)/ - ((yh1)**5.5*(R + Sqrt(yh1))**7) - - (19070251200.*(-S + x)**4)/ - ((yh1)**4.5*(R + Sqrt(yh1))**7) + - (681080400.*(-S + x)**2)/ - ((yh1)**3.5*(R + Sqrt(yh1))**7) + - (188561973600.*(-S + x)**14)/ - ((yh1)**10*(R + Sqrt(yh1))**6) - - (702420919200.*(-S + x)**12)/ - ((yh1)**9*(R + Sqrt(yh1))**6) + - (1022642220600.*(-S + x)**10)/ - ((yh1)**8*(R + Sqrt(yh1))**6) - - (730458729000.*(-S + x)**8)/ - ((yh1)**7*(R + Sqrt(yh1))**6) + - (262215954000.*(-S + x)**6)/ - ((yh1)**6*(R + Sqrt(yh1))**6) - - (42908065200.*(-S + x)**4)/ - ((yh1)**5*(R + Sqrt(yh1))**6) + - (2383781400.*(-S + x)**2)/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - 16216200./ - ((yh1)**3*(R + Sqrt(yh1))**6) + - (235702467000.*(-S + x)**14)/ - ((yh1)**10.5*(R + Sqrt(yh1))**5) - - (903112610400.*(-S + x)**12)/ - ((yh1)**9.5*(R + Sqrt(yh1))**5) + - (1363522960800.*(-S + x)**10)/ - ((yh1)**8.5*(R + Sqrt(yh1))**5) - - (1022642220600.*(-S + x)**8)/ - ((yh1)**7.5*(R + Sqrt(yh1))**5) + - (393323931000.*(-S + x)**6)/ - ((yh1)**6.5*(R + Sqrt(yh1))**5) - - (71513442000.*(-S + x)**4)/ - ((yh1)**5.5*(R + Sqrt(yh1))**5) + - (4767562800.*(-S + x)**2)/ - ((yh1)**4.5*(R + Sqrt(yh1))**5) - - 48648600/ - ((yh1)**3.5*(R + Sqrt(yh1))**5) + - (274986211500.*(-S + x)**14)/ - ((yh1)**11*(R + Sqrt(yh1))**4) - - (1072446224850.*(-S + x)**12)/ - ((yh1)**10*(R + Sqrt(yh1))**4) + - (1655706452400.*(-S + x)**10)/ - ((yh1)**9*(R + Sqrt(yh1))**4) - - (1278302775750.*(-S + x)**8)/ - ((yh1)**8*(R + Sqrt(yh1))**4) + - (511321110300.*(-S + x)**6)/ - ((yh1)**7*(R + Sqrt(yh1))**4) - - (98330982750.*(-S + x)**4)/ - ((yh1)**6*(R + Sqrt(yh1))**4) + - (7151344200.*(-S + x)**2)/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - 85135050/ - ((yh1)**4*(R + Sqrt(yh1))**4) + - (302484832650.*(-S + x)**14)/ - ((yh1)**11.5*(R + Sqrt(yh1))**3) - - (1191606916500.*(-S + x)**12)/ - ((yh1)**10.5*(R + Sqrt(yh1))**3) + - (1862669758950.*(-S + x)**10)/ - ((yh1)**9.5*(R + Sqrt(yh1))**3) - - (1460917458000.*(-S + x)**8)/ - ((yh1)**8.5*(R + Sqrt(yh1))**3) + - (596541295350.*(-S + x)**6)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3) - - (117997179300.*(-S + x)**4)/ - ((yh1)**6.5*(R + Sqrt(yh1))**3) + - (8939180250.*(-S + x)**2)/ - ((yh1)**5.5*(R + Sqrt(yh1))**3) - - 113513400./ - ((yh1)**4.5*(R + Sqrt(yh1))**3) + - (316234143225.*(-S + x)**14)/ - ((yh1)**12*(R + Sqrt(yh1))**2) - - (1251187262325.*(-S + x)**12)/ - ((yh1)**11*(R + Sqrt(yh1))**2) + - (1966151412225.*(-S + x)**10)/ - ((yh1)**10*(R + Sqrt(yh1))**2) - - (1552224799125.*(-S + x)**8)/ - ((yh1)**9*(R + Sqrt(yh1))**2) + - (639151387875.*(-S + x)**6)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - - (127830277575.*(-S + x)**4)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + - (9833098275.*(-S + x)**2)/ - ((yh1)**6*(R + Sqrt(yh1))**2) - - 127702575./ - ((yh1)**5*(R + Sqrt(yh1))**2) + - (316234143225.*(-S + x)**14)/ - ((yh1)**12.5*(R + Sqrt(yh1))) - - (1251187262325.*(-S + x)**12)/ - ((yh1)**11.5*(R + Sqrt(yh1))) + - (1966151412225.*(-S + x)**10)/ - ((yh1)**10.5*(R + Sqrt(yh1))) - - (1552224799125.*(-S + x)**8)/ - ((yh1)**9.5*(R + Sqrt(yh1))) + - (639151387875.*(-S + x)**6)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - - (127830277575.*(-S + x)**4)/ - ((yh1)**7.5*(R + Sqrt(yh1))) + - (9833098275.*(-S + x)**2)/ - ((yh1)**6.5*(R + Sqrt(yh1))) - - 127702575./((yh1)**5.5*(R + Sqrt(R**2 + & (-S + x)**2))) bl_c13=fun13 end ! ********************************************************* real*8 function BL_C14(x,R,S) ! 12th derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + x)**2 fun14=(-6227020800.*(-S + x)**15)/ - ((yh1)**7*(R + Sqrt(yh1))**14) - - (43589145600.*(-S + x)**15)/ - ((yh1)**7.5*(R + Sqrt(yh1))**13) + - (50295168000.*(-S + x)**13)/ - ((yh1)**6.5*(R + Sqrt(yh1))**13) - - (163459296000.*(-S + x)**15)/ - ((yh1)**8*(R + Sqrt(yh1))**12) + - (326918592000.*(-S + x)**13)/ - ((yh1)**7*(R + Sqrt(yh1))**12) - - (163459296000.*(-S + x)**11)/ - ((yh1)**6*(R + Sqrt(yh1))**12) - - (435891456000.*(-S + x)**15)/ - ((yh1)**8.5*(R + Sqrt(yh1))**11) + - (1144215072000.*(-S + x)**13)/ - ((yh1)**7.5*(R + Sqrt(yh1))**11) - - (980755776000.*(-S + x)**11)/ - ((yh1)**6.5*(R + Sqrt(yh1))**11) + - (272432160000.*(-S + x)**9)/ - ((yh1)**5.5*(R + Sqrt(yh1))**11) - - (926269344000.*(-S + x)**15)/ - ((yh1)**9*(R + Sqrt(yh1))**10) + - (2860537680000.*(-S + x)**13)/ - ((yh1)**8*(R + Sqrt(yh1))**10) - - (3187456272000.*(-S + x)**11)/ - ((yh1)**7*(R + Sqrt(yh1))**10) + - (1498376880000.*(-S + x)**9)/ - ((yh1)**6*(R + Sqrt(yh1))**10) - - (245188944000.*(-S + x)**7)/ - ((yh1)**5*(R + Sqrt(yh1))**10) - - (1667284819200.*(-S + x)**15)/ - ((yh1)**9.5*(R + Sqrt(yh1))**9) + - (5721075360000.*(-S + x)**13)/ - ((yh1)**8.5*(R + Sqrt(yh1))**9) - - (7437397968000.*(-S + x)**11)/ - ((yh1)**7.5*(R + Sqrt(yh1))**9) + - (4495130640000.*(-S + x)**9)/ - ((yh1)**6.5*(R + Sqrt(yh1))**9) - - (1225944720000.*(-S + x)**7)/ - ((yh1)**5.5*(R + Sqrt(yh1))**9) + - (114421507200.*(-S + x)**5)/ - ((yh1)**4.5*(R + Sqrt(yh1))**9) - - (2639867630400.*(-S + x)**15)/ - ((yh1)**10*(R + Sqrt(yh1))**8) + - (9725828112000.*(-S + x)**13)/ - ((yh1)**9*(R + Sqrt(yh1))**8) - - (13945121190000.*(-S + x)**11)/ - ((yh1)**8*(R + Sqrt(yh1))**8) + - (9739449720000.*(-S + x)**9)/ - ((yh1)**7*(R + Sqrt(yh1))**8) - - (3371347980000.*(-S + x)**7)/ - ((yh1)**6*(R + Sqrt(yh1))**8) + - (514896782400.*(-S + x)**5)/ - ((yh1)**5*(R + Sqrt(yh1))**8) - - (23837814000.*(-S + x)**3)/ - ((yh1)**4*(R + Sqrt(yh1))**8) - - (3771239472000.*(-S + x)**15)/ - ((yh1)**10.5*(R + Sqrt(yh1))**7) + - (14588742168000.*(-S + x)**13)/ - ((yh1)**9.5*(R + Sqrt(yh1))**7) - - (22312193904000.*(-S + x)**11)/ - ((yh1)**8.5*(R + Sqrt(yh1))**7) + - (17044037010000.*(-S + x)**9)/ - ((yh1)**7.5*(R + Sqrt(yh1))**7) - - (6742695960000.*(-S + x)**7)/ - ((yh1)**6.5*(R + Sqrt(yh1))**7) + - (1287241956000.*(-S + x)**5)/ - ((yh1)**5.5*(R + Sqrt(yh1))**7) - - (95351256000.*(-S + x)**3)/ - ((yh1)**4.5*(R + Sqrt(yh1))**7) + - (1459458000.*(-S + x))/ - ((yh1)**3.5*(R + Sqrt(yh1))**7) - - (4949751807000.*(-S + x)**15)/ - ((yh1)**11*(R + Sqrt(yh1))**6) + - (19799007228000.*(-S + x)**13)/ - ((yh1)**10*(R + Sqrt(yh1))**6) - - (31608941364000.*(-S + x)**11)/ - ((yh1)**9*(R + Sqrt(yh1))**6) + - (25566055515000.*(-S + x)**9)/ - ((yh1)**8*(R + Sqrt(yh1))**6) - - (10956880935000.*(-S + x)**7)/ - ((yh1)**7*(R + Sqrt(yh1))**6) + - (2359943586000.*(-S + x)**5)/ - ((yh1)**6*(R + Sqrt(yh1))**6) - - (214540326000.*(-S + x)**3)/ - ((yh1)**5*(R + Sqrt(yh1))**6) + - (5108103000.*(-S + x))/ - ((yh1)**4*(R + Sqrt(yh1))**6) - - (6049696653000.*(-S + x)**15)/ - ((yh1)**11.5*(R + Sqrt(yh1))**5) + - (24748759035000.*(-S + x)**13)/ - ((yh1)**10.5*(R + Sqrt(yh1))**5) - - (40640067468000.*(-S + x)**11)/ - ((yh1)**9.5*(R + Sqrt(yh1))**5)+ - (34088074020000.*(-S + x)**9)/ - ((yh1)**8.5*(R + Sqrt(yh1))**5)- - (15339633309000.*(-S + x)**7)/ - ((yh1)**7.5*(R + Sqrt(yh1))**5)+ - (3539915379000.*(-S + x)**5)/ - ((yh1)**6.5*(R + Sqrt(yh1))**5)- - (357567210000.*(-S + x)**3)/ - ((yh1)**5.5*(R + Sqrt(yh1))**5)+ - (10216206000.*(-S + x))/ - ((yh1)**4.5*(R + Sqrt(yh1))**5)- - (6957151150950.*(-S + x)**15)/ - ((yh1)**12*(R + Sqrt(yh1))**4) + - (28873552207500.*(-S + x)**13)/ - ((yh1)**11*(R + Sqrt(yh1))**4) - - (48260080118250.*(-S + x)**11)/ - ((yh1)**10*(R + Sqrt(yh1))**4) + - (41392661310000.*(-S + x)**9)/ - ((yh1)**9*(R + Sqrt(yh1))**4) - - (19174541636250.*(-S + x)**7)/ - ((yh1)**8*(R + Sqrt(yh1))**4) + - (4601889992700.*(-S + x)**5)/ - ((yh1)**7*(R + Sqrt(yh1))**4) - - (491654913750.*(-S + x)**3)/ - ((yh1)**6*(R + Sqrt(yh1))**4) + - (15324309000.*(-S + x))/ - ((yh1)**5*(R + Sqrt(yh1))**4) - - (7589619437400.*(-S + x)**15)/ - ((yh1)**12.5*(R + Sqrt(yh1))**3)+ - (31760907428250.*(-S + x)**13)/ - ((yh1)**11.5*(R + Sqrt(yh1))**3)- - (53622311242500.*(-S + x)**11)/ - ((yh1)**10.5*(R + Sqrt(yh1))**3)+ - (46566743973750.*(-S + x)**9)/ - ((yh1)**9.5*(R + Sqrt(yh1))**3) - - (21913761870000.*(-S + x)**7)/ - ((yh1)**8.5*(R + Sqrt(yh1))**3) + - (5368871658150.*(-S + x)**5)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3) - - (589985896500.*(-S + x)**3)/ - ((yh1)**6.5*(R + Sqrt(yh1))**3) + - (19155386250.*(-S + x))/ - ((yh1)**5.5*(R + Sqrt(yh1))**3) - - (7905853580625.*(-S + x)**15)/ - ((yh1)**13*(R + Sqrt(yh1))**2) + - (33204585038625.*(-S + x)**13)/ - ((yh1)**12*(R + Sqrt(yh1))**2) - - (56303426804625.*(-S + x)**11)/ - ((yh1)**11*(R + Sqrt(yh1))**2) + - (49153785305625.*(-S + x)**9)/ - ((yh1)**10*(R + Sqrt(yh1))**2) - - (23283371986875.*(-S + x)**7)/ - ((yh1)**9*(R + Sqrt(yh1))**2) + - (5752362490875.*(-S + x)**5)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - - (639151387875.*(-S + x)**3)/ - ((yh1)**7*(R + Sqrt(yh1))**2) + - (21070924875.*(-S + x))/ - ((yh1)**6*(R + Sqrt(yh1))**2) - - (7905853580625.*(-S + x)**15)/ - ((yh1)**13.5*(R + Sqrt(yh1))) + - (33204585038625.*(-S + x)**13)/ - ((yh1)**12.5*(R + Sqrt(yh1))) - - (56303426804625.*(-S + x)**11)/ - ((yh1)**11.5*(R + Sqrt(yh1))) + - (49153785305625.*(-S + x)**9)/ - ((yh1)**10.5*(R + Sqrt(yh1))) - - (23283371986875.*(-S + x)**7)/ - ((yh1)**9.5*(R + Sqrt(yh1))) + - (5752362490875.*(-S + x)**5)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - - (639151387875.*(-S + x)**3)/ - ((yh1)**7.5*(R + Sqrt(yh1))) + - (21070924875.*(-S + x))/ - ((yh1)**6.5*(R + Sqrt(yh1))) bl_c14=fun14 end ! ********************************************************* real*8 function BL_C15(x,R,S) ! 15th derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + x)**2 fun15=(87178291200.*(-S + x)**16)/ - ((yh1)**7.5*(R + Sqrt(yh1))**15) + - (653837184000.*(-S + x)**16)/ - ((yh1)**8*(R + Sqrt(yh1))**14) - - (747242496000.*(-S + x)**14)/ - ((yh1)**7*(R + Sqrt(yh1))**14) + - (2615348736000.*(-S + x)**16)/ - ((yh1)**8.5*(R + Sqrt(yh1))**13) - - (5230697472000.*(-S + x)**14)/ - ((yh1)**7.5*(R + Sqrt(yh1))**13) + - (2615348736000.*(-S + x)**12)/ -((yh1)**6.5*(R + Sqrt(yh1))**13) + - (7410154752000.*(-S + x)**16)/ -((yh1)**9*(R + Sqrt(yh1))**12) - - (19615115520000.*(-S + x)**14)/ -((yh1)**8*(R + Sqrt(yh1))**12) + - (16999766784000.*(-S + x)**12)/ -((yh1)**7*(R + Sqrt(yh1))**12) - - (4794806016000.*(-S + x)**10)/ -((yh1)**6*(R + Sqrt(yh1))**12) + - (16672848192000.*(-S + x)**16)/ -((yh1)**9.5*(R + Sqrt(yh1))**11) - - (52306974720000.*(-S + x)**14)/ -((yh1)**8.5*(R + Sqrt(yh1))**11) + - (59499183744000.*(-S + x)**12)/ -((yh1)**7.5*(R + Sqrt(yh1))**11) - - (28768836096000.*(-S + x)**10)/ -((yh1)**6.5*(R + Sqrt(yh1))**11) + - (4903778880000.*(-S + x)**8)/ -((yh1)**5.5*(R + Sqrt(yh1))**11) + - (31678411564800.*(-S + x)**16)/ -((yh1)**10*(R + Sqrt(yh1))**10) - - (111152321280000.*(-S + x)**14)/ -((yh1)**9*(R + Sqrt(yh1))**10) + - (148747959360000.*(-S + x)**12)/ -((yh1)**8*(R + Sqrt(yh1))**10) - - (93498717312000.*(-S + x)**10)/ -((yh1)**7*(R + Sqrt(yh1))**10) + - (26970783840000.*(-S + x)**8)/ -((yh1)**6*(R + Sqrt(yh1))**10) - - (2746116172800.*(-S + x)**6)/ -((yh1)**5*(R + Sqrt(yh1))**10) + - (52797352608000.*(-S + x)**16)/ -((yh1)**10.5*(R + Sqrt(yh1))**9) - - (200074178304000.*(-S + x)**14)/ -((yh1)**9.5*(R + Sqrt(yh1))**9) + - (297495918720000.*(-S + x)**12)/ -((yh1)**8.5*(R + Sqrt(yh1))**9) - - (218163673728000.*(-S + x)**10)/ -((yh1)**7.5*(R + Sqrt(yh1))**9) + - (80912351520000.*(-S + x)**8)/ -((yh1)**6.5*(R + Sqrt(yh1))**9) - - (13730580864000.*(-S + x)**6)/ -((yh1)**5.5*(R + Sqrt(yh1))**9) + - (762810048000.*(-S + x)**4)/ -((yh1)**4.5*(R + Sqrt(yh1))**9) + - (79196028912000.*(-S + x)**16)/ -((yh1)**11*(R + Sqrt(yh1))**8) - - (316784115648000.*(-S + x)**14)/ -((yh1)**10*(R + Sqrt(yh1))**8) + - (505743061824000.*(-S + x)**12)/ -((yh1)**9*(R + Sqrt(yh1))**8) - - (409056888240000.*(-S + x)**10)/ -((yh1)**8*(R + Sqrt(yh1))**8) + - (175310094960000.*(-S + x)**8)/ -((yh1)**7*(R + Sqrt(yh1))**8) - - (37759097376000.*(-S + x)**6)/ -((yh1)**6*(R + Sqrt(yh1))**8) + - (3432645216000.*(-S + x)**4)/ -((yh1)**5*(R + Sqrt(yh1))**8) - - (81729648000.*(-S + x)**2)/ -((yh1)**4*(R + Sqrt(yh1))**8) + - (108894539754000.*(-S + x)**16)/ -((yh1)**11.5*(R + Sqrt(yh1))**7) - - (452548736640000.*(-S + x)**14)/ -((yh1)**10.5*(R + Sqrt(yh1))**7) + - (758614592736000.*(-S + x)**12)/ -((yh1)**9.5*(R + Sqrt(yh1))**7) - - (654491021184000.*(-S + x)**10)/ -((yh1)**8.5*(R + Sqrt(yh1))**7) + - (306792666180000.*(-S + x)**8)/ -((yh1)**7.5*(R + Sqrt(yh1))**7) - - (75518194752000.*(-S + x)**6)/ -((yh1)**6.5*(R + Sqrt(yh1))**7) + - (8581613040000.*(-S + x)**4)/ -((yh1)**5.5*(R + Sqrt(yh1))**7) - - (326918592000.*(-S + x)**2)/ -((yh1)**4.5*(R + Sqrt(yh1))**7) + - 1459458000./ -((yh1)**3.5*(R + Sqrt(yh1))**7) + - (139143023019000.*(-S + x)**16)/ -((yh1)**12*(R + Sqrt(yh1))**6) - - (593970216840000.*(-S + x)**14)/ -((yh1)**11*(R + Sqrt(yh1))**6) + - (1029548375856000.*(-S + x)**12)/ -((yh1)**10*(R + Sqrt(yh1))**6) - - (927195613344000.*(-S + x)**10)/ -((yh1)**9*(R + Sqrt(yh1))**6) + - (460188999270000.*(-S + x)**8)/ -((yh1)**8*(R + Sqrt(yh1))**6) - - (122717066472000.*(-S + x)**6)/ -((yh1)**7*(R + Sqrt(yh1))**6) + - (15732957240000.*(-S + x)**4)/ -((yh1)**6*(R + Sqrt(yh1))**6) - - (735566832000.*(-S + x)**2)/ -((yh1)**5*(R + Sqrt(yh1))**6) + - 5108103000./ -((yh1)**4*(R + Sqrt(yh1))**6) + - (166971627622800.*(-S + x)**16)/ -((yh1)**12.5*(R + Sqrt(yh1))**5) - - (725963598360000.*(-S + x)**14)/ -((yh1)**11.5*(R + Sqrt(yh1))**5) + - (1286935469820000.*(-S + x)**12)/ -((yh1)**10.5*(R + Sqrt(yh1))**5) - - (1192108645728000.*(-S + x)**10)/ -((yh1)**9.5*(R + Sqrt(yh1))**5) + - (613585332360000.*(-S + x)**8)/ -((yh1)**8.5*(R + Sqrt(yh1))**5) - - (171803893060800.*(-S + x)**6)/ -((yh1)**7.5*(R + Sqrt(yh1))**5) + - (23599435860000.*(-S + x)**4)/ -((yh1)**6.5*(R + Sqrt(yh1))**5) - - (1225944720000.*(-S + x)**2)/ -((yh1)**5.5*(R + Sqrt(yh1))**5) + - 10216206000./ -((yh1)**4.5*(R + Sqrt(yh1))**5) + - (189740485935000.*(-S + x)**16)/ -((yh1)**13*(R + Sqrt(yh1))**4) - - (834858138114000.*(-S + x)**14)/ -((yh1)**12*(R + Sqrt(yh1))**4) + - (1501424714790000.*(-S + x)**12)/ -((yh1)**11*(R + Sqrt(yh1))**4) - - (1415629016802000.*(-S + x)**10)/ -((yh1)**10*(R + Sqrt(yh1))**4) + - (745067903580000.*(-S + x)**8)/ -((yh1)**9*(R + Sqrt(yh1))**4) - - (214754866326000.*(-S + x)**6)/ -((yh1)**8*(R + Sqrt(yh1))**4) + - (30679266618000.*(-S + x)**4)/ -((yh1)**7*(R + Sqrt(yh1))**4) - - (1685673990000.*(-S + x)**2)/ -((yh1)**6*(R + Sqrt(yh1))**4) + - 15324309000./ -((yh1)**5*(R + Sqrt(yh1))**4) + - (205552193096250.*(-S + x)**16)/ -((yh1)**13.5*(R + Sqrt(yh1))**3) - - (910754332488000.*(-S + x)**14)/ -((yh1)**12.5*(R + Sqrt(yh1))**3) + - (1651567186269000.*(-S + x)**12)/ -((yh1)**11.5*(R + Sqrt(yh1))**3) - - (1572921129780000.*(-S + x)**10)/ -((yh1)**10.5*(R + Sqrt(yh1))**3) + - (838201391527500.*(-S + x)**8)/ -((yh1)**9.5*(R + Sqrt(yh1))**3) - - (245434132944000.*(-S + x)**6)/ -((yh1)**8.5*(R + Sqrt(yh1))**3) + - (35792477721000.*(-S + x)**4)/ -((yh1)**7.5*(R + Sqrt(yh1))**3) - - (2022808788000.*(-S + x)**2)/ -((yh1)**6.5*(R + Sqrt(yh1))**3) + - 19155386250./ -((yh1)**5.5*(R + Sqrt(yh1))**3) + - (213458046676875.*(-S + x)**16)/ -((yh1)**14*(R + Sqrt(yh1))**2) - - (948702429675000.*(-S + x)**14)/ -((yh1)**13*(R + Sqrt(yh1))**2) + - (1726638422008500.*(-S + x)**12)/ -((yh1)**12*(R + Sqrt(yh1))**2) - - (1651567186269000.*(-S + x)**10)/ -((yh1)**11*(R + Sqrt(yh1))**2) + - (884768135501250.*(-S + x)**8)/ -((yh1)**10*(R + Sqrt(yh1))**2) - - (260773766253000.*(-S + x)**6)/ -((yh1)**9*(R + Sqrt(yh1))**2) + - (38349083272500.*(-S + x)**4)/ -((yh1)**8*(R + Sqrt(yh1))**2) - - (2191376187000.*(-S + x)**2)/ -((yh1)**7*(R + Sqrt(yh1))**2) + - 21070924875./ -((yh1)**6*(R + Sqrt(yh1))**2) + - (213458046676875.*(-S + x)**16)/ -((yh1)**14.5*(R + Sqrt(yh1))) - - (948702429675000.*(-S + x)**14)/ -((yh1)**13.5*(R + Sqrt(yh1))) + - (1726638422008500.*(-S + x)**12)/ -((yh1)**12.5*(R + Sqrt(yh1))) - - (1651567186269000.*(-S + x)**10)/ -((yh1)**11.5*(R + Sqrt(yh1))) + - (884768135501250.*(-S + x)**8)/ -((yh1)**10.5*(R + Sqrt(yh1))) - - (260773766253000.*(-S + x)**6)/ -((yh1)**9.5*(R + Sqrt(yh1))) + - (38349083272500.*(-S + x)**4)/ -((yh1)**8.5*(R + Sqrt(yh1))) - - (2191376187000.*(-S + x)**2)/ -((yh1)**7.5*(R + Sqrt(yh1))) + - 21070924875./ - ((yh1)**6.5*(R + Sqrt(yh1))) bl_c15=fun15 end ! ********************************************************* real*8 function BL_C16(x,R,S) ! 16th derv implicit real*8 (a-h,o-z) ! yh1=R**2 + (-S + x)**2 fun16=(-1307674368000.*(-S + x)**17)/ - ((yh1)**8*(R + Sqrt(yh1))**16) - - (10461394944000.*(-S + x)**17)/ - ((yh1)**8.5*(R + Sqrt(yh1))**15) + - (11856247603200.*(-S + x)**15)/ - ((yh1)**7.5*(R + Sqrt(yh1))**15) - - (44460928512000.*(-S + x)**17)/ - ((yh1)**9*(R + Sqrt(yh1))**14) + - (88921857024000.*(-S + x)**15)/ - ((yh1)**8*(R + Sqrt(yh1))**14) - - (44460928512000.*(-S + x)**13)/ - ((yh1)**7*(R + Sqrt(yh1))**14) - - (133382785536000.*(-S + x)**17)/ - ((yh1)**9.5*(R + Sqrt(yh1))**13) + - (355687428096000.*(-S + x)**15)/ - ((yh1)**8.5*(R + Sqrt(yh1))**13) - - (311226499584000.*(-S + x)**13)/ - ((yh1)**7.5*(R + Sqrt(yh1))**13) + - (88921857024000.*(-S + x)**11)/ - ((yh1)**6.5*(R + Sqrt(yh1))**13) - - (316784115648000.*(-S + x)**17)/ - ((yh1)**10*(R + Sqrt(yh1))**12) + - (1007781046272000.*(-S + x)**15)/ - ((yh1)**9*(R + Sqrt(yh1))**12) - - (1167099373440000.*(-S + x)**13)/ - ((yh1)**8*(R + Sqrt(yh1))**12) + - (577992070656000.*(-S + x)**11)/ - ((yh1)**7*(R + Sqrt(yh1))**12) - - (101889627840000.*(-S + x)**9)/ - ((yh1)**6*(R + Sqrt(yh1))**12) - - (633568231296000.*(-S + x)**17)/ - ((yh1)**10.5*(R + Sqrt(yh1))**11) + - (2267507354112000.*(-S + x)**15)/ - ((yh1)**9.5*(R + Sqrt(yh1))**11) - - (3112264995840000.*(-S + x)**13)/ - ((yh1)**8.5*(R + Sqrt(yh1))**11) + - (2022972247296000.*(-S + x)**11)/ - ((yh1)**7.5*(R + Sqrt(yh1))**11) - - (611337767040000.*(-S + x)**9)/ - ((yh1)**6.5*(R + Sqrt(yh1))**11) + - (66691392768000.*(-S + x)**7)/ - ((yh1)**5.5*(R + Sqrt(yh1))**11) - - (1108744404768000.*(-S + x)**17)/ - ((yh1)**11*(R + Sqrt(yh1))**10) + - (4308263972812800.*(-S + x)**15)/ - ((yh1)**10*(R + Sqrt(yh1))**10) - - (6613563116160000.*(-S + x)**13)/ - ((yh1)**9*(R + Sqrt(yh1))**10) + - (5057430618240000.*(-S + x)**11)/ - ((yh1)**8*(R + Sqrt(yh1))**10) - - (1986847742880000.*(-S + x)**9)/ - ((yh1)**7*(R + Sqrt(yh1))**10) + - (366802660224000.*(-S + x)**7)/ -((yh1)**6*(R + Sqrt(yh1))**10) - - (23341987468800.*(-S + x)**5)/ -((yh1)**5*(R + Sqrt(yh1))**10) - - (1742312636064000.*(-S + x)**17)/ -((yh1)**11.5*(R + Sqrt(yh1))**9) + - (7180439954688000.*(-S + x)**15)/ - ((yh1)**10.5*(R + Sqrt(yh1))**9) - - (11904413609088000.*(-S + x)**13)/ - ((yh1)**9.5*(R + Sqrt(yh1))**9) + - (10114861236480000.*(-S + x)**11)/ - ((yh1)**8.5*(R + Sqrt(yh1))**9) - - (4635978066720000.*(-S + x)**9)/ -((yh1)**7.5*(R + Sqrt(yh1))**9) + - (1100407980672000.*(-S + x)**7)/ - ((yh1)**6.5*(R + Sqrt(yh1))**9) - - (116709937344000.*(-S + x)**5)/ - ((yh1)**5.5*(R + Sqrt(yh1))**9) + - (3705077376000.*(-S + x)**3)/ -((yh1)**4.5*(R + Sqrt(yh1))**9) - - (2504574414342000.*(-S + x)**17)/ -((yh1)**12*(R + Sqrt(yh1))**8) + - (10770659932032000.*(-S + x)**15)/ -((yh1)**11*(R + Sqrt(yh1))**8) - - (18848654881056000.*(-S + x)**13)/ - ((yh1)**10*(R + Sqrt(yh1))**8) + - (17195264102016000.*(-S + x)**11)/ - ((yh1)**9*(R + Sqrt(yh1))**8) - - (8692458875100000.*(-S + x)**9)/ - ((yh1)**8*(R + Sqrt(yh1))**8) + - (2384217291456000.*(-S + x)**7)/ - ((yh1)**7*(R + Sqrt(yh1))**8) - - (320952327696000.*(-S + x)**5)/ - ((yh1)**6*(R + Sqrt(yh1))**8) + - (16672848192000.*(-S + x)**3)/ - ((yh1)**5*(R + Sqrt(yh1))**8) - - (173675502000.*(-S + x))/ - ((yh1)**4*(R + Sqrt(yh1))**8) - - (3339432552456000.*(-S + x)**17)/ -((yh1)**12.5*(R + Sqrt(yh1))**7) + - (14809657406544000.*(-S + x)**15)/ -((yh1)**11.5*(R + Sqrt(yh1))**7) - - (26926649830080000.*(-S + x)**13)/ -((yh1)**10.5*(R + Sqrt(yh1))**7) + - (25792896153024000.*(-S + x)**11)/ - ((yh1)**9.5*(R + Sqrt(yh1))**7) - - (13907934200160000.*(-S + x)**9)/ -((yh1)**8.5*(R + Sqrt(yh1))**7) + - (4172380260048000.*(-S + x)**7)/ - ((yh1)**7.5*(R + Sqrt(yh1))**7) - - (641904655392000.*(-S + x)**5)/ - ((yh1)**6.5*(R + Sqrt(yh1))**7) + - (41682120480000.*(-S + x)**3)/ - ((yh1)**5.5*(R + Sqrt(yh1))**7) - - (694702008000.*(-S + x))/ -((yh1)**4.5*(R + Sqrt(yh1))**7) - - (4174290690570000.*(-S + x)**17)/ - ((yh1)**13*(R + Sqrt(yh1))**6) + - (18923451130584000.*(-S + x)**15)/ - ((yh1)**12*(R + Sqrt(yh1))**6) - - (35341227901980000.*(-S + x)**13)/ - ((yh1)**11*(R + Sqrt(yh1))**6) + - (35004644779104000.*(-S + x)**11)/ - ((yh1)**10*(R + Sqrt(yh1))**6) - - (19702906783560000.*(-S + x)**9)/ - ((yh1)**9*(R + Sqrt(yh1))**6) + - (6258570390072000.*(-S + x)**7)/ - ((yh1)**8*(R + Sqrt(yh1))**6) - - (1043095065012000.*(-S + x)**5)/ - ((yh1)**7*(R + Sqrt(yh1))**6) + - (76417220880000.*(-S + x)**3)/ - ((yh1)**6*(R + Sqrt(yh1))**6) - - (1563079518000.*(-S + x))/ - ((yh1)**5*(R + Sqrt(yh1))**6) - - (4933252634310000.*(-S + x)**17)/ -((yh1)**13.5*(R + Sqrt(yh1))**5) + - (22708141356700800.*(-S + x)**15)/ -((yh1)**12.5*(R + Sqrt(yh1))**5) - - (43194834102420000.*(-S + x)**13)/ - ((yh1)**11.5*(R + Sqrt(yh1))**5) + - (43755805973880000.*(-S + x)**11)/ - ((yh1)**10.5*(R + Sqrt(yh1))**5) - - (25332308721720000.*(-S + x)**9)/ - ((yh1)**9.5*(R + Sqrt(yh1))**5) + - (8344760520096000.*(-S + x)**7)/ - ((yh1)**8.5*(R + Sqrt(yh1))**5) - - (1460333091016800.*(-S + x)**5)/ - ((yh1)**7.5*(R + Sqrt(yh1))**5) + - (114625831320000.*(-S + x)**3)/ - ((yh1)**6.5*(R + Sqrt(yh1))**5) - - (2605132530000.*(-S + x))/ - ((yh1)**5.5*(R + Sqrt(yh1))**5) - - (5549909213598750.*(-S + x)**17)/ - ((yh1)**14*(R + Sqrt(yh1))**4) + - (25804706087160000.*(-S + x)**15)/ - ((yh1)**13*(R + Sqrt(yh1))**4) - - (49674059217783000.*(-S + x)**13)/ - ((yh1)**12*(R + Sqrt(yh1))**4) + - (51048440302860000.*(-S + x)**11)/ -((yh1)**11*(R + Sqrt(yh1))**4) - - (30082116607042500.*(-S + x)**9)/ -((yh1)**10*(R + Sqrt(yh1))**4) + - (10132923488688000.*(-S + x)**7)/ - ((yh1)**9*(R + Sqrt(yh1))**4) - - (1825416363771000.*(-S + x)**5)/ - ((yh1)**8*(R + Sqrt(yh1))**4) + - (149013580716000.*(-S + x)**3)/ - ((yh1)**7*(R + Sqrt(yh1))**4) - - (3582057228750.*(-S + x))/ - ((yh1)**6*(R + Sqrt(yh1))**4) - - (5976825306952500.*(-S + x)**17)/ -((yh1)**14.5*(R + Sqrt(yh1))**3) + - (27955098261090000.*(-S + x)**15)/ - ((yh1)**13.5*(R + Sqrt(yh1))**3) - - (54189882783036000.*(-S + x)**13)/ -((yh1)**12.5*(R + Sqrt(yh1))**3) + - (56153284333146000.*(-S + x)**11)/ - ((yh1)**11.5*(R + Sqrt(yh1))**3) - - (33424574007825000.*(-S + x)**9)/ - ((yh1)**10.5*(R + Sqrt(yh1))**3) + - (11399538924774000.*(-S + x)**7)/ - ((yh1)**9.5*(R + Sqrt(yh1))**3) - - (2086190130024000.*(-S + x)**5)/ - ((yh1)**8.5*(R + Sqrt(yh1))**3) + - (173849177502000.*(-S + x)**3)/ - ((yh1)**7.5*(R + Sqrt(yh1))**3) - - (4298468674500.*(-S + x))/ - ((yh1)**6.5*(R + Sqrt(yh1))**3) - - (6190283353629375.*(-S + x)**17)/ - ((yh1)**15*(R + Sqrt(yh1))**2) + - (29030294348055000.*(-S + x)**15)/ - ((yh1)**14*(R + Sqrt(yh1))**2) - - (56447794565662500.*(-S + x)**13)/ - ((yh1)**13*(R + Sqrt(yh1))**2) + - (58705706348289000.*(-S + x)**11)/ - ((yh1)**12*(R + Sqrt(yh1))**2) - - (35095802708216250.*(-S + x)**9)/ - ((yh1)**11*(R + Sqrt(yh1))**2) + - (12032846642817000.*(-S + x)**7)/ - ((yh1)**10*(R + Sqrt(yh1))**2) - - (2216577013150500.*(-S + x)**5)/ - ((yh1)**9*(R + Sqrt(yh1))**2) + - (186266975895000.*(-S + x)**3)/ - ((yh1)**8*(R + Sqrt(yh1))**2) - - (4656674397375.*(-S + x))/ - ((yh1)**7*(R + Sqrt(yh1))**2) - - (6190283353629375.*(-S + x)**17)/ - ((yh1)**15.5*(R + Sqrt(yh1))) + - (29030294348055000.*(-S + x)**15)/ - ((yh1)**14.5*(R + Sqrt(yh1))) - - (56447794565662500.*(-S + x)**13)/ - ((yh1)**13.5*(R + Sqrt(yh1))) + - (58705706348289000.*(-S + x)**11)/ - ((yh1)**12.5*(R + Sqrt(yh1))) - - (35095802708216250.*(-S + x)**9)/ - ((yh1)**11.5*(R + Sqrt(yh1))) + - (12032846642817000.*(-S + x)**7)/ - ((yh1)**10.5*(R + Sqrt(yh1))) - - (2216577013150500.*(-S + x)**5)/ - ((yh1)**9.5*(R + Sqrt(yh1))) + - (186266975895000.*(-S + x)**3)/ - ((yh1)**8.5*(R + Sqrt(yh1))) - - (4656674397375.*(-S + x))/ - ((yh1)**7.5*(R + Sqrt(yh1))) bl_c16=fun16 end ! ********************************************************* real*8 function BL_Q1(r,ss,S) ! 1th order implicit real*8 (a-h,o-z) ! aux0=(-S + ss)/Sqrt(r**2 + (-S + ss)**2) bl_q1=aux0 end ! ********************************************************* real*8 function BL_Q2(r,ss,S) ! 2nd order implicit real*8 (a-h,o-z) ! aux0= -(((S - ss)*(-S + ss))/ - ((r**2 + (-S + ss)**2)**1.5)) - - 1/(Sqrt(r**2 + (-S + ss)**2)) bl_q2=aux0 end ! ********************************************************* real*8 function BL_Q3(r,ss,S) ! 3rd order implicit real*8 (a-h,o-z) ! aux0= (3*(S - ss)*(-S + ss)**2)/ - ((r**2 + (-S + ss)**2)**2.5) - - (S - ss)/((r**2 + (-S + ss)**2)**1.5) + - (2*(-S + ss))/((r**2 + (-S + ss)**2)**1.5) bl_q3=aux0 end ! ********************************************************* real*8 function BL_Q4(r,ss,S) ! 4th order implicit real*8 (a-h,o-z) ! aux0= (-15*(S - ss)*(-S + ss)**3)/ - (r*(r**2 + (-S + ss)**2)**3.5) + - (9*(S - ss)*(-S + ss))/(r*(r**2 + (-S + ss)**2)**2.5) - - (9*(-S + ss)**2)/(r*(r**2 + (-S + ss)**2)**2.5) + - 3/(r*(r**2 + (-S + ss)**2)**1.5) bl_q4=aux0 end ! ********************************************************* real*8 function BL_Q5(r,ss,S) ! 5th order implicit real*8 (a-h,o-z) ! aux0= (105*(S - ss)*(-S + ss)**4)/ - (r*(r**2 + (-S + ss)**2)**4.5) - - (90*(S - ss)*(-S + ss)**2)/ - (r*(r**2 + (-S + ss)**2)**3.5) + - (60*(-S + ss)**3)/(r*(r**2 + (-S + ss)**2)**3.5) + - (9*(S - ss))/(r*(r**2 + (-S + ss)**2)**2.5) - - (36*(-S + ss))/(r*(r**2 + (-S + ss)**2)**2.5) bl_q5=aux0 end ! ********************************************************* real*8 function BL_Q6(r,ss,S) ! 6th order implicit real*8 (a-h,o-z) ! aux0= (-945*(S - ss)*(-S + ss)**5)/ - (r*(r**2 + (-S + ss)**2)**5.5) + - (1050*(S - ss)*(-S + ss)**3)/ - (r*(r**2 + (-S + ss)**2)**4.5) - - (525*(-S + ss)**4)/(r*(r**2 + (-S + ss)**2)**4.5) - - (225*(S - ss)*(-S + ss))/ - (r*(r**2 + (-S + ss)**2)**3.5) + - (450*(-S + ss)**2)/(r*(r**2 + (-S + ss)**2)**3.5) - - 45/(r*(r**2 + (-S + ss)**2)**2.5) bl_q6=aux0 end ! ********************************************************* real*8 function BL_Q7(r,ss,S) ! 7th order implicit real*8 (a-h,o-z) ! aux0= (10395*(S - ss)*(-S + ss)**6)/ - (r*(r**2 + (-S + ss)**2)**6.5) - - (14175*(S - ss)*(-S + ss)**4)/ - (r*(r**2 + (-S + ss)**2)**5.5) + - (5670*(-S + ss)**5)/(r*(r**2 + (-S + ss)**2)**5.5) + - (4725*(S - ss)*(-S + ss)**2)/ - (r*(r**2 + (-S + ss)**2)**4.5) - - (6300*(-S + ss)**3)/(r*(r**2 + (-S + ss)**2)**4.5) - - (225*(S - ss))/(r*(r**2 + (-S + ss)**2)**3.5) + - (1350*(-S + ss))/(r*(r**2 + (-S + ss)**2)**3.5) bl_q7=aux0 end ! ********************************************************* real*8 function BL_Q8(r,ss,S) ! 8th order implicit real*8 (a-h,o-z) ! aux0= (-135135*(S - ss)*(-S + ss)**7)/ - (r*(r**2 + (-S + ss)**2)**7.5) + - (218295*(S - ss)*(-S + ss)**5)/ - (r*(r**2 + (-S + ss)**2)**6.5) - - (72765*(-S + ss)**6)/(r*(r**2 + (-S + ss)**2)**6.5) - - (99225*(S - ss)*(-S + ss)**3)/ - (r*(r**2 + (-S + ss)**2)**5.5) + - (99225*(-S + ss)**4)/(r*(r**2 + (-S + ss)**2)**5.5) + - (11025*(S - ss)*(-S + ss))/ - (r*(r**2 + (-S + ss)**2)**4.5) - - (33075*(-S + ss)**2)/(r*(r**2 + (-S + ss)**2)**4.5) + - 1575/(r*(r**2 + (-S + ss)**2)**3.5) bl_q8=aux0 end ! ********************************************************* real*8 function BL_Q9(r,ss,s) ! 9th order implicit real*8 (a-h,o-z) ! aux0= (2027025*(S - ss)*(-S + ss)**8)/ - (r*(r**2 + (-S + ss)**2)**8.5) - - (3783780*(S - ss)*(-S + ss)**6)/ - (r*(r**2 + (-S + ss)**2)**7.5) + - (1081080*(-S + ss)**7)/(r*(r**2 + (-S + ss)**2)**7.5) + - (2182950*(S - ss)*(-S + ss)**4)/ - (r*(r**2 + (-S + ss)**2)**6.5) - - (1746360*(-S + ss)**5)/(r*(r**2 + (-S + ss)**2)**6.5) - - (396900*(S - ss)*(-S + ss)**2)/ - (r*(r**2 + (-S + ss)**2)**5.5) + - (793800*(-S + ss)**3)/(r*(r**2 + (-S + ss)**2)**5.5) + - (11025*(S - ss))/(r*(r**2 + (-S + ss)**2)**4.5) - - (88200*(-S + ss))/(r*(r**2 + (-S + ss)**2)**4.5) bl_q9=aux0 end ! ********************************************************* real*8 function BL_Z1(r,ss,s) ! 1th order implicit real*8 (a-h,o-z) ! aux0=(-S + ss)/(r**2 + (-S + ss)**2)**1.5 bl_z1=aux0 end ! ********************************************************* real*8 function BL_Z2(r,ss,s) ! 2nd order implicit real*8 (a-h,o-z) ! aux0= (-3*(-S + ss)**2)/(r**2 + (-S + ss)**2)**2.5 + - (r**2 + (-S + ss)**2) bl_z2=aux0 end ! ********************************************************* real*8 function BL_Z3(r,ss,s) ! 3rd order implicit real*8 (a-h,o-z) ! aux0= (15*(-S + ss)**3)/(r**2 + (-S + ss)**2)**3.5 - - (9*(-S + ss))/(r**2 + (-S + ss)**2)**2.5 bl_z3=aux0 end ! ********************************************************* real*8 function BL_Z4(r,ss,s) ! 4th order implicit real*8 (a-h,o-z) ! aux0= (-105*(-S + ss)**4)/(r**2 + (-S + ss)**2)**4.5 + - (90*(-S + ss)**2)/(r**2 + (-S + ss)**2)**3.5 - - 9/(r**2 + (-S + ss)**2)**2.5 bl_z4=aux0 end ! ********************************************************* real*8 function BL_Z5(r,ss,s) ! 5th order implicit real*8 (a-h,o-z) ! aux0= (945*(-S + ss)**5)/(r**2 + (-S + ss)**2)**5.5 - - (1050*(-S + ss)**3)/(r**2 + (-S + ss)**2)**4.5 + - (225*(-S + ss))/(r**2 + (-S + ss)**2)**3.5 bl_z5=aux0 end ! ********************************************************* real*8 function BL_Z6(r,ss,s) ! 6th order implicit real*8 (a-h,o-z) ! aux0= (-10395*(-S + ss)**6)/(r**2 + (-S + ss)**2)**6.5 + - (14175*(-S + ss)**4)/(r**2 + (-S + ss)**2)**5.5 - - (4725*(-S + ss)**2)/(r**2 + (-S + ss)**2)**4.5 + - 225/(r**2 + (-S + ss)**2)**3.5 bl_z6=aux0 end ! ********************************************************* real*8 function BL_Z7(r,ss,s) ! 7th order implicit real*8 (a-h,o-z) ! aux0= (135135*(-S + ss)**7)/(r**2 + (-S + ss)**2)**7.5 - - (218295*(-S + ss)**5)/(r**2 + (-S + ss)**2)**6.5 + - (99225*(-S + ss)**3)/(r**2 + (-S + ss)**2)**5.5 - - (11025*(-S + ss))/(r**2 + (-S + ss)**2)**4.5 bl_z7=aux0 end ! ********************************************************* real*8 function BL_Z8(r,ss,s) ! 8th order implicit real*8 (a-h,o-z) ! aux0= (-2027025*(-S + ss)**8)/(r**2 + (-S + ss)**2)**8.5 + - (3783780*(-S + ss)**6)/(r**2 + (-S + ss)**2)**7.5 - - (2182950*(-S + ss)**4)/(r**2 + (-S + ss)**2)**6.5 + - (396900*(-S + ss)**2)/(r**2 + (-S + ss)**2)**5.5 - - 11025/(r**2 + (-S + ss)**2)**4.5 bl_z8=aux0 end ! ********************************************************* real*8 function BL_Z9(r,ss,s) ! 9th order implicit real*8 (a-h,o-z) ! aux0= (34459425*(-S + ss)**9)/(r**2 + (-S + ss)**2)**9.5 - - (72972900*(-S + ss)**7)/(r**2 + (-S + ss)**2)**8.5 + - (51081030*(-S + ss)**5)/(r**2 + (-S + ss)**2)**7.5 - - (13097700*(-S + ss)**3)/(r**2 + (-S + ss)**2)**6.5 + - (893025*(-S + ss))/(r**2 + (-S + ss)**2)**5.5 bl_z9=aux0 end c ********************************************************* subroutine BROD(fpar,xyz,e,b) ! ! Get field from bent axial current rod ! Use interpolated x-y field from predefined grid ! implicit double precision (a-h,o-z) include 'icommon.inc' c save jxlo,jylo real*8 xyz(3),e(3),b(3),fpar(15) c kk = fpar(2) ! grid number interp = fpar(10) ! interpolation model rring = fpar(12) ! do i=1,3 e(i) = 0d0 b(i) = 0d0 end do htrk = 1d0 / rring ! y = xyz(2) call HUNT(zgr(1,kk),nzgr(kk),y,jylo) ! Watch out for initial point exactly on left edge of grid if( jylo .eq. 0 ) jylo = 1 if( y .ge. zgr(nzgr(kk),kk) ) go to 900 ! B=0 beyond grid x = xyz(1) call HUNT(rgr(1,kk),nrgr(kk),x,jxlo) if( jxlo .eq. 0 ) jxlo = 1 if( x .ge. rgr(nrgr(kk),kk) ) go to 900 ! if( interp .le. 1 ) then ! bi-linear interpolation call BILINEAR(jylo,jxlo,y,x,kk,bz,br) ! else if( interp.eq.2 .or. interp.eq.3 ) then call BFIELDINT(jylo,jxlo,interp,y,x,kk,bz,br) ! end if ! end of test on interpolation ! b(1) = br b(2) = bz c 900 continue end ! ********************************************************* subroutine BSHEET(z,clen,rho,a,bz,br) ! ! Get field from a single solenoidal current sheet (unnormalized) ! ! Ref. M. Garrett, J. Appl. Phys. 34:2567, 1963 ! ! z distance of observation point to center of sheet ! clen length of the sheet ! rho radius of observation point ! a radius of sheet ! bz longitudinal field component ! br radial field component ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 k ! data p2/1.570796327d0/ ! ! Dont allow observation points in the sheet if( ABS(z).le.clen/2. .and. rho.eq.a ) then if(termout) then write(*,*) 'Asking for B on a sheet, terminate' write(*,'(a,3f10.4)')'Zo,Lsh,Rsh: ',z,clen,a write(2,*) 'Asking for B on a sheet, terminate' write(2,'(a,3f10.4)')'Zo,Lsh,Rsh: ',z,clen,a end if stop ! bz = 0.d0 ! br = 0.d0 ! iflag = -59 ! go to 900 end if ! c2 = -4.d0*a*rho / (a+rho)**2 f2 = (a-rho) / (2.d0*a) ! ! Left-end of sheet u = z - clen/2.d0 ! distance of particle to left end of sheet r1 = SQRT( (a + rho)**2 + u*u ) k = 2.d0 / r1 * SQRT( a*rho ) f1 = u*a / ((a+rho)*r1) ele = ELLE(piov2,k) elk = ELLF(piov2,k) if( rho .ne. a ) then elp = ELLPI(piov2,c2,k) bz1 = f1 * ( elk + f2*(elp-elk) ) if( rho .ne. 0. ) then f3 = r1 / (4.d0*rho) br1 = f3 * ( 2.d0*(elk-ele) - k*k*elk ) else br1 = 0.d0 end if else bz1 = 0.d0 if( rho .ne. 0. ) then f3 = r1 / (4.d0*rho) br1 = f3 * ( 2.d0*(elk-ele) - k*k*elk ) else br1 = 0.d0 end if end if ! ! Right-end of sheet u = z + clen/2.d0 ! distance of particle to right end of sheet r1 = SQRT( (a + rho)**2 + u*u ) k = 2.d0 / r1 * SQRT( a*rho ) f1 = u*a / ((a+rho)*r1) ele = ELLE(piov2,k) elk = ELLF(piov2,k) if( rho .ne. a ) then elp = ELLPI(piov2,c2,k) bz2 = f1 * ( elk + f2*(elp-elk) ) if( rho .ne. 0. ) then f3 = r1 / (4.d0*rho) br2 = f3 * ( 2.d0*(elk-ele) - k*k*elk ) else br2 = 0.d0 end if else bz2 = 0.d0 if( rho .ne. 0. ) then f3 = r1 / (4.d0*rho) br2 = f3 * ( 2.d0*(elk-ele) - k*k*elk ) else br2 = 0.d0 end if end if ! bz = -bz1 + bz2 br = br1 - br2 ! Sign convention ! Bz is +ve if field line parallel to dipole moment of solenoid ! Br is +ve if field line is diverging from the solenoid axis ! 900 continue end c ********************************************************* subroutine BS1ORD(x,y,h,g,b1) ! ! Compute 1st order bent solenoid fields ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b1(3) ! t02 = -a1p0 - h*a0p0 - g*b0p0 ! b1(1) = a0p0 + (a1p0-bsp1/2.)*x + b1p0*y ! b1(2) = b0p0 + b1p0*x + (t02-bsp1/2.)*y ! b1(3) = bsp0 + (a0p1-h*bsp0)*x + (b0p1-g*bsp0)*y ! end c ********************************************************* subroutine BS2ORD(x,y,h,hp,g,gp,b2) ! ! Compute 2nd order bent solenoid fields ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b2(3) ! a1 = a1p0 - bsp1/2. a1p = a1p1 - bsp2/2. ! t02 = -a1 - bsp1 - h*a0p0 - g*b0p0 t02p = -a1p - bsp2 - h*a0p1 - hp*a0p0 - g*b0p1 - gp*b0p0 t12 = -h*a1 - 2.*a2p0 + h*h*a0p0 + h*g*b0p0 - g*b1p0 & + 2.*h*bsp1 - a0p2 + hp*bsp0 t03 = g*a1 -2.*b2p0 - h*b1p0 + 2.*g*h*a0p0 +2.*g*g*b0p0 & - b0p2 + 3.*g*bsp1 + gp*bsp0 ! b2(1) = a0p0 + a1*x + b1p0*y + a2p0*x*x + 2.*b2p0*x*y + t12/2.*y*y ! b2(2) = b0p0 + b1p0*x + t02*y + b2p0*x*x + t12*x*y + t03/2.*y*y ! b2(3) = bsp0 + (a0p1-h*bsp0)*x + (b0p1-g*bsp0)*y & + (a1p/2.-h*a0p1+h*h*bsp0)*x*x & + (b1p1-h*b0p1-g*a0p1+2.*h*g*bsp0)*x*y & + (t02p/2.-g*b0p1+g*g*bsp0)*y*y ! end c ********************************************************* subroutine BS2ORDH(x,y,h,b2) ! ! Compute 2nd order bent solenoid fields ! Constant curvature h ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b2(3) ! t02 = -a1p0 - h*a0p0 t02p = -a1p1 - h*a0p1 t12 = -a2p0 - h*a1p0 + h*bsp1 - h*t02/2. - a0p2/2. t03 = -b2p0 - h*b1p0/2. - b0p2/2. ! b2(1) = a0p0 + (a1p0-bsp1/2.)*x + b1p0*y + a2p0*x*x + 2.*b2p0*x*y & + t12*y*y ! b2(2) = b0p0 + b1p0*x + (t02-bsp1/2.)*y + b2p0*x*x + 2.*t12*x*y & + t03*y*y ! b2(3) = bsp0 + (a0p1-h*bsp0)*x + b0p1*y & + (a1p1/4.-h*a0p1+h*h*bsp0)*x*x & + (b1p1/2.-h*b0p1)*x*y & + (t02p/4.-bsp2/4.)*y*y ! end c ********************************************************* subroutine BS3ORD(x,y,h,hp,hpp,g,gp,gpp,b3) ! ! Compute 3rd order bent solenoid fields ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b3(3) ! a1 = a1p0 - bsp1/2. a1p = a1p1 - bsp2/2. a1pp = a1p2 - bsp3/2. a3 = a3p0 + bsp3/16. ! ! Auxiliary quantities h3 = 1.d0 + h*x + g*y t02 = -h*a0p0 - a1 - g*b0p0 - bsp1 t02p = -a1p - bsp2 - h*a0p1 - hp*a0p0 - g*b0p1 - gp*b0p0 t02pp = -a1pp-bsp3-h*a0p2-2*hp*a0p1-hpp*a0p0-g*b0p2-2*gp*b0p1 & -gpp*b0p0 t03 = -g*a1 -2.*b2p0 - h*b1p0 - b0p2 + g*bsp1 + gp*bsp0 & - 2*g*t02 t03p = -g*a1p - gp*a1 - 2*b2p1 - h*b1p1 - hp*b1p0 - b0p3 & + g*bsp2 + gp*bsp1 + gp*bsp1 + gpp*bsp0 - 2*g*t02p - 2*gp*t02 t12 = -2*h*a1 - 2.*a2p0 - g*b1p0 + h*bsp1 - a0p2 + hp*bsp0 & - h*t02 t12p = -2*h*a1p - 2*hp*a1 - 2*a2p1 - g*b1p1 - gp*b1p0 & + h*bsp2 + hp*bsp1 - a0p3 + hp*bsp1 + hpp*bsp0 & - h*t02p - hp*t02 t13 = -2*g*t12 - 6*b3p0 - 2*g*a2p0 - 4*h*b2p0 + g*a0p2 & - b1p2 - h*t03 - 2*g*h*bsp1 + h*b0p2 & + gp*a0p1 + hp*b0p1 - 2*bsp0*(h*gp+g*hp) t22 = 2*(-h*h*bsp1 - 3*h*a2p0 - g*b2p0 - h*t12 + h*a0p2 - 3*a3 & - a1pp/2 - 2*h*hp*bsp0 + hp*a0p1 ) t04 = 2*(-3*g*t03/2 - bsp1*g*g + b0p2*g - t02pp/2 - t22/2 & - 2*g*b2p0 - h*t12/2 + gp*b0p1 - 2*g*gp*bsp0 ) ! ! Field components b3(1) = a0p0 + a1*x + b1p0*y & + a2p0*x*x + 2.*b2p0*x*y + t12/2.*y*y & + a3*x**3 + 3.*b3p0*x*x*y + t22/2.*x*y*y + t13/6.*y**3 ! b3(2) = b0p0 + b1p0*x + t02*y & + b2p0*x*x + t12*x*y + t03/2.*y*y & + b3p0*x**3 + t22/2.*x*x*y + t13/2.*x*y*y + t04/6.*y**3 ! b3(3) = bsp0+a0p1*x+b0p1*y+a1p/2*x*x+b1p1*x*y+t02p/2*y*y & +a2p1/3*x**3+b2p1*x*x*y+t12p/2*x*y*y+t03p/6*y**3 b3(3) = b3(3) / h3 ! end c ********************************************************* subroutine BS3ORDH(x,y,h,b3) ! ! Compute 3rd order bent solenoid fields ! Constant curvature h ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b3(3) ! ! Auxiliary quantities t02 = -a1p0 - h*a0p0 t02p = -a1p1 - h*a0p1 t02pp = -a1p2 - h*a0p2 t12 = -a2p0 - h*a1p0 + h*bsp1 - h*t02/2. - a0p2/2. t12p = -a2p1 - h*a1p1 + h*bsp2 - h*t02p/2 - a0p3/2 t03 = -b2p0 - h*b1p0/2. - b0p2/2. t03p = -b2p1 - h*b1p1/2 - b0p3/2 t22 = -a3p0 - h*a2p0 -2*h*t12/3 - a1p2/12 + h*a0p2/3 - h*h*bsp1/3 t13 = -b3p0 - 2*h*b2p0/3 - b1p2/12 + h*b0p2/6 t04 = -t22 - h*t12/3 - t02pp/12 ! ! Field components b3(1) = a0p0 + (a1p0-bsp1/2.)*x + b1p0*y & + a2p0*x*x + 2.*b2p0*x*y + t12*y*y & + (a3p0+bsp3/16.)*x**3 + 3.*b3p0*x*x*y + (3.*t22+bsp3/16.)*x*y*y & + t13*y**3 ! b3(2) = b0p0 + b1p0*x + (t02-bsp1/2.)*y & + b2p0*x*x + 2.*t12*x*y + t03*y*y & + b3p0*x**3 + (3.*t22+bsp3/16)*x*x*y + 3.*t13*x*y*y & + (t04+bsp3/16.)*y**3 ! b3(3) = bsp0 + (a0p1-h*bsp0)*x + b0p1*y & + (a1p1/4.-h*a0p1-bsp2/4.+h*h*bsp0)*x*x & + (b1p1/2.-h*b0p1)*x*y + (t02p/4.-bsp2/4.)*y*y & + (a2p1/6.-h*a1p1/4.+h*h*a0p1-h*bsp2/4.-h*h*h*bsp0)*x**3 & + (b2p1/2.-h*b1p1/2.+h*h*b0p1)*x*x*y & + (t12p/2.-h*t02p/4.-h*bsp2/4.)*x*y*y + t03p*y**3 ! end c ********************************************************* subroutine BS4ORD(x,y,h,hp,hpp,hppp,g,gp,gpp,gppp,b4) ! ! Compute 4th order bent solenoid fields ! Uses J. Gallardo Maxima results, 23 Dec 02 ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b4(3) ! a1 = a1p0 - bsp1/2. a1p = a1p1 - bsp2/2. a1pp = a1p2 - bsp3/2. a1ppp = a1p3 - bsp4/2. a3 = a3p0 + bsp3/16. a3p = a3p1 + bsp4/16. ! c Independent quantities: bs a0 b0 a1 b1 a2 b2 a3 b3 a4 b4 c auxiliary quantities aux0=(1.+h*x+g*y) aux1=(h*a0p0+g*b0p0) aux1p=(hp*a0p0+h*a0p1+gp*b0p0+g*b0p1) aux1pp=(hpp*a0p0+2.*hp*a0p1+h*a0p2+gpp*b0p0+2.*gp*b0p1 & +g*b0p2) c !jcg aux1ppp=(hppp*a0p0+hpp*a0p1+2.*hpp*a0p1+2.*hp*a0p2+hp*a0p2+ & h*a0p3+gppp*b0p0+gpp*b0p1+2.*gpp*b0p1+2.*gp*b0p2+gp*b0p2+ & g*b0p3) a02=-a1-aux1-bsp1 a02p=-a1p-aux1p-bsp2 a02pp=-a1pp-aux1pp-bsp3 c !jcg a02ppp=-a1ppp-aux1ppp-bsp4 aux2=(h*a1+g*b1p0) aux2p=(hp*a1+h*a1p+gp*b1p0+g*b1p1) aux2pp=(hpp*a1+2.*hp*a1p+h*a1pp+gpp*b1p0+2.*gp*b1p1+g*b1p2) c !jcg aux2ppp=(hpp*a1p+hppp*a1+2.*hp*a1pp+2.*hpp*a1p+h*a1ppp+hp*a1pp+ & gpp*b1p1+gppp*b1p0+2.*gp*b1p2+2.*gpp*b1p1+g*b1p3+gp*b1p2) a12=-2*a2p0-aux2+h*aux1+hp*bsp0-a0p2+2*h*bsp1 a12p=-2*a2p1-aux2p+hp*aux1+h*aux1p+hpp*bsp0+hp*bsp1-a0p3 * +2*hp*bsp1+2*h*bsp2 a12pp=-2*a2p2-aux2pp+hpp*aux1+2*hp*aux1p+h*aux1pp+hppp*bsp0 & +2*hpp*bsp1+hp*bsp2-a0p4+2*hpp*bsp1+4*hp*bsp2+2*h*bsp3 aux3=(h*b1p0+g*a02) aux3p=(hp*b1p0+h*b1p1+gp*a02+g*a02p) aux3pp=(hpp*b1p0+2.*hp*b1p1+h*b1p2+gpp*a02+ * 2.*gp*a02p+g*a02pp) a03=-2*b2p0-aux3+g*aux1+gp*bsp0-b0p2+2*g*bsp1 a03p=-2*b2p1-aux3p+gp*aux1+g*aux1p+gpp*bsp0+gp*bsp1-b0p3 & +2*gp*bsp1+2*g*bsp2 a03pp=-2*b2p2-aux3pp+gpp*aux1+2*gp*aux1p+g*aux1pp & +gppp*bsp0+2.*gpp*bsp1+gp*bsp2-b0p4 & +2*gpp*bsp1+4*gp*bsp2+2*g*bsp3 aux4=(h*2*a2p0+g*2*b2p0) aux4p = 2.*(h*a2p1+hp*a2p0+g*b2p1+gp*b2p0) a22=-6*a3-aux4+2*h*aux2-2*h**2*aux1 * +hp*a0p1-3*h*hp*bsp0-a1pp+4*h*a0p2-6*h**2*bsp1 a22p = -6*a3p-aux4p+2*hp*aux2+2*h*aux2p-4*h*hp*aux1 & -2*h*h*aux1p+hpp*a0p1+hp*a0p2-3*hp*hp*bsp0-3*h*hpp*bsp0 & -3*h*hp*bsp1-a1ppp+4*hp*a0p2+4*h*a0p3-12*h*hp*bsp1-6*h*h*bsp2 aux5=(2*h*b2p0+g*a12) aux5p=(2*hp*b2p0+2*h*b2p1+gp*a12+g*a12p) a13=-6*b3p0-aux5-h*aux3-g*aux2+hp*b0p1+gp*a0p1 * -3*h*gp*bsp0-b1p2+2*h*b0p2+2*g*a0p2 a13p = -6*b3p1-aux5p-hp*aux3-h*aux3p-gp*aux2-g*aux2p & +hpp*b0p1+hp*b0p2+gpp*a0p1+gp*a0p2-3*hp*gp*bsp0-3*h*gpp*bsp0 & -3*h*gp*bsp1-b1p3+2*hp*b0p2+2*h*b0p3+2*gp*a0p2+2*g*a0p3 aux6=(h*a12+g*a03) aux6p = h*a12p+hp*a12+g*a03p+gp*a03 ! a04=-a22-aux6+2*g*aux3+2*g*g*aux1+2*gp*b0p1+a1p+aux1p+bsp3 & +4*g*b0p2-6*g*g*bsp1 a04p = -a22p-aux6p+2*g*aux3p+2*gp*aux3+2*g*g*aux1p+4*g*aux1 & +2*gp*b0p2+2*gpp*b0p1+a1pp+aux1pp+bsp4+4*g*b0p3+4*gp*b0p2 & -6*g*g*bsp2-12*g*bsp1 aux7 = 6*(h*a3+g*b3p0) aux7p = 6*(h*a3p+hp*a3+g*b3p1+gp*b3p0) a32=-24*a4p0-aux7+3*h*aux4-6*h**2*aux3+6*h**3*aux1+ * 3*hp*a1p-18*h*hp*a0p1+36*h**2*hp*bsp0-2*a2p2+6*h*a1pp- * 18*h**2*a0p2+24*h**3*bsp1 aux8 = 6*h*b3p0+g*a22 aux8p = 6*hp*b3p0+6*h*b3p1+gp*a22+g*a22p a23=-24*b4p0-aux8+2*h*aux5 +g*aux4-2*h**2*aux3- * 2*h*g*aux4+2*(hp*b1p1+gp*a1p)-6*h*(hp*b0p1+gp*a0p1) & -6*h*hp*a0p1+12*h*h*gp*bsp0+24*h*h*hp*bsp0-2*b2p2 & +4*h*b1p2+2*g*a1pp-6*h*h*b0p2-12*h*g*a0p2 aux9=(h*a22+g*a13) aux9p = h*a22p+hp*a22+g*a13p+gp*a13 a14=-a32-aux9+h*aux6+2*g*aux5+hp*a02p+2*gp*b0p1- * 6*(h*gp*b0p1+g*gp*a0p1)+24*h*g*gp*bsp0-a12pp & +2*h*a02pp+4*g*b1p2-12*h*g*b0p2-6*g**2*a0p2 aux10=(h*a13+g*a04) a05=-a23-aux10+3*g*aux6-6*g**2*aux3+6*g**3*aux1+3*gp*a02p- * 18*g*gp*b0p1+36*g**2*gp*bsp0-a03pp+6*g*a02pp-18*g**2*b0p2+ * 24*g**3*bsp1 c define the field components b4(1) = a0p0+x*(a1+x*(a2p0+x*(a3+x*(a4p0))))+ * y*(b1p0+y*(a12/2.+y*(a13/6.+y*(a14/24.))))+ * x*y*(2*b1p0+x*(3*b3p0+x*(4*b4p0))+y*(a22/2.+y*(a23/6.))+ * x*y*a32/4.) ! b4(2) = b0p0+x*(b1p0+x*(b2p0+x*(b3p0+x*(b4p0))))+ * y*(a02+y*(a03/2.+y*(a04/6.+y*(a05/24.))))+ * x*y*(a12+x*(a22/2.+x*(a32/6.))+y*(a13/2.+y*(a14/6.))+ * x*y*a23/4.) ! fld3 = bsp0+x*(a0p1+x*(a1p/2.+x*(a2p1/3.+x*a3p/4.)))+ * y*(b0p1+y*(a02p/2.+y*(a03p/6.+y*a04p/24.)))+ * x*y*( b1p1+x*(b2p1+x*b3p1)+y*(a12p/2.+y*a13p/6.)+ * x*y*a22p/4. ) b4(3) = fld3/aux0 ! end c ********************************************************* subroutine BS4ORDH(x,y,h,b4) ! ! Compute 4th order bent solenoid fields ! Constant curvature h ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b4(3) ! a1 = a1p0 - bsp1/2. a1p = a1p1 - bsp2/2. a1pp = a1p2 - bsp3/2. a1ppp = a1p3 - bsp4/2. a3 = a3p0 + bsp3/16. a3p = a3p1 + bsp4/16. ! ! Auxiliary quantities h3 = 1.d0 + h*x t02 = -h*a0p0 - a1 - bsp1 t02p = -a1p - bsp2 - h*a0p1 t02pp = -a1pp - bsp3 - h*a0p2 t02ppp = -a1ppp - bsp4 - h*a0p3 t03 = -2.*b2p0 - h*b1p0 - b0p2 t03p = -2*b2p1 - h*b1p1 - b0p3 t03pp = -2*b2p2 - h*b1p2 - b0p4 t12 = -2*h*a1 - 2.*a2p0 + h*bsp1 - a0p2 - h*t02 t12p = -2*h*a1p - 2*a2p1 + h*bsp2 - a0p3 - h*t02p t12pp = -2*h*a1pp - 2*a2p2 + h*bsp3 - a0p4 - h*t02pp t13 = -6*b3p0 - 4*h*b2p0 - b1p2 - h*t03 + h*b0p2 t13p = -6*b3p1 - 4*h*b2p1 - b1p3 - h*t03p + h*b0p3 t22 = 2*(-h*h*bsp1 - 3*h*a2p0 - h*t12 + h*a0p2 - 3*a3 - a1pp/2) t22p = 2*(-h*h*bsp2 - 3*h*a2p1 - h*t12p + h*a0p3 - 3*a3p & - a1ppp/2 ) t23 = 2*(-h*h*b0p2 + h*b1p2 - 9*h*b3p0 - 12*b4p0 - b2p2 & - h*t13) t32 = 6*(-h*t22/2 + h/2*a1pp - 4*a4p0 - 4*h*a3 - a2p2/3 & + h**3*bsp1 - h*h*a0p2) t04 = -t02pp - t22 - h*t12 t04p = -t02ppp - t22p - h*t12p t14 = 2*(-h*t22 - t32/2 + h/2*t02pp - t12pp/2 - h/2*t04) t05 = -t23 - t03pp - h*t13 ! ! Field components b4(1) = a0p0 + a1*x + b1p0*y & + a2p0*x*x + 2.*b2p0*x*y + t12/2.*y*y & + a3*x**3 + 3.*b3p0*x*x*y + t22/2.*x*y*y + t13/6.*y**3 & + a4p0*x**4 + 4*b4p0*x**3*y + t32/4*x*x*y*y + t23/6*x*y**3 & + t14/24*y**4 ! b4(2) = b0p0 + b1p0*x + t02*y & + b2p0*x*x + t12*x*y + t03/2.*y*y & + b3p0*x**3 + t22/2.*x*x*y + t13/2.*x*y*y + t04/6.*y**3 & + b4p0*x**4 + t32/6*x**3*y + t23/4*x*x*y*y & + t14/6*x*y**3 + t05/24*y**4 ! b4(3) = bsp0 + a0p1*x + b0p1*y & + a1p/2*x*x + b1p1*x*y + t02p/2*y*y & + a2p1/3*x**3 + b2p1*x*x*y + t12p/2*x*y*y + t03p/6*y**3 & + a3p/4*x**4 + b3p1*x**3*y + t22p/4*x*x*y*y & + t13p/6*x*y**3 + t04p/24*y**4 ! b4(3) = b4(3) / h3 ! end c ********************************************************* subroutine BS5ORD(x,y,kp0,kp1,kp2,kp3,kp4,b5) ! ! Compute 5th order bent solenoid fields ! Uses C. Wang Mathmatica results, 14 October 2001 ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 b5(3) real*8 kp0,kp1,kp2,kp3,kp4 ! ! Note that odd skew harmonics (OSH) are explicitly broken into true transverse harmonics ! and solenoid parts in the following expansion. ! a1 = a1p0 a1p = a1p1 a1pp = a1p2 a1ppp = a1p3 a1pppp = a1p4 a3 = a3p0 a3p = a3p1 a3pp = a3p2 a5 = a5p0 ! bx1 = a0p0 + (a1 - bsp1/2)*x + a2p0*x**2 + (a3 + bsp3/16)*x**3 2 + a4p0*x**4 + 3 (a5 - bsp5/384)*x**5 + b1p0*y + 2*b2p0*x*y + 3*b3p0*x**2*y + 4 4*b4p0*x**3*y + 5*b5p0*x**4*y + 5 ((-2*a0p2 - 4*a2p0 - 2*a1*kp0 + 5*bsp1*kp0 + 2*a0p0*kp0**2 + 6 2*bsp0*kp1)*y**2)/4 + (-a1pp/2 - 3*a3 + bsp3/16 + 2*a0p2*kp0 - 7 a2p0*kp0 + a1*kp0**2 - (7*bsp1*kp0**2)/2 - a0p0*kp0**3 8 + a0p1*kp1 - 9 3*bsp0*kp0*kp1)*x*y**2 + ((-16*a2p2 - 192*a4p0 + 48*a1pp*kp0 - A 48*a3*kp0 - 27*bsp3*kp0 - 144*a0p2*kp0**2 + 48*a2p0*kp0**2 - B 48*a1*kp0**3 + 216*bsp1*kp0**3 + 48*a0p0*kp0**4 + 24*a1p*kp1 - C 12*bsp2*kp1 - 144*a0p1*kp0*kp1 D + 288*bsp0*kp0**2*kp1)*x**2*y**2)/32 ! bx2 = (-a3pp/2 - 10*a5 - bsp5/192 + (4*a2p2*kp0)/3 2 - 2*a4p0*kp0 - 3 3*a1pp*kp0**2 + 2*a3*kp0**2 + (13*bsp3*kp0**2)/8 4 + 8*a0p2*kp0**3 - 5 2*a2p0*kp0**3 + 2*a1*kp0**4 - 11*bsp1*kp0**4 - 2*a0p0*kp0**5 + 6 (2*a2p1*kp1)/3 - 3*a1p*kp0*kp1 + (3*bsp2*kp0*kp1)/2 + 7 12*a0p1*kp0**2*kp1 - 20*bsp0*kp0**3*kp1)*x**3*y**2 + 8 ((-b1p2 - 6*b3p0 + 2*b0p2*kp0 - 2*b2p0*kp0 + b1p0*kp0**2 9 + b0p1*kp1)* A y**3)/6 + ((-b2p2 - 12*b4p0 + 2*b1p2*kp0 - 3*b3p0*kp0 B - 3*b0p2*kp0**2 + C 2*b2p0*kp0**2 - b1p0*kp0**3 + b1p1*kp1 D - 3*b0p1*kp0*kp1)*x*y**3)/3 ! bx3 = ((-b3p2 - 20*b5p0 + 2*b2p2*kp0 - 4*b4p0*kp0 - 3*b1p2*kp0**2 2 + 3*b3p0*kp0**2 + 4*b0p2*kp0**3 - 2*b2p0*kp0**3 + b1p0*kp0**4 3 + b2p1*kp1 - 4 3*b1p1*kp0*kp1 + 6*b0p1*kp0**2*kp1)*x**2*y**3)/2 + 5 ((4*a0p4 + 16*a2p2 + 96*a4p0 - 24*a1pp*kp0 + 48*a3*kp0 6 - bsp3*kp0 + 7 40*a0p2*kp0**2 - 24*a2p0*kp0**2 + 12*a1*kp0**3 - 70*bsp1*kp0**3 8 - 12*a0p0*kp0**4 - 8*a1p*kp1 - 20*bsp2*kp1 + 28*a0p1*kp0*kp1 - 9 116*bsp0*kp0**2*kp1 - 12*a0p0*kp1**2 + 4*a1*kp2 - 18*bsp1*kp2 - A 16*a0p0*kp0*kp2 - 4*bsp0*kp3)*y**4)/96 ! bx4 = ((16*a1pppp + 192*a3pp + 1920*a5 - bsp5 - 128*a0p4*kp0 2 - 320*a2p2*kp0 + 3 768*a4p0*kp0 + 448*a1pp*kp0**2 - 480*a3*kp0**2 4 + 66*bsp3*kp0**2 - 5 896*a0p2*kp0**3 + 288*a2p0*kp0**3 - 192*a1*kp0**4 6 + 1376*bsp1*kp0**4 + 7 192*a0p0*kp0**5 - 192*a0p3*kp1 - 128*a2p1*kp1 8 + 336*a1p*kp0*kp1 + 9 792*bsp2*kp0*kp1 - 1184*a0p1*kp0**2*kp1 + 3040*bsp0*kp0**3*kp1 - A 96*a1*kp1**2 + 528*bsp1*kp1**2 + 480*a0p0*kp0*kp1**2 B - 128*a0p2*kp2 + C 32*a2p0*kp2 - 128*a1*kp0*kp2 + 704*bsp1*kp0*kp2 + D 320*a0p0*kp0**2*kp2 + 320*bsp0*kp1*kp2 - 32*a0p1*kp3 + E 160*bsp0*kp0*kp3)*x*y**4)/384 ! bx5 = ((b1p4 + 12*b3p2 + 120*b5p0 - 4*b0p4*kp0 - 12*b2p2*kp0 2 + 48*b4p0*kp0 + 3 10*b1p2*kp0**2 - 18*b3p0*kp0**2 - 16*b0p2*kp0**3 4 + 6*b2p0*kp0**3 - 5 3*b1p0*kp0**4 - 6*b0p3*kp1 - 4*b2p1*kp1 + 7*b1p1*kp0*kp1 - 6 29*b0p1*kp0**2*kp1 - 3*b1p0*kp1**2 - 4*b0p2*kp2 + 2*b2p0*kp2 - 7 4*b1p0*kp0*kp2 - b0p1*kp3)*y**5)/120 by1 = b0p0 + b1p0*x + b2p0*x**2 + 2 b3p0*x**3 + b4p0*x**4 + b5p0*x**5 + ((-2*a1 - bsp1 3 - 2*a0p0*kp0)*y)/2 + 4 (-a0p2 - 2*a2p0 - a1*kp0 + (5*bsp1*kp0)/2 + a0p0*kp0**2 5 + bsp0*kp1)*x* 6 y + (-a1pp/2 - 3*a3 + bsp3/16 + 2*a0p2*kp0 - a2p0*kp0 7 + a1*kp0**2 - 8 (7*bsp1*kp0**2)/2 - a0p0*kp0**3 + a0p1*kp1 9 - 3*bsp0*kp0*kp1)*x**2*y + A (-a2p2/3 - 4*a4p0 + a1pp*kp0 - a3*kp0 - (9*bsp3*kp0)/16 - B 3*a0p2*kp0**2 + a2p0*kp0**2 - a1*kp0**3 + (9*bsp1*kp0**3)/2 + C a0p0*kp0**4 + (a1p*kp1)/2 - (bsp2*kp1)/4 - 3*a0p1*kp0*kp1 + D 6*bsp0*kp0**2*kp1)*x**3*y ! by2 = (-a3pp/4 - 5*a5 - bsp5/384 + 2 (2*a2p2*kp0)/3 - a4p0*kp0 - (3*a1pp*kp0**2)/2 + a3*kp0**2 + 3 (13*bsp3*kp0**2)/16 + 4*a0p2*kp0**3 - a2p0*kp0**3 + a1*kp0**4 - 4 (11*bsp1*kp0**4)/2 - a0p0*kp0**5 + (a2p1*kp1)/3 5 - (3*a1p*kp0*kp1)/2 + 6 (3*bsp2*kp0*kp1)/4 + 6*a0p1*kp0**2*kp1 7 - 10*bsp0*kp0**3*kp1)*x**4*y + 8 ((-b0p2 - 2*b2p0 - b1p0*kp0)*y**2)/2 + 9 ((-b1p2 - 6*b3p0 + 2*b0p2*kp0 - 2*b2p0*kp0 + b1p0*kp0**2 A + b0p1*kp1)*x*y**2)/2 + ((-b2p2 - 12*b4p0 + 2*b1p2*kp0 B - 3*b3p0*kp0 - 3*b0p2*kp0**2 + C 2*b2p0*kp0**2 - b1p0*kp0**3 + b1p1*kp1 D - 3*b0p1*kp0*kp1)*x**2*y**2)/2 ! by3 = ((-b3p2 - 20*b5p0 + 2*b2p2*kp0 - 4*b4p0*kp0 - 3*b1p2*kp0**2 2 + 3*b3p0*kp0**2 + 4*b0p2*kp0**3 - 2*b2p0*kp0**3 + b1p0*kp0**4 3 + b2p1*kp1 - 4 3*b1p1*kp0*kp1 + 6*b0p1*kp0**2*kp1)*x**3*y**2)/2 + 5 ((16*a1pp + 48*a3 + 3*bsp3 - 16*a0p2*kp0 + 32*a2p0*kp0 - 6 8*a1*kp0**2 + 36*bsp1*kp0**2 + 8*a0p0*kp0**3 + 40*bsp0*kp0*kp1 7 + 8*a0p0*kp2)*y**3)/48 + ((4*a0p4 + 16*a2p2 + 96*a4p0 8 - 24*a1pp*kp0 + 9 48*a3*kp0 - bsp3*kp0 + 40*a0p2*kp0**2 - 24*a2p0*kp0**2 + A 12*a1*kp0**3 - 70*bsp1*kp0**3 - 12*a0p0*kp0**4 - 8*a1p*kp1 - B 20*bsp2*kp1 + 28*a0p1*kp0*kp1 - 116*bsp0*kp0**2*kp1 C - 12*a0p0*kp1**2 + D 4*a1*kp2 - 18*bsp1*kp2 - 16*a0p0*kp0*kp2 E - 4*bsp0*kp3)*x*y**3)/24 ! by4 = (a1pppp/12 + a3pp + 10*a5 - bsp5/192 - (2*a0p4*kp0)/3 - 2 (5*a2p2*kp0)/3 + 4*a4p0*kp0 + (7*a1pp*kp0**2)/3 3 - (5*a3*kp0**2)/2 + 4 (11*bsp3*kp0**2)/32 - (14*a0p2*kp0**3)/3 + (3*a2p0*kp0**3)/2 - 5 a1*kp0**4 + (43*bsp1*kp0**4)/6 + a0p0*kp0**5 - a0p3*kp1 - 6 (2*a2p1*kp1)/3 + (7*a1p*kp0*kp1)/4 + (33*bsp2*kp0*kp1)/8 - 7 (37*a0p1*kp0**2*kp1)/6 + (95*bsp0*kp0**3*kp1)/6 - (a1*kp1**2)/2 8 + (11*bsp1*kp1**2)/4 + (5*a0p0*kp0*kp1**2)/2 - (2*a0p2*kp2)/3 + 9 (a2p0*kp2)/6 - (2*a1*kp0*kp2)/3 + (11*bsp1*kp0*kp2)/3 + A (5*a0p0*kp0**2*kp2)/3 + (5*bsp0*kp1*kp2)/3 - (a0p1*kp3)/6 + B (5*bsp0*kp0*kp3)/6)*x**2*y**3 ! by5 = ((b0p4 + 4*b2p2 + 24*b4p0 - 2*b1p2*kp0 + 12*b3p0*kp0 2 + 4*b0p2*kp0**2 - 3 2*b2p0*kp0**2 + b1p0*kp0**3 + 5*b0p1*kp0*kp1 + b1p0*kp2)*y**4)/24 4 + ((b1p4 + 12*b3p2 + 120*b5p0 - 4*b0p4*kp0 - 12*b2p2*kp0 5 + 48*b4p0*kp0 + 6 10*b1p2*kp0**2 - 18*b3p0*kp0**2 - 16*b0p2*kp0**3 + 6*b2p0*kp0**3 7 - 3*b1p0*kp0**4 - 6*b0p3*kp1 - 4*b2p1*kp1 + 7*b1p1*kp0*kp1 - 8 29*b0p1*kp0**2*kp1 - 3*b1p0*kp1**2 - 4*b0p2*kp2 + 2*b2p0*kp2 - 9 4*b1p0*kp0*kp2 - b0p1*kp3)*x*y**4)/24 + A ((-48*a1pppp - 288*a3pp - 1920*a5 - 5*bsp5 + 144*a0p4*kp0 + B 192*a2p2*kp0 - 1152*a4p0*kp0 - 336*a1pp*kp0**2 + 288*a3*kp0**2 C - 134*bsp3*kp0**2 + 720*a0p2*kp0**3 - 192*a2p0*kp0**3 D + 144*a1*kp0**4 - E 1096*bsp1*kp0**4 - 144*a0p0*kp0**5 + 256*a0p3*kp1 F - 240*a1p*kp0*kp1 - G 1080*bsp2*kp0*kp1 + 976*a0p1*kp0**2*kp1 - 2576*bsp0*kp0**3*kp1 + H 128*a1*kp1**2 - 832*bsp1*kp1**2 - 528*a0p0*kp0*kp1**2 + I 144*a0p2*kp2 - 96*a2p0*kp2 + 144*a1*kp0*kp2 - 936*bsp1*kp0*kp2 J - 304*a0p0*kp0**2*kp2 - 560*bsp0*kp1*kp2 - 224*bsp0*kp0*kp3 - K 16*a0p0*kp4)*y**5)/1920 bz1 = bsp0 + (a0p1 - bsp0*kp0)*x + 2 (a1p/2 - bsp2/4 - a0p1*kp0 + bsp0*kp0**2)*x**2 + 3 (a2p1/3 - ((2*a1p - bsp2)*kp0)/4 + a0p1*kp0**2 4 - bsp0*kp0**3)*x**3 + 5 ((16*a3p + bsp4)/64 - (a2p1*kp0)/3 + ((2*a1p - bsp2)*kp0**2)/4 6 - a0p1*kp0**3 + bsp0*kp0**4)*x**4 + (a4p1/5 7 - ((16*a3p + bsp4)*kp0)/64 + 8 (a2p1*kp0**2)/3 - ((2*a1p - bsp2)*kp0**3)/4 + a0p1*kp0**4 9 - bsp0*kp0**5)* A x**5 + b0p1*y + (b1p1 - b0p1*kp0)*x*y + (b2p1 - b1p1*kp0 B + b0p1*kp0**2)* C x**2*y + (b3p1 - b2p1*kp0 + b1p1*kp0**2 - b0p1*kp0**3)*x**3*y + D (b4p1 - b3p1*kp0 + b2p1*kp0**2 - b1p1*kp0**3 E + b0p1*kp0**4)*x**4*y ! bz2 = ((-2*a1p - bsp2 - 2*a0p1*kp0 - 2*a0p0*kp1)*y**2)/4 + 2 ((-2*a0p3 - 4*a2p1 + 4*a0p1*kp0**2 - 2*a1*kp1 + 7*bsp1*kp1 + 3 6*kp0*(bsp2 + a0p0*kp1) + 2*bsp0*kp2)*x*y**2)/4 + 4 ((-8*a1ppp - 48*a3p + bsp4 - 48*a0p1*kp0**3 + 48*a0p2*kp1 - 5 16*a2p0*kp1 - 48*bsp0*kp1**2 + 8*kp0**2*(2*a1p - 13*bsp2 - 6 12*a0p0*kp1) + 16*a0p1*kp2 + 8*kp0*(6*a0p3 + 2*a2p1 + 6*a1*kp1 7 - 27*bsp1*kp1 - 8*bsp0*kp2))*x**2*y**2)/32 + 8 (2*a0p1*kp0**4 + kp0**3*(-a1p + (11*bsp2)/2 + 5*a0p0*kp1) + 9 kp0*((3*a1ppp)/4 + a3p - (5*bsp4)/16 - 6*a0p2*kp1 A + (3*a2p0*kp1)/2 + B (15*bsp0*kp1**2)/2 - 2*a0p1*kp2) + kp0**2*(-3*a0p3 - 3*a1*kp1 + C (33*bsp1*kp1)/2 + 5*bsp0*kp2) + (-16*a2p3 - 192*a4p1 + D 72*a1pp*kp1 - 48*a3*kp1 - 39*bsp3*kp1 - 144*a0p1*kp1**2 + E 24*a1p*kp2 - 12*bsp2*kp2)/96)*x**3*y**2 ! bz3 = ((-b0p3 - 2*b2p1 - b1p1*kp0 - b1p0*kp1)*y**3)/6 + 2 ((-b1p3 - 6*b3p1 + 2*b1p1*kp0**2 + 3*b0p2*kp1 - 2*b2p0*kp1 + 3 3*kp0*(b0p3 + b1p0*kp1) + b0p1*kp2)*x*y**3)/6 + 4 ((-b2p3 - 12*b4p1 - 3*b1p1*kp0**3 + 3*b1p2*kp1 - 3*b3p0*kp1 - 5 3*b0p1*kp1**2 + 2*kp0**2*(b2p1 - 3*(b0p3 + b1p0*kp1)) + b1p1*kp2 6 + kp0*(3*b1p3 + 3*b3p1 - 12*b0p2*kp1 + 6*b2p0*kp1 7 - 4*b0p1*kp2))*x**2* 8 y**3)/6 + ((16*a1ppp + 48*a3p + 3*bsp4 + 8*a0p1*kp0**3 9 - 16*a0p2*kp1 + A 32*a2p0*kp1 + 40*bsp0*kp1**2 + 4*kp0**2*(-2*a1p + 9*bsp2 + B 6*a0p0*kp1) + 8*a0p1*kp2 + 8*kp0*(-2*a0p3 + 4*a2p1 - 2*a1*kp1 + C 14*bsp1*kp1 + 5*bsp0*kp2) + 8*a0p0*kp3)*y**4)/192 ! bz4 = ((8*a0p5 + 32*a2p3 + 192*a4p1 - 32*a0p1*kp0**4 - 64*a1pp*kp1 2 + 96*a3*kp1 - 42*bsp3*kp1 + 32*a0p1*kp1**2 + 3 8*kp0**3*(4*a1p - 22*bsp2 - 15*a0p0*kp1) - 8*a1p*kp2 - 4 76*bsp2*kp2 - 80*a0p0*kp1*kp2 - 4*kp0**2*(-24*a0p3 + 20*a2p1 - 5 22*a1*kp1 + 191*bsp1*kp1 + 68*bsp0*kp2) + 8*a1*kp3 - 6 44*bsp1*kp3 + kp0*(-64*a1ppp + 48*a3p - 5*bsp4 + 232*a0p2*kp1 - 7 128*a2p0*kp1 - 504*bsp0*kp1**2 + 16*a0p1*kp2 - 40*a0p0*kp3) - 8 8*bsp0*kp4)*x*y**4)/192 + ((b0p5 + 4*b2p3 + 24*b4p1 + b1p1*kp0**3 9 - 2*b1p2*kp1 + 12*b3p0*kp1 + 5*b0p1*kp1**2 + A kp0**2*(4*b0p3 - 2*b2p1 + 3*b1p0*kp1) + b1p1*kp2 + B kp0*(-2*b1p3 + 12*b3p1 + 13*b0p2*kp1 - 4*b2p0*kp1 + 5*b0p1*kp2) + C b1p0*kp3)*y**5)/120 ! b5(1) = bx1 + bx2 + bx3 + bx4 + bx5 b5(2) = by1 + by2 + by3 + by4 + by5 b5(3) = bz1 + bz2 + bz3 + bz4 ! end c ********************************************************* subroutine CAV4_GET(xyz,tt,rp,sphi,cphi,e,b) c ! Interpolate model=4 cavity fields ! implicit double precision (a-h,o-z) include 'icommon.inc' c save izlo,jrlo real*8 xyz(3),e(3),b(3) integer*2 sym(2) data epss/1e-6/ ! z = xyz(3) ! relative z into cavity c if( rp .gt. rdm4 ) then ! lost it radially ifail = -37 go to 900 end if c c Leave the RF field zero in the drift tubes if((z.lt.(zcavm4-gm4)/2.).or.(z.gt.(zcavm4+gm4)/2.)) go to 900 c c Find field components in the gap by interpolating zg = ABS(z - zcavm4/2.) ! relative z from center of gap call HUNT(zcrf,nzom4,zg+epss,izlo) ! find appropriate grid element call HUNT(rcrf,nrom4,rp+epss,jrlo) if( (izlo.lt.1) .or. (izlo.gt.nzom4-1) .or. (jrlo.lt.1) .or. + jrlo.gt.nrom4-1) then ifail = -38 go to 900 end if c c Electric fields u = ( zg - zcrf(izlo) ) / delzm4 ! bilinear interpolation v = ( rp - rcrf(jrlo) ) / delrm4 er = (1.-u)*(1.-v)*ercrf(jrlo,izlo) + u*(1.-v)*ercrf(jrlo,izlo+1) + + u*v*ercrf(jrlo+1,izlo+1) + (1.-u)*v*ercrf(jrlo+1,izlo) ez = (1.-u)*(1.-v)*ezcrf(jrlo,izlo) + u*(1.-v)*ezcrf(jrlo,izlo+1) + + u*v*ezcrf(jrlo+1,izlo+1) + (1.-u)*v*ezcrf(jrlo+1,izlo) c c Magnetic fields c B is defined over a grid del/2 smaller on all sides than E c We can interpolate for small rp by assuming B=0 on axis c We have to extrapolate B into the other border regions c dr2 = delrm4 / 2. dz2 = delzm4 / 2. g2 = gm4 / 2. ! if( rp .lt. dr2 ) then bp = bpcrf(1,izlo) / dr2 * rp c else if( rp .gt. rdm4-dr2) then slp = (bpcrf(nrom4-1,izlo) - bpcrf(nrom4-2,izlo)) / delrm4 bp = slp * (rp - (rdm4-dr2)) + bpcrf(nrom4-1,izlo) c else if( zg.lt.dz2 .and. rp.ge.dr2 .and. rp.le.rdm4-dr2 ) then slp = (bpcrf(jrlo,2) - bpcrf(jrlo,1)) / delzm4 bp = -slp * (dz2-zg) + bpcrf(jrlo,1) c else if( zg.gt.g2-dz2 .and. rp.ge.dr2 .and. rp.le.rdm4-dr2) then slp = (bpcrf(jrlo,nzom4-1) - bpcrf(jrlo,nzom4-2)) / delzm4 bp = slp * (zg-(g2-dz2)) + bpcrf(jrlo,nzom4-1) c else zb = zg - dz2 ! use shifted grid for B rb = rp - dr2 u = ( zb - zcrf(izlo) ) / delzm4 v = ( rb - rcrf(jrlo) ) / delrm4 bp = (1.-u)*(1.-v)*bpcrf(jrlo,izlo) & + u*(1.-v)*bpcrf(jrlo,izlo+1) & + u*v*bpcrf(jrlo+1,izlo+1) + (1.-u)*v*bpcrf(jrlo+1,izlo) end if c c Determine field symmetry, depending on which quadrant particle is in do i=1,2 sym(i) = 1 end do c if( z .lt. zcavm4/2. )then ! 2nd and 3rd quad sym(1) = -1 sym(2) = -1 end if c arg = omrf * tt + phase carg = COS(arg) sarg = SIN(arg) c e(1) = er * cphi * carg * sym(1) e(2) = er * sphi * carg * sym(2) e(3) = ez * carg b(1) = -bp * sphi * sarg b(2) = bp * cphi * sarg ! 900 continue return end c ********************************************************* subroutine CAV4_SET(fpar) c ! Set up model=4 cavity fields ! implicit double precision (a-h,o-z) include 'icommon.inc' ! parameter ( nsrc = 15 ) real*8 fpar(15),rsrc(nsrc),zsrc(nsrc) c zcavm4 = fpar(8) ! total length of cavity gm4 = fpar(9) ! gap rdm4 = fpar(10) ! radius of drift tube rn = fpar(11) ! radius of nose piece ! c Generate field grid ----------------------------------- g2 = gm4 / 2. rcn = rdm4 + rn delzm4 = g2 / nzom4 delrm4 = rdm4 / nrom4 dr2 = delrm4 / 2. dz2 = delzm4 / 2. ! c Generate relative source point locations on nose th = 0. dth = pi / 2. / nsrc c do i=1,nsrc rsrc(i) = rn * cos(th) zsrc(i) = rn * sin(th) th = th + dth end do c c Calculate electric field components over grid of observation c points in the gap c do 200 iz=1,nzom4 zcrf(iz) = (iz-1) * delzm4 c do 150 ir=1,nrom4 rcrf(ir) = (ir-1) * delrm4 ez = 0. er = 0. c do 100 is=1,nsrc dr = rcrf(ir) - ( rcn - rsrc(is) ) ! upper right nose dz = zcrf(iz) - ( g2 - zsrc(is) ) r3 = ( dr*dr + dz*dz )**(1.5) ez = ez - dz / r3 er = er - dr / r3 c dr = rcrf(ir) - ( rcn - rsrc(is) ) ! upper left nose dz = zcrf(iz) - ( -g2 + zsrc(is) ) r3 = ( dr*dr + dz*dz )**(1.5) ez = ez + dz / r3 er = er + dr / r3 c dr = rcrf(ir) - ( -rcn + rsrc(is) ) ! lower left nose dz = zcrf(iz) - ( -g2 + zsrc(is) ) r3 = ( dr*dr + dz*dz )**(1.5) ez = ez + dz / r3 er = er + dr / r3 c dr = rcrf(ir) - ( -rcn + rsrc(is) ) ! lower right nose dz = zcrf(iz) - ( g2 - zsrc(is) ) r3 = ( dr*dr + dz*dz )**(1.5) ez = ez - dz / r3 er = er - dr / r3 100 continue c if( ir .eq. 1 .and. iz .eq. 1 ) escal = emagrf / ez ezcrf(ir,iz) = ez * escal ercrf(ir,iz) = er * escal 150 continue c 200 continue c c Calculate magnetic field from Maxwell equations c do 300 iz=1,nzom4-1 c do 250 ir=1,nrom4-1 derdz = ercrf(ir,iz+1)-ercrf(ir,iz)+ercrf(ir+1,iz+1) & -ercrf(ir+1,iz) derdz = derdz / 2. / delzm4 dezdr = ezcrf(ir+1,iz)-ezcrf(ir,iz)+ezcrf(ir+1,iz+1) & -ezcrf(ir,iz+1) dezdr = dezdr / 2. / delrm4 bpcrf(ir,iz) = (dezdr - derdz) / omrf 250 continue c 300 continue c 900 continue return end c ********************************************************* subroutine COIL(fpar,xyz,e,b) c c Compute field of circular coil segment c implicit double precision (a-h,o-z) include 'icommon.inc' c save jsold,jcold,jzlo,jrlo real*8 xyz(3),e(3),b(3),fpar(15) character*80 title character*10 fname data bnorm/2.e-7/ data jcold/-1/,jsold/-1/ c model = fpar(1) c c Compute field from sum of coil segments ! do i=1,3 e(i) = 0. b(i) = 0. end do ! rp = sqrt( xyz(1)**2 + xyz(2)**2 ) ! if( xyz(1)*xyz(2) .ne. 0. ) then if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) else phi = 0. end if sphi = sin(phi) cphi = cos(phi) ! --------------------------------------------------------- if( model .eq. 1 ) then ! exact field from single loop a = fpar(3) ! coil radius curr = fpar(4) ! current z = xyz(3) - fpar(2) ! call BLOOP(a,z,rp,bz,br) ! b(1) = br * cphi* curr b(2) = br * sphi * curr b(3) = bz * curr ! do i=1,3 b(i) = b(i) * bnorm end do go to 900 end if ! -------------------------------------------------------- ! Initialize for model 2 ! if((kbcell.eq.0 .and. jsrg.ne.jsold ) .or. & (kbcell.eq.1 .and. jcel.ne.jcold ) ) then ! if( model.eq.2 ) then ! need data file mfile = fpar(2) cscale = fpar(11) if( mfile .lt. 0 ) then ! allow short cut for polarity flip mfile = -mfile fcsign = -1. else fcsign = 1. end if write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open COIL data file' iflag = -1 go to 900 end if write(2,*)' COIL SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr if( ioc .ne. 0 ) then write(2,*)'Error reading COIL data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NCOILS = ',nelgr ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),relgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading COIL data file' iflag = -1 go to 900 end if curelgr(i) = curelgr(i) * fcsign * cscale write(2,50) idm,dzelgr(i),relgr(i),curelgr(i) 50 format(' #,DZ,R,I: ',i5,3e12.4) end do ! close(mfile) ! end if ! end of test on model = 2 ! jcold = jcel jsold = jsrg end if ! end of test on kbcell ... ! End of initialization ! -------------------------------------------------------- if( model .eq. 2 ) then ! Exact field from sum of coils in data file ! do i=1,nelgr a = relgr(i) z = xyz(3) - dzelgr(i) ! call BLOOP(a,z,rp,bz,br) ! b(1) = b(1) + br * cphi * curelgr(i) b(2) = b(2) + br * sphi * curelgr(i) b(3) = b(3) + bz * curelgr(i) end do c do i=1,3 b(i) = b(i) * bnorm end do ! go to 900 end if ! -------------------------------------------------------- if( model .eq. 3 ) then ! Interpolated field from predefined grid kk = fpar(2) ! grid number interp = fpar(3) ! interpolation model z = xyz(3) call HUNT(zgr(1,kk),nzgr(kk),z,jzlo) ! Watch out for initial point exactly on left edge of grid if( jzlo .eq. 0 ) jzlo = 1 call HUNT(rgr(1,kk),nrgr(kk),rp,jrlo) if( jrlo .eq. 0 ) jrlo = 1 if( z .ge. zgr(nzgr(kk),kk) ) go to 900 ! B=0 beyond grid if( rp .ge. rgr(nrgr(kk),kk) ) go to 900 ! if( interp .le. 1 ) then ! bi-linear interpolation call BILINEAR(jzlo,jrlo,z,rp,kk,bz,br) ! else if( interp.eq.2 .or. interp.eq.3 ) then ! Polynomial interpolation ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,bzgr,mxzgr,mxrgr,interp, ! & z,rp,bz,df) ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,brgr,mxzgr,mxrgr,interp, ! & z,rp,br,df) ! call BFIELDINT(jzlo,jrlo,interp,z,rp,kk,bz,br) ! end if ! end of test on interpolation ! b(1) = br * cphi b(2) = br * sphi b(3) = bz end if ! end of test on model c 900 continue return end c ********************************************************* subroutine DIPOLE(fpar,xx,e,b) c c Compute field of vertical sector dipole c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 xx(3),e(3),b(3),fpar(15) real*8 lamd,lamq,lams,k ! do i=1,3 e(i) = 0. b(i) = 0.d0 end do ! model = fpar(1) dip = fpar(2) ! dipole strength pref = fpar(4) ! reference momentum ! ! -------------------------------------------------------- if( model .eq. 1 ) then ! constant field (hard edge), 2nd order cf = fpar(3) ! curvature factor if( cf .eq. 0.d0 ) then htrk = pcharge * dip / pref else htrk = cf end if hptrk = 0.d0 b1 = fpar(5) b2 = fpar(6) x = xx(1) y = xx(2) ! b(1) = b1*y + 2.*b2*x*y b(2) = dip + b1*x + b2*x*x - (b2+htrk*b1/2.)*y*y ! ! -------------------------------------------------------- else if( model .eq. 2 ) then ! 3rd order dTANH(s) dipole plus multipoles ! Use multipole model of C. Wang MC#207 cf = fpar(3) ! curvature factor clend = fpar(5) ! length elend = fpar(6) ! dipole end length lamd = fpar(7) ! dipole attenuation length quad = fpar(8) ! quadrupole component strength sex = fpar(9) ! sextupole component strength clenq = fpar(10) elenq = fpar(11) lamq = fpar(12) clens = fpar(13) elens = fpar(14) lams = fpar(15) ! x = xx(1) y = xx(2) z = xx(3) - elend b0 = dip * FTANH(z,clend,lamd,0) b0p = dip * FTANH(z,clend,lamd,1) b0pp = dip * FTANH(z,clend,lamd,2) b0ppp = dip * FTANH(z,clend,lamd,3) ! if( cf .eq. 0.d0 ) then htrk = pcharge * b0 / pref hptrk = pcharge * b0p / pref else htrk = cf hptrk = 0. end if ! h = htrk hp = hptrk ! if( quad .ne. 0. ) then z = xx(3) - elenq b1 = quad * FTANH(z,clenq,lamq,0) b1p = quad * FTANH(z,clenq,lamq,1) b1pp = quad * FTANH(z,clenq,lamq,2) else b1 = 0. b1p = 0. b1pp = 0. end if ! if( sex .ne. 0. ) then z = xx(3) - elens b2 = sex * FTANH(z,clens,lams,0) b2p = sex * FTANH(z,clens,lams,1) else b2 = 0. b2p = 0. end if ! t1 = b1pp+2.*h*(b2-b0pp)-h*h*b1-hp*b0p ! b(1) = b1*y + 2.*b2*x*y - t1/6.*y**3 b(2) = b0 + b1*x + b2*x*x & - (2.*b2+b0pp+h*b1)/2.*y*y - t1/2.*x*y*y b(3) = b0p*y + (b1p-h*b0p)*x*y & + (b2p-h*b1p+h*h*b0p)*x*x*y & - (2.*b2p+b0ppp+h*b1p+hp*b1)/6.*y**3 ! ! -------------------------------------------------------- else if( model .eq. 3 ) then ! combined function dipole k = fpar(3) ! radial exponent r0 = fpar(5) ! distance from machine center to ICOOL ref circle cf = fpar(6) ! curvature factor nbffag = fpar(7) ! flag for radial-sector negative-bend FFAG x = xx(1) r = r0 + x ! radius of particle in LAB if( x .le. -r0) then iflag = -82 go to 900 end if by0 = dip * (r0/r)**k ! By on the x-z midplane a2 = -k*k/2. a4 = -a2/12.*(4. + 4.*k + k*k) a6 = -a4/30.*(16. + 8.*k + k*k) a8 = -a6/56.*(36. + 12.*k + k*k) w = xx(2) / r b(2) = by0*( 1. + a2*w**2 + a4*w**4 + a6*w**6 + a8*w**8 ) b(1) = -by0*( k*w + a2*(k+2.)/3.*w**3 + a4*(k+4.)/5.*w**5 & + a6*(k+6.)/7.*w**7 + a8*(k+8.)/9.*w**9 ) ! if( cf .eq. 0.d0 ) then htrk = pcharge * dip / pref else htrk = cf end if ! ! -------------------------------------------------------- else if( model .eq. 4 ) then ! hard edge with rotated pole faces zc = fpar(7) ! axial center d = fpar(8) ! half length bt1 = fpar(9)*pi/180. ! entrance face rotation bt2 = fpar(10)*pi/180. ! exit face rotation w = fpar(11) ! half horizontal aperture ht = fpar(12) ! half vertical aperture ! htrk = 0. ! initialize locally x = xx(1) y = xx(2) z = xx(3) ! if( ABS(x).gt.w .or. ABS(y).gt.ht ) then iflag = -10 go to 900 end if ! ! n.b. x < 0 on side closest the center of curvature zl = zc - d - x*TAN(bt1) zr = zc + d + x*TAN(bt2) ! if( z.ge.zl .and. z.le.zr ) then ! inside dipole htrk = pcharge / pref * dip b(2) = dip end if ! ! Enable vertical bending when crossing pole faces if( .not.edge1 .and. z.ge.zl ) then b(2) = dip dpyedge = pcharge*b(2)*y*TAN(bt1) edge1 = .true. edgekick = .true. end if ! if( .not.edge2 .and. z.ge.zr ) then b(2) = dip dpyedge = pcharge*b(2)*y*TAN(bt2) edge2 = .true. edgekick = .true. end if ! ELSE IF (MODEL.EQ.5) THEN C C PARAMETERS: C 2 : DIPOLE FIELD (T) C 3 : (SIGNED) CURVATURE OF REFERENCE ORBIT (M^{-1}) C 4 : ENTRANCE POLE FACE ANGLE (DEG) C 5 : EXIT POLE FACE ANGLE (DEG) C 6 : FLAG FOR END FIELD EDGE KICKS C 0 : KICKS AT BOTH ENTRANCE AND EXIT C 1 : KICK AT EXIT ONLY C 2 : KICK AT ENTRANCE ONLY C 3 : NO KICKS C C NOTE THAT THE HEAVY LIFTING HERE REALLY HAPPENS IN BEG_OF_REGION C AND END_OF_REGION C HTRK = FPAR(3) HPTRK = 0D0 B(1) = 0D0 B(2) = DIP end if ! end test on model ! -------------------------------------------------------- ! 900 continue return end c ********************************************************* subroutine DIPOLE_H(fpar,x,e,b) c c Compute field of horizontal sector dipole ! cf. field expansion in routine BS3ORD c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 x(3),e(3),b(3),fpar(15) real*8 lamd,lamq,lams ! do i=1,3 e(i) = 0. b(i) = 0. end do htrk = 0. hptrk = 0. ! model = fpar(1) dip = fpar(2) ! dipole strength pref = fpar(4) ! reference momentum ! ! -------------------------------------------------------- if( model .eq. 1 ) then ! hard edge, 2nd order a1 = fpar(5) a2 = fpar(6) gtrk = pcharge / pref * dip g = gtrk t02 = -a1 t03 = -2*g*t02 - g*a1 t12 = -a2 ! b(1) = dip + a1*xx + a2*xx*xx + t12*y*y b(2) = t02*y + 2.*t12*xx*y + t03*y*y ! -------------------------------------------------------- else if( model .eq. 2 ) then ! dTANH(s) dipole plus multipoles ! clend = fpar(5) ! length elend = fpar(6) ! dipole end length lamd = fpar(7) ! dipole attenuation length g0 = fpar(8) ! quadrupole component strength s0 = fpar(9) ! setupole component strength clenq = fpar(10) elenq = fpar(11) lamq = fpar(12) clens = fpar(13) elens = fpar(14) lams = fpar(15) ! xx = x(1) y = x(2) z = x(3) - elend ! a0 = dip * FTANH(z,clend,lamd,0) a0p = dip * FTANH(z,clend,lamd,1) a0pp = dip * FTANH(z,clend,lamd,2) a0ppp = dip * FTANH(z,clend,lamd,3) gtrk = pcharge / pref * a0 gptrk = pcharge / pref * a0p g = gtrk gp = gptrk ! if( g0 .ne. 0. ) then z = x(3) - elenq a1 = g0 * FTANH(z,clenq,lamq,0) a1p = g0 * FTANH(z,clenq,lamq,1) a1pp = g0 * FTANH(z,clenq,lamq,2) else a1 = 0. a1p = 0. a1pp = 0. end if ! if( s0 .ne. 0. ) then z = x(3) - elens a2 = s0 * FTANH(z,clens,lams,0) a2p = s0 * FTANH(z,clens,lams,1) else a2 = 0. a2p = 0. end if ! t02 = -a1 t02p = -a1p t02pp = -a1pp t03 = -2*g*t02 - g*a1 t03p = -2*g*t02p - 2*gp*t02 - g*a1p - gp*a1 t12 = -a0pp - 2*a2 t12p = -a0ppp - 2*a2p t13 = -2.*g*t12 - 2*g*a2 + g*a0pp t22 = - a1pp t04 = -3*g*t03 - t02pp - t22 ! b(1) = a0 + a1*xx + a2*xx*xx + t12/2.*y*y & + t22/2.*xx*y*y + t13/6.*y**3 b(2) = t02*y + t12*xx*y + t03/2*y*y & + t22/2*xx*xx*y + t13/2*xx*y*y + t04/6*y**3 b(3) = a0p*xx + a1p/2*xx*xx + t02p/2*y*y & + a2p/3*xx**3 + t12p/2*xx*y*y + t03p/6*y**3 ! end if ! end test on model ! -------------------------------------------------------- ! 900 continue return end c ********************************************************* integer function EDGE_ACCEL_13(x,p,t,q,m,om,ph,g) ! ! Hard edge focus for open cell accelerator model 13 ! Ref. J.S. Berg 22 Nov 2004 double precision g,m,om,p(3),ph,q,t,x(3) double precision xps(6) integer rc,impmid,rfendi,rfendv external rfendv rc = rfendi(g,om,ph,q,m) if (rc.ne.0) then EDGE_ACCEL_13 = rc return end if xps(1) = x(1) xps(2) = p(1) xps(3) = x(2) xps(4) = p(2) xps(5) = t xps(6) = sqrt(p(1)*p(1)+p(2)*p(2)+p(3)*p(3)+m*m) rc = impmid(xps,rfendv) if (rc.ne.0) then EDGE_ACCEL_13 = rc return end if if (rc.eq.0) then x(1) = xps(1) x(2) = xps(3) p(1) = xps(2) p(2) = xps(4) t = xps(5) p(3) = sqrt(xps(6)*xps(6)-p(1)*p(1)-p(2)*p(2)-m*m) end if EDGE_ACCEL_13 = 0 return end c ********************************************************* subroutine EFIELD(fpar,xyz,e,b) ! ! Make a static electric field ! Ref. J.S. Berg 7 Aug 2006 ! ! Model 1: Given field in the box (xmin,ymin,zmin),(xmax,ymax,zmax) ! Parameters: ! 2 : Ex ! 3 : Ey ! 4 : Ez ! 5 : xmin ! 6 : xmax ! 7 : ymin ! 8 : ymax ! 9 : zmin ! 10 : zmax ! 11 : Respect curvature=0, force rectangular=1 ! Specify curvature=2, spec curv+force rect=3 ! 12 : Curvature if 11 is 2 or 3 ! implicit double precision (a-h,o-z) include 'icommon.inc' ! double precision fpar(15),xyz(3),e(3),b(3) double precision eps,thcyl,rcyl,xrect,zrect,zmid integer curvflag parameter (eps=2d0**(-27)) ! b(1) = 0d0 b(2) = 0d0 b(3) = 0d0 e(1) = 0d0 e(2) = 0d0 e(3) = 0d0 curvflag = nint(fpar(11)) ! if (curvflag/2.eq.1) htrk = fpar(12) if (mod(curvflag,2).eq.1.and.abs(htrk*slen).lt.eps) $ curvflag = curvflag-1 ! if (fpar(7).le.xyz(2).and.xyz(2).le.fpar(8)) then if (mod(curvflag,2).eq.0.and. $ fpar(5).le.xyz(1).and.xyz(1).le.fpar(6).and. $ fpar(9).le.xyz(3).and.xyz(3).le.fpar(10)) then e(1) = fpar(2)/(1.0+htrk*xyz(1)) e(2) = fpar(3) e(3) = fpar(4) else if (mod(curvflag,2).eq.1) then zmid = 0.5d0*(fpar(9)+fpar(10)) thcyl = htrk*(xyz(3)-zmid) rcyl = xyz(1)+1d0/htrk xrect = xyz(1) - 2d0*rcyl*sin(0.5d0*thcyl)**2 zrect = rcyl*sin(thcyl) if (fpar(5).le.xrect.and.xrect.le.fpar(6).and. $ fpar(9)-zmid.le.zrect.and.zrect.le.fpar(10)-zmid) then e(1) = fpar(2) e(2) = fpar(3) e(3) = fpar(4) end if end if end if ! end c ********************************************************* subroutine ELLTHREE(phi,en,ak,ele,elk,elp) ! ! combines elle,ellk,and ellpi calculations, which call ! rd,rf, and rj, into one loop; ! assumes positive arguments cc,q,1+enss to rd,rf,rj ! implicit double precision (a-h,o-z) PARAMETER (ERRTOLrj=.05,ERRTOLrf=0.08,ERRTOLrd=0.05, *TINY=1.d-20,BIG=1.d20,THIRD=1./3.,C1=3./14.,C2=1./3.,C3=3./22., *C4=3./26.,C5=.75*C3,C6=1.5*C4,C7=.5*C2,C8=C3+C3, *A1=1./24.,A2=.1,A3=3./44.,A4=1./14., *B1=3./14.,B2=1./6.,B3=9./22.,B4=3./26.,B5=.25*C3,B6=1.5*C4) CU USES rc c s=sin(phi) enss=en*s*s cc=cos(phi)**2 q=(1.-s*ak)*(1.+s*ak) x=cc y=q z=1. p=1.+enss ! Require that x,y,z,p are positive if(min(x,y,z,p).lt.0..or.min(x+y,x+z,y+z,p,z).lt.TINY.or. *max(x,y,z,p).gt.BIG) then stop 'invalid arguments in rj' end if ! sumrd=0. sumrj=0. fac=1. 1 continue sqrtx=sqrt(x) sqrty=sqrt(y) sqrtz=sqrt(z) alamb=sqrtx*(sqrty+sqrtz)+sqrty*sqrtz alpha=(p*(sqrtx+sqrty+sqrtz)+sqrtx*sqrty*sqrtz)**2 beta=p*(p+alamb)**2 sumrj=sumrj+fac*rc(alpha,beta) sumrd=sumrd+fac/(sqrtz*(z+alamb)) fac=.25*fac x=.25*(x+alamb) y=.25*(y+alamb) z=.25*(z+alamb) p=.25*(p+alamb) averj=.2*(x+y+z+p+p) averd=.2*(x+y+3.*z) averf=THIRD*(x+y+z) dxrj=(averj-x)/averj dyrj=(averj-y)/averj dzrj=(averj-z)/averj dprj=(averj-p)/averj dxrd=(averd-x)/averd dyrd=(averd-y)/averd dzrd=(averd-z)/averd dxrf=(averf-x)/averf dyrf=(averf-y)/averf dzrf=(averf-z)/averf if(max(abs(dxrj),abs(dyrj),abs(dzrj),abs(dprj)).gt.ERRTOLrj)goto 1 if(max(abs(dxrf),abs(dyrf),abs(dzrf)).gt.ERRTOLrf)goto 1 if(max(abs(dxrd),abs(dyrd),abs(dzrd)).gt.ERRTOLrd)goto 1 c ea=dxrj*(dyrj+dzrj)+dyrj*dzrj eb=dxrj*dyrj*dzrj ec=dprj**2 ed=ea-3.*ec ee=eb+2.*dprj*(ea-ec) rj=3.*sumrj+fac*(1.+ed*(-C1+C5*ed-C6*ee)+ *eb*(C7+dprj*(-C8+dprj*C4))+ *dprj*ea*(C2-dprj*C3)-C2*dprj*ec)/(averj*sqrt(averj)) c eard=dxrd*dyrd ebrd=dzrd*dzrd ecrd=eard-ebrd edrd=eard-6.*ebrd eerd=edrd+ecrd+ecrd rd=3.*sumrd+fac*(1.+edrd*(-B1+B5*edrd-B6*dzrd*eerd)+ *dzrd*(B2*eerd+dzrd*(-B3*ecrd+dzrd*B4*eard)))/(averd*sqrt(averd)) c e2=dxrf*dyrf-dzrf**2 e3=dxrf*dyrf*dzrf rf=(1.+(A1*e2-A2-A3*e3)*e2+A4*e3)/sqrt(averf) c ele=s*(rf-((s*ak)**2)*rd/3.) elk=s*rf elp=s*(rf-enss*rj/3.) c return end c ********************************************************* double precision function FINDZC(omega,grad,slen,zstep,pz,pm) ! ! Used in phasemodel 2 ! Ref. J.S. Berg 22 Nov 2004 ! implicit none double precision grad,omega,pm,pz,slen,zstep double precision de0,de1,dde,e,phi0,phi1 double precision c data c/2.99792458d8/ ! e = sqrt(pz*pz+pm*pm) phi0 = -0.5d0*omega*slen*e/(c*pz) de0 = grad*slen 100 de1 = de0 phi1 = phi0 call zcfcn(phi0,omega,grad,slen,zstep,e,pm,de0,dde) if (de0.eq.0d0 .or. abs(de0).ge.abs(de1)) goto 200 phi0 = phi0 - de0/dde goto 100 ! 200 if (de0.eq.0d0) then findzc = phi0/omega else findzc = phi1/omega end if ! end c ********************************************************* subroutine FOFO(fpar,xyz,e,b) c c Get field from solenoidal FOFO lattice c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 xyz(3),e(3),b(3) real*8 fpar(15) c do i=1,3 b(i) = 0. e(i) = 0. end do ! model = fpar(1) ! --------------------------------------------------------- if( model .eq. 1 ) then ! linear ramp bmag = fpar(2) bcen = fpar(3) period = fpar(4) rp = SQRT( xyz(1)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) else phi = 0. end if ! d = xyz(3) - fpar(5) d = MOD( d+period, period ) p4 = period / 4. g = bmag / p4 ! if( d .le. p4 ) then b(3) = bcen + g * d br = -0.5 * g * rp else if( (p4.lt.d) .and. (d.le.2*p4) )then b(3) = bcen + g * ( 2.*p4 - d ) br = 0.5 * g * rp else if( (2*p4.lt.d) .and. (d.le.3*p4) )then b(3) = bcen - g * ( d - 2.*p4 ) br = 0.5 * g * rp else b(3) = bcen - g * ( period - d ) br = -0.5 * g * rp end if ! b(1) = br * cos(phi) b(2) = br * sin(phi) ! --------------------------------------------------------- else if( model .eq. 2 ) then ! sinusoidal ramp with Bessel r bmag = fpar(2) bcen = fpar(3) period = fpar(4) rp = sqrt( xyz(1)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) else phi = 0. end if ! d = xyz(3) - fpar(5) xk = 2. * pi / period b(3) = bcen + bmag * SIN( xk*d ) * BESSI0( xk*rp ) br= -bmag * COS( xk*d ) * BESSI1( xk*rp ) b(1) = br * cos(phi) b(2) = br * sin(phi) end if ! end of test on model c 900 continue return end c ********************************************************* function FOURIER(cs,sn,am,bm,mxord,fn,per,idrv) ! ! Calculate derivatives of function specified by its Fourier coefficients ! Function defined up to 7th derivative ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 am(0:mxord),bm(0:mxord),cs(0:mxord),sn(0:mxord) ! fourier = 0.d0 ! if( idrv .eq. 0 ) then do n=0,mxord fourier = fourier + am(n)*cs(n) + bm(n)*sn(n) end do fourier = fn * fourier ! else if( idrv .eq. 1 ) then do n=0,mxord fourier = fourier + n*( -am(n)*sn(n) + bm(n)*cs(n) ) end do fourier = fn * twopi/per * fourier ! else if( idrv .eq. 2 ) then do n=0,mxord fourier = fourier + n**2 * (-am(n)*cs(n) - bm(n)*sn(n)) end do fourier = fn * (twopi/per)**2 * fourier ! else if( idrv .eq. 3 ) then do n=0,mxord fourier = fourier + n**3 * (am(n)*sn(n) - bm(n)*cs(n)) end do fourier = fn * (twopi/per)**3 * fourier ! else if( idrv .eq. 4 ) then do n=0,mxord fourier = fourier + n**4 * (am(n)*cs(n) + bm(n)*sn(n)) end do fourier = fn * (twopi/per)**4 * fourier ! else if( idrv .eq. 5 ) then do n=0,mxord fourier = fourier + n**5 * (-am(n)*sn(n) + bm(n)*cs(n)) end do fourier = fn * (twopi/per)**5 * fourier ! else if( idrv .eq. 6 ) then do n=0,mxord fourier = fourier + n**6 * (-am(n)*cs(n) - bm(n)*sn(n)) end do fourier = fn * (twopi/per)**6 * fourier ! else if( idrv .eq. 7 ) then do n=0,mxord fourier = fourier + n**7 * (am(n)*sn(n) - bm(n)*cs(n)) end do fourier = fn * (twopi/per)**7 * fourier end if ! end c ********************************************************* function FTANH(z,clen,lamb,order) ! ! Calculates difference of TANH functions and derivatives ! for z dependence of field strengths ! implicit double precision (a-h,o-z) real*8 lamb integer order ! d = z / lamb g = (z - clen) / lamb scd = 1. / COSH(d) scg = 1. / COSH(g) thd = TANH(d) thg = TANH(g) ! if( order .eq. 0 ) then ftanh = 0.5 * ( thd - thg ) ! else if( order .eq. 1 ) then ftanh = 0.5/lamb * ( scd**2 - scg**2 ) ! else if( order .eq. 2 ) then ftanh = 1./lamb**2 * ( -thd * scd**2 + thg * scg**2 ) ! else if( order .eq. 3 ) then ftanh = 1./lamb**3 * ( -scd**4 + 2.*thd**2 * scd**2 & + scg**4 - 2.*thg**2 * scg**2 ) ! else if( order .eq. 4 ) then ftanh = 4./lamb**4 * ( 2.*scd**4 * thd & - thd**3 * scd**2 - 2.*scg**4 * thg + thg**3 * scg**2 ) ! else if( order .eq. 5 ) then ftanh = 4./lamb**5 * ( -11.* scd**4 * thd**2 & + 2.*scd**6 + 2.*thd**4 * scd**2 & + 11.*scg**4 * thg**2 - 2.*scg**6 - 2.*thg**4 * scg**2 ) ! else if( order .eq. 6 ) then ftanh = 8./lamb**6 * ( 26.*scd**4 * thd**3 & - 17.*thd * scd**6 - 2.*scd**2 * thd**5 & - 26.*thg**3 * scg**4 + 17.*scg**6 * thg & + 2.*thg**5 * scg**2 ) ! else if( order .eq. 7 ) then ftanh = 8./lamb**7 * ( -114.* scd**4 * thd**4 & + 180.*scd**6 * thd**2 - 17.* scd**8 & + 4.*scd**2 * thd**6 + 114.*scg**4 * thg**4 & - 180.*thg**2 * scg**6 + 17.*scg**8 - 4.*thg**6 * scg**2 ) ! end if ! return end c ********************************************************* subroutine GRID(kk,grtype,fpar) ! ! Set up a field grid ! implicit double precision (a-h,o-z) dimension fpar(15) character*4 grtype ! if( grtype.eq.'STUS' .or. grtype.eq.'G43D' ) then ! read in 3D field map call GRID3D(grtype,fpar) ! else if( grtype .eq. 'BROD' ) then ! make new x-y field map on grid call GRIDXY(kk,grtype,fpar) ! else ! make new r-z field map on a grid call GRIDRZ(kk,grtype,fpar) end if ! end c ********************************************************* subroutine GRIDRZ(kk,grtype,fpar) ! ! Set up r-z field grid ! implicit double precision (a-h,o-z) include 'icommon.inc' dimension fpar(15) character*80 title,line character*10 fname character*4 grtype data bnorm/4.e-7/ ! mfile = fpar(2) delzgr(kk) = fpar(3) delrgr(kk) = fpar(4) zglen = fpar(5) rglen = fpar(6) dzcut = fpar(7) kfile = fpar(9) cscale = fpar(11) ! if( mfile .lt. 0 ) then ! allow short cut for polarity flip mfile = -mfile fcsign = -1. else fcsign = 1. end if ! write(fname,10) mfile 10 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open GRID data file' iflag = -1 go to 900 end if ! ! --------------------------------------------------------- ! Set up grid parameters if( grtype .eq. 'BLOC' ) then write(2,*)' BLOCK SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NBLOCKS = ',nelgr ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),zlelgr(i), & relgr(i),r2elgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading BLOCK data file' iflag = -1 go to 900 end if curelgr(i) = curelgr(i) * fcsign * cscale write(2,20) idm,dzelgr(i),zlelgr(i), & relgr(i),r2elgr(i),curelgr(i) 20 format(' #,DZ,L,R1,R2,I:',i4,5e12.4) end do ! ! ......................................................... else if( grtype .eq. 'COIL' ) then write(2,*)' COIL SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr if( ioc .ne. 0 ) then write(2,*)'Error reading COIL data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NCOILS = ',nelgr ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),relgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading COIL data file' iflag = -1 go to 900 end if curelgr(i) = curelgr(i) * fcsign * cscale write(2,30) idm,dzelgr(i),relgr(i),curelgr(i) 30 format(' #,DZ,R,I: ',i5,3e12.4) end do ! ! ......................................................... else if( grtype .eq. 'SHEE' ) then if( cscale .eq. 0. ) cscale = 1. ! write(2,*)' SHEET SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr,fcelgr if( ioc .ne. 0 ) then write(2,*)'Error reading SHEET data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NSHEETS = ',nelgr,' fcur = ', & fcelgr ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),zlelgr(i), & relgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading SHEET data file' iflag = -1 go to 900 end if curelgr(i) = curelgr(i) * fcelgr * fcsign * cscale write(2,40) idm,dzelgr(i),zlelgr(i), & relgr(i),curelgr(i) 40 format(' #,DZ,L,R,I: ',i5,4e12.4) end do ! ! ......................................................... else if( grtype .eq. 'SOL' ) then ! user-supplied field grid if( fpar(10) .eq. 0.d0 ) then ! longitudinal shift ksh = 0 else ksh = fpar(10) end if ! write(2,*)' R-Z MAP SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nzgr(kk),nrgr(kk) if( ioc .ne. 0 ) then write(2,*)'Error reading SOL data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,i5)')' R-Z MAP: NZGR,NRGR',nzgr(kk),nrgr(kk) if( nzgr(kk).gt.mxzgr .or. nrgr(kk).gt.mxrgr ) then write(2,*)'Number of SOL z or r grid points exceeds max' iflag = -1 go to 900 end if ! do i=1,nzgr(kk) do j=1,nrgr(kk) read(mfile,*,iostat=ioc)idm,jdm,zgr(i,kk),rgr(j,kk), & bz,br if( ioc .ne. 0 ) then write(2,*)'Error reading SOL field data file' iflag = -1 go to 900 end if kz = KSHIF(i,ksh,nzgr(kk)) bzgr(kz,j,kk) = bz * fcsign * cscale brgr(kz,j,kk) = br * fcsign * cscale end do end do ! delzgr(kk) = zgr(2,kk) - zgr(1,kk) delrgr(kk) = rgr(2,kk) - rgr(1,kk) ! go to 800 ! ! ......................................................... else if( grtype .eq. 'G4RZ' ) then ! G4BL cylinder grid if( fpar(10) .eq. 0.d0 ) then ! longitudinal shift ksh = 0 else ksh = fpar(10) end if ! write(2,*)' R-Z MAP SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) line ! title write(2,'(a80)') line read(mfile,'(a80)',iostat=ioc) line ! param read(mfile,'(a80)',iostat=ioc) line ! cylinder call PARSE_G4RZ(line,nr,nz) if( nr*nz .le. 0 ) then write(2,*)'Error parsing G4RZ data file' iflag = -1 go to 900 else nzgr(kk) = nz nrgr(kk) = nr end if read(mfile,'(a80)',iostat=ioc) line ! data write(2,'(a,i5,5x,i5)')' R-Z MAP: NZGR,NRGR',nz,nr if( nz.gt.mxzgr .or. nr.gt.mxrgr ) then write(2,*)'Number of G4RZ z or r grid points exceeds max' iflag = -1 go to 900 end if ! do i=1,nz do j=1,nr read(mfile,*,iostat=ioc) r,z,br,bz if( ioc .ne. 0 ) then write(2,*)'Error reading G4RZ field data file' iflag = -1 go to 900 end if rgr(j,kk) = r / 1d3 ! [m] zgr(j,kk) = z / 1d3 kz = KSHIF(i,ksh,nzgr(kk)) bzgr(kz,j,kk) = bz * fcsign * cscale brgr(kz,j,kk) = br * fcsign * cscale end do end do ! delzgr(kk) = zgr(2,kk) - zgr(1,kk) delrgr(kk) = rgr(2,kk) - rgr(1,kk) ! go to 800 ! end if ! end of test on grtype ! close(mfile) ! ! --------------------------------------------------------- ! Compute B values on grid for COIL, SHEET, and BLOCK ! nzgr(kk) = ANINT( zglen / delzgr(kk) + 1. ) nrgr(kk) = ANINT( rglen / delrgr(kk) + 1. ) if( fpar(10) .eq. 0.d0 ) then ! longitudinal shift ksh = 0 else ksh = fpar(10) end if ! do iz=1,nzgr(kk) if( termout .and. MOD(iz,20).eq.1 ) then write(*,80)' Computing grid:',iz,' /',nzgr(kk),' done' 80 format(a,i5,a,i5,a) end if kz = KSHIF(iz,ksh,nzgr(kk)) zgr(kz,kk) = (kz-1) * delzgr(kk) ! do ir=1,nrgr(kk) rgr(ir,kk) = (ir-1) * delrgr(kk) ! bzgr(kz,ir,kk) = 0. brgr(kz,ir,kk) = 0. ! --------------------------------------------------------- if( grtype .eq. 'BLOC' ) then do 100 ic=1,nelgr clen = zlelgr(ic) a1 = relgr(ic) a2 = r2elgr(ic) cur = curelgr(ic) z = zgr(kz,kk) - dzelgr(ic) - clen/2. ! if( fpar(12) .eq. 1.d0 ) then call BBLOCS(z,rgr(ir,kk),clen,a1,a2,cur,bz,br) else call BBLOCI(z,rgr(ir,kk),clen,a1,a2,bu,bv) br = bv * bnorm * cur * 1.d6 bz = bu * bnorm * cur * 1.d6 end if if( iflag .ne. 0 ) go to 900 ! bzgr(kz,ir,kk) = bzgr(kz,ir,kk) + bz brgr(kz,ir,kk) = brgr(kz,ir,kk) + br 100 continue ! ! ......................................................... else if(grtype .eq. 'COIL' ) then ! Check that the bozo isn't asking for a grid point on a coil do i=1,nelgr if( zgr(kz,kk).eq.dzelgr(i) .and. rgr(ir,kk).eq.relgr(i)) then write(2,*) ' Grid point requested on a coil' iflag = -1 go to 900 end if end do ! do 120 ic=1,nelgr a = relgr(ic) z = zgr(kz,kk) - dzelgr(ic) if( ABS(z) .gt. dzcut ) go to 120 ! call BLOOP( a,z,rgr(ir,kk),bz,br ) ! bzgr(kz,ir,kk) = bzgr(kz,ir,kk) + bz * curelgr(ic) brgr(kz,ir,kk) = brgr(kz,ir,kk) + br * curelgr(ic) 120 continue ! sum on coils ! bzgr(kz,ir,kk) = bzgr(kz,ir,kk) * bnorm * 0.5 brgr(kz,ir,kk) = brgr(kz,ir,kk) * bnorm * 0.5 ! ! ......................................................... else if(grtype .eq. 'SHEE' ) then ! Check that the bozo isn't asking for a grid point on a sheet do i=1,nelgr zend = dzelgr(i) + zlelgr(i) if( zgr(kz,kk).ge.dzelgr(i) .and. zgr(kz,kk).le.zend & .and. rgr(ir,kk).eq.relgr(i) ) then write(2,*) ' Grid point requested on a current sheet' iflag = -1 go to 900 end if end do ! do 140 ic=1,nelgr clen = zlelgr(ic) a = relgr(ic) z = zgr(kz,kk) - dzelgr(ic) - clen/2. if( zgr(kz,kk) .lt. dzelgr(ic)-dzcut ) go to 140 if( zgr(kz,kk) .gt. dzelgr(ic)+clen+dzcut ) go to 140 ! call BSHEET( z,clen,rgr(ir,kk),a,bz,br ) if( iflag .ne. 0 ) go to 900 ! bzgr(kz,ir,kk) = bzgr(kz,ir,kk) + bz * curelgr(ic) brgr(kz,ir,kk) = brgr(kz,ir,kk) + br * curelgr(ic) 140 continue ! sum on sheets ! bzgr(kz,ir,kk) = bzgr(kz,ir,kk) * bnorm brgr(kz,ir,kk) = brgr(kz,ir,kk) * bnorm end if ! end of test on grtype ! --------------------------------------------------------- ! end do ! r grid ! end do ! z grid ! ! --------------------------------------------------------- 800 continue if( kfile.gt.19 .and. kfile.lt.100 ) then ! Save field grid to file write(fname,10) kfile open(unit=kfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open field grid output file' iflag = -1 go to 900 end if write(kfile,'(a80)') title write(kfile,'(2i8)') nzgr(kk),nrgr(kk) do i=1,nzgr(kk) do j=1,nrgr(kk) write(kfile,'(2i6,4e13.5)') i,j,zgr(i,kk),rgr(j,kk), & bzgr(i,j,kk),brgr(i,j,kk) end do end do close(kfile) end if ! end of test on kfile ! 900 continue end c ********************************************************* subroutine GRIDXY(kk,grtype,fpar) ! ! Set up x-y field grid ! Save data in the r-z arrays: x <--> r and y <--> z ! implicit double precision (a-h,o-z) include 'icommon.inc' dimension fpar(15) character*80 title character*10 fname character*4 grtype,dum data bnorm/4d-7/ ! delrgr(kk) = fpar(3) delzgr(kk) = fpar(4) xlo = fpar(5) ylo = fpar(6) nx = fpar(7) ny = fpar(8) kfile = fpar(9) rring = fpar(12) rrod = fpar(13) cur = fpar(14) bdip = fpar(15) ! nzgr(kk) = ny nrgr(kk) = nx dum = grtype ! ! Compute B values on grid do iz=1,ny if( termout .and. MOD(iz,20).eq.1 ) then write(*,80)' Computing grid:',iz,' /',nzgr(kk),' done' 80 format(a,i5,a,i5,a) end if zgr(iz,kk) = ylo + (iz-1) * delzgr(kk) u = zgr(iz,kk) ! do ir=1,nx rgr(ir,kk) = xlo + (ir-1) * delrgr(kk) nsh = nrgr(kk) bz = 0d0 br = 0d0 ! do 150 is=1,nsh ! loop over current sheets xsheet = rgr(1,kk) + 0.5d0*delrgr(kk) + (is-1)*delrgr(kk) if( ABS(xsheet) .gt. rrod ) go to 150 arg = rrod**2-xsheet**2 ys = SQRT(arg) ! height of rod surface at this x clen = 2d0 * ys robs = rring + rgr(ir,kk) rsheet = rring + xsheet ! call BSHEET(u,clen,robs,rsheet,bu,bv) ! bz = bz + bu * cur / rsheet ! J ~ 1/r br = br + bv * cur / rsheet 150 continue ! bzgr(iz,ir,kk) = bz * bnorm / nsh + bdip brgr(iz,ir,kk) = br * bnorm / nsh ! ! --------------------------------------------------------- ! end do ! r grid ! end do ! z grid ! ! --------------------------------------------------------- 800 continue if( kfile.gt.19 .and. kfile.lt.100 ) then ! Save field grid to file write(fname,840) kfile 840 format('for0',i2,'.dat') open(unit=kfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open field grid output file' iflag = -1 go to 900 end if write(kfile,'(a80)') title write(kfile,'(2i8)') nzgr(kk),nrgr(kk) do i=1,nzgr(kk) do j=1,nrgr(kk) write(kfile,'(2i6,4e13.5)') i,j,zgr(i,kk),rgr(j,kk), & bzgr(i,j,kk),brgr(i,j,kk) end do end do close(kfile) end if ! end of test on kfile ! 900 continue end ! ********************************************************* subroutine GRID3D(grtype,fpar) ! ! Set up 3D field grid ! implicit double precision (a-h,o-z) include 'icommon.inc' dimension fpar(15) character*80 title,line character*30 utitle character*10 fname character*4 grtype ! mfile = fpar(2) icf = fpar(3) ! curvature flag pref = fpar(4) fac = fpar(5) if( fpar(10) .eq. 0.d0 ) then ! longitudinal shift ksh = 1 else ksh = fpar(10) end if ! write(fname,40) mfile 40 format('for0',i2,'.dat') ! --------------------------------------------------------- ! Set up grid parameters if( grtype .eq. 'STUS' ) then ! if( fpar(9) .eq. 0. ) then ! formatted B grid file open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open user static field' iflag = -1 go to 900 end if ! read(mfile,'(a80)') title read(mfile,*) nxgr,nygr,nsgr,hmap if( nxgr.gt.mxxg .or. nygr.gt.mxyg .or. nsgr.gt.mxsg ) then write(2,*)'Dimensions exceed max limit: user static field' iflag = -1 go to 900 end if ! read(mfile,*) (xgr(j),j=1,nxgr) read(mfile,*) (ygr(j),j=1,nygr) read(mfile,*) (sgr(j),j=1,nsgr) ! nrecords = nxgr*nygr*nsgr ! do i=1,nrecords read(mfile,*)ix,iy,iz,bx,by,bz kz = KSHIF(iz,ksh,nsgr) bxgr(ix,iy,kz) = bx * fac bygr(ix,iy,kz) = by * fac bsgr(ix,iy,kz) = bz * fac end do ! ! ......................................................... else if( fpar(9) .gt. 0. ) then ! unformatted grid or splines open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & form='unformatted',status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open user static field' iflag = -1 go to 900 end if ! read(mfile,iostat=ioc) utitle read(mfile) nxgr,nygr,nsgr,hmap if( nxgr.gt.mxxg .or. nygr.gt.mxyg .or. nsgr.gt.mxsg ) then write(2,*)'Dimensions exceed max limit: user static field' iflag = -1 go to 900 end if ! read(mfile) (xgr(j),j=1,nxgr) read(mfile) (ygr(j),j=1,nygr) read(mfile) (sgr(j),j=1,nsgr) ! if( fpar(9) .eq. 1 ) then ! unformatted B grid nrecords = nxgr*nygr*nsgr do i=1,nrecords read(mfile)ix,iy,iz,bx,by,bz kz = KSHIF(iz,ksh,nsgr) bxgr(ix,iy,kz) = bx * fac bygr(ix,iy,kz) = by * fac bsgr(ix,iy,kz) = bz * fac end do ! else ! unformatted spline coefficients read(mfile) nbfsp,nordsp read(mfile) msksp,mcbsp read(mfile) (sknotsp(j),j=1,msksp) read(mfile) (cbxsp(j),j=1,mcbsp) read(mfile) (cbysp(j),j=1,mcbsp) read(mfile) (cbssp(j),j=1,mcbsp) end if ! end of test on fpar(9)=1,2 ! end if ! end of test on fpar(9) close(mfile) ! if( fpar(8) .eq. 1. ) hmap = -hmap ! if( termout ) then write(2,'(a,i2)')'Reading STUS file ',mfile if( fpar(9) .eq. 0 ) then write(2,'(a80)') title else write(2,*) utitle end if write(2,'(a,3i5,f13.5)') 'NXGR,NYGR,NSGR,HMAP: ', & nxgr,nygr,nsgr,hmap if( fpar(9) .eq. 2 ) then write(2,'(a,3i5,3x,3i5)')' NBF,ORD: ',nbfsp,nordsp end if end if ! if( icf.eq.1 .and. hmap.eq.0.d0 ) then call HUNT(xgr,nxgr,0.d0,jx) if( jx .eq. 0 ) jx = 1 call HUNT(ygr,nygr,0.d0,jy) if( jy .eq. 0 ) jy = 1 do iz=1,nsgr htrkgr(iz) = pcharge / pref * bygr(jx,jy,iz) gtrkgr(iz) = pcharge / pref * bxgr(jx,jy,iz) end do end if ! ! --------------------------------------------------------- else if( grtype .eq. 'G43D' ) then open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open G4BL 3D field' iflag = -1 go to 900 end if ! read(mfile,'(a80)') title read(mfile,'(a80)') line ! param read(mfile,'(a80)') line ! grid call PARSE_G43D(line,nxgr,nygr,nsgr) if( nxgr*nygr*nsgr .le. 0 ) then write(2,*) 'Error parsing G4BL data file' iflag = -1 go to 900 end if read(mfile,'(a80)') line ! data write(2,'(a,3i5)') 'NXGR,NYGR,NSGR: ',nxgr,nygr,nsgr if( nxgr.gt.mxxg .or. nygr.gt.mxyg .or. nsgr.gt.mxsg ) then write(2,*)'Dimensions exceed max limit: G4BL 3D field' iflag = -1 go to 900 end if ! do ix=1,nxgr do iy=1,nygr do iz=1,nsgr ! read(mfile,*) x,y,z,bx,by,bz xgr(ix) = x / 1d3 ! [m] ygr(iy) = y / 1d3 ! [m] sgr(iz) = z / 1d3 ! [m] kz = KSHIF(iz,ksh,nsgr) bxgr(ix,iy,kz) = bx * fac bygr(ix,iy,kz) = by * fac bsgr(ix,iy,kz) = bz * fac end do end do end do ! if( icf .eq. 1 ) then call HUNT(xgr,nxgr,0.d0,jx) if( jx .eq. 0 ) jx = 1 call HUNT(ygr,nygr,0.d0,jy) if( jy .eq. 0 ) jy = 1 do iz=1,nsgr htrkgr(iz) = pcharge / pref * bygr(jx,jy,iz) gtrkgr(iz) = pcharge / pref * bxgr(jx,jy,iz) end do end if ! end if ! end of test on grtype ! dxgr = xgr(2) - xgr(1) dygr = ygr(2) - ygr(1) dsgr = sgr(2) - sgr(1) ! 900 continue end c ********************************************************* subroutine HELFLD(rp,phi,phi0,z,pern,bref,r0,an,bn,n,bx,by,bz) ! ! Compute helical field for multipole order n ! ! Ref: Tominaka, NIM A 459:398, 2001, eq. 40. ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 k,kr,kz,an(20),bn(20),pern(20) ! k = ABS( twopi / pern(n) ) kr = k * rp ! Allow sign of kz => twist orientation if( pern(n) .gt. 0. ) then kz = k * z else kz = -k * z end if cphi = COS(phi) sphi = SIN(phi) ! bnorm = bref * r0 br = 0. bphi = 0. bz = 0. ! ar = n * kr az = n * (phi - phi0 - kz) caz = COS(az) saz = SIN(az) fac = FACTRL(n) * (2./(n*k*r0))**n ! if( n .eq. 1 ) then birp = 0.5 * ( BESSI0(ar) + BESSI(2,ar) ) bz = -fac*k*BESSI1(ar) * (an(n)*saz + bn(n)*caz) ! else if( n .eq. 2 ) then birp = 0.5 * ( BESSI1(ar) + BESSI(3,ar) ) bz = -fac*k*BESSI(2,ar) * (an(n)*saz + bn(n)*caz) ! else birp = 0.5 * ( BESSI(n-1,ar) + BESSI(n+1,ar) ) bz = -fac*k*BESSI(n,ar) * (an(n)*saz + bn(n)*caz) end if ! br = fac*k*birp * (-an(n)*caz + bn(n)*saz) if( rp .gt. 0. ) then bphi = -bz / kr else bphi = 0. end if ! br = bnorm * br bphi = bnorm * bphi bz = bnorm * bz ! bx = br*cphi - bphi*sphi by = bphi*cphi + br*sphi ! end c ********************************************************* subroutine HELIX(fpar,xyz,e,b) c c Get field from helical magnet ! implicit double precision (a-h,o-z) include 'icommon.inc' c external HX_BR,HX_BPHI,HX_BS,HX_BX0,HX_BY0 real*8 xyz(3),e(3),b(3) real*8 fpar(15) real*8 an(20),bn(20),pern(20) dimension ds0(5),c(5) real*8 kr,k,ka,kz character*80 title character*10 fname save jsold,jcold,a,per save bs1,bs2,bs3,elens1,clens,elens2,zhel0,hellen save nmult,an,bn,pern,npts ! data jcold/-1/,jsold/-1/,zhel0/0./ ! model = fpar(1) ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! Initialization ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ if( (kbcell.eq.0 .and. jsrg.ne.jsold) .or. & (kbcell.eq.1 .and. jcel.ne.jcold) ) then bs1 = fpar(8) ! solenoid field at start of region bs2 = fpar(9) ! central solenoid field bs3 = fpar(10) ! solenoid field at end of region elens1 = fpar(11) ! solenoid entrance end length clens = fpar(12) ! solenoidal central field length elens2 = fpar(13) ! solenoid exit end length if( fpar(14) .eq. 1. ) zhel0 = sr0 ! reset parameter hellen = elens1 + clens + elens2 ! total helix length ! ! --------------------------------------------------------- if( model .eq. 4 ) then mfile = fpar(2) write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open multipole file in HELIX' iflag = -1 go to 900 end if write(2,*)' Opened file',mfile, ' for helical multipoles' read(mfile,'(a80)') title read(mfile,*) nmult do i=1,nmult read(mfile,*) bn(i),an(i) end do close(mfile) ! ! --------------------------------------------------------- else if( model .eq. 5 ) then mfile = fpar(2) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open helix multipole input file' iflag = -1 go to 900 end if write(2,*) ' Opened file ',mfile,' for helix multipoles' ! read(mfile,'(a80)') title read(mfile,*) npts ! do i=1,npts read(mfile,*) soa(i),bsoa(i),a0oa(i),b0oa(i), & a1oa(i),b1oa(i),a2oa(i),b2oa(i),a3oa(i),b3oa(i) end do ! close( mfile ) ! end if ! end of test on model ! jsold = jsrg jcold = jcel end if ! end of test on initialization ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! Enter here for tracking ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ do i=1,3 e(i) = 0. b(i) = 0. end do ! z = xyz(3) dz = z - zhel0 rp = SQRT(xyz(1)**2 + xyz(2)**2) if( rp .gt. 10. ) then ! protect against I0 blow up iflag = -72 go to 900 end if c if (rp .eq. 0.) then phi = 0. sphi = 0. cphi = 1. else phi = ATAN2(xyz(2),xyz(1)) sphi = SIN(phi) cphi = COS(phi) end if ! ! Put external solenoid field in {b(1),b(2),b(3)} ! Put helix field in {bx,by,bz} ! ! --------------------------------------------------------- if (model .eq. 1) then ! simple rotating dipole btm = fpar(2) ! maximum dipole field per = fpar(3) ! helix period if( per .eq. 0. ) then write(2,*)'Period equal 0 in HELIX' iflag = -1 go to 900 end if elent = fpar(4) ! taper length for dipole field ! ! Solenoidal field contribution if( dz .lt. elens1 ) then ! entrance region gs = (bs2 - bs1) / elens1 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs1 + gs * dz else if( dz .gt. elens1+clens ) then ! exit region gs = (bs3 - bs2) / elens2 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs2 + gs * (dz-elens1-clens) else ! central region b(1) = 0. b(2) = 0. b(3) = bs2 end if ! ! Dipole field taper clent = hellen - 2.*elent if( dz .lt. elent ) then br = btm * dz / elent else if( dz .gt. elent+clent ) then br = btm * (hellen-dz) / elent else br = btm end if ! k = ABS( twopi / per ) ! Allow sign of kz => twist orientation if( per .gt. 0. ) then kz = k * dz else kz = -k * dz end if ! bx = br * SIN(kz) by = br * COS(kz) bz = 0. ! ! --------------------------------------------------------- else if (model .eq. 2) then ! helical current sheet ! Ref. T. Tominaka et al, NIM A459:398, 2001, eq. 39 cur0 = fpar(2) ! current a = fpar(3) ! helical sheet radius per = fpar(4) ! helix period if( per .eq. 0. ) then write(2,*)'Period equal 0 in HELIX' iflag = -1 go to 900 end if phi0 = fpar(5) * pi/180. ! starting azimuth elent = fpar(7) ! taper length for transverse field nterms = fpar(15) ! # of terms in series ! ! Solenoidal field contribution if( dz .lt. elens1 ) then ! entrance region gs = (bs2 - bs1) / elens1 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs1 + gs * dz else if( dz .gt. elens1+clens ) then ! exit region gs = (bs3 - bs2) / elens2 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs2 + gs * (dz-elens1-clens) else ! central region b(1) = 0. b(2) = 0. b(3) = bs2 end if ! ! Transverse field taper clent = hellen - 2.*elent if( dz .lt. elent ) then cur = cur0 * dz / elent else if( dz .gt. elent+clent ) then cur = cur0 * (hellen-dz) / elent else cur = cur0 end if ! k = ABS( twopi / per ) ka = k * a kr = k * rp ! Allow sign of kz => twist orientation if( per .gt. 0. ) then kz = k * dz else kz = -k * dz end if ! bnorm = permfp/twopi*cur*k br = 0. bphi = 0. bz = 0. ! do m=1,nterms ! do n=1,m aa = n * ka ar = n * kr az = n * (phi - phi0 - kz) ! if( n .eq. 1 ) then bknp = -(BESSK0(aa) + BESSK1(aa)/ka ) birp = 0.5*( BESSI0(ar) + BESSI(2,ar) ) termz = n * bknp * BESSI1(ar) * COS(az) else if( n .eq. 2 ) then bknp = -(BESSK1(aa) + BESSK(2,aa)/ka ) birp = 0.5*( BESSI1(ar) + BESSI(3,ar) ) termz = n * bknp * BESSI(2,ar) * COS(az) else ! n.b. BESSK and BESSI require n > 1 bknp = -(BESSK(n-1,aa) + BESSK(n,aa)/ka ) birp = 0.5*( BESSI(n-1,ar) + BESSI(n+1,ar) ) termz = n * bknp * BESSI(n,ar) * COS(az) end if ! termr = n * bknp * birp * SIN(az) if( rp .gt. 0. ) then termp = termz / rp else termp = 0. end if ! br = br + termr bphi = bphi + termp bz = bz + termz end do ! end do ! br = 2.*bnorm * ka/nterms * br bphi = 2.*bnorm * a/nterms * bphi bz = 2.*ka/nterms * bz bz = bnorm*(1. - bz) ! bx = br*cphi - bphi*sphi by = bphi*cphi + br*sphi ! ! --------------------------------------------------------- else if( model .eq. 3 ) then ! multifilar helical multipoles ! Ref. J. Gallardo, note 21 January 2005 ! For this model the solenoid end field distances must be symmetric a = fpar(2) ! helical radius per = fpar(3) ! helix period if( per .eq. 0. ) then write(2,*)'Period equal 0 in HELIX' iflag = -1 go to 900 end if g1 = fpar(4) ! dipole strength g2 = fpar(5) ! quadrupole strength g3 = fpar(6) ! sextupole strength elent = fpar(7) ! taper length for transverse field g4 = fpar(13) ! octupole strength nterms = fpar(15) ! number of terms in series ! ! Solenoidal field contribution elens2 = elens1 if( dz .lt. elens1 ) then ! entrance region gs = (bs2 - bs1) / elens1 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs1 + gs * dz else if( dz .gt. elens1+clens ) then ! exit region gs = (bs3 - bs2) / elens2 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs2 + gs * (dz-elens1-clens) else ! central region b(1) = 0. b(2) = 0. b(3) = bs2 end if ! ! Transverse field taper clent = hellen - 2.*elent if( dz .lt. elent ) then ftft = dz / elent else if( dz .gt. elent+clent ) then ftft = (hellen-dz) / elent else ftft = 1.d0 end if ! k = ABS( twopi / per ) ka = k * a kr = k * rp ! Allow sign of kz => twist orientation if( per .gt. 0. ) then kz = k * dz else kz = -k * dz end if br = 0.d0 bphi = 0.d0 bz = 0.d0 ! ........................................................ if( g1 .ne. 0.d0 ) then ! bifilar helical dipole bnorm = -g1 / (ka * BESSK0(ka) + BESSK1(ka)) ! do j=0,nterms n = 2*j+1 ! n = 1,3,... aa = n * ka ar = n * kr az = n * (phi - kz) ! if( n .eq. 1 ) then bknp = -(BESSK0(aa) + BESSK1(aa)/ka ) birp = 0.5*( BESSI0(ar) + BESSI(2,ar) ) termz = n * bknp * BESSI1(ar) * COS(az) else ! n.b. BESSK and BESSI require n > 1 bknp = -(BESSK(n-1,aa) + BESSK(n,aa)/ka ) birp = 0.5*( BESSI(n-1,ar) + BESSI(n+1,ar) ) termz = n * bknp * BESSI(n,ar) * COS(az) end if ! termr = n * bknp * birp * SIN(az) if( rp .gt. 0. ) then termp = termz / rp else termp = 0. end if ! br = br + termr*bnorm bphi = bphi + termp*bnorm bz = bz + termz*bnorm end do ! ! ........................................................ else if( g2 .ne. 0.d0 ) then ! 4-filar helical quadrupole arg = 2.d0 * ka denom = ka*BESSK1(arg) + BESSK(2,arg) bnorm = -g2 / (2d0*k*denom) ! do j=0,nterms n = 2*(2*j+1) ! n = 2,6,... aa = n * ka ar = n * kr az = n * (phi - kz) ! if( n .eq. 2 ) then bknp = -(BESSK1(aa) + BESSK(2,aa)/ka ) birp = 0.5*( BESSI1(ar) + BESSI(3,ar) ) termz = n * bknp * BESSI(2,ar) * COS(az) else ! n.b. BESSK and BESSI require n > 1 bknp = -(BESSK(n-1,aa) + BESSK(n,aa)/ka ) birp = 0.5*( BESSI(n-1,ar) + BESSI(n+1,ar) ) termz = n * bknp * BESSI(n,ar) * COS(az) end if ! termr = n * bknp * birp * SIN(az) if( rp .gt. 0. ) then termp = termz / rp else termp = 0. end if ! br = br + termr*bnorm bphi = bphi + termp*bnorm bz = bz + termz*bnorm end do ! ! ........................................................ else if( g3 .ne. 0.d0 ) then ! 6-filar helical sextupole arg = 3.d0 * ka denom = ka*BESSK(2,arg) + BESSK(3,arg) bnorm = -g3 / (27d0/8d0*k*k*denom) ! do j=0,nterms n = 3*(2*j+1) ! n = 3,9,... aa = n * ka ar = n * kr az = n * (phi - kz) ! bknp = -(BESSK(n-1,aa) + BESSK(n,aa)/ka ) birp = 0.5*( BESSI(n-1,ar) + BESSI(n+1,ar) ) termz = n * bknp * BESSI(n,ar) * COS(az) ! termr = n * bknp * birp * SIN(az) if( rp .gt. 0. ) then termp = termz / rp else termp = 0. end if ! br = br + termr*bnorm bphi = bphi + termp*bnorm bz = bz + termz*bnorm end do ! ! ........................................................ else if( g4 .ne. 0.d0 ) then ! 8-filar helical octupole arg = 4.d0 * ka denom = ka*BESSK(3,arg) + BESSK(4,arg) bnorm = -g4 / (16d0/3d0*k**3*denom) ! do j=0,nterms n = 4*(2*j+1) ! n = 4,12,... aa = n * ka ar = n * kr az = n * (phi - kz) ! bknp = -(BESSK(n-1,aa) + BESSK(n,aa)/ka ) birp = 0.5*( BESSI(n-1,ar) + BESSI(n+1,ar) ) termz = n * bknp * BESSI(n,ar) * COS(az) ! termr = n * bknp * birp * SIN(az) if( rp .gt. 0. ) then termp = termz / rp else termp = 0. end if ! br = br + termr*bnorm bphi = bphi + termp*bnorm bz = bz + termz*bnorm end do end if ! end of test on multipole strength ! br = 2d0 * ftft * ka * br bphi = 2d0 * ftft * a * bphi bz = -2d0 * ka * ftft * bz bx = br*cphi - bphi*sphi by = bphi*cphi + br*sphi ! ! --------------------------------------------------------- else if( model .eq. 4 ) then ! user supplied multipoles ! Ref. T. Tominaka et al, NIM A459:398, 2001, eq. 40 mfile = fpar(2) r0 = fpar(3) ! helical reference radius hper = fpar(4) ! helix period if( hper .eq. 0. ) then write(2,*)'Period equal 0 in HELIX' iflag = -1 go to 900 end if phi0 = fpar(5) * pi/180. ! starting azimuth bref0 = fpar(6) ! elent = fpar(7) ! taper length for transverse field ! ! Solenoidal field contribution if( dz .lt. elens1 ) then ! entrance region gs = (bs2 - bs1) / elens1 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs1 + gs * dz else if( dz .gt. elens1+clens ) then ! exit region gs = (bs3 - bs2) / elens2 br = 0.5 * gs * rp b(1) = br * cphi b(2) = br * sphi b(3) = bs2 + gs * (dz-elens1-clens) else ! central region b(1) = 0. b(2) = 0. b(3) = bs2 end if ! ! Transverse field taper clent = hellen - 2.*elent if( dz .lt. elent ) then bref = bref0 * dz / elent else if( dz .gt. elent+clent ) then bref = bref0 * (hellen-dz) / elent else bref = bref0 end if ! bx = 0d0 by = 0d0 bz = 0d0 ! do n=1,nmult pern(n) = hper ! call HELFLD(rp,phi,phi0,dz,pern,bref,r0,an,bn,n,bhx,bhy,bhz) ! bx = bx + bhx by = by + bhy bz = bz + bhz end do ! ! --------------------------------------------------------- else if( model .eq. 5 ) then ! user multipoles including torsion ! Ref. J. Gallardo, MC note 522, 3 March 2008 ! iord = fpar(3) hrad = fpar(4) hpitch = fpar(5) ! call HUNT(soa,npts,dz,jt) if( jt.lt.0 .or. jt.gt.npts ) then ! outside range go to 900 else if( jt .eq. 0 ) then jt = 1 else if( jt .eq. npts ) then jt = npts - 1 end if ! nord = 5 kk = MIN( MAX(jt-nord/2,1),npts-nord ) do i=1,5 ds0(i) = soa(kk+i-1) - soa(jt) end do ! ds = s - soa(jt) ds = dz - soa(jt) ! rcf 6/10/09 ! call POLCOF(ds0,bsoa(kk),nord,c) bsp0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 bsp1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 bsp2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 bsp3 = 6.*c(4) + 24.*c(5)*ds ! call POLCOF(ds0,a0oa(kk),nord,c) a0p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 a0p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 a0p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 a0p3 = 6.*c(4) + 24.*c(5)*ds ! call POLCOF(ds0,b0oa(kk),nord,c) b0p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b0p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 b0p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 b0p3 = 6.*c(4) + 24.*c(5)*ds ! call POLCOF(ds0,a1oa(kk),nord,c) a1p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 a1p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 a1p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 ! call POLCOF(ds0,b1oa(kk),nord,c) b1p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b1p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 b1p2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 ! call POLCOF(ds0,a2oa(kk),nord,c) a2p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 a2p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 ! call POLCOF(ds0,b2oa(kk),nord,c) b2p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 b2p1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 ! call POLCOF(ds0,a3oa(kk),nord,c) a3p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 ! call POLCOF(ds0,b3oa(kk),nord,c) b3p0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 ! htrk = hrad / (hrad**2 + hpitch**2) tortion = hpitch / (hrad**2 + hpitch**2) ! x = xyz(1) y = xyz(2) h = htrk t = tortion ! if( iord .eq. 1 ) then a02 = -a1p0 - bsp1 -h*a0p0 c a11 = b1p0 + t*bsp0 ! bx = a0p0 + x*a1p0 + y*(b1p0+t*bsp0) by = b0p0 + x*(b1p0-t*bsp0) + y*a02 bz = bsp0 + x*(-t*b0p0+a0p1-2.*h*bsp0) + y*(t*a0p0+b0p1) ! else if( iord .eq. 2 ) then a02 = -a1p0 - bsp1 -h*a0p0 c a11 = b1p0 + t*bsp0 a12 = -2.*a2p0 -2.*h*a1p0 + t*t*a0p0 - a0p2 + h*bsp1 - h*a02 & + 2.*t*b0p1 c a21 = 2.*b2p0 - 2.*t*t*b0p0 + 2.*t*a0p1 - 4.*t*h*bsp0 a03 = -2.*b2p0 - h*b1p0 + t*h*bsp0 + t*t*b0p0 & - 2.*t*a0p1 - b0p2 c a11p = b1p1 a02p = -a1p1 - bsp2 - h*a0p1 ! bx = a0p0 + x*a1p0 + y*(b1p0+t*bsp0) + x**2*a2p0 & + x*y*(2.*b2p0-t*t*b0p0+t*a0p1-2.*t*h*bsp0) & + y**2*(0.5*a12+t*t*a0p0+t*b0p1) by = b0p0 + x*(b1p0-t*bsp0)+y*a02 + x**2*(b2p0+t*t*b0p0- & t*a0p1+2*t*h*bsp0) & + x*y*(a12-t*t*a0p0-t*b0p1) + y**2*(0.5*a03) bz = bsp0 + x*(-t*b0p0+a0p1-2.*h*bsp0) + y*(t*a0p0+b0p1) & + x**2*(-t*b1p0+2.*t*h*b0p0+0.5*a1p1-2.*h*a0p1+3.*h*h*bsp0) & + x*y*(t*a1p0-2.*t*h*a0p0-t*a02+b1p1-2.*h*b0p1) & + y**2*(t*b1p0+0.5*a02p) ! else if( iord .eq. 3 ) then a02 = -a1p0 - bsp1 -h*a0p0 c a11 = b1p0 + t*bsp0 a12 = -2.*a2p0 -2.*h*a1p0 + t*t*a0p0 - a0p2 + h*bsp1 - h*a02 & + 2.*t*b0p1 cjcg a21 = 2.*b2p0 - 2.*t*t*b0p0 + 2.*t*a0p1 - 4.*t*h*bsp0 a03 = -2.*b2p0 - h*b1p0 + t*h*bsp0 + t*t*b0p0 & - 2.*t*a0p1 - b0p2 a22 = -6.*a3p0-6.*h*a2p0-2.*t*t*h*a0p0-2.*h*a12-4.*t*h*b0p1 & +2.*t*t*a1p0-2.*t*t*a02-4.*t*b1p1+2.*h*a0p2 & -2.*h*h*bsp1-a1p2 c a11p = b1p1 + t*bsp1 c a11pp = b1p2 + t*bsp2 a02p = -a1p1 - bsp2 - h*a0p1 a02pp = -a1p2 - bsp3 - h*a0p2 a04 = -a22-h*a12+2.*t*t*h*a0p0-2.*t*t*a1p0+2.*t*t*a02 & +2.*t*h*b0p1-4.*t*b1p1-a02pp c jcg a31 = 6.*b3p0 - 6.*t*t*a11 + 12.*t*t*h*b0p0 + 3.*t*a1p1 c & -12.*t*h*a0p1+18.*t*h*h*bsp0 a13 = -6.*b3p0-4.*h*b2p0-2.*t*t*h*b0p0+3.*t*h*a0p1 & +2.*t*h*h*bsp0 & +4.*t*t*b1p0-2.*t*a1p1+2.*t*a02p-h*a03+h*b0p2-b1p2 a12p = -2.*a2p1 -2.*h*a1p1 + t*t*a0p1 - a0p3 + h*bsp2 & - h*a02p + 2.*t*b0p2 cjcg a21p = 2.*b2p1 - 2.*t*t*b0p1 + 2.*t*a0p2 - 4.*t*h*bsp1 a03p = -2.*b2p1 - h*b1p1 + t*h*bsp1 + t*t*b0p1 & - 2.*t*a0p2 - b0p3 ! bx = a0p0 + x*a1p0 + y*(b1p0+t*bsp0) + x**2*a2p0 c jcg & + x*y*(2.*b2p0-t*t*b0p0+t*a0p1-2.*t*h*bsp0) & + x*y*(2.*b2p0-t*t*b0p0+t*a0p1-2.*t*h*bsp0) & + y**2*(0.5*a12+t*t*a0p0+t*b0p1) + x**3*a3p0 & + x*x*y*(3.*b3p0+3.*t*h*h*bsp0-t*t*b1p0+2.*t*t*h*b0p0 & +0.5*t*a1p1-2.*t*h*a0p1) & + x*y*y*(0.5*a22+t*t*a1p0-2.*t*t*h*a0p0-t*t*a02 & +t*b1p1-2.*t*h*b0p1) & + y**3*(a13/6.+t*t*b1p0+0.5*t*a02p) ! by = b0p0 + x*(b1p0-t*bsp0) + y*a02 & + x**2*(b2p0+ t*t*b0p0-t*a0p1+2.*t*h*bsp0) ! jcg & + x*y*(a12-t*t*a0p0-t*b0p1) + y**2*(0.5*a03) & + x**3*(b3p0+t*t*b1p0-2.*t*t*h*b0p0-0.5*t*a1p1+2.*t*h*a0p1 ! jcg & - 3.*t*h*h*bsp0) ! jcg & + x*x*y*(0.5*a22-t*t*a1p0+2.*t*t*h*a0p0+t*t*a02 & -t*b1p1+2.*t*h*b0p1) & + x*y*y*(0.5*a13-t*t*b1p0-0.5*t*a02p) + y**3*(a04/6.) ! bz = bsp0 + x*(-t*b0p0+a0p1-2.*h*bsp0) + y*(t*a0p0+b0p1) & + x**2*(-t*b1p0+2.*t*h*b0p0+0.5*a1p1-2.*h*a0p1+3.*h*h*bsp0) & + x*y*(t*a1p0-2.*t*h*a0p0-t*a02+b1p1-2.*h*b0p1) & + y**2*(t*b1p0+0.5*a02p) & + x**3*(-t*b2p0+2.*t*h*b1p0-3.*t*h*h*b0p0+a2p1/3. & -h*a1p1+3.*h*h*a0p1-4.*h*h*h*bsp0) & + x*x*y*(t*a2p0-2.*t*h*a1p0+3.*t*h*h*a0p0-t*a12 & +2.*t*h*a02+b2p1-2.*h*b1p1+3.*h*h*b0p1) & + x*y*y*(-0.5*t*a03+2.*t*b2p0-2.*t*h*b1p0+0.5*a12p-h*a02p) & + y**3*(0.5*t*a12+a03p/6.) ! end if ! end if ! end of test on models ! ------------------------------------------------------- ! b(1) = b(1) + bx ! add helix to solenoid b(2) = b(2) + by b(3) = b(3) + bz ! 900 continue end ! ********************************************************* subroutine HORN(fpar,xyz,e,b) ! ! Get magnetic field from horn ! implicit double precision (a-h,o-z) include 'icommon.inc' c save jsold,jcold,jzlo real*8 xyz(3),e(3),b(3) real*8 fpar(15) character*80 title character*10 fname data bnorm/2.e-7/ data jcold/-1/,jsold/-1/ c model = fpar(1) ! do i=1,3 e(i) = 0. b(i) = 0. end do ! rp = SQRT( xyz(1)**2 + xyz(2)**2 ) if( rp .eq. 0. ) go to 900 if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) else phi = 0. end if sphi = SIN(phi) cphi = COS(phi) ! ! -------------------------------------------------------- ! <<<<< Horn Initialization >>>>> ! -------------------------------------------------------- ! ! Initialize for model 2 ! if( (kbcell.eq.0 .and. jsrg.ne.jsold ) .or. & (kbcell.eq.1 .and. jcel.ne.jcold ) ) then ! if( model.eq.2 ) then ! need data file(s) ngr = fpar(2) ! do igr=1,ngr mfile = fpar(igr+3) write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open HORN data input file' iflag = -1 go to 900 end if ! write(2,*)' HORN SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) curhorn(igr) if( ioc .ne. 0 ) then write(2,*)'Error reading CURRENT in data input file' iflag = -1 go to 900 end if curhorn(igr) = curhorn(igr) * fpar(3) write(2,'(a,e14.6)') 'Scaled current = ',curhorn(igr) read(mfile,*,iostat=ioc) nptshorn(igr) if( ioc .ne. 0 ) then write(2,*)'Error reading NPTSHORN in data input file' iflag = -1 go to 900 end if write(2,'(a,i5)') ' NPOINTS = ',nptshorn(igr) if( nptshorn(igr).le.0 .or. nptshorn(igr).gt.mxptshorn )then write(2,*) 'Error in number of HORN data points' iflag = -1 go to 900 end if ! do i=1,nptshorn(igr) read(mfile,*,iostat=ioc)idm,zhorn(i,igr), & rihorn(i,igr),rohorn(i,igr) if( ioc .ne. 0 ) then write(2,*)'Error reading HORN data input file' iflag = -1 go to 900 end if zhorn(i,igr) = zhorn(i,igr) / 100. ! [m] rihorn(i,igr) = rihorn(i,igr) / 100. ! [m] rohorn(i,igr) = rohorn(i,igr) / 100. ! [m] write(2,50) idm,zhorn(i,igr),rihorn(i,igr),rohorn(i,igr) 50 format(' #,DZ,Ri,Ro [m]: ',i5,3e12.4) end do ! close(mfile) end do ! end of loop on igr ! end if ! end of test on model = 2 ! ........................................................ ! jcold = jcel jsold = jsrg end if ! end of test on kbcell ! ! Initialization finished ! ! -------------------------------------------------------- ! <<<<< Get horn field value at current location >>>>> ! -------------------------------------------------------- ! if( model .eq. 1 ) then ! analytic horn model thet1 = fpar(2)*pi/180. ! starting angle thet2 = fpar(3)*pi/180. ! ending angle cur = fpar(4) ! current in horn ! z = xyz(3) if( z .eq. 0. ) go to 900 thet = ATAN2(rp,z) if( thet.le.thet1 .or. thet.ge.thet2 ) go to 900 ! bphi = bnorm * cur / rp b(1) = -bphi * sphi b(2) = bphi * cphi ! -------------------------------------------------------- else if( model .eq. 2 ) then ! Field from user-specified shape in data file ! igr = 1 if( fpar(2).eq.2 .and. rp.gt.fpar(6) ) igr = 2 ! cur = curhorn(igr) ! current in horn z = xyz(3) call HUNT(zhorn(1,igr),nptshorn(igr),z,jzlo) ! Watch out for initial point exactly on left edge of grid if( jzlo .eq. 0 ) jzlo = 1 if( jzlo .gt. mxptshorn ) then iflag = -81 go to 900 end if if( z .ge. zhorn(nptshorn(igr),igr) ) go to 900 ! B=0 beyond grid ! Local radii of horn fac = (z-zhorn(jzlo,igr)) /(zhorn(jzlo+1,igr)-zhorn(jzlo,igr)) rh1 = rihorn(jzlo,igr) & + fac * (rihorn(jzlo+1,igr)-rihorn(jzlo,igr)) rh2 = rohorn(jzlo,igr) & + fac * (rohorn(jzlo+1,igr)-rohorn(jzlo,igr)) ! if( rp.le.rh1 .or. rp.ge.rh2 ) go to 900 ! bphi = bnorm * cur / rp b(1) = -bphi * sphi b(2) = bphi * cphi end if ! end of test on model ! -------------------------------------------------------- 900 continue end c ********************************************************* function HX_BPHI(t) ! ! Bphi integrand for helix model 2 ! ! Ref. J. Gallardo, BNL-CAP-MUON note, 4 Sept. 2000 ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 k,num ! a = radiushx k = alphahx p0 = phi0hx r = rrhx p = phihx s = sshx ! az = p - k*t + p0 caz = COS(az) saz = SIN(az) om = SQRT(r*r+a*a-2.*r*a*caz) num = r - a*k*(s-t)*saz - a*caz denom = (om**2+(s-t)**2)**(1.5) hx_bphi = num / denom end c ********************************************************* function HX_BR(t) ! ! Br integrand for helix model 2 ! ! Ref. J. Gallardo, BNL-CAP-MUON note, 4 Sept. 2000 ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 k,num ! a = radiushx k = alphahx p0 = phi0hx r = rrhx p = phihx s = sshx ! az = p - k*t + p0 caz = COS(az) saz = SIN(az) om = SQRT(r*r+a*a-2.*r*a*caz) num = saz - k*(s-t)*caz denom = (om**2+(s-t)**2)**(1.5) hx_br = num / denom end c ********************************************************* function HX_BS(t) ! ! Bs integrand for helix model 2 ! ! Ref. J. Gallardo, BNL-CAP-MUON note, 4 Sept. 2000 ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 k,num ! a = radiushx k = alphahx p0 = phi0hx r = rrhx p = phihx s = sshx ! az = p - k*t + p0 caz = COS(az) saz = SIN(az) om = SQRT(r*r+a*a-2.*r*a*caz) num = a - r*caz denom = (om**2+(s-t)**2)**(1.5) hx_bs = num / denom end c ********************************************************* function HX_BX0(t) ! ! Bx on-axis integrand for helix model 2 ! ! Ref. J. Gallardo, BNL-CAP-MUON note, 4 Sept. 2000 ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 k,num ! a = radiushx k = alphahx p0 = phi0hx s = sshx ! az = - k*t + p0 num = SIN(az) denom = (a**2+(s-t)**2)**(1.5) hx_bx0 = num / denom end c ********************************************************* function HX_BY0(t) ! ! By on-axis integrand for helix model 2 ! ! Ref. J. Gallardo, BNL-CAP-MUON note, 4 Sept. 2000 ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 k,num ! a = radiushx k = alphahx p0 = phi0hx s = sshx ! az = - k*t + p0 num = (s-t)*COS(az) denom = (a**2+(s-t)**2)**(1.5) hx_by0 = num / denom end c ********************************************************* subroutine KICKER(fpar,xyz,tt,e,b) c c Compute fields for kickers and deflection cavities c implicit double precision (a-h,o-z) include 'icommon.inc' save jsold,tstart,tk,bmag real*8 xyz(3),e(3),b(3),x(3) real*8 fpar(15),tk(mxnbx),bmag(mxnbx) character*10 fname data jsold/-1/ c do i=1,3 e(i) = 0.d0 b(i) = 0.d0 end do ! model = fpar(1) ! --------------------------------------------------------- if( model .eq. 1 ) then ! TM210 rectangular cavity om = fpar(2) * twopi * 1.e6 ! RF angular frequency emagrf = fpar(3) * 1.e6 ! electric field strength phase = fpar(4) * pi / 180. ! initial phase a = fpar(5) ! width ! if( ABS( xp(1) ) .gt. a ) then iflag = -47 go to 900 end if ! wl = twopi * cc / om d = slen targ = om * tt + phase ! TM210 mode arg = (2./wl)**2 - (2./a)**2 if( arg .gt. 0. ) then h = 1. / SQRT( arg ) else write(2,*)'Inconsistent parameters for TM210 kicker' iflag = -1 go to 900 end if ! if( ABS( xp(2) ) .gt. h ) then iflag = -47 go to 900 end if ! ! Translate to coordinate system with origin at corner of box x(1) = xyz(1) + a / 2. x(2) = xyz(2) + h / 2. ! warg = 2. * pi * x(1) / a harg = pi * x(2) / h ! e(3) = emagrf * SIN(warg) * SIN(harg) * COS(targ) b(1) = -emagrf/cc * wl/2./h*SIN(warg)* COS(harg)*SIN(targ) b(2) = emagrf/cc *wl/a * COS(warg)* SIN(harg) * SIN(targ) ! --------------------------------------------------------- else if( model .eq. 2 ) then ! time-dependent transverse kicker if( jsrg .ne. jsold ) then ! initialize at start of region mfile = fpar(2) offset = fpar(3) ! ! tpref is the time the ref particle EXITED the region ! We want tstart to refer to the time the ref particle ENTERED xm = mass( ABS(iptypref) ) y = ppref(3) / xm bet = y / SQRT(1. + y*y) dt = slen / (bet*cc) tstart = tpref - dt + offset ! if( mfile .gt. 19 ) then ! input waveform file write(fname,80) mfile 80 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open KICKER input waveform file' iflag = -1 go to 900 end if read(mfile,*) npts if( npts.le.0 .or. npts.gt.mxnbx ) then write(2,*)'NPTS error in KICKER input waveform file' iflag = -1 go to 900 end if do i=1,npts read(mfile,*) tk(i),bmag(i) end do close(mfile) end if ! end of test on mfile jsold = jsrg ! flag initialization done end if ! end of test on jsrg ! ! ......................................................... nord = fpar(4) phi = fpar(5) * pi / 180. tr = tp - tstart if( mfile .gt. 19 ) then call HUNT(tk,npts,tr,jt) if( jt.eq.0 .or. jt.eq.npts ) go to 900 ! outside range ! if( nord .eq. 1 ) then ! linear interpolation bm = bmag(jt) + (bmag(jt+1)-bmag(jt)) / & (tk(jt+1)-tk(jt)) * (tr-tk(jt)) else ! quadratic interpolation kk = MIN( MAX(jt-nord/2,1),npts-nord ) call POLINT(tk(kk),bmag(kk),nord+1,tr,bm,df) end if ! else ! compute Bmag analytically b0 = fpar(6) slp = fpar(7) dt = fpar(8) boff = fpar(9) bm = boff + slp*b0*tr/dt end if ! b(1) = bm * COS(phi) b(2) = bm * SIN(phi) ! --------------------------------------------------------- else if( model .eq. 3 ) then ! TE011 rectangular cavity om = fpar(2) * twopi * 1.d6 ! RF angular frequency emagrf = fpar(3) * 1.d6 ! electric field strength phase = fpar(4) * pi / 180.d0 ! initial phase ! wl = twopi * cc / om d = slen targ = om * tt + phase ! TE011 mode arg = (2.d0/wl)**2 - (1.d0/d)**2 if( arg .gt. 0. ) then h = 1.d0 / SQRT( arg ) else write(2,*)'Inconsistent parameters for TE011 kicker' iflag = -1 go to 900 end if ! if( ABS( xp(2) ) .gt. h ) then iflag = -47 go to 900 end if ! ! Translate to coordinate system with origin at corner of box x(2) = xyz(2) + h / 2.d0 x(3) = xyz(3) ! harg = pi * x(2) / h zarg = pi * x(3) / d ! e(1) = -emagrf * SIN(harg)*SIN(zarg) * SIN(targ) b(2) = -emagrf/cc *wl/2d0/d * SIN(harg)*COS(zarg) * COS(targ) b(3) = emagrf/cc *wl/2d0/h * COS(harg)*SIN(zarg) * COS(targ) ! end if ! end of test on model ! --------------------------------------------------------- 900 continue return end ! ********************************************************* function KSHIF(iz,ksh,nz) ! ! Return shifted, wrapped index for an array ! ! IZ: index {1, 2, ... NZ} of input array ! KSH: shift {0 ... NZ-1} ! NZ: size of arrays ! ! Properties ! 1) KSHIF(iz,0,nz) = iz ! 2) if KSH=k => shift output by k ! 3) KSHIF(ksh,ksh,nz) = 1 ! 4) wrap output index around NZ ! if( ksh .eq. 0 ) then kshif = iz else if( iz .lt. ksh ) then kshif = iz + nz - ksh + 1 else kshif = iz - ksh + 1 end if end if ! end c ********************************************************* subroutine PARSE_G4RZ(line,nr,nz) ! ! Parse input line for G4BL cylinder field grid ! character*80 line character*5 cn ! if( line(1:8) .ne. 'cylinder' ) then nr = 0 go to 900 end if ! do j=9,80 if( line(j:j+2) .eq. 'nR=' ) then jlo = j+3 do k=1,8 if( line(jlo+k:jlo+k) .eq. ' ' ) then jhi = jlo + k - 1 exit end if end do cn = line(jlo:jhi) read(cn,'(i5)') nr go to 50 end if ! end do ! 50 continue do j=9,80 if( line(j:j+2) .eq. 'nZ=' ) then jlo = j+3 do k=1,8 if( line(jlo+k:jlo+k) .eq. ' ' ) then jhi = jlo + k - 1 exit end if end do cn = line(jlo:jhi) read(cn,'(i5)') nz go to 900 end if ! end do ! 900 continue end c ********************************************************* subroutine PARSE_G43D(line,nx,ny,nz) ! ! Parse input line for G4BL grid field grid ! character*80 line character*5 cn ! if( line(1:4) .ne. 'grid' ) then nx = 0 go to 900 end if ! do j=5,80 if( line(j:j+2) .eq. 'nX=' ) then jlo = j+3 do k=1,8 if( line(jlo+k:jlo+k) .eq. ' ' ) then jhi = jlo + k - 1 exit end if end do cn = line(jlo:jhi) read(cn,'(i5)') nx go to 50 end if ! end do ! 50 continue ! do j=5,80 if( line(j:j+2) .eq. 'nY=' ) then jlo = j+3 do k=1,8 if( line(jlo+k:jlo+k) .eq. ' ' ) then jhi = jlo + k - 1 exit end if end do cn = line(jlo:jhi) read(cn,'(i5)') ny go to 70 end if ! end do ! 70 continue ! do j=5,80 if( line(j:j+2) .eq. 'nZ=' ) then jlo = j+3 do k=1,8 if( line(jlo+k:jlo+k) .eq. ' ' ) then jhi = jlo + k - 1 exit end if end do cn = line(jlo:jhi) read(cn,'(i5)') nz go to 900 end if ! end do ! 900 continue end c ********************************************************* subroutine PILLBOX(rp,sphi,cphi,tt,e,b) ! ! Get field of TM010 pillbox RF cavity ! implicit double precision (a-h,o-z) include 'icommon.inc' real*8 kc dimension e(3),b(3) ! kc = omrf / cc rcav = j01 / kc if( rp .gt. rcav ) then iflag = -39 go to 900 end if rarg = kc * rp targ = omrf * tt + phase ! e(3) = emagrf * BESSJ0(rarg) * COS(targ) bphi = -emagrf / cc * BESSJ1(rarg) * SIN(targ) ! b(1) = -bphi * sphi b(2) = bphi * cphi ! 900 continue end c ********************************************************* subroutine QUADRUPOLE(fpar,x,e,b) c c Get field of magnetic quadrupole at position x ! Use multipole model of C. Wang MC#207 c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 x(3),e(3),b(3) real*8 lamq,lams,fpar(15),k,kp,kpp,m,mp c do i=1,3 b(i) = 0. e(i) = 0. end do ! model = fpar(1) b1 = fpar(2) ! gradient strength [T/m] xx = x(1) y = x(2) ! --------------------------------------------------------- if( model .eq. 1) then ! constant gradient, 2nd order ! b2 = fpar(3) ! sextupole strength [T/m^2] ! b(1) = b1*y + 2.*b2*xx*y b(2) = b1*xx + b2*xx*xx - b2*y*y ! ! --------------------------------------------------------- else if( model .eq. 2 ) then ! dTANH(s) quad plus sextupole, 3rd order clenq = fpar(3) ! length elenq = fpar(5) ! end length lamq = fpar(6) ! end attenuation length s0 = fpar(7) ! setupole component strength clens = fpar(8) ! sextupole length elens = fpar(9) ! sextupole end length lams = fpar(10) ! sextupole attenuation length ! s = x(3) - elenq k = b1 * FTANH(s,clenq,lamq,0) kp = b1 * FTANH(s,clenq,lamq,1) kpp = b1 * FTANH(s,clenq,lamq,2) ! if( s0 .ne. 0. ) then s = x(3) - elens m = s0 * FTANH(s,clens,lams,0) mp = s0 * FTANH(s,clens,lams,1) else m = 0. mp = 0. end if ! b(1) = k*y + 2*m*xx*y - kpp/6.*y**3 b(2) = k*xx + m*(xx*xx-y*y) - kpp/2.*xx*y*y b(3) = kp*xx*y + mp*xx*xx*y - mp/3*y**3 end if ! end test on model ! --------------------------------------------------------- 900 continue return end c ********************************************************* subroutine QUAD_SKEW(fpar,x,e,b) c c Field from skew quadrupole ! Use multipole model of C. Wang MC#207 c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 x(3),e(3),b(3) real*8 lamq,fpar(15),k,kp,kpp c do i=1,3 b(i) = 0. e(i) = 0. end do ! model = fpar(1) a1 = fpar(2) ! gradient strength [T/m] xx = x(1) y = x(2) ! --------------------------------------------------------- if( model .eq. 1) then ! constant gradient, 2nd order ! a2 = fpar(3) ! skew sextupole strength [T/m^2] ! b(1) = a1*xx + a2*xx*xx - a2*y*y b(2) = -a1*y - 2.*a2*xx*y ! --------------------------------------------------------- else if( model .eq. 2 ) then ! dTANH(s) quad only, 3rd order clenq = fpar(3) ! length elenq = fpar(5) ! end length lamq = fpar(6) ! end attenuation length ! s = x(3) - elenq k = a1 * FTANH(s,clenq,lamq,0) kp = a1 * FTANH(s,clenq,lamq,1) kpp = a1 * FTANH(s,clenq,lamq,2) ! b(1) = k*xx - kpp/2.*xx*y*y b(2) = -k*y - kpp/2.*xx*xx*y + kpp/3.*y**3 b(3) = kp/2.*(xx*xx - y*y) end if ! end test on model ! --------------------------------------------------------- 900 continue end c ********************************************************* subroutine RFDIP ! ! Get field from pillbox or SuperFish cavity in constant dipole field ! Cavity is centered axially in sector shaped region with angle 2*alpha ! Cavity can be offset radially by d ! implicit double precision (a-h,o-z) include 'icommon.inc' dimension fpar(15) real*8 kc dimension xyz(3),e(3),b(3) save d,cavl,rcav,rho,ctyp ! --------------------------------------------------------- entry RFDIP_INIT(fpar) ! set up and check consistency ! if( htrk .eq. 0.d0 ) then write(2,*) 'Called ACCEL model 11 with 0 curvature' stop end if ! d = fpar(5) ! radial offset ! d is positive when center of cavity is to left of reference trajectory cavl = fpar(6) ! cavity length ctyp = fpar(7) if( ctyp .eq. 1.d0 ) then ! SuperFish cavity call SFISH_SET(fpar) rcav = rccrf else ! pillbox kc = omrf / cc rcav = j01 / kc end if ! rho = 1.d0 / ABS(htrk) ! if( htrk .lt. 0.d0 ) then uf = rho - d - rcav else uf = rho + d - rcav end if ! if( uf .le. 0.d0 ) then write(2,*) 'ACCEL model 11', & ' WARNING: cavity radius is bigger than space available' ! stop end if ! zf = uf * TAN(slen/2.d0/rho) if( cavl .gt. 2.d0*zf ) then write(2,*) 'ACCEL model 11', & ' WARNING: cavity length is bigger than space available' ! stop end if if( cavl .eq. 0.d0 ) cavl = 2.d0 * zf return ! --------------------------------------------------------- entry RFDIP_TRK(xyz,rp,sphi,cphi,tt,e,b) ! track particle ! ! Find location of particle (up,vp,zp) in cavity coordinate system bet = (xyz(3) - slen/2.d0) / rho ! particle angle from center of region ! BET is negative coming into the cavity and positive going out sbet = SIN(bet) if( htrk .lt. 0.d0 ) then zp = sbet * (rho - rp*cphi ) else zp = sbet * (rho + rp*cphi ) end if ! if( ABS(zp) .le. cavl/2.d0 ) then cbet = COS(bet) ! if( htrk .lt. 0.d0 ) then up = rho - d - rho*cbet + rp*cphi*cbet else up = -rho - d + rho*cbet + rp*cphi*cbet end if ! vp = rp * sphi ! ruv = SQRT( up*up + vp*vp ) if( ruv .ge. rcav ) then iflag = -39 go to 900 end if if( ruv .gt. 0.d0 ) then cphuv = up / ruv sphuv = vp / ruv else cphuv = 1.d0 sphuv = 0.d0 end if ! if( ctyp .eq. 1.d0 ) then ! SuperFish cavity call SFISH_GET(zp,tt,ruv,sphuv,cphuv,e,b) else ! pillbox call PILLBOX(ruv,sphuv,cphuv,tt,e,b) end if ! end if ! 900 continue end c ********************************************************* integer function RFENDI(es,omega,phi,q,m) C C rfendi : initialize the RF ends vector field C es : gradient (V/m) C omega : 2 pi frequency (s^{-1}) C phi : phase offset (radians), 0 is crest C q : Charge, in units of the electron charge (proton : +1, electron: -1) C m : mc^2, in eV C The vector field is rfendv, with order x,px,y,py,t,E C C If the momentuma, energy, es, and m are all scaled by the same factor, C the answers will be consistent C ! Ref. J.S. Berg 22 Nov 2004 ! double precision BESSJ0,BESSJ1 double precision es,m,omega,phi,q double precision x(6),y(6) double precision k,m2,om,ph,v double precision b0,b1,carg,et,kx,kp,ps,sarg double precision c integer rfendv save k,m2,om,ph,v data c/2.99792458d8/ k = omega/c m2 = m*m om = omega ph = phi v = 0.5d0*q*es rfendi = 0 return entry RFENDV(x,y) et = x(2)*x(2)+x(4)*x(4)+m2 ps = x(6)*x(6)-et if (ps.lt.0d0) then rfendv = 1 return end if ps = sqrt(ps) b1 = k*sqrt(x(1)*x(1)+x(3)*x(3)) b0 = bessj0(b1) b1 = b1*bessj1(b1) sarg = om*x(5)+ph carg = cos(sarg) sarg = sin(sarg) kx = v/(om*k*ps*ps*ps)*carg*b1 kp = -v*x(6)/ps*carg*b0 y(1) = kx*c*x(2)*x(6) y(2) = kp*x(1) y(3) = kx*c*x(4)*x(6) y(4) = kp*x(3) y(5) = kx*et y(6) = -v*x(6)/(k*ps)*sarg*b1 rfendv = 0 end c ********************************************************* subroutine ROD(fpar,xyz,e,b) c c Field from current-carrying axial rod c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 xyz(3),e(3),b(3) real*8 fpar(15) c do i=1,3 b(i) = 0. e(i) = 0. end do ! model = fpar(1) brr = fpar(2) ! Bphi at outer radius rr = fpar(3) ! radius of rod [m] clen = fpar(4) ! length of rod [m] rp = SQRT( xyz(1)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then phi = ATAN2(xyz(2),xyz(1)) else phi = 0. end if ! if( rp .le. rr ) then bm = brr * rp / rr else bm = brr * rr / rp end if c if( model .eq. 1 ) then ! central field only ! if( xyz(3).lt.0. .or. xyz(3).gt.clen ) go to 900 b(1) = bm * SIN(phi) b(2) = -bm * COS(phi) end if c ........................................................ if( (model .eq. 2) .or. (model .eq. 3) ) then c Correct for fringe field at ends z = xyz(3) ! distance into rod zend = fpar(5) ! end length ! if( z.lt.0. .or. z.gt.clen+2*zend ) go to 900 c b(1) = bm * sin(phi) b(2) = -bm * cos(phi) ! if( z .lt. zend ) then b(1) = b(1) * z / zend b(2) = b(2) * z / zend c else if( z .gt. clen-zend ) then f = (clen - z) / zend b(1) = b(1) * f b(2) = b(2) * f end if c ........................................................ if( model .eq. 3 ) then c Correct field for minimal nonlinearity c Use A. Lennox, FNAL pbar note #269, Fig. 5 x = rp / rr if( x .lt. 0.3 ) then a2 = 0.13753 a1 = 0.66775 a0 = -0.00173 else if( x .gt. 0.8 ) then a2 = -1.00000 a1 = 2.37000 a0 = -0.65000 else a2 = 0. a1 = 0.78571 a0 = -0.02548 end if c f = (a2*x*x + a1*x + a0) b(1) = b(1) * f b(2) = b(2) * f end if end if ! end of test on model=2,3 c if( model .eq. 4 ) then ! dTANH(z) elen = fpar(5) xlam = fpar(6) z = xyz(3) if( z.le.0. .or. z.gt.clen+2.*elen ) go to 900 z = z - elen bphi = bm * FTANH(z,clen,xlam,0) b(1) = bphi * SIN(phi) b(2) = -bphi * COS(phi) end if ! if( model .eq. 5 ) then ! tapered rod bc = fpar(2) ! central field [T] rc = fpar(3) ! central radius [m] zc = fpar(4) ! central length r1 = fpar(5) ! entrance radius [m] z1 = fpar(6) ! entrance length r2 = fpar(7) ! exit radius [m] z2 = fpar(8) ! exit length ! z = xyz(3) if( z .lt. z1 ) then ! entrance taper rr = r1 + (rc-r1)/z1 * z else if( z .gt. (z1+zc) ) then ! exit taper rr = rc + (r2-rc)/z2 * (z - z1 - zc) else ! central region rr = rc end if ! if( rp .ge. rr ) then iflag = -40 go to 900 end if ! if( rp .gt. 0. ) then sphi = SIN(phi) cphi = COS(phi) bm = bc * (rc / rr) * (rp / rr) b(1) = bm * sphi b(2) = -bm * cphi else b(1) = 0. b(2) = 0. end if ! end if ! end of test on model ! 900 continue return end ! ********************************************************* subroutine scx(x,ez,by) c c Return x dependence of ez and by at coordinate x c Must be initialized first at iniscx c Entry points scxezs and scxbys to get only one field component c ! Auxiliary routine for SECCAV ! implicit none double precision c,j1c,y1c,yc,zc save c,j1c,y1c,yc,zc double precision x,ez,by,y,z,arg double precision ezeval,byeval double precision bessj0,bessj1,bessy0,bessy1 ezeval(arg) = c*(bessy1(arg)*j1c - bessj1(arg)*y1c) byeval(arg) = c*(bessy0(arg)*j1c - bessj0(arg)*y1c) arg = zc*(yc+x) ez = ezeval(arg) by = byeval(arg) return entry scxezs(x,ez) arg = zc*(yc+x) ez = ezeval(arg) return entry scxbys(x,by) arg = zc*(yc+x) by = byeval(arg) return c c Initialize scx routine; argument of bessel functions are c z(y+x) c entry iniscx(y,z) c = 2d0*atan(1d0)*y*z j1c = bessj1(z*(y-0.5d0)) y1c = bessy1(z*(y-0.5d0)) yc = y zc = z end ! ********************************************************* double precision function scxby(x) c c Function call form for the scxbys entry point above c ! Auxiliary routine for SECCAV ! implicit none double precision x,by call scxbys(x,by) scxby = by end ! ********************************************************* double precision function scxez(x) c c Function call form for the scxezs entry point above c ! Auxiliary routine for SECCAV ! implicit none double precision x,ez call scxezs(x,ez) scxez = ez end ! ********************************************************* subroutine SCXfmatr(w,h,f,rho) c c Update width w and height h, maintaining aspect ratio, so that c cavity frequency is f. Radius of curvature rho. c ! Auxiliary routine for SECCAV ! implicit none double precision w,h,f,rho double precision wc2,pi,k0,k1,k2,f0,f1,err0,err1 double precision scz1 pi = 4d0*atan(1d0) wc2 = (pi*f/1.49896229d8)**2 k1 = sqrt(((pi/h)**2 + (pi/w)**2)/wc2) k0 = 0.99d0*k1 f0 = (pi/(k0*h))**2+(scz1(rho/(k0*w))/(k0*w))**2-wc2 err0 = 1d0 c Secant iteration 100 f1 = (pi/(k1*h))**2+(scz1(rho/(k1*w))/(k1*w))**2-wc2 c Quit if found solution or solution stops improving if (f1.eq.0d0.or.abs(f1).ge.abs(f0)) goto 200 k2 = (f1*k0-f0*k1)/(f1-f0) err1 = abs(k2-k1) c Quit if solution stops improving if (err1.eq.0d0.or.err1.ge.err0) goto 200 err0 = err1 k0 = k1 k1 = k2 f0 = f1 goto 100 200 w = w*k1 h = h*k1 end ! ********************************************************* double precision function scxopt(w,h,f,rho,xoff) c c Find zero in x of by. Arguments are as for fmatr, with the c exception that rho is the radius of curvature of the reference c orbit where by=0; the rho that fmatr is called with is the c radius of curvature at the cavity center. Also, the radius of c curvature that this function is called with has a sign. See c the top of this file for a description of that. c Additionally, xoff is the displacement of the orbit relative to c the reference orbit we compute here. Positive means the reference c orbit should be on the side of by=0 toward the center of c curvature when rho>0, signs opposite when rho<0. xoff+ c is the displacement of the reference orbit toward the center c of curvature when rho>0, signs flip when rho<0. c w and h are updated c ! Auxiliary routine for SECCAV ! implicit none double precision w,h,f,rho,xoff double precision rho0,rho1,y,err0,err1,rhosgn double precision scxzer,scz1 rhosgn = sign(1d0,rho) rho0 = abs(rho)+rhosgn*xoff err0 = abs(rho) 100 call SCXfmatr(w,h,f,rho0) y = rho0/w rho1 = abs(rho) + rhosgn*xoff - w*scxzer(y,scz1(y)) err1 = abs(rho1-rho0) if (err1.eq.0.0.or.err1.ge.err0) goto 200 rho0 = rho1 err0 = err1 goto 100 200 scxopt = sign(rho1,rho)-rho-xoff end ! ********************************************************* double precision function scxzer(y,z) c c Find zero in x of by. c ! Auxiliary routine for SECCAV ! implicit none double precision y,z double precision x0,x1,x2,by0,by1,err0,err1 double precision scxby call iniscx(y,z) by0 = scxby(0d0) x0 = 1d-2 by1 = scxby(x0) x1 = -by0*x0/(by1-by0) by0 = by1 err0 = 5d-1 c Secant method 100 by1 = scxby(x1) c Quit if found solution or solution stops improving if (by1.eq.0d0.or.abs(by1).ge.abs(by0)) goto 200 x2 = (by1*x0-by0*x1)/(by1-by0) err1 = abs(x2-x1) c Quit if solution stops improving if (err1.eq.0d0.or.err1.ge.err0) goto 200 err0 = err1 x0 = x1 x1 = x2 by0 = by1 goto 100 200 scxzer = x1 end ! ********************************************************* double precision function scz1(y) ! c Find lowest zero of c Y1(z(y+1/2))J1(z(y-1/2))-J1(z(y+1/2))Y1(z(y-1/2)) c ! Auxiliary routine for SECCAV ! implicit none double precision y,yphalf,ymhalf double precision z,err,z1,err1,f,df double precision j0p,j0m,j1p,j1m,y0p,y0m,y1p,y1m double precision bessj0,bessj1,bessy0,bessy1 yphalf = y+0.5d0 ymhalf = y-0.5d0 z = 4d0*atan(1d0) err = z c Newton iteration 100 j0p = bessj0(z*yphalf) j0m = bessj0(z*ymhalf) y0p = bessy0(z*yphalf) y0m = bessy0(z*ymhalf) j1p = bessj1(z*yphalf) j1m = bessj1(z*ymhalf) y1p = bessy1(z*yphalf) y1m = bessy1(z*ymhalf) f = y1p*j1m-j1p*y1m c Quit if found the zero if (f.eq.0d0) goto 200 df = yphalf*(y0p*j1m-j0p*y1m)+ymhalf*(y1p*j0m-j1p*y0m)-2d0*f/z z1 = z-f/df err1 = abs(z-z1) c Quit if no longer improving if (err1.eq.0d0.or.err1.ge.err) goto 200 z = z1 err = err1 goto 100 200 scz1 = z return end ! ********************************************************* subroutine seccav(x,t,b,e) c c Routines for handling rectangular cross-section sector c cavities. c c Version 020429b c c To use, first call c c initsc(p,phi,rho) c c for each cavity, where c c p(2) = cavity frequency (MHz) c p(3) = cavity gradient (MV/m) c p(5) = horizontal cavity displacement (m) c p(6) = trial width c p(7) = trial height c phi = cavity phase (rad) c rho = radius of curvature of reference orbit c c The cavity height and width are adjusted to make the frequency c correct in such a way that the aspect ratio remains equal to the c aspect ratio defined by the trial width and height. Positive c rho means that the center of curvature is to the right if you are c looking in the direction of particle motion. Positive x is always c to the left, positive y is always up (again, when looking in the c direction of particle motion). Positive p(5) means that the c cavity is displaced to the left that distance from the point where c B_y=0. c c To find the fields, call c c seccav(x,t,b,e) c c where c c x(1) is the horizontal coordinate (m), in a curved coordinate c system c x(2) is the vertical coordinate (m) c t is the time (s) c b(3) are the (x,y,z) components of the magnetic field c e(3) are the (x,y,z) components of the electric field c c Report all suspected bugs to J. Scott Berg c or c c $Date: 2002/04/29 22:08:59 $ c $Revision: 1.3 $ c $Source: /home/jsberg/cvsroot/icool/sectorcav/seccav.f,v $ c c c Evaluate fields at transverse coordinates x (in m), c time t (in s), putting fields in b (T) and e (V/m). c Must call initsc first with current cavity parameters. c See initsc (below) for details. c implicit none c Parameters to seccav double precision b(3),e(3),t,x(2) c Working variables for seccav double precision by,c,cy,ez,phase,s,sy c Parameters for initsc double precision p(15),phiin,rho c Working variables for initsc double precision clight,eig,ht,y c Saved variables set by initsc double precision bxc,byc,ezc,ky,omega,phi,rhosgn,wd,xdisp double precision scxopt,scxez,scz1 save bxc,byc,ezc,ky,omega,phi,rhosgn,wd,xdisp phase = phi + omega*t c = cos(phase) s = sin(phase) cy = cos(ky*x(2)) sy = -sin(ky*x(2)) call scx(rhosgn*(x(1)-xdisp)/wd,ez,by) e(1) = 0d0 e(2) = 0d0 e(3) = ezc*ez*cy*c b(1) = bxc*ez*sy*s b(2) = byc*by*cy*s b(3) = 0d0 return c c p are the input parameters from the icool input deck c p(2) = frequency (MHz) c p(3) = gradient (MV/m) c p(4) = phase offset (deg) IGNORE! c p(5) = horizontal displacement (m) c p(6) = width (m) c p(7) = height (m) c entry initsc(p,phiin,rho) omega = p(2)*8d6*atan(1d0) phi = phiin wd = p(6) ht = p(7) if (wd.eq.0d0.and.ht.eq.0d0) then wd = 1d0 ht = 1d0 endif rhosgn = sign(1d0,rho) xdisp = p(5) + scxopt(wd,ht,p(2)*1d6,rho,p(5)) y = abs(rho+xdisp)/wd eig = scz1(y) call iniscx(y,eig) ezc = p(3)*1d6/scxez(-rhosgn*xdisp/wd) clight = 2.99792458d8 ky = 4d0*atan(1d0)/ht bxc = -ezc*ky/omega byc = ezc*eig*rhosgn/(omega*wd) return end ! ********************************************************* subroutine SC_INIT(spcharge) ! ! Initialize space-charge for this region ! ! If we want space charge, set spcharge = true ! implicit double precision (a-h,o-z) include 'icommon.inc' ! logical spcharge ! if( (spacelev.eq.1 .and. ftag(1).eq.'ACCE') & .or. (spacelev.ge.2) ) then spcharge = .true. end if ! return end ! ********************************************************* subroutine SC_KICK(loc) ! ! Give particles in bunch a space charge kick ! implicit double precision (a-h,o-z) include 'icommon.inc' ! dimension stem(2,2),esc(3),bsc(3) dimension pba(10),pzmn(10),pzmx(10) ! ! Break bunch into momentum bins in LAB pzmin = 1e9 ! find momentum limits pzmax = -1e9 do 50 ip=1,npart call READ_EVT(ip) if( ipflg .ne. 0 ) go to 50 if( pp(3) .lt. pzmin ) pzmin = pp(3) if( pp(3) .gt. pzmax ) pzmax = pp(3) pol(1) = tp ! save LAB t in pol array pol(3) = pp(3) ! save LAB pz in pol array ! => you cant do spin and space charge simultaneously call WRITE_EVT(ip) 50 continue ! ! Initialize binning nbins = 10 dpbin = (pzmax-pzmin) / nbins do ib=1,nbins pzmn(ib) = pzmin + (ib-1)*dpbin pzmx(ib) = pzmn(ib) + dpbin pba(ib) = 0d0 end do good = 0d0 ! ! Separate particles into bins binmin = 15. ! minimum number of particle per bin 60 continue ! do 70 ip=1,npart call READ_EVT(ip) if( ipflg .ne. 0 ) go to 70 good = good + 1d0 do ib=1,nbins if( pp(3).ge.pzmn(ib) .and. pp(3).le.pzmx(ib) ) then pba(ib) = pba(ib) + 1d0 exit end if end do 70 continue ! ! Check if too few particles are in a bin do ib=1,nbins if( pba(ib).lt.binmin .and. nbins.gt.1 ) then if( ib .eq. 1 ) then ! merge bin 1 with 2 do i=1,nbins-1 pzmx(i) = pzmx(i+1) end do do i=2,nbins-1 pzmn(i) = pzmn(i+1) end do ! else if( ib .eq. nbins ) then ! merge bin N with N-1 pzmx(nbins-1) = pzmx(nbins) ! else ! interior bin ibl = ib - 1 ibh = ib + 1 if( pba(ibl) .le. pba(ibh) ) then ! merge bin ib with ib-1 do i=ib,nbins-1 pzmn(i) = pzmn(i+1) end do do i=ib-1,nbins-1 pzmx(i) = pzmx(i+1) end do else ! merge bin ib with ib+1 do i=ib+1,nbins-1 pzmn(i) = pzmn(i+1) end do do i=ib,nbins-1 pzmx(i) = pzmx(i+1) end do end if end if ! end of test on ib ! nbins = nbins - 1 good = 0d0 do jb=1,nbins pba(jb) = 0d0 end do go to 60 end if end do ! end of loop on momentum bins ! ! Distribute the bunch charge do ib=1,nbins pba(ib) = pba(ib) / good * parbunsc end do ! ! --------------------------------------------------------- do ip=1,npart ! remember initial z position call READ_EVT(ip) if( ipflg .ne. 0 ) cycle zinit = xp(3) exit end do ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! do ib=1,nbins ! loop over momentum bins parb = pba(ib) ! ! Get average time and pz in LAB for this bin do i=1,2 do j=1,2 stem(i,j) = 0d0 end do end do stev = 0d0 ! do 120 ip=1,npart call READ_EVT(ip) if( ipflg .ne. 0 ) go to 120 if( pol(3).lt.pzmn(ib) .or. pol(3).gt.pzmx(ib) ) goto 120 stev = stev + 1.d0 stem(1,1) = stem(1,1) + pol(1) stem(2,1) = stem(2,1) + pol(3) 120 continue ! if( stev .gt. 0. ) then tavg = stem(1,1) / stev pzavg = stem(2,1) / stev gambun = SQRT( 1d0 + (pzavg/pmass)**2 ) betbun = pzavg / pmass / gambun zregbun = slen / 2d0 ! length of kick in LAB zregbun = zregbun / gambun ! length of kick in BUNCH else write(2,'(a,2i5)') 'No particles in SC_KICK: ',jsrg,loc if(termout) then write(*,'(a,2i5)') 'No particles in SC_KICK',jsrg,loc end if go to 800 end if ! --------------------------------------------------------- ! Get bunch size in BUNCH frame stem(1,1) = 0d0 stem(2,1) = 0d0 ! do 140 ip=1,npart call READ_EVT(ip) if( ipflg .ne. 0 ) go to 140 if( pol(3).lt.pzmn(ib) .or. pol(3).gt.pzmx(ib) ) goto 140 dt = pol(1) - tavg ! in LAB if( trainsc ) then ! collapse bunch train count = 0.5 + dt * frfbunsc * 1.e6 if( count .lt. 0. ) count = count - 1 dt = dt - INT(count) * 1.e-6 / frfbunsc end if pm = SQRT( pp(1)**2 + pp(2)**2 + pol(3)**2 ) ! in LAB e = SQRT(pm**2 + pmass**2) gamma = e / pmass ! Transform z-plane values into a spatial bunch in LAB at fixed time term = cc * dt / pmass / gamma xp(1) = xp(1) + pp(1)*term xp(2) = xp(2) + pp(2)*term xp(3) = pol(3)*term ! ! Lorentz transform position and momentum of particle into BUNCH frame xp(3) = xp(3) * gambun call LOREN4(betbun,pp,e) ! call WRITE_EVT(ip) ! save BUNCH frame values for this bin ! ! Compute rms length and radius of bunch in BUNCH frame stem(1,1) = stem(1,1) + xp(3) stem(2,1) = stem(2,1) + xp(3) * xp(3) r = SQRT( xp(1)**2 + xp(2)**2 ) stem(1,2) = stem(1,2) + r stem(2,2) = stem(2,2) + r * r 140 continue ! zavg = stem(1,1) / stev arg = stem(2,1)/stev - zavg*zavg if( arg .gt. 0d0 ) then szbsc = SQRT(arg) else write(2,'(a,2i5)') 'Cant find sigma length in SC_KICK',jsrg,loc if(termout) then write(*,'(a,2i5)') 'Cant find sigma length in SC_KICK',jsrg,loc end if go to 800 end if ! ravg = stem(1,2) / stev arg = stem(2,2)/stev - ravg*ravg if( arg .gt. 0d0 ) then srbsc = SQRT(arg) else write(2,'(a,2i5)') 'Cant find sigma radius in SC_KICK',jsrg,loc if(termout) then write(*,'(a,2i5)') 'Cant find sigma radius in SC_KICK',jsrg,loc end if go to 800 end if ! --------------------------------------------------------- ! Redo bunch statistics after throwing out 4-sigma tail particles do i=1,2 do j=1,2 stem(i,j) = 0d0 end do end do stev = 0d0 ! do 240 ip=1,npart call READ_EVT(ip) if( ipflg .ne. 0 ) go to 240 if( pol(3).lt.pzmn(ib) .or. pol(3).gt.pzmx(ib) ) goto 240 r = SQRT( xp(1)**2 + xp(2)**2 ) if( r .gt. 4.*srbsc ) go to 240 if( ABS(xp(3)-zavg) .gt. 4.*szbsc ) go to 240 ! stev = stev + 1d0 stem(1,1) = stem(1,1) + xp(3) stem(2,1) = stem(2,1) + xp(3) * xp(3) r = SQRT( xp(1)**2 + xp(2)**2 ) stem(1,2) = stem(1,2) + r stem(2,2) = stem(2,2) + r * r 240 continue ! if( stev .gt. 0d0 ) then zavg = stem(1,1) / stev arg = stem(2,1)/stev - zavg*zavg if( arg .gt. 0d0 ) then szbsc = SQRT(arg) else write(2,'(a,2i5)') 'Cant find sigma length in SC_KICK', & jsrg,loc if(termout) then write(*,'(a,2i5)') 'Cant find sigma length in SC_KICK', & jsrg,loc end if go to 800 endif ! ravg = stem(1,2) / stev arg = stem(2,2)/stev - ravg*ravg if( arg .gt. 0d0 ) then srbsc = SQRT(arg) else write(2,'(a,2i5)') 'Cant find sigma radius in SC_KICK', & jsrg,loc if(termout) then write(*,'(a,2i5)') 'Cant find sigma radius in SC_KICK', & jsrg,loc end if go to 800 endif end if ! --------------------------------------------------------- ! Apply space charge kick to all the particles do 560 ip=1,npart call READ_EVT(ip) if( ipflg .ne. 0 ) go to 560 if( pol(3).lt.pzmn(ib) .or. pol(3).gt.pzmx(ib) ) goto 560 slongsc = xp(3) radsc = SQRT( xp(1)**2 + xp(2)**2 ) ! ! Get field due to space charge in the BUNCH frame call SPACE_CHARGE(parb,esc,bsc) ! ! Kick particle momentum in BUNCH frame call VMAG(pp,pm) e = SQRT(pm**2 + pmass**2) gamma = e / pmass call PKICK(gamma,zregbun,esc,bsc) ! update momentum array PP ! ! Transform particle back to LAB frame at fixed time call LOREN4(-betbun,pp,e) ! ! Transform particles in LAB bunch back to their initial z plane call VMAG(pp,pm) gamma = SQRT( 1. + (pm/pmass)**2 ) dt = pol(1) - tavg term = cc * dt / pmass / gamma do i=1,2 xp(i) = xp(i) - pp(i)*term end do tp = pol(1) xp(3) = zinit ! call WRITE_EVT(ip) ! back in LAB frame for this bin ! if( ip .le. nprnt ) then ! print out current data write(2,550) loc,jsrg,ip,iptyp,xp,pp,tp 550 format(' Space charge loc:',i4,' region:',i4, & ' , particle:',i4,' , type:',i4,/,7e11.4) end if ! 560 continue end do ! end of loop on momentum bins ! --------------------------------------------------------- go to 900 ! 800 continue ! iflag = -34 ! 900 continue return end c ********************************************************* subroutine SEXTUPOLE(fpar,x,e,b) c c Get field of sextupole magnet ! Use multipole expansion of C. Wang, MC#207 c implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 x(3),e(3),b(3) real*8 lamb,fpar(15),m,mp c do i=1,3 b(i) = 0. e(i) = 0. end do ! model = fpar(1) b2 = fpar(2) ! sextupole strength [T/m^2] xx = x(1) y = x(2) ! --------------------------------------------------------- if( model .eq. 1) then ! constant strength, 2nd order ! b(1) = 2.* b2 * xx * y b(2) = b2 * (xx**2 - y**2) ! --------------------------------------------------------- else if( model .eq. 2 ) then ! dTANH(s) model, 3rd order clen = fpar(3) ! central length elen = fpar(5) ! end length lamb = fpar(6) ! end attenuation length ! s = x(3) - elen m = b2 * FTANH(s,clen,lamb,0) mp = b2 * FTANH(s,clen,lamb,1) ! b(1) = 2.*m*xx*y b(2) = m*(xx*xx - y*y) b(3) = mp*xx*xx*y - mp/3.*y**3 ! end if ! end test on model ! --------------------------------------------------------- 900 continue return end c ********************************************************* subroutine SFISH_GET(zz,tt,rp,sphi,cphi,e,b) c ! Interpolate Superfish cavity fields ! implicit double precision (a-h,o-z) include 'icommon.inc' c save izlo,jrlo real*8 e(3),b(3) ! if( rp .gt. rccrf ) then ! lost it radially ifail = -57 go to 900 end if ! ! Get z in coordinates of cavity grid if( ascrf .lt. 0.5 ) then ! axial symmetric case zp = ABS(zz - z2crf) ! relative z from center of gap else ! no axial symmetry zp = zz - z2crf end if ! call LOCATE2(zcrf,nzcrf,zp,izlo) call LOCATE2(rcrf,nrcrf,rp,jrlo) if( zp .eq. zcrf(nzcrf) ) izlo = nzcrf - 1 if( rp .eq. rcrf(nrcrf) ) jrlo = nrcrf - 1 ! if( (izlo.ge.1) .and. (izlo.le.nzcrf-1) .and. & (jrlo.ge.1) .and. (jrlo.le.nrcrf-1) ) then u = ( zp - zcrf(izlo) ) / dzcrf ! bilinear interpolation v = ( rp - rcrf(jrlo) ) / drcrf ! ez = (1.-u)*(1.-v)*ezcrf(jrlo,izlo) + u*(1.-v)*ezcrf(jrlo,izlo+1) & + u*v*ezcrf(jrlo+1,izlo+1) + (1.-u)*v*ezcrf(jrlo+1,izlo) er = (1.-u)*(1.-v)*ercrf(jrlo,izlo) + u*(1.-v)*ercrf(jrlo,izlo+1) & + u*v*ercrf(jrlo+1,izlo+1) + (1.-u)*v*ercrf(jrlo+1,izlo) bp = (1.-u)*(1.-v)*bpcrf(jrlo,izlo) + u*(1.-v)*bpcrf(jrlo,izlo+1) & + u*v*bpcrf(jrlo+1,izlo+1) + (1.-u)*v*bpcrf(jrlo+1,izlo) else ez = 0d0 er = 0d0 bp = 0d0 end if ! if( ascrf .lt. 0.5 ) then ! axial symmetry if( zz .lt. z2crf ) then ! cavity half outside the SF grid ersym = -1. ! Er symmetry factor else ersym = 1. end if er = er * ersym end if ! arg = omrf * tt + phase carg = COS(arg) sarg = SIN(arg) ! e(1) = er * sfnrm * cphi * carg e(2) = er * sfnrm * sphi * carg e(3) = ez * sfnrm * carg b(1) = -bp * sfnrm * sphi * sarg b(2) = bp * sfnrm * cphi * sarg ! 900 continue end c ********************************************************* subroutine SFISH_SET(fpar) ! ! Load RF cavity fields from Superfish SF7 Parmela file ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 fpar(15) character*10 fname ! sfnrm = fpar(9) emagrf = sfnrm * 1d6 mfile = fpar(8) if( mfile .eq. isffile ) then go to 900 else isffile = mfile end if ! rccrf = fpar(10) z2crf = fpar(11) ascrf = fpar(12) ! write(fname,105) mfile 105 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open Superfish rf cavity field' iflag = -1 go to 900 end if write(2,107) mfile,j7rec 107 format(' Opened Superfish RF file: ',i3,' for region: ',i5) ! ! ......................................................... read(mfile,*,iostat=ioc) zmin,zmax,nzint if( ioc .ne. 0 ) then write(2,*)'Error reading Superfish data file' iflag = -1 go to 900 end if if( nzint .le. 0 ) then write(2,*)'Error in Superfish: nzint' iflag = -1 go to 900 end if if( nzint .ge. mxzcrf ) then write(2,*)'Error in Superfish: nzint' iflag = -1 go to 900 end if ! read(mfile,*,iostat=ioc) dum if( ioc .ne. 0 ) then write(2,*)'Error reading Superfish data file' iflag = -1 go to 900 end if ! read(mfile,*,iostat=ioc) rmin,rmax,nrint if( ioc .ne. 0 ) then write(2,*)'Error reading Superfish data file' iflag = -1 go to 900 end if if( nrint .le. 0 ) then write(2,*)'Error in Superfish: nrint' iflag = -1 go to 900 end if if( nrint .ge. mxrcrf ) then write(2,*)'Error in Superfish: nrint' iflag = -1 go to 900 end if ! dzcrf = (zmax-zmin)*1d-2 / nzint ! [m] drcrf = (rmax-rmin)*1d-2 / nrint ! [m] nzcrf = nzint + 1 nrcrf = nrint + 1 ! do i=1,nzcrf zcrf(i) = zmin*1d-2 + (i-1)*dzcrf end do ! do i=1,nrcrf rcrf(i) = (i-1)*drcrf end do nrec = (nzint+1) * (nrint+1) ! do ir=1,nrcrf do iz=1,nzcrf read(mfile,*,iostat=ioc) ezcrf(ir,iz),ercrf(ir,iz), & dum,bpcrf(ir,iz) if( ioc .ne. 0 ) then write(2,*)'Error reading Superfish field' iflag = -1 go to 900 end if ezcrf(ir,iz) = ezcrf(ir,iz) * 1d6 ! [V/m] ercrf(ir,iz) = ercrf(ir,iz) * 1d6 ! [V/m] bpcrf(ir,iz) = -bpcrf(ir,iz) * 4d-7*pi ! [T] end do end do ! close(mfile) ! ......................................................... ! 900 continue end ! ********************************************************* subroutine SHEET(fpar,xyz,e,b) ! ! Get field from solenoidal current sheet ! implicit double precision (a-h,o-z) include 'icommon.inc' c save jsold,jcold,jrlo,jzlo,jgmdlo,zmapnxt real*8 xyz(3),e(3),b(3),fpar(15) character*80 title character*10 fname data bnorm/4.e-7/ data zmapnxt/-1.e9/ data eps/1.d-9/ data jcold/-1/,jsold/-1/ c model = fpar(1) ! do i=1,3 e(i) = 0. b(i) = 0. end do ! rp = SQRT( xyz(1)**2 + xyz(2)**2 ) if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) else phi = 0. end if sphi = sin(phi) cphi = cos(phi) ! ! -------------------------------------------------------- ! <<<<< Sheet Initialization >>>>> ! -------------------------------------------------------- ! ! Initialize for models 2 and 3 ! if( model.eq.2 .or. model.eq.3 ) then ! need data file jrd = fpar(10) if( (kbcell.eq.0 .and. jsrg.ne.jsold .and. jrd.eq.0) .or. & (kbcell.eq.1 .and. jcel.ne.jcold .and. jrd.eq.0) ) then ! mfile = fpar(2) cscale = fpar(11) if( cscale .eq. 0. ) cscale = 1. if( mfile .lt. 0 ) then ! allow short cut for polarity flip mfile = -mfile fcsign = -1. else fcsign = 1. end if write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open SHEET data file' iflag = -1 go to 900 end if ! write(2,*)' SHEET SUMMARY for file',mfile read(mfile,'(a80)',iostat=ioc) title write(2,'(a80)') title read(mfile,*,iostat=ioc) nelgr,fcelgr if( ioc .ne. 0 ) then write(2,*)'Error reading SHEET data file' iflag = -1 go to 900 end if write(2,'(a,i5,5x,a,e13.5)') ' NSHEETS = ',nelgr,' fcur = ', & fcelgr ! do i=1,nelgr read(mfile,*,iostat=ioc)idm,dzelgr(i),zlelgr(i), & relgr(i),curelgr(i) if( ioc .ne. 0 ) then write(2,*)'Error reading SHEET data file' iflag = -1 go to 900 end if curelgr(i) = curelgr(i) * fcelgr * fcsign * cscale write(2,50) idm,dzelgr(i),zlelgr(i), & relgr(i),curelgr(i) 50 format(' #,DZ,L,R,I: ',i5,4e12.4) end do ! close(mfile) ! ! ........................................................ if( model .eq. 3 ) then ! need grid set up kk = 1 delzgr(kk) = fpar(3) delrgr(kk) = fpar(4) zglen = fpar(5) rglen = fpar(6) dzcut = fpar(7) mfile = fpar(9) if( delzgr(kk).le.0. .or. delrgr(kk).le.0. ) then write(2,*) 'Grid dz or dr LE 0 in SHEET' iflag = -1 go to 900 end if nzgr(kk) = ANINT( zglen / delzgr(kk) + 1. ) nrgr(kk) = ANINT( rglen / delrgr(kk) + 1. ) if( nzgr(kk).gt.mxzgr .or. nrgr(kk).gt.mxrgr .or. & krad.ne.1) then write(2,*)'Number of SHEET z or r grid points exceeds max' iflag = -1 go to 900 end if ! do iz=1,nzgr(kk) if( termout .and. MOD(iz,20).eq.1 ) then write(*,80)' Computing grid:',iz,' /',nzgr(kk),' done' 80 format(a,i5,a,i5,a) end if zgr(iz,kk) = (iz-1) * delzgr(kk) ! do ir=1,nrgr(kk) rgr(ir,kk) = (ir-1) * delrgr(kk) ! ! Check that the bozo isn't asking for a grid point on a sheet do i=1,nelgr zend = dzelgr(i) + zlelgr(i) if( zgr(iz,kk).ge.dzelgr(i) .and. zgr(iz,kk).le.zend & .and. rgr(ir,kk).eq.relgr(i) ) then write(2,*) ' Grid point requested on a current sheet' jcold = jcel jsold = jsrg iflag = -1 go to 900 end if end do ! bzgr(iz,ir,kk) = 0. brgr(iz,ir,kk) = 0. ! do 100 ic=1,nelgr clen = zlelgr(ic) if( zgr(iz,kk) .lt. dzelgr(ic)-dzcut ) go to 100 if( zgr(iz,kk) .gt. dzelgr(ic)+clen+dzcut ) go to 100 a = relgr(ic) z = zgr(iz,kk) - dzelgr(ic) - clen/2. ! call BSHEET( z,clen,rgr(ir,kk),a,bz,br ) if( iflag .ne. 0 ) go to 900 ! bzgr(iz,ir,kk) = bzgr(iz,ir,kk) + bz * curelgr(ic) brgr(iz,ir,kk) = brgr(iz,ir,kk) + br * curelgr(ic) 100 continue ! sum on sheets ! bzgr(iz,ir,kk) = bzgr(iz,ir,kk) * bnorm brgr(iz,ir,kk) = brgr(iz,ir,kk) * bnorm end do ! r grid ! end do ! z grid end if ! end of test on model 3 ! ......................................................... if( model .eq. 3 ) then if( mfile.gt.19 .and. mfile.lt.100 ) then ! Save field grid to file write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open field grid output file from SHEET' iflag = -1 go to 900 end if write(mfile,'(a80)') title write(mfile,'(2i8)') nzgr(kk),nrgr(kk) do i=1,nzgr(kk) do j=1,nrgr(kk) write(mfile,'(2i6,4e13.5)') i,j,zgr(i,kk),rgr(j,kk), & bzgr(i,j,kk),brgr(i,j,kk) end do end do close(mfile) end if ! end of test on mfile end if ! end of test on model=3 ! jcold = jcel jsold = jsrg end if ! end of test on kbcell end if ! end of test on model = 2 or 3 ! --------------------------------------------------------- ! ! Initialize for model 4 ! Use absolute z coordinates here ! This model must be used with a defined reference particle ! The first particle to get here will be a reference particle ! In phasemodel=3 the z location will be set to the end of the region ! without tracking, so we must correct for that ! if( jpar.le.0 .and. phmodref.eq.3 ) then zz = xp(3) - slen else zz = xp(3) end if dzz = zz - zmapnxt + eps ! ! if( model.eq.4 .and. zz.ge.zmapnxt .and. jpar.le. 0 )then if( model.eq.4 .and. dzz.ge.0.d0 .and. jpar.le. 0 )then ! kk = 1 ! Get grid parameters corresponding to current location ! call HUNT(z1md,nmapmd,zz,jgmdlo) call HUNT(z1md,nmapmd,zz+eps,jgmdlo) if( jgmdlo .eq. 0 ) jgmdlo = 1 z1g = z1md(jgmdlo) dzg = dzmd(jgmdlo) delzgr(kk) = dzg nzgr(kk) = nzmd(jgmdlo) drg = drmd(jgmdlo) delrgr(kk) = drg nrgr(kk) = nrmd(jgmdlo) dzcut = dzcutmd(jgmdlo) mfile = fpar(2) write(2,'(a,f12.4,a)')'Computing grid at z = ',z1g,' m' ! do iz=1,nzgr(kk) if( termout .and. MOD(iz,20).eq.1 ) then write(*,80)' Computing grid:',iz,' /',nzgr(kk),' done' end if zgr(iz,kk) = z1g + (iz-1)*dzg ! Get first useful magnet for finding field at this grid location ! z = zgr(iz,kk) - dzcut ! call HUNT(z1mc,nmagmc,z,jcmdlo) ! if( jcmdlo .eq. 0 ) jcmdlo = 1 ! do ir=1,nrgr(kk) rgr(ir,kk) = (ir-1) * drg bzgr(iz,ir,kk) = 0. brgr(iz,ir,kk) = 0. ! do 200 im=1,nmagmc if( zgr(iz,kk) .lt. z1mc(im)-dzcut ) go to 200 if( zgr(iz,kk) .gt. z1mc(im)+dzmc(im)+dzcut ) go to 200 z = zgr(iz,kk) - z1mc(im) - dzmc(im)/2. drsh = drmc(im) / nshmc(im) ! radial distance between sheets curs = 1.e6 * curdmc(im) * drsh ! do is=1,nshmc(im) a = r1mc(im) + 0.5*drsh + (is-1)*drsh call BSHEET(z,dzmc(im),rgr(ir,kk),a,bz,br) bzgr(iz,ir,kk) = bzgr(iz,ir,kk) + bz*curs brgr(iz,ir,kk) = brgr(iz,ir,kk) + br*curs end do ! end of loop on sheets ! 200 continue ! end of loop on magnets ! bzgr(iz,ir,kk) = bzgr(iz,ir,kk) * bnorm brgr(iz,ir,kk) = brgr(iz,ir,kk) * bnorm end do ! end of loop on r ! continue end do ! end of loop on z ! if( jgmdlo+1 .le. nmapmd ) then zmapnxt = z1md(jgmdlo+1) else zmapnxt = 1.e9 end if ! ......................................................... mfile = fpar(2) if( mfile.gt.19 .and. mfile.lt.100 ) then ! Save field grid to file write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open field grid output file from SHEET' iflag = -1 go to 900 end if write(mfile,'(a)') 'Magnetic field grid' write(mfile,'(2i8)') nzgr(kk),nrgr(kk) do i=1,nzgr(kk) do j=1,nrgr(kk) write(mfile,'(2i6,4e13.5)') i,j,zgr(i,kk),rgr(j,kk), & bzgr(i,j,kk),brgr(i,j,kk) end do end do close(mfile) end if ! end of test on mfile ! end if ! end of test on model = 4 ! ! Initialization finished ! ! -------------------------------------------------------- ! <<<<< Get sheet field value at current location >>>>> ! -------------------------------------------------------- ! if( model .eq. 1 ) then ! exact field from single sheet a = fpar(3) ! radius of sheet clen = fpar(4) ! length of sheet curr = fpar(5) ! current z = xyz(3) - fpar(2) ! call BSHEET(z,clen,rp,a,bz,br) if( iflag .ne. 0 ) go to 900 ! b(1) = br * cphi* curr b(2) = br * sphi * curr b(3) = bz * curr ! do i=1,3 b(i) = b(i) * bnorm end do ! ! -------------------------------------------------------- else if( model .eq. 2 ) then ! Exact field from sum of sheets in data file ! do i=1,nelgr clen = zlelgr(i) a = relgr(i) z = xyz(3) - dzelgr(i) - clen/2. ! call BSHEET( z,clen,rp,a,bz,br ) if( iflag .ne. 0 ) go to 900 ! b(1) = b(1) + br * cphi* curelgr(i) b(2) = b(2) + br * sphi * curelgr(i) b(3) = b(3) + bz * curelgr(i) end do c do i=1,3 b(i) = b(i) * bnorm end do ! ! -------------------------------------------------------- else if( model.ge.3 .or. model.le.5 ) then ! take B from grid ! Interpolated field from an r-z grid if( model .eq. 3 ) then kk = 1 interp = fpar(8) z = xyz(3) else if( model .eq. 4 ) then kk = 1 interp = fpar(3) z = xp(3) else if( model .eq. 5 ) then kk = fpar(2) ! grid number interp = fpar(3) ! interpolation model if( fpar(4) .eq. 1. ) then z = xp(3) ! absolute mode else z = xyz(3) ! relative mode end if end if ! call HUNT(zgr(1,kk),nzgr(kk),z,jzlo) ! Watch out for initial point exactly on left edge of grid if( jzlo .eq. 0 ) jzlo = 1 call HUNT(rgr(1,kk),nrgr(kk),rp,jrlo) if( jrlo .eq. 0 ) jrlo = 1 if( z .ge. zgr(nzgr(kk),kk) ) go to 900 ! B=0 beyond grid if( rp .ge. rgr(nrgr(kk),kk) ) go to 900 ! if( interp .le. 1 ) then ! bi-linear interpolation call BILINEAR(jzlo,jrlo,z,rp,kk,bz,br) ! else if( interp.eq.2 .or. interp.eq.3 ) then ! Polynomial interpolation ! ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,bzgr,mxzgr,mxrgr,interp, ! & z,rp,bz,df) ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,brgr,mxzgr,mxrgr,interp, ! & z,rp,br,df) ! call BFIELDINT(jzlo,jrlo,interp,z,rp,kk,bz,br) ! end if ! end of test on interpolation ! b(1) = br * cphi b(2) = br * sphi b(3) = bz end if ! end of test on model ! -------------------------------------------------------- 900 continue return end c ********************************************************* subroutine SOLENOID(fpar,xyz,e,b) c c Get field from solenoid c implicit double precision (a-h,o-z) include 'icommon.inc' c dimension zsol(1000),bsol(1000) dimension ds0(8),c(8) save bsc,bsd,per,fns,cs,sn,iord save zsol,bsol,npts save sold save z1,b1,pow,a1,a2,a3 ! real*8 cs(0:199),sn(0:199),bsc(0:199),bsd(0:199) real*8 xyz(3),e(3),b(3) real*8 lamb,fpar(15) character*10 fname data bnorm/2.e-7/ data i11/0/ ! to test for first entry to model 11 (KTM) ! do i=1,3 e(i) = 0. b(i) = 0. end do ! model = fpar(1) b0 = fpar(2) ! field at center [T] clen = fpar(3) ! central length [m] elen = fpar(4) ! end length [m] ! rp = SQRT( xyz(1)**2 + xyz(2)**2 ) ! if( rp .gt. 0. ) then phi = ATAN2( xyz(2),xyz(1) ) sphi = SIN(phi) cphi = COS(phi) else phi = 0. sphi = 0. cphi = 1. end if ! ! --------------------------------------------------------- if( model .eq. 1) then ! linear ends + constant center Bz bstart = fpar(5) elen2 = fpar(6) bend = fpar(7) z = xyz(3) ! if( elen .gt. 0. )then g = (b0-bstart) / elen else g = 0. end if ! if( elen2 .gt. 0. )then g2 = (bend-b0) / elen2 else g2 = 0. end if ! if( z .lt. elen ) then ! entrance fringe br = -0.5 * g * rp b(1) = br * cphi b(2) = br * sphi b(3) = g * z + bstart else if( z .gt. (elen+clen) ) then ! exit fringe br = 0.5 * g2 * rp b(1) = br * cphi b(2) = br * sphi b(3) = b0 + g2 * (z - elen - clen) else ! central region b(1) = 0. b(2) = 0. b(3) = b0 end if ! --------------------------------------------------------- else if( model .eq. 2 ) then ! hyperbolic tangent Bz order = fpar(5) lamb = fpar(6) bcen = fpar(7) z = xyz(3) - elen bz0 = bcen + b0 * FTANH(z,clen,lamb,0) dbz1 = b0 * FTANH(z,clen,lamb,1) br = -rp /2. * dbz1 bz = bz0 ! if( order .gt. 1 ) then ! 3rd order dbz2 = b0 * FTANH(z,clen,lamb,2) dbz3 = b0 * FTANH(z,clen,lamb,3) br = br + rp**3 / 16. * dbz3 bz = bz - rp**2 / 4. * dbz2 end if ! if( order .gt. 3 ) then ! 5th order dbz4 = b0 * FTANH(z,clen,lamb,4) dbz5 = b0 * FTANH(z,clen,lamb,5) br = br - rp**5 / 384. * dbz5 bz = bz + rp**4 / 64. * dbz4 end if ! if( order .gt. 5 ) then ! 7th order dbz6 = b0 * FTANH(z,clen,lamb,6) dbz7 = b0 * FTANH(z,clen,lamb,7) br = br + rp**7 / 18432. * dbz7 bz = bz - rp**6 / 2304. * dbz6 end if ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! --------------------------------------------------------- else if( model .eq. 3 ) then ! sum of circular coils ! ncoils = fpar(5) a = fpar(6) ! ! Determine current that gives B0 in the center of solenoid ! if( ncoils .eq. 1 ) then cur = 2. * a * b0 / permfp ! else dzc = clen / (ncoils-1) off = -clen/2.- dzc ! do i=1,ncoils z = off + dzc*i bz = a*a / (a*a + z*z)**(1.5) b(3) = b(3) + bz end do c b(3) = b(3) * permfp / 2. cur = b0 / b(3) b(3) = 0. end if ! ! Compute field at particle location ! if( ncoils .eq. 1 ) then z = xyz(3) - elen ! call BLOOP(a,z,rp,bz,br) ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! else off = -elen + dzc ! do i=1,ncoils z = xyz(3) + off - dzc*i ! call BLOOP( a,z,rp,bz,br ) ! b(1) = b(1) + br * cphi b(2) = b(2) + br * sphi b(3) = b(3) + bz end do end if c do i=1,3 b(i) = b(i) * cur * bnorm end do ! ! --------------------------------------------------------- else if( model .eq. 4 ) then ! exact field from single sheet a = fpar(5) ! radius of sheet z = xyz(3) - elen c2 = clen / 2. ! ! Determine current density that gives B0 in center of solenoid cur = b0/permfp/c2 * SQRT(c2*c2+a*a) ! ! Compute field at particle location ! call BSHEET(z,clen,rp,a,bz,br) if( iflag .ne. 0 ) go to 900 ! b(1) = br * cphi* cur * 2.*bnorm b(2) = br * sphi * cur * 2.*bnorm b(3) = bz * cur * 2.*bnorm ! ! --------------------------------------------------------- else if( model .eq. 5 ) then ! exact field from single block r1 = fpar(5) ! inner radius of block r2 = fpar(6) ! outer radius of block c2 = clen / 2. z = xyz(3) - elen ! ! Determine current density that gives B0 in center of solenoid cur = (r2+SQRT(r2**2+c2**2))/(r1+SQRT(r1**2+c2**2)) cur = b0/permfp/c2 / LOG(cur) / 1e6 ! BBLOCK wants current in [A/mm^2] ! ! Compute field at particle location ! call BBLOCI(z,rp,clen,r1,r2,bu,bv) br = bv * bnorm * cur * 2.d6 bz = bu * bnorm * cur * 2.d6 if( iflag .ne. 0 ) go to 900 ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! ! --------------------------------------------------------- else if( model .eq. 6 ) then ! field from user-supplied map ! Interpolated field from predefined grid kk = fpar(2) ! grid number interp = fpar(3) ! interpolation model z = xyz(3) call HUNT(zgr(1,kk),nzgr(kk),z,jzlo) ! Watch out for initial point exactly on left edge of grid if( jzlo .eq. 0 ) jzlo = 1 call HUNT(rgr(1,kk),nrgr(kk),rp,jrlo) if( jrlo .eq. 0 ) jrlo = 1 if( z .ge. zgr(nzgr(kk),kk) ) go to 900 ! B=0 beyond grid if( rp .ge. rgr(nrgr(kk),kk) ) go to 900 ! if( interp .le. 1 ) then ! bi-linear interpolation call BILINEAR(jzlo,jrlo,z,rp,kk,bz,br) ! else if( interp.eq.2 .or. interp.eq.3 ) then ! Polynomial interpolation ! ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,bzgr,mxzgr,mxrgr,interp, ! & z,rp,bz,df) ! call BIPOLYINT(jzlo,jrlo,zgr,rgr,brgr,mxzgr,mxrgr,interp, ! & z,rp,br,df) ! call BFIELDINT(jzlo,jrlo,interp,z,rp,kk,bz,br) ! end if ! end of test on interpolation ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! --------------------------------------------------------- else if( model .eq. 7 ) then ! tapered solenoid bc = fpar(2) ! central field [T] rc = fpar(3) ! central radius [m] zc = fpar(4) ! central length b1 = fpar(5) ! entrance field [T] r1 = fpar(6) ! entrance radius [m] z1 = fpar(7) ! entrance length b2 = fpar(8) ! exit field [T] r2 = fpar(9) ! exit radius [m] z2 = fpar(10) ! exit length ! if( z1 .gt. 0. ) then g1 = (bc-b1) / z1 end if if( z2 .gt. 0. ) then g2 = (b2-bc) / z2 end if ! z = xyz(3) ! if( z .lt. z1 ) then ! entrance fringe br = -0.5 * g1 * rp b(1) = br * cphi b(2) = br * sphi b(3) = b1 + g1 * z rs = r1 + (rc-r1)/z1 * z c else if( z .gt. (z1+zc) ) then ! exit fringe br = -0.5 * g2 * rp b(1) = br * cphi b(2) = br * sphi b(3) = bc + g2 * (z - z1 - zc) rs = rc + (r2-rc)/z2 * (z - z1 - zc) c else ! central region b(1) = 0. b(2) = 0. b(3) = bc rs = rc end if c if( rp .ge. rs ) then iflag = -31 end if ! ! --------------------------------------------------------- else if( model .eq. 8 ) then ! hard edge ! The logic for this model is in routines BEG_OF_REGION and END_OF_REGION b(3) = fpar(2) ! ! --------------------------------------------------------- else if( model .eq. 9 ) then ! from file of Fourier coefficients ! if( (kbcell.eq.0.and.jsrg.ne.jsold) .or. & (kbcell.eq.1.and.jcel.ne.jcold) ) then mfile = fpar(2) iord = fpar(3) sclbs = fpar(4) write(fname,40) mfile 40 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open SOL Fourier input file' iflag = -1 go to 900 end if write(2,*) ' Opened file ',mfile,' for solenoid' ! read(mfile,'(a80)') title read(mfile,*) per,fns read(mfile,*) mxord ! do i=0,mxord read(mfile,*) idum,bsc(i),bsd(i) end do ! fns = fns * sclbs ! close( mfile ) ! jsold = jsrg jcold = jcel end if ! end of test on kbcell ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! Enter here for tracking for model 9 ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ s = xyz(3) ! ! Compute trig moments for this s location if( s .ne. sold ) then do n=0,mxord arg = twopi * n * s / per cs(n) = COS(arg) sn(n) = SIN(arg) end do ! ! Define the longitudinal functions and their derivatives bz0 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,0) dbz1 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,1) dbz2 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,2) dbz3 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,3) dbz4 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,4) dbz5 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,5) dbz6 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,6) dbz7 = FOURIER(cs,sn,bsc,bsd,mxord,fns,per,7) ! sold = s end if ! end of test on S .NE. SOLD ! br = -rp /2. * dbz1 bz = bz0 ! if( iord .gt. 1 ) then ! 3rd order br = br + rp**3 / 16. * dbz3 bz = bz - rp**2 / 4. * dbz2 end if ! if( iord .gt. 3 ) then ! 5th order br = br - rp**5 / 384. * dbz5 bz = bz + rp**4 / 64. * dbz4 end if ! if( iord .gt. 5 ) then ! 7th order br = br + rp**7 / 18432. * dbz7 bz = bz - rp**6 / 2304. * dbz6 end if ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! ! --------------------------------------------------------- else if( model .eq. 10 ) then ! from file of on-axis field ! if( (kbcell.eq.0.and.jsrg.ne.jsold) .or. & (kbcell.eq.1.and.jcel.ne.jcold) ) then mfile = fpar(2) iord = fpar(3) sclbs = fpar(4) write(fname,40) mfile open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Couldnt open SOL on-axis input file' iflag = -1 go to 900 end if write(2,*) ' Opened file ',mfile,' for solenoid' ! read(mfile,'(a80)') title read(mfile,*) npts ! do i=1,npts read(mfile,*) zsol(i),bsol(i) bsol(i) = bsol(i) * sclbs end do ! close( mfile ) ! jsold = jsrg jcold = jcel end if ! end of test on kbcell ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! Enter here for tracking for model 10 ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ s = xyz(3) ! call HUNT(zsol,npts,s,jt) if( jt.lt.0 .or. jt.gt.npts ) then ! outside range go to 900 else if( jt .eq. 0 ) then jt = 1 else if( jt .eq. npts ) then jt = npts - 1 end if ! nord = 8 kk = MIN( MAX(jt-nord/2,1),npts-nord ) do i=1,8 ds0(i) = zsol(kk+i-1) - zsol(jt) end do ds = s - zsol(jt) ! call POLCOF(ds0,bsol(kk),nord,c) ! bz0 = c(1) + c(2)*ds + c(3)*ds**2 + c(4)*ds**3 + c(5)*ds**4 & + c(6)*ds**5 + c(7)*ds**6 + c(8)*ds**7 dbz1 = c(2) + 2.*c(3)*ds + 3.*c(4)*ds**2 + 4.*c(5)*ds**3 & + 5.*c(6)*ds**4 + 6.*c(7)*ds**5 + 7.*c(8)*ds**6 dbz2 = 2.*c(3) + 6.*c(4)*ds + 12.*c(5)*ds**2 & + 20.*c(6)*ds**3 + 30.*c(7)*ds**4 + 42.*c(8)*ds**5 dbz3 = 6.*c(4) + 24.*c(5)*ds + 60.*c(6)*ds**2 & + 120.*c(7)*ds**3 + 210.*c(8)*ds**4 dbz4 = 24.*c(5) + 120.*c(6)*ds + 360.*c(7)*ds**2 & + 840.*c(8)*ds**3 dbz5 = 120.*c(6) + 720.*c(7)*ds + 2520.*c(8)*ds**2 dbz6 = 720.*c(7) + 5040.*c(8)*ds dbz7 = 5040.*c(8) ! br = -rp /2. * dbz1 bz = bz0 ! if( iord .gt. 1 ) then ! 3rd order br = br + rp**3 / 16. * dbz3 bz = bz - rp**2 / 4. * dbz2 end if ! if( iord .gt. 3 ) then ! 5th order br = br - rp**5 / 384. * dbz5 bz = bz + rp**4 / 64. * dbz4 end if ! if( iord .gt. 5 ) then ! 7th order br = br + rp**7 / 18432. * dbz7 bz = bz - rp**6 / 2304. * dbz6 end if ! b(1) = br * cphi b(2) = br * sphi b(3) = bz ! ! --------------------------------------------------------- else if( model .eq. 11) then ! field from inverse cubic ! if (i11 .eq. 0) then ! initialize i11 = 1 pow = fpar(2) if (pow .eq. 0.) then write(2,*) 'SOL model 11 power cannot be 0' stop end if iord = fpar(3) z1 = fpar(4) ! z of beginning of region b1 = fpar(5) ! b_z(0,z1) = axial field at z1 db1 = fpar(6) ! = d b_z(0,z1) / dz at z1 z2 = fpar(7) b2 = fpar(8) if (b2 .le. 0.) then write(2,*) 'SOL model 11 field at z2 must be positive' stop end if db2 = fpar(9) a1 = - db1/b1/pow z21 = z2 - z1 temp1 = ((b1/b2)**(1./pow) - 1.) / z21**2 temp2 = a1 / z21 a2 = 3.*temp1 - 2.*temp2 a3 = (-2.*temp1 + temp2) / z21 end if ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! Enter here for tracking for model 11 ! ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++ ! ! The inverse cubic fit ! dz = xyz(3) - z1 cubic = 1. + a1*dz + a2*dz**2 + a3*dz**3 if (cubic .le. 0.) then write(2,*) 'SOL model 11: cubic went negative' stop end if bz = b1 / cubic**pow temp3 = - pow*bz / cubic temp4 = a1 + 2.*a2*dz + 3.*a3*dz**2 dbdz1 = temp3*temp4 br = -rp / 2. * dbz1 ! if( iord .gt. 1 ) then ! 3rd order temp5 = 2.*a2 + 6.*a3*dz temp6 = -(pow + 1.)*temp3 / cubic dbz2 = temp3*temp5 + temp6*temp4**2 bz = bz - rp**2 / 4. * dbz2 temp7 = -(pow + 2.)*temp6 / cubic dbz3 = 6.*a3*temp3 + 2.*temp6*temp4*temp5 + temp7*temp4**3 br = br + rp**3 / 16. * dbz3 end if b(1) = br * cphi b(2) = br * sphi b(3) = bz ! end if ! end of test on models ! -------------------------------------------------------- 900 continue end ! ********************************************************* subroutine SPACE_CHARGE(parb,e,b) ! ! Compute space-charge field at given location in a bunch ! ! SZBSC,SRBSC: sigma length and radius of bunch ! SLONGSC,RADSC: axial and radial position of this particle ! implicit double precision (a-h,o-z) include 'icommon.inc' ! real*8 e(3),b(3),me ! logical first ! data epsi_0/8.8542d-12/ ! data first/.true./ ! ! if( first ) then ! open(unit=106,file='ellipin.dat',status='unknown') ! open(unit=107,file='ellipout.dat',status='unknown') ! first = .false. ! end if ! if( .not.(xp(1).eq.0. .and. xp(2).eq.0.) ) then phi = ATAN2( xp(2),xp(1) ) else phi = 0. end if cphi = COS( phi ) sphi = SIN( phi ) ! spcon = 1 / (4 pi eps0) * e * (# mu/bunch) spcon = 1.438d-9 * parb ! ......................................................... if( spacelev .eq. 2 ) then ! Uniformly-charged ellipsoidal bunch in free space ! in rest frame of bunch ! Ref. Reiser, p.404-5; Fubiani et al, PRSTAB 9, 064402 (2006) ! a = 1.36*srbsc ! ellipsoid radius zm = 1.36*szbsc ! ellipsoid half-length r = radsc ! radius of particle z = slongsc ! axial distance of particle from center ! aoz = a / zm if( aoz .lt. 0.01 ) aoz = 0.01d0 ! constrain aspect ratio if( aoz .gt. 0.99 ) aoz = 0.99d0 rho0 = 3d0*parb*1.602d-19 / (4d0*pi*a*a*zm) ! if( ABS(z) .gt. zm ) go to 100 ! outside the bunch rell = a * SQRT(1d0-(z/zm)**2) ! bunch radius at this z if( r .gt. rell ) go to 100 ! outside the bunch ! ! Get field INSIDE the bunch psi = SQRT( 1d0 - aoz**2 ) ! from Reiser me = 0.5d0/psi * LOG((1d0+psi)/(1d0-psi)) - 1d0 me = me * (1d0-psi**2) / psi**2 er = rho0 * (1d0-me) * r / (2d0*epsi_0) e(3) = rho0 * me * z / epsi_0 ! ! errr1=1.d0/(psi*psi) ! from solving Eq. A13 in Fubiani (J. Gallardo) ! errr1=errr1-aoz**2/psi**3*LOG((1.d0+psi)/aoz) ! errr=rho0/(2d0*epsi_0)*r*errr1 ! ezzz=rho0/(2d0*epsi_0)*z*2d0*(1d0-errr1) ! ! write(106,*)jsrg,r,z,(er-errr),(e(3)-ezzz) go to 110 ! 100 continue ! get field OUTSIDE the bunch r2=r*r ! from solving Eq. A13 in Fubiani (J. Gallardo) z2=z*z a2=a*a zm2=zm*zm aux=r2+z2-a2-zm2 xlambda=aux+sqrt(aux**2+4.d0*(zm2*r2+a2*(z2-zm2))) xlambda=0.5d0*xlambda aux1=a2/xlambda aux2=zm2/xlambda aux3=aoz**2 err1=1.d0/(SQRT(xlambda)*zm2*(1d0+aux1)*(1d0-aux3)) err2=SQRT(1d0+aux2) err3=SQRT(xlambda)*(1d0+aux1)/(zm*SQRT(1d0-aux3))* $ LOG((zm*SQRT((1d0-aux3))+ $ SQRT(xlambda)*(1d0+aux2))/SQRT(xlambda*(1d0+aux1))) err4=err1*(err2-err3) er=rho0/(2d0*epsi_0)*a2*zm*r*err4 e(3)=rho0/epsi_0*a2*zm*z* $ (1d0/((xlambda)**1.5*(1d0+aux1)*SQRT(1d0+aux2))-err4) ! ! dist = SQRT(r*r + z*z) ! from point charge approximation ! emag = spcon / dist**2 ! erpc = emag * r / dist ! ezpc = emag * z / dist ! ! write(107,'(i4,2f12.4,2e12.4)')jsrg,r,z,(er-erpc),(e(3)-ezpc) ! 110 continue e(1) = er * cphi e(2) = er * sphi b(1) = 0d0 b(2) = 0d0 b(3) = 0d0 end if ! end of test on spacelev ! --------------------------------------------------------- ! return end c ********************************************************* subroutine STUS(fpar,xx,e,b) ! ! Get field from user-supplied 3D static field grid ! implicit double precision (a-h,o-z) include 'icommon.inc' dimension xx(3),e(3),b(3),fpar(15) dimension scr(kyor*ksor+kxor) ! ! model = fpar(1) ! do i=1,3 e(i) = 0. b(i) = 0. end do ! x = xx(1) y = xx(2) s = xx(3) ! inter = fpar(2) ! interpolation flag icf = fpar(3) ! curvature flag ! call LOCATE2(xgr,nxgr,x,jx) ! returns lower grid corner call LOCATE2(ygr,nygr,y,jy) call LOCATE2(sgr,nsgr,s,jz) if( jx*jy*jz .eq. 0 ) go to 900 ! if( icf .eq. 1 ) then ! curved grid if( hmap .ne. 0.d0 ) then ! curvature is in grid file htrk = hmap else ! calculate it from local p and B htrk = htrkgr(jz) + (s-sgr(jz))*(htrkgr(jz+1)-htrkgr(jz)) gtrk = gtrkgr(jz) + (s-sgr(jz))*(gtrkgr(jz+1)-gtrkgr(jz)) end if end if ! if( inter .eq. 0 ) then ! simple linear u = (x-xgr(jx)) / dxgr v = (y-ygr(jy)) / dygr w = (s-sgr(jz)) / (dsgr*(1.+htrk*x)) ! call TRILINEAR(u,v,w,jx,jy,jz,bxgr,mxxg,mxyg,mxsg,b(1)) call TRILINEAR(u,v,w,jx,jy,jz,bygr,mxxg,mxyg,mxsg,b(2)) call TRILINEAR(u,v,w,jx,jy,jz,bsgr,mxxg,mxyg,mxsg,b(3)) ! else if( inter .eq. 1 ) then ! use B-spline ! b(1) = BS3EV(xx,cbxsp,nbfsp,nordsp,sknotsp,scr) b(2) = BS3EV(xx,cbysp,nbfsp,nordsp,sknotsp,scr) b(3) = BS3EV(xx,cbssp,nbfsp,nordsp,sknotsp,scr) ! else if( inter .ge. 2 ) then ! Find offset into grid arrays kx = MIN( MAX(jx-inter/2,1), nxgr-inter ) ky = MIN( MAX(jy-inter/2,1), nygr-inter ) kz = MIN( MAX(jz-inter/2,1), nsgr-inter ) ! call TRIPOLY(kx,ky,kz,xgr,ygr,sgr,bxgr,mxxg,mxyg,mxsg, & inter,x,y,s,htrk,b(1),df) call TRIPOLY(kx,ky,kz,xgr,ygr,sgr,bygr,mxxg,mxyg,mxsg, & inter,x,y,s,htrk,b(2),df) call TRIPOLY(kx,ky,kz,xgr,ygr,sgr,bsgr,mxxg,mxyg,mxsg, & inter,x,y,s,htrk,b(3),df) ! end if ! end of test on inter ! 900 continue end c ********************************************************* subroutine WAVEFORM(fpar) ! ! Compute induction linac voltage waveform ! based on current value of E-t phase space distribution ! ! Adapted from a routine written by Charles Kim, 1999. ! ! subroutine WAVEFORM(XMEAN,XSTAND,ii7) ! PARAMETER (NBOX=51) implicit double precision (a-h,o-z) include 'icommon.inc' ! ! real xmean(8),xstand(8) ! DIMENSION XBOX(mxnbx),YBOX(mxnbx),XBOXSAV(mxnbx),YBOXSAV(mxnbx), ! & XBOXSAV2(mxnbx),YBOXSAV2(mxnbx) ! DIMENSION WBOX(mxnbx),WBOXSAV(mxnbx),WBOXSAV2(mxnbx) ! DIMENSION PZBOX(mxnbx),PZBOXSAV(mxnbx),PZBOXSAV2(mxnbx), ! & deltPz(mxnbx) real*8 fpar(15) real*8 XBOX(mxnbx),YBOX(mxnbx),WBOX(mxnbx),PZBOX(mxnbx), & deltPz(mxnbx) DIMENSION Vmovingave(mxnbx) logical dopen character*10 fname ! real*8 xyz(3),ttt ! SAVE BOXSAV,YBOXSAV,WBOXSAV,XBOXSAV2,YBOXSAV2,WBOXSAV2,tbeg,tend SAVE tbeg,tend save ngap_left,vinduc0,vinduc1,pz_target ! ! user defined parameters ! fparm(1,1)=7 : for induction linac ! fparm(2,1)=100 : number of gaps ! fparm(3,1)=-0.0005 : start voltage MV ! fparm(4,1)= 0.002 : max swing MV ! automatically determined parameters ! fparm(5,1)=0.10 : =slen gap length meters determined by sregion parameters ! fparm(6,1)=dx : time step = (timmax - timmin ) / 50 ! fparm(7,1)=pulse duration ! ! data nplt/1/, ifnduc/0/ data dopen/.false./ cyf ! data tgtkin/0.121/ ! ! User specified parameters ngaps = fpar(2) vstart = fpar(3) vswing = fpar(4) toffset = fpar(5) tgtkin = fpar(6) tscale = fpar(7) fp8 = fpar(8) nbox = fpar(9) mfile = fpar(11) killflag = fpar(12) ! if( fpar(13) .eq. 1 ) then ! restart ifnduc = 0 else ifnduc = 1 end if ! ! n.b. the following code uses calls to READ_EVT ! We assume that the current particle data has not been changed since ! READ_EVT was originally called in SIMULATE. ! ! xyz(1)=0.0d0 ! xyz(2)=0.0d0 ! xyz(3)=0.0d0 ! ! INITIALIZING HISTOGRAM FOR TIME DISTRIBUTION AND WAVEFORM ! ! Get mean time of bunch ! bunchmean = 0. wtsum = 0. ! do 50 i=1,npart call READ_EVT(i) if( ipflg .ne. 0 ) go to 50 cyf c select mu+ with kin.Enetgy < 250 MeV c if( iabs(iptyp) .ne. 2.or.pp(3).gt.0.3397) go to 50 cyfend bunchmean = bunchmean + evtwt*tp wtsum = wtsum + evtwt 50 continue ! if( wtsum .gt. 0. ) then bunchmean = bunchmean / wtsum else bunchmean = 0. end if ! DX=tscale/NBOX ! timmax = timmin + tscale ! RF ! timmean=0.5*(timmin+timmax) ! DO N=1,NBOX ! XBOX(N)=xmean(7)+timmean+DX*(N-0.5*NBOX) ! time XBOX(N)= bunchmean + toffset + DX*(N-0.5*NBOX) ! time YBOX(N)=0.0 ! surviving number WBOX(N)=0.0 ! weight PZBOX(N)=0.0 ! average momentum deltPz(N)=0.0 ! VINDUC(N)=FPARM(3,1) VINDUC(N) = vstart END DO ! ! Calculate histograms ! do ip=1,npart do 100 ip=1,npart call READ_EVT(ip) ! if( ipflg .ne. 0) cycle cyf select muons if(iabs(iptyp) .ne. 2.or.pp(3).gt.0.3397) go to 100 cyfend if( ipflg .ne. 0) go to 100 ! if( iptypestat.ne.0 .and. iptyp.ne.iptypestat) cycle ekine=sqrt(mass(2)**2+pp(1)**2+pp(2)**2+pp(3)**2)-mass(2) DO N=1,NBOX-1 IF(Tp.GE.XBOX(N) .AND. Tp.LT.XBOX(N+1)) then WBOX(N)=WBOX(N)+evtwt cyf PZBOX(N)=PZBOX(N)+evtwt*PP(3) cyf deltPz(N)=deltPz(N)+evtwt*(PP(3)**2) PZBOX(N)=PZBOX(N)+evtwt*ekine deltPz(N)=deltPz(N)+evtwt*ekine**2 END IF END DO ! end do 100 continue ! ! averaging 1 c do n=3,nbox-2 c pzbox(n)=0.2*(pzbox(n-2)+pzbox(n-1)+pzbox(n)+pzbox(n+1) c & +pzbox(n+2)) c wbox(n)=0.2*(wbox(n-2)+wbox(n-1)+wbox(n)+wbox(n+1)+wbox(n+2)) c deltpz(n)=0.2*(deltpz(n-2)+deltpz(n-1)+deltpz(n)+deltpz(n+1) c & +deltpz(n+2)) c end do ! averaging 2 c do n=3,nbox-2 c pzbox(n)=0.2*(pzbox(n-2)+pzbox(n-1)+pzbox(n)+pzbox(n+1) c & +pzbox(n+2)) c wbox(n)=0.2*(wbox(n-2)+wbox(n-1)+wbox(n)+wbox(n+1)+wbox(n+2)) c deltpz(n)=0.2*(deltpz(n-2)+deltpz(n-1)+deltpz(n)+deltpz(n+1) c & +deltpz(n+2)) c end do ! normalizing by weights DO N=1,NBOX IF(wbox(n).GT.0.0) then pzbox(n)=PZBOX(n)/wbox(n) else PZBOX(N)=0.0 END IF END DO ! DO N=1,NBOX IF(wbox(n).GT.0.0) then ARG=deltPz(N)/wbox(n)-pzbox(n)**2 if(ARG.GT.0.0) THEN deltPz(N)=SQRT(ARG) else deltPz(N)=0.0 end if else deltPz(N)=0.0 END IF END DO ! ! averaging 3 c do n=3,nbox-2 c pzbox(n)=0.2*(pzbox(n-2)+pzbox(n-1)+pzbox(n)+pzbox(n+1) c & +pzbox(n+2)) c wbox(n)=0.2*(wbox(n-2)+wbox(n-1)+wbox(n)+wbox(n+1)+wbox(n+2)) c deltpz(n)=0.2*(deltpz(n-2)+deltpz(n-1)+deltpz(n)+deltpz(n+1) c & +deltpz(n+2)) c end do ! averaging 4 c do n=3,nbox-2 c pzbox(n)=0.2*(pzbox(n-2)+pzbox(n-1)+pzbox(n)+pzbox(n+1) c & +pzbox(n+2)) c wbox(n)=0.2*(wbox(n-2)+wbox(n-1)+wbox(n)+wbox(n+1)+wbox(n+2)) c deltpz(n)=0.2*(deltpz(n-2)+deltpz(n-1)+deltpz(n)+deltpz(n+1) c & +deltpz(n+2)) c end do ! ! find the maximum pz pzmax=-1.0e-20 do n=1,nbox if(pzbox(n).gt.pzmax) then pzmax=pzbox(n) nmax=n end if end do ! ! optimizing for the pz/time range for given accel rate and length of induction linac if(ifnduc.eq.0) then ! initialize waveform parameters at the biginning of the induction linac ! NGAP_LEFT=fparm(2,1) NGAP_LEFT= ngaps ! vinduc0=fparm(3,1)! deceleration at the beginning of the bunch due to core reset vinduc0= vstart ! deceleration at the beginning of the bunch due to core reset ! vinduc1=fparm(3,1) + fparm(4,1) vinduc1= vstart + vswing ! v0=NGAP_LEFT*fparm(3,1) v0=NGAP_LEFT* vstart end if ! pzrange=NGAP_LEFT*fparm(4,1) pzrange=NGAP_LEFT* vswing ! ! now determining optimum pz/time window ymax=-1.0d20 klast=0 ! ! do k=nmax,nbox do 150 k=nmax,nbox nbeg=0 nend=0 ! if(pzbox(k).eq.0) cycle if(pzbox(k).eq.0) go to 150 pz1=pzbox(k) pz2=pzbox(k)-pzrange if(pz2.lt.0.0) then ! pz1=pz1-pz2 ! RF pz2=0.0 klast=1 end if ! do n=k,nbox if(pzbox(n).gt.pz2 .and. pzbox(n) .le. pz1) then ybox(k)=ybox(k)+wbox(n) if(nbeg.eq.0) then nbeg=n end if nend=n end if end do ! if(ymax.lt.ybox(k)) then ymax=ybox(k) kbeg=nbeg kend=nend pzbeg=pz1 pzend=pz2 tbeg=xbox(kbeg) end if nbeg=0 nend=0 ! if(klast.eq.1) exit if(klast.eq.1) go to 160 ! end do 150 continue ! 160 continue cyf if(ifnduc.eq.0) then cyf pz_target=pzbeg+v0 c c put the target kin. energy pz_target=tgtkin cyfend cyf end if ! do k=kbeg,nbox dpz= (pzbox(kbeg)-pzbox(k)) dt=xbox(k)-xbox(kbeg) ! if(dpz.le.pzrange .and. dt.lt.fparm(7,1)) then if(dpz.le.pzrange .and. dt.lt.tscale ) then kend=k pzend=pzbox(k) tend=xbox(kend) end if end do ! pz window and target pz are determined ! end of range optimization ! ! creating waveform with the voltage limit ! do n=1,nbox do 180 n=1,nbox tinduc(n)=xbox(n) if( n .lt. nmax ) go to 180 ! RF IF(WBOX(N).GT.0) THEN vinduc(n) = (pz_target - PZBOX(n))/ngap_left if (vinduc(n).lt.vinduc0) vinduc(n) = vinduc0 END IF ! tinduc(n)=xbox(n) ! end do 180 continue ! ! MOVING AVERAGING do n=kbeg+1,NBOX-1 Vmovingave(n)=0.0 m=0 ! do k=-1,1,1 do 200 k=-1,1,1 ! if(n-k.lt.1) cycle if(n-k.lt.1) go to 200 ! if(n-k.gt.nbox) cycle if(n-k.gt.nbox) go to 200 Vmovingave(n)=Vmovingave(n)+vinduc(n-k) m=m+1 ! end do 200 continue ! if(m.gt.0) Vmovingave(n)=Vmovingave(n)/m end do ! do n=kbeg+1,NBOX-1 if(Vmovingave(n).lt.1.5e-3) then vinduc(n)=Vmovingave(n) else vinduc(n)=1.5e-3 end if end do ! c! FINAL AVERAGING c do n=kbeg+2,NBOX-2 c Vmovingave(n)=0.0 c m=0 c! do k=-2,2,1 c do 250 k=-2,2,1 c! if(n-k.lt.1) cycle c if(n-k.lt.1) go to 250 c! if(n-k.gt.nbox) cycle c if(n-k.gt.nbox) go to 250 c Vmovingave(n)=Vmovingave(n)+vinduc(n-k) c m=m+1 c! end do c250 continue c if(m.gt.0) Vmovingave(n)=Vmovingave(n)/m c end do c! c do n=kbeg+2,NBOX-2 c vinduc(n)=Vmovingave(n) c end do ! ! fparm(6,1)=dx tstepnd = dx ! if(ngap_left.gt.1) then ngap_left=ngap_left-1 ! else ! ifnduc=0 ! end of an induction linac. ! ready for the next induction linac end if ! ! time limit ! tend=tbeg+fparm(7,1) tend=tbeg+tscale ! do n=1,nbox if(xbox(n) .ge. tend) then ! vinduc(n) = (vinduc(n-1)+fparm(3,1))*fparm(8,1)-fparm(3,1) vinduc(n) = (vinduc(n-1)+vstart)*fp8-vstart else if (xbox(n) .le. tbeg) then vinduc(n) = vinduc0 end if end do ! ! Is full swing not achieved ? Vinducmin= 1.0e20 Vinducmax=-1.0e20 do n=1,nbox if(vinduc(n).gt.Vinducmax) Vinducmax=vinduc(n) if(vinduc(n).lt.Vinducmin) Vinducmin=vinduc(n) end do Vinducdif=Vinducmax-Vinducmin ! ! if(vinducdif.gt.0.0 .and. vinducdif.lt.fparm(4,1)) then if(vinducdif.gt.0.0 .and. vinducdif.lt.vswing) then ! scaling = fparm(4,1)/Vinducdif scaling = vswing/Vinducdif do n=1,nbox ! vinduc(n)=(vinduc(n)-fparm(3,1))*scaling + fparm(3,1) vinduc(n)=(vinduc(n)-vstart)*scaling + vstart end do end if ! do n=1,nbox if(vinducdif.le.0.0) then vinduc(n)=vinduc1 end if end do ! ! nplt=nplt+1 ! write(335,*) ngap_left,tbeg,tend ! ! kill particles out of the energy limit just once at the beginning ! if(ifnduc.eq.0) then if( killflag.eq.1 .and. ifnduc.eq.0 ) then pzmax=pzbeg pzmin=pzend ! do ip=1,npart do 280 ip=1,npart call read_evt(ip) if( ipflg .ne. 0 ) go to 280 ! cyf if(pp(3).lt.pzmin .or. pp(3).gt. pzmax) then if(ekine.lt.pzmin .or. ekine.gt. pzmax) then ! ipflg=-99 ipflg = -61 call write_evt(ip) write(2,260) jsrg,ip,ipflg 260 format(' FAILURE:JSRG,IP,IPFLG: ',3i8) nfail = nfail + 1 end if ! end do 280 continue ! ifnduc=1 end if ! ! saving if( mfile .gt. 19 ) then ! save diagnostics to file if( .not.dopen ) then write(fname,300) mfile 300 format('for0',i2,'.dat') open(unit=mfile, & file=pathname(1:LASTNB(pathname)) // fname, & status='unknown',iostat=ioc) if( ioc .ne. 0 ) then write(2,*)'Coldnt open waveform diagnostics file' iflag = -1 go to 900 end if dopen = .true. end if ! do n=1,NBOX write(mfile,'(2i6,8e15.7)') jsrg,n,tinduc(n),vinduc(n),pzbox(n), & wbox(n),ybox(n),deltPz(n),vinduc0,vinduc1 end do end if ! 900 continue if( jpar .eq. 0 ) then ! restore reference particle values call READ_EVT(0) end if ! return end c ********************************************************* subroutine WIGGLER(fpar,xyz,e,b) c c Get field from wiggler ! implicit double precision (a-h,o-z) include 'icommon.inc' c real*8 xyz(3),e(3),b(3) real*8 scale(3) real*8 fpar(15),lamb,Lw,kappa,kr,kx,ky,kz dimension bess(0:4),tderiv(0:4),rfact(0:4),sfact(0:4) dimension c(3) integer n,nmax ! do i=1,3 e(i) = 0. b(i) = 0. end do ! model = fpar(1) z = xyz(3) rp = SQRT(xyz(1)**2 + xyz(2)**2) c if (rp .eq. 0.) then phi = 0. sphi = 0. cphi = 1. else phi = ATAN2(xyz(2),xyz(1)) sphi = SIN(phi) cphi = COS(phi) end if c ! --------------------------------------------------------- if( model .eq. 1 ) then ! simple rotating dipole b0 = fpar(2) per = fpar(4) if( per .eq. 0d0 ) stop 'Zero wiggler period' kz = twopi * z / per b(1) = b0 * COS(kz) b(2) = b0 * SIN(kz) b(3) = fpar(3) ! ! ......................................................... else if (model .eq. 2) then ! Linear ends + constant center wiggler fields ! Adapted from G. Penn (20 June 2000) ! modified 15 November 2000 clen = fpar(2) ! central length (m) elen = fpar(3) ! end length (m) Lw = fpar(4) ! wiggler period (m) c Lw>0, right-hand twist; Lw<0, left-hand twist; don't allow Lw=0 if (Lw .eq. 0.) stop 'Zero wiggler period' c kappa = 2.*pi/Lw kr = kappa * rp c b0 = fpar(5) ! ramp in solenoid field (T) bcen = fpar(6) ! uniform solenoid field (T) c solenoid field if (elen .gt. 0.) then g = b0 / elen else g = 0. end if c if( z .lt. elen ) then ! entrance fringe br = -0.5 * g * rp b(1) = br * cphi b(2) = br * sphi b(3) = g * z + bcen else if( z .gt. (elen+clen) ) then ! exit fringe br = 0.5 * g * rp b(1) = br * cphi b(2) = br * sphi b(3) = b0 - g * (z - elen - clen) + bcen else ! central region b(1) = 0. b(2) = 0. b(3) = b0 + bcen end if c c wiggler field by harmonic scale(1) = 2./kappa scale(2) = 1./kappa**2 scale(3) = (16./27.)/kappa**3 c nmax = 3 ! up to third order - dipole, quadrupole, sextupole c parameters for each order given in triplets as below n = 0 do while (n .lt. nmax) n=n+1 bw0=fpar(7+3*(n-1)) ! amount wiggler field ramps up by bwcen=fpar(8+3*(n-1)) ! uniform wiggler field phi0=fpar(9+3*(n-1))*pi/180. c initial rotation (given in degrees) c units for b0 and bcen depend on mode order: c n=1, dipole, T; n=2, quadrupole, T/m; n=3, sextupole, T/m^2 if ((bw0 .ne. 0.) .or. (bwcen .ne. 0.)) then c this mode order gets used if (elen .gt. 0.) then g = bw0 / elen else g = 0. end if c if (z .lt. elen) then ! entrance fringe bw = bwcen + g * z dbdz = g else if (z .gt. (elen+clen)) then ! exit fringe bw = bwcen + bw0 - g * (z - elen - clen) dbdz = -g else ! central region bw = bwcen + bw0 dbdz = 0. end if c bess(m)=I_m(n*kr), modified Bessel function bess(0)=BESSI0(n*kr) bess(1)=BESSI1(n*kr) bess(2) = BESSI(2,n*kr) bess(3) = BESSI(3,n*kr) bess(4) = BESSI(4,n*kr) cw = bw*scale(n) dcdz = dbdz*scale(n) psi = n*(phi - phi0 - kappa*z) cpsi = COS(psi) spsi = SIN(psi) bphi = cw*(n*kappa/2.)*(bess(n-1)-bess(n+1))*cpsi & -0.5*dcdz*(bess(n-1)+(2*n-1)*bess(n+1))*spsi br = cw*(n*kappa/2.)*(bess(n-1)+bess(n+1))*spsi & +0.5*dcdz*(bess(n-1)-(2*n-1)*bess(n+1)+2.*n*kr*bess(n))*cpsi bz = - kr * bphi + dcdz * bess(n) * spsi b(1) = b(1) + br * cphi - bphi * sphi b(2) = b(2) + bphi * cphi + br * sphi b(3) = b(3) + bz end if ! end of test for use mode order end do ! loop on mode order c ......................................................... else if (model .eq. 3) then ! hyperbolic tangent ramp c this model only allows solenoid + dipole wiggler field component ! Adapted from G. Penn (20 June 2000) ! modified 15 November 2000 clen = fpar(2) ! central length (m) elen = fpar(3) ! end length (m) Lw = fpar(4) ! wiggler period (m) c Lw>0, right-hand twist; Lw<0, left-hand twist; don't allow Lw=0 if (Lw .eq. 0.) stop 'Zero wiggler period' c kappa = 2.*pi/Lw kr = kappa * rp c lamb = fpar(5) ! scale length for TANH b0 = fpar(6) ! ramp in solenoid field (T) bcen = fpar(7) ! uniform solenoid field (T) bw0=fpar(8) ! amount wiggler field ramps up by bwcen=fpar(9) ! uniform wiggler field phi0=fpar(10)*pi/180. c initial rotation (given in degrees) zshift = z - elen c bess(m)=I_m(kr), modified Bessel function bess(0)=BESSI0(kr) bess(1)=BESSI1(kr) bess(2) = BESSI(2,kr) bess(3) = BESSI(3,kr) bess(4) = BESSI(4,kr) c tderiv(m)=mth derivative of tanh function, c 0.5(tanh(zshift/lamb) - tanh((clen-zshift)/lamb)) tderiv(0) = FTANH(zshift,clen,lamb,0) tderiv(1) = FTANH(zshift,clen,lamb,1) tderiv(2) = FTANH(zshift,clen,lamb,2) tderiv(3) = FTANH(zshift,clen,lamb,3) tderiv(4) = FTANH(zshift,clen,lamb,4) bz0 = bcen + b0 * tderiv(0) dbz1 = b0 * tderiv(1) dbz2 = b0 * tderiv(2) dbz3 = b0 * tderiv(3) dbz4 = b0 * tderiv(4) ! br = - rp /2. * dbz1 + rp**3 / 16. * dbz3 bz = bz0 - rp**2 / 4. * dbz2 + rp**4 / 64. * dbz4 b(1) = br * cphi b(2) = br * sphi b(3) = bz ! psi = phi - phi0 - kappa*z cpsi = COS(psi) spsi = SIN(psi) ! rfact(0) = bess(0)-bess(2) rfact(1) = -2.*bess(2) rfact(2) = -kr*bess(1)+3.*bess(2) rfact(3) = kr**2*bess(2)/3.-kr*bess(3) rfact(4) = (kr**3*bess(1)-6.*kr**2*bess(2) & +15.*kr*bess(3))/12. sfact(0) = bess(0)+bess(2) sfact(1) = -2.*(kr*bess(1)-bess(2)) sfact(2) = -kr**2*bess(0)+2.*kr*bess(1)-3.*bess(2) sfact(3) = kr**3*bess(1)/3.-2.*kr**2*bess(2)/3. & +kr*bess(3) sfact(4) = (kr**4*bess(0)-3.*kr**3*bess(1) & +9.*kr**2*bess(2)-15.*kr*bess(3))/12. c bphi=(bwcen+bw0*tderiv(0))*rfact(0)*cpsi bphi=bphi+bw0*tderiv(1)*(rfact(1)-rfact(0))*spsi/kappa bphi=bphi+bw0*tderiv(2)*(rfact(2)+rfact(1))*cpsi/kappa**2 bphi=bphi+bw0*tderiv(3)*(rfact(3)-rfact(2))*spsi/kappa**3 bphi=bphi+bw0*tderiv(4)*(rfact(4)+rfact(3))*cpsi/kappa**4 bz=-kr*bphi bz=bz+kr*bw0*tderiv(1)*rfact(0)*spsi/kappa bz=bz+kr*bw0*tderiv(2)*(rfact(0)-rfact(1))*cpsi/kappa**2 bz=bz+kr*bw0*tderiv(3)*(rfact(1)+rfact(2))*spsi/kappa**3 bz=bz+kr*bw0*tderiv(4)*(rfact(2)-rfact(3))*cpsi/kappa**4 br=(bwcen+bw0*tderiv(0))*sfact(0)*spsi br=br+bw0*tderiv(1)*(sfact(0)-sfact(1))*cpsi/kappa br=br+bw0*tderiv(2)*(sfact(1)+sfact(2))*spsi/kappa**2 br=br+bw0*tderiv(3)*(sfact(2)-sfact(3))*cpsi/kappa**3 br=br+bw0*tderiv(4)*(sfact(3)+sfact(4))*spsi/kappa**4 b(1) = b(1) + br * cphi - bphi * sphi b(2) = b(2) + bphi * cphi + br * sphi b(3) = b(3) + bz ! ......................................................... else if( model .eq. 4 ) then ! planar wiggler ! J. Gallardo 27 March 2006; R. Walker, CERN 95-06, p.814 per = fpar(2) ky = fpar(3) bsol = fpar(5) nterms = fpar(6) c(1) = fpar(7) c(2) = fpar(8) c(3) = fpar(9) phi = fpar(10) ! x = xyz(1) y = xyz(2) if( per .le. 0d0 ) stop 'Bad wiggler period' kz = twopi / per if( ky .le. 0d0 ) stop 'Bad wiggler ky' ! arg = ky*ky - kz*kz if( arg .ge. 0d0 ) then kx = SQRT(arg) do i=1,nterms bx = -c(i)*kx/ky*SIN(i*kx*x)*SINH(i*ky*y)*COS(i*kz*z+phi) by = c(i)*COS(i*kx*x)*COSH(i*ky*y)*COS(i*kz*z+phi) bz = -c(i)*kz/ky*COS(i*kx*x)*SINH(i*ky*y)*SIN(i*kz*z+phi) b(1) = b(1) + bx b(2) = b(2) + by b(3) = b(3) + bz end do ! else if( arg .lt. 0d0 ) then kx = SQRT(-arg) do i=1,nterms bx = c(i)*kx/ky*SINH(i*kx*x)*SINH(i*ky*y)*COS(i*kz*z+phi) by = c(i)*COSH(i*kx*x)*COSH(i*ky*y)*COS(i*kz*z+phi) bz = -c(i)*kz/ky*COSH(i*kx*x)*SINH(i*ky*y)*SIN(i*kz*z+phi) b(1) = b(1) + bx b(2) = b(2) + by b(3) = b(3) + bz end do end if b(3) = b(3) + bsol ! ! ......................................................... else if( model .eq. 5 ) then ! helical wiggler per = fpar(2) b0 = fpar(3) bsol = fpar(4) if( per .le. 0d0 ) stop 'Bad wiggler period' kz = twopi / per rarg = kz * rp zarg = phi - kz*z ! br = b0*(BESSI0(rarg)+BESSI(2,rarg))*SIN(zarg) bz = -2.*b0*BESSI1(rarg)*COS(zarg) if( rp .gt. 0d0 ) then bp = -bz / rarg else bp = 0d0 end if ! b(1) = br*cphi - bp*sphi b(2) = br*sphi + bp*cphi b(3) = bz + bsol ! end if ! end of test on models ! ------------------------------------- 900 continue return end ! ********************************************************* subroutine ZCFCN(phi,omega,grad,slen,zstep,e,pm,de,dde) ! ! Used in phasemodel 2 ! Ref. J.S. Berg 22 Nov 2004 ! implicit none double precision de,dde,e,grad,omega,phi,pm,slen,zstep double precision ds,x(3) integer i,ns logical varstep data varstep/.false./ ! if (varstep) then write (*,*) 'varstep not yet implemented for phase model 2' else if ( zstep .ge. slen ) then ds = slen ns = 1 else ns = NINT(slen/zstep) ds = slen/DBLE(ns) end if x(1) = 0d0 x(2) = e x(3) = 0d0 do 100 i=1,ns call zcrk4(x,ds,phi,omega,grad,pm) 100 continue end if ! de = x(2)-e dde = x(3) ! end ! ********************************************************* subroutine ZCRK4(x,ds,phi,omega,grad,pm) ! ! Used in phasemodel 2 ! Ref. J.S. Berg 25 Oct 2004 implicit none double precision ds,grad,omega,phi,pm,x(3) double precision v(4,3),xx(3,3) double precision c data c/2.99792458d8/ ! C This assumes a constant gradient C One can use with other models having sinusoidal RF by replacing C grad with the gradient as a function of position v(1,1) = x(2)/(c*sqrt(x(2)**2-pm**2)) v(1,2) = grad*sin(omega*x(1)+phi) v(1,3) = grad*cos(omega*x(1)+phi) xx(1,1) = x(1) + 0.5d0*ds*v(1,1) xx(1,2) = x(2) + 0.5d0*ds*v(1,2) xx(1,3) = x(3) + 0.5d0*ds*v(1,3) v(2,1) = xx(1,2)/(c*sqrt(xx(1,2)**2-pm**2)) v(2,2) = grad*sin(omega*xx(1,1)+phi) v(2,3) = grad*cos(omega*xx(1,1)+phi) xx(2,1) = x(1) + 0.5d0*ds*v(2,1) xx(2,2) = x(2) + 0.5d0*ds*v(2,2) xx(2,3) = x(3) + 0.5d0*ds*v(2,3) v(3,1) = xx(2,2)/(c*sqrt(xx(2,2)**2-pm**2)) v(3,2) = grad*sin(omega*xx(2,1)+phi) v(3,3) = grad*cos(omega*xx(2,1)+phi) xx(3,1) = x(1) + ds*v(3,1) xx(3,2) = x(2) + ds*v(3,2) xx(3,3) = x(3) + ds*v(3,3) v(4,1) = xx(3,2)/(c*sqrt(xx(3,2)**2-pm**2)) v(4,2) = grad*sin(omega*xx(3,1)+phi) v(4,3) = grad*cos(omega*xx(3,1)+phi) x(1) = x(1) + ds/6d0*(v(1,1)+2d0*(v(2,1)+v(3,1))+v(4,1)) x(2) = x(2) + ds/6d0*(v(1,2)+2d0*(v(2,2)+v(3,2))+v(4,2)) x(3) = x(3) + ds/6d0*(v(1,3)+2d0*(v(2,3)+v(3,3))+v(4,3)) ! end