-
Notifications
You must be signed in to change notification settings - Fork 39
/
Copy pathFS_trrun.f90
7229 lines (6563 loc) · 243 KB
/
FS_trrun.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
LOGICAL FUNCTION is_thin()
double precision node_value, el
el = node_value('l ')
is_thin = el.eq.0d0;
END FUNCTION is_thin
LOGICAL FUNCTION is_drift()
double precision node_value, code
code = node_value('mad8_type ')
is_drift = code.eq.1;
END FUNCTION is_drift
LOGICAL FUNCTION is_dipole()
double precision node_value, code
code = node_value('mad8_type ')
is_dipole = code.eq.3;
END FUNCTION is_dipole
LOGICAL FUNCTION is_matrix()
double precision node_value, code
code = node_value('mad8_type ')
is_matrix = code.eq.4;
END FUNCTION is_matrix
LOGICAL FUNCTION is_quad()
double precision node_value, code
code = node_value('mad8_type ')
is_quad = code.eq.5;
END FUNCTION is_quad
subroutine trrun(switch,turns,orbit0,rt,part_id,last_turn, &
last_pos,z,dxt,dyt,last_orbit,eigen,coords,e_flag,code_buf,l_buf)
use twtrrfi
use bbfi
use time_varfi
use spch_bbfi
use twiss0fi
use name_lenfi
use trackfi
use fasterror
use SCdat
implicit none
!----------------------------------------------------------------------*
! Purpose: *
! Interface RUN and DYNAP command to tracking routine *
! *
!-- Input: *
! switch (int) 1: RUN, 2: DYNAP fastune *
! turns (int) number of turns to track *
! orbit0 (dp. array) start of closed orbit *
! rt (dp. matrix) one-turn matrix *
!-- Output: *
! eigen dp(6,6) eigenvector matrix *
! coords dp(6,0:turns,npart) (only switch > 1) particle coords. *
! e_flag (int) (only switch > 1) 0: OK, else part. lost *
!-- Buffers: *
! part_id int(npart) *
! last_turn int(npart) *
! last_pos dp(npart) *
! z dp(6,npart) *
! last_orbit dp(6,npart) *
! code_buf int(nelem) local mad-8 code storage *
! l_buf dp(nelem) local length storage *
!----------------------------------------------------------------------*
logical onepass,onetable,last_out,info,aperflag,doupdate,first, &
bb_sxy_update,virgin_state,emittance_update, &
checkpnt_restart,fast_error_func,exit_loss_turn,I_OPEN
integer j,code,restart_sequ,advance_node, &
node_al_errors,n_align,nlm,jmax,j_tot,turn,turns,i,k,get_option, &
ffile,SWITCH,nint,ndble,nchar,part_id(*),last_turn(*),char_l, &
segment, e_flag, nobs,lobs,int_arr(1),tot_segm,code_buf(*), &
tot_turn,max_part,I_NUMBER
parameter(max_part=20000)
integer part_id_keep(max_part),last_turn_keep(max_part)
double precision tmp_d,orbit0(6),orbit(6),el,rt(6,6),re(6,6), &
al_errors(align_max),z(6,*),zz(6),dxt(*),dyt(*),eigen(6,6),sum, &
node_value,get_variable,last_pos(*),last_orbit(6,*), &
maxaper(6),get_value,obs_orb(6),coords(6,0:turns,*),l_buf(*), &
betx_start,bety_start,alfx_start,alfy_start,gamx_start,gamy_start,&
dx_start,dpx_start,dy_start,dpy_start,deltap, &
N_ions_in_beam, Npart_gain, t_rms, pt_rms, &
N_ions_ini, n_ions_macro, sigma_z_ini, z_factor, &
N_ions_for_bb,z_keep(6,max_part),ex_rms0,ey_rms0,sigma_p0,sigma_z0
double precision sgpr(6,6),npart,sc3d_emit(5)
data ex_rms0,ey_rms0,sigma_p0,sigma_z0 / 0d0, 0d0, 0d0, 0d0 /
character(12) tol_a, char_a
character(20) text
integer, external :: get_file_unit
!VVK 20100321 -------------------------------------------------
integer i_part ! local counter
double precision Summ_t_mean ! local for mean value
double precision Summ_t_square ! local for rms value
!-------------------------------------------------------------------
integer iiii,jjjj
!hbu
double precision spos
!hbu
character(4) vec_names(7)
!hbu
character(name_len) el_name
character(120) msg
data tol_a,char_a / 'maxaper ', ' ' /
!hbu
data vec_names / 'x', 'px', 'y', 'py', 't', 'pt','s' /
save first
data first / .true. /
save tot_turn
data tot_turn / 0 /
save betx_start, bety_start, &
alfx_start, alfy_start, &
gamx_start, gamy_start, &
dx_start, dpx_start, &
dy_start, dpy_start, &
N_ions_in_beam, Npart_gain, t_rms, pt_rms, &
N_ions_ini, n_ions_macro, sigma_z_ini, z_factor, &
N_ions_for_bb,z_keep,part_id_keep, &
last_turn_keep,jmax,segment, &
ex_rms0,ey_rms0,sigma_p0,sigma_z0
data betx_start, bety_start, &
alfx_start, alfy_start, &
gamx_start, gamy_start, &
dx_start, dpx_start, &
dy_start, dpy_start / 1d0, 1d0, 0d0, 0d0, &
0d0, 0d0, 0d0, 0d0, 0d0, 0d0 /
logical is_drift, is_thin, is_quad, is_matrix, is_dipole
!***************************************************************************
sr = reshape(&
(/0d0,1d0,0d0,0d0,0d0,0d0, -1d0,0d0, 0d0,0d0, 0d0,0d0,&
0d0,0d0,0d0,1d0,0d0,0d0, 0d0,0d0,-1d0,0d0, 0d0,0d0,&
0d0,0d0,0d0,0d0,0d0,1d0, 0d0,0d0, 0d0,0d0,-1d0,0d0/), shape(sr))
sr = TRANSPOSE(sr)
!frs 07.06.18: Needed for sc_3d_kick
arad = get_value('probe ','arad ')
!frs 02.05.16: Needed for sc_3d_beamsize
!-------added by Yipeng SUN 01-12-2008--------------
deltap = get_value('probe ','deltap ')
if(deltap.eq.0d0) then
onepass = get_option('onepass ') .ne. 0
if(onepass) then
else
call trclor(switch,orbit0)
endif
else
endif
!---- Set fast_error_func flag to use faster error function
!---- including tables. Thanks to late G. Erskine
fast_error_func = get_option('fast_error_func ') .ne. 0
if(fast_error_func.and..not.fasterror_on) then
call wzset
fasterror_on = .true.
endif
!---- get options for space charge variables
exit_loss_turn = get_option('exit_loss_turn ') .ne. 0
bb_sxy_update = get_option('bb_sxy_update ') .ne. 0
sc_3d_kick = get_option('sc_3d_kick ') .ne. 0
sc_3d_beamsize = get_option('sc_3d_beamsize ') .ne. 0
if(sc_3d_beamsize) then
i_spch = 0 !a special spch-update counter
j = restart_sequ()
11 continue
code = node_value('mad8_type ')
if(code.eq.22) then
i_spch = i_spch+1
sc_intstr=2d0*0.5d0**1.5d0*arad * get_value('probe ', 'npart ') / get_value('probe ','gamma ')
sc_charge(i_spch)=node_value('charge ')
do j=1,6
sc_map(i_spch,j,j)=1d0
enddo
endif
if (advance_node().ne.0) then
j=j+1
go to 11
endif
i_spch = 0 !a special spch-update counter
endif
sc_3d_periodic = get_option('sc_3d_periodic ') .ne. 0
! sc_3d_damp = get_option('sc_3d_damp ') .ne. 0 ! replaced by periodic approach
! sc_3d_damp_amp = get_value('run ','sc_3d_damp_amp ') ! replaced by periodic approach
checkpnt_restart = get_value('run ', 'checkpnt_restart ') .ne. 0d0
emittance_update = get_option('emittance_update ') .ne. 0
if(sc_3d_kick) emittance_update=.FALSE.
virgin_state = get_value('run ', 'virgin_state ') .ne. 0d0
if(switch.eq.1) then
! 2015-Feb-23 16:20:19 ghislain: open file only when necessary
if (bb_sxy_update) then
inquire (file='checkpoint_restart.dat', OPENED=I_OPEN, NUMBER=I_NUMBER)
if (I_OPEN) close(I_NUMBER)
unit_chpt=get_file_unit(lu_max)
if(checkpnt_restart) then
open(unit_chpt,file='checkpoint_restart.dat',form='unformat&
&ted',status='old')
else
open(unit_chpt,file='checkpoint_restart.dat',form='unformatted')
endif
endif
else
bb_sxy_update = .false.
checkpnt_restart = .false.
endif
if(bb_sxy_update) then
if(virgin_state) first=.true.
call table_input( &
betx_start, bety_start, &
alfx_start, alfy_start, &
gamx_start, gamy_start, &
dx_start, dpx_start, &
dy_start, dpy_start)
if(sc_3d_kick.and.first) call mymap()
if(first) call make_bb6d_ixy(turns)
endif
if(fsecarb) then
write (msg,*) 'Second order terms of arbitrary Matrix not '// &
'allowed for tracking.'
call aawarn('trrun: ',msg)
return
endif
aperflag = get_option('aperture ') .ne. 0
e_flag = 0
! flag to avoid double entry of last line
last_out = .false.
onepass = get_option('onepass ') .ne. 0
onetable = get_option('onetable ') .ne. 0
info = get_option('info ') * get_option('warn ') .gt. 0
if(onepass) call m66one(rt)
call trbegn(rt,eigen)
if (switch .eq. 1) then
ffile = get_value('run ', 'ffile ')
else
ffile = 1
endif
! for one_table case
if(first) then
segment = 0
endif
tot_segm = turns / ffile + 1
if (mod(turns, ffile) .ne. 0) tot_segm = tot_segm + 1
if(first) then
if(checkpnt_restart) then
read(unit_chpt,END=100) jmax
read(unit_chpt,END=100) Ex_rms0
read(unit_chpt,END=100) Ey_rms0
do i = 1, jmax
do j=1,6
read(unit_chpt,END=100) z(j,i)
enddo
enddo
read(unit_chpt,END=100) sigma_t
read(unit_chpt,END=100) mean_t
read(unit_chpt,END=100) N_ini
else
call trinicmd(switch,orbit0,eigen,jmax,z,turns,coords)
N_ini = jmax
endif
R_part = dble(jmax)/dble(N_ini)
!--- set particle id
!$OMP PARALLEL PRIVATE(k)
!$OMP DO
do k=1,jmax
part_id(k) = k
enddo
!$OMP END DO
!$OMP END PARALLEL
else
if(jmax.eq.0) then
call aawarn('trrun: ',&
'Particle Number equals zero: exit from trrun')
return
endif
if(jmax.gt.max_part) then
write(text, '(1p,d20.12)') max_part
call aafail('TRRUN: ','Fatal: Maximum Particle Number exceeded =' // text)
endif
!$OMP PARALLEL PRIVATE(i,j)
!$OMP DO
do i = 1, jmax
last_turn(i)=last_turn_keep(i)
part_id(i)=part_id_keep(i)
do j=1,6
z(j,i)=z_keep(j,i)
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
endif
!--- jmax may be reduced by particle loss - keep number in j_tot
j_tot = jmax
!--- get vector of six coordinate maxapers (both RUN and DYNAP)
call comm_para(tol_a, nint, ndble, nchar, int_arr, maxaper, char_a, char_l)
!hbu--- init info for tables initial s position is 0
!hbu initial s position is 0
spos=0
!hbu start of line, element 0
nlm=0
!hbu
el_name='start '
!---- Initialize kinematics and orbit
bet0 = get_value('beam ','beta ')
betas = get_value('probe ','beta ')
gammas= get_value('probe ','gamma ')
bet0i = 1d0 / bet0
beti = 1d0 / betas
dtbyds = get_value('probe ','dtbyds ')
deltas = get_variable('track_deltap ')
dorad = get_value('probe ','radiate ') .ne. 0d0
dodamp = get_option('damp ') .ne. 0
dorand = get_option('quantum ') .ne. 0
if(first) then
!--- enter start coordinates in summary table
if (switch .eq. 1) then
t_max=get_value('probe ','circ ')/get_value('run ', 'track_harmon ')/betas
pt_max = get_value('run ', 'deltap_max ')
pt_max=(sqrt((betas*(pt_max+1d0))**2+1d0/gammas**2)-1d0)*beti
else
t_max=1d20
pt_max=1d20
endif
do i = 1,j_tot
if(abs(z(5,i)).gt.t_max) then
write(text, '(1p,d13.5,a1,i6)') t_max,"p",i
call aafail('TRACK_INITIAL: ','Fatal: T-Coordinate larger than' // text)
endif
if(abs(z(6,i)).gt.pt_max) then
write(text, '(1p,d13.5,a1,i6)') pt_max,"p",i
call aafail('TRACK_INITIAL: ','Fatal: PT-Coordinate larger than' // text)
endif
tmp_d = i
call double_to_table_curr('tracksumm ', 'number ', tmp_d)
tmp_d = tot_turn
call double_to_table_curr('tracksumm ', 'turn ', tmp_d)
do j = 1, 6
tmp_d = z(j,i) - orbit0(j)
call double_to_table_curr('tracksumm ', vec_names(j), tmp_d)
enddo
!hbu add s
call double_to_table_curr('tracksumm ',vec_names(7),spos)
call augment_count('tracksumm ')
enddo
endif
!--- enter first turn, and possibly eigen in tables
if (switch .eq. 1) then
if (onetable) then
if(first) then
call track_pteigen(eigen)
!hbu add s, node id and name
call tt_putone(jmax, tot_turn, tot_segm, segment, part_id, &
z, orbit0,spos,nlm,el_name)
endif
else
do i = 1, jmax
!hbu
call tt_puttab(part_id(i), 0, 1, z(1,i), orbit0,spos)
enddo
endif
endif
call dcopy(orbit0,orbit,6)
doupdate = get_option('update ') .ne. 0
!--- loop over turns
nobs = 0
if(bb_sxy_update) then
trrun_nt=0
if(first) then
time_var_m_cnt=0
time_var_p_cnt=0
time_var_c_cnt=0
! <<N_macro_part_ini=N_macro_surv+N_macro_lost>>
N_ions_in_beam=get_value('probe ', 'npart ') !BEAM->NPART
if(N_ions_in_beam .LT. 0d0) call aafail('TRRUN: ','N_ions_in_beam .LE. 0.0')
Npart_gain = get_value('run ', 'n_part_gain ')
N_ions_ini=Npart_gain*N_ions_in_beam
N_macro_surv=jmax ! = number of START lines submitted
n_ions_macro=N_ions_ini/N_macro_surv
N_for_I=N_macro_surv ! at start (to be redefined in Ixy)
if(N_macro_surv .GT. N_macro_max) call aafail('TRRUN: ',&
'Number N_macro_surv exceeds N_macro_max (array size)')
if(N_macro_surv .GT. N_macro_surv) call aafail('TRRUN: ',&
'Number START-lines exceeds the initial number of macroparticles N_macro_surv')
t_rms = get_value('run ', 'sigma_z ')*beti
pt_rms = get_value('run ', 'deltap_rms ')
pt_rms=(sqrt((betas*(pt_rms+1d0))**2+1d0/gammas**2)-1d0)*beti
sigma_z_ini=t_rms !betas: BEAM->BETA
sigma_z=sigma_z_ini !at start (to be redefined in Ixy)
sigma_p=pt_rms !default
z_factor=1d0 !at start sigma_z_ini/sigma_z
Ex_rms=get_value('probe ', 'ex ') !BEAM->Ex
Ey_rms=get_value('probe ', 'ey ') !BEAM->Ey
if(checkpnt_restart.and.emittance_update) then
Ex_rms=Ex_rms0
Ey_rms=Ey_rms0
endif
first=.false.
endif
endif
do turn = 1, turns
if (doupdate) call trupdate(turn)
j = restart_sequ()
if(bb_sxy_update) then
trrun_nt=turn
time_var_m_lnt=0
time_var_p_lnt=0
time_var_c_lnt=0
N_macro_surv=jmax
i_spch = 0 !a special spch-update counter
ex_rms0=ex_rms
ey_rms0=ey_rms
sigma_z0=sigma_z
sigma_p0=sigma_p
!sigma_p0 = sigma_p !CM, 3/11/14
!fill, table=Ixy_unsorted; column=i_macro_part, Ix, Iy, dpi, z_part;
!new on 3/31/14:
!frs on 04.06.2016 - fixing
!a) bug concerning sigma_p
!b) Filling data in file bb6d_ixy.txt even for "emittance_update = !.false.",
! obviously without update!
!c) Fixing checkpnt_restart for "emittance_update = !.false." which
! worked for ".true." alright.
! if (emittance_update) then
if(sc_3d_kick) then
call SigmaMatrixFit6Dv7(z,jmax,sgpr,sc3d_emit)
if(sc_3d_beamsize) then
call SigmaTransport2(sgpr)
else
call SigmaTransport(sgpr)
endif
call double_to_table_curr('bb6d_ixy ', 'turn ', dble(tot_turn+turn))
call double_to_table_curr('bb6d_ixy ', 'n_macro_surv ', dble(N_ini))
call double_to_table_curr('bb6d_ixy ', 'n_for_i ', dble(jmax))
call double_to_table_curr('bb6d_ixy ', 'ex_rms ', sc3d_emit(1))
call double_to_table_curr('bb6d_ixy ', 'ey_rms ', sc3d_emit(2))
call double_to_table_curr('bb6d_ixy ', 'ez_rms ', sc3d_emit(3))
call double_to_table_curr('bb6d_ixy ', 'sigma_ct ', sc3d_emit(4))
call double_to_table_curr('bb6d_ixy ', 'sigma_pt ', sc3d_emit(5))
else
call ixy_calcs(betas, orbit0, z, &
betx_start, bety_start, &
alfx_start, alfy_start, &
gamx_start, gamy_start, &
dx_start, dpx_start, &
dy_start, dpy_start)
call ixy_fitting()
call double_to_table_curr('bb6d_ixy ', 'turn ', dble(tot_turn+turn))
call double_to_table_curr('bb6d_ixy ', 'n_macro_surv ', dble(n_macro_surv))
call double_to_table_curr('bb6d_ixy ', 'n_for_i ', dble(n_for_i))
call double_to_table_curr('bb6d_ixy ', 'ex_rms ', ex_rms)
call double_to_table_curr('bb6d_ixy ', 'ey_rms ', ey_rms)
call double_to_table_curr('bb6d_ixy ', 'ez_rms ', 0d0)
call double_to_table_curr('bb6d_ixy ', 'sigma_ct ', sigma_z0)
call double_to_table_curr('bb6d_ixy ', 'sigma_pt ', sigma_p0)
endif
call augment_count('bb6d_ixy ')
if(sigma_p0.eq.0d0) sigma_p0=sigma_p
!frs on 04.06.2016 - fixing
!a) bug concerning sigma_p
!b) Filling data in file bb6d_ixy.txt even for "emittance_update = !.false.",
! obviously without update!
!c) Fixing checkpnt_restart for "emittance_update = !.false." which
! worked for ".true." alright.
!new on 3/31/14:
! endif
sigma_z=sigma_z0
sigma_p=sigma_p0
!frs on 04.06.2016 - fixing
!a) bug concerning sigma_p
!b) Filling data in file bb6d_ixy.txt even for "emittance_update = !.false.",
! obviously without update!
!c) Fixing checkpnt_restart for "emittance_update = !.false." which
! worked for ".true." alright.
! sigma_p=sigma_p0
if((sigma_z.GT.0d0).AND.(sigma_z_ini.GT.0d0)) then
z_factor=sigma_z_ini/sigma_z
else
z_factor=1D0
endif
N_ions_for_bb=n_ions_macro*N_for_I*z_factor
if(N_ions_in_beam.le.0d0) then
rat_bb_n_ions=0d0
else
rat_bb_n_ions=N_ions_for_bb/N_ions_in_beam
endif
if(.not.emittance_update) then
ex_rms=ex_rms0
ey_rms=ey_rms0
rat_bb_n_ions=1d0
endif
if(idnint(time_var_m_nt(time_var_m_cnt+1)).eq.tot_turn+turn) then
time_var_m=.true.
else
time_var_m=.false.
endif
if(idnint(time_var_p_nt(time_var_p_cnt+1)).eq.tot_turn+turn) then
time_var_p=.true.
else
time_var_p=.false.
endif
if(idnint(time_var_c_nt(time_var_c_cnt+1)).eq.tot_turn+turn) then
time_var_c=.true.
else
time_var_c=.false.
endif
endif
nlm = 0
sum=0d0
!frs on 04.06.2016 - fixing
!a) bug concerning sigma_p
!b) Filling data in file bb6d_ixy.txt even for "emittance_update = !.false.",
! obviously without update!
!c) Fixing checkpnt_restart for "emittance_update = !.false." which
! worked for ".true." alright.
if(bb_sxy_update) then
if(emittance_update.or.(.not.emittance_update.and.mean_t.eq.0d0.and.sigma_t.eq.0d0)) then
! if((sigma_t.eq.0d0.or.emittance_update).and.bb_sxy_update) then
!VVK 20100321 -------- Find RMS-value of t ----------------------
! if we do 1-turn tracking, orbit0(5)=0 always
Summ_t_mean=0d0
do i_part=1, jmax
if (abs(z(5,i_part)) .GE. 0d0) then
Summ_t_mean=Summ_t_mean+z(5,i_part)
else
Print *, 'NaN z(5,i) ? :', i_part, z(5,i_part)
endif
enddo
mean_t=Summ_t_mean/dble(jmax)
Summ_t_square=0d0
do i_part=1, jmax
if (abs(z(5,i_part)) .GE. 0d0) &
Summ_t_square=Summ_t_square+(z(5,i_part)-mean_t)**2
enddo
sigma_t=sqrt(Summ_t_square/dble(jmax))
if(abs(sigma_t).eq.0d0) then
sigma_t=t_max/2d0
call aawarn('TTRUN Frozen SC: sigma_t = zero: ','sigma_t set to L/track_harmon/betas/2')
endif
endif
sigma_t=sigma_z
!-----------------------------------------------------------------
endif
!--- loop over nodes
10 continue
bbd_pos=j
if (turn .eq. 1) then
code = node_value('mad8_type ')
if(code.eq.39) code=15
if(code.eq.38) code=24
el = node_value('l ')
code_buf(nlm+1) = code
l_buf(nlm+1) = el
!hbu get current node name
call element_name(el_name,len(el_name))
if (.not.(is_drift() .or. is_thin() .or. is_quad() .or. is_dipole() .or. is_matrix())) then
print*," "
print*,"code: ",code," el: ",el," THICK ELEMENT FOUND"
sum = node_value('name ')
print*," "
print*,"Track dies nicely"
print*,"Thick lenses will get nowhere"
print*,"MAKETHIN will save you"
print*," "
print*," "
! Better to use stop. If use return, irrelevant tracking output
! is printed
stop
! call aafail('TRRUN: Fatal ', &
! '----element with length found : CONVERT STRUCTURE WITH '// &
! 'MAKETHIN')
endif
else
el = l_buf(nlm+1)
code = code_buf(nlm+1)
endif
if (switch .eq. 1) then
nobs = node_value('obs_point ')
endif
!-------- Misalignment at beginning of element (from twissfs.f)
if (code .ne. 1) then
call dzero(al_errors, align_max)
n_align = node_al_errors(al_errors)
!print *, "n_align = ", n_align
!print *, "align_max = ", align_max
if (n_align .ne. 0) then
do i = 1, jmax
call dcopy(z(1,i),zz(1),6)
call tmali1(zz(1),al_errors, betas, gammas,z(1,i), re(1,1))
enddo
endif
endif
!-------- Track through element // suppress dxt 13.12.04
call ttmap(switch,code,el,z,jmax,dxt,dyt,sum,tot_turn+turn,part_id, &
last_turn,last_pos, last_orbit,aperflag,maxaper,al_errors,onepass)
!-------- Misalignment at end of element (from twissfs.f)
!frs on 04.06.2016 - fixing
!a) bug concerning sigma_p
!b) Filling data in file bb6d_ixy.txt even for "emittance_update = !.false.",
! obviously without update!
!c) Fixing checkpnt_restart for "emittance_update = !.false." which
! worked for ".true." alright.
ex_rms0=ex_rms
ey_rms0=ey_rms
sigma_z0=sigma_z
sigma_p0=sigma_p
if (bb_sxy_update .and. is_lost) then
call ixy_calcs(betas, orbit0, z, &
betx_start, bety_start, &
alfx_start, alfy_start, &
gamx_start, gamy_start, &
dx_start, dpx_start, &
dy_start, dpy_start)
call ixy_fitting()
is_lost = .false.
endif
if( .not. emittance_update) then
ex_rms=ex_rms0
ey_rms=ey_rms0
rat_bb_n_ions=1d0
endif
sigma_z=sigma_z0
sigma_p=sigma_p0
!frs on 04.06.2016 - fixing
!a) bug concerning sigma_p
!b) Filling data in file bb6d_ixy.txt even for "emittance_update = !.false.",
! obviously without update!
!c) Fixing checkpnt_restart for "emittance_update = !.false." which
! worked for ".true." alright.
if (code .ne. 1) then
if (n_align .ne. 0) then
do i = 1, jmax
call dcopy(z(1,i),zz(1),6)
call tmali2(el,zz(1),al_errors, betas, gammas,z(1,i),re(1,1))
enddo
endif
endif
nlm = nlm+1
if (nobs .gt. 0) then
call dzero(obs_orb,6)
call get_node_vector('obs_orbit ', lobs, obs_orb)
if (lobs .lt. 6) &
call aafail('TRACK: Fatal ', 'obs. point orbit not found')
if (onetable) then
!hbu
spos=sum
!hbu get current node name
call element_name(el_name,len(el_name))
!hbu
call tt_putone(jmax, tot_turn+turn, tot_segm, segment, part_id, &
z, obs_orb,spos,nlm,el_name)
else
if (mod(turn, ffile) .eq. 0) then
do i = 1, jmax
!hbu add spos
call tt_puttab(part_id(i), turn, nobs, z(1,i), obs_orb, &
spos)
enddo
endif
endif
endif
if (advance_node().ne.0) then
j=j+1
go to 10
endif
!--- end of loop over nodes
if (switch .eq. 1) then
if (mod(turn, ffile) .eq. 0) then
if (turn .eq. turns) last_out = .true.
if (onetable) then
!hbu
spos=sum
!hbu spos added
!2013-Nov-05 15:05:33 ghislain: get current node name
call element_name(el_name,len(el_name))
call tt_putone(jmax, tot_turn+turn, tot_segm, segment, part_id, &
z, orbit0,spos,nlm,el_name)
else
do i = 1, jmax
!hbu
call tt_puttab(part_id(i), turn, 1, z(1,i), orbit0,spos)
enddo
endif
endif
else
!$OMP PARALLEL PRIVATE(i,j)
!$OMP DO
do i = 1, jmax
do j = 1, 6
coords(j,turn,i) = z(j,i) - orbit0(j)
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
endif
if (jmax .eq. 0 .or. (switch .gt. 1 .and. jmax .lt. j_tot)) &
goto 20
if (switch .eq. 2 .and. info) then
if (mod(turn,100) .eq. 0) print *, 'turn :', turn
endif
if(exit_loss_turn .and. lost_in_turn) then
lost_in_turn = .false.
write(text, '(i6)') turn
call aawarn('TRRUN Frozen SC: ',&
'Stop in Loss Turn: '//text)
goto 20
endif
enddo
!--- end of loop over turns
20 continue
if (switch .gt. 1 .and. jmax .lt. j_tot) then
e_flag = 1
return
endif
do i = 1, jmax
last_turn(part_id(i)) = min(turns+tot_turn, turn+tot_turn)
last_pos(part_id(i)) = sum
do j = 1, 6
last_orbit(j,part_id(i)) = z(j,i)
enddo
enddo
turn = min(turn, turns)
if(bb_sxy_update) then
tot_turn=tot_turn+turn
!$OMP PARALLEL PRIVATE(i,j)
!$OMP DO
do i = 1, jmax
part_id_keep(i)=part_id(i)
last_turn_keep(i)=last_turn(i)
do j=1,6
z_keep(j,i)=z(j,i)
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
endif
!--- enter last turn in tables if not done already
if (.not. last_out) then
if (switch .eq. 1) then
if (onetable) then
!hbu
spos=sum
!hbu get current node name
call element_name(el_name,len(el_name))
!hbu spos added
call tt_putone(jmax, tot_turn+turn, tot_segm, segment, part_id, &
z, orbit0,spos,nlm,el_name)
else
do i = 1, jmax
!hbu
call tt_puttab(part_id(i), turn, 1, z(1,i), orbit0,spos)
enddo
endif
endif
endif
!--- Write checkpoint_restart data
if (bb_sxy_update) then
rewind unit_chpt
write(unit_chpt) jmax
write(unit_chpt) Ex_rms
write(unit_chpt) Ey_rms
do i = 1, jmax
do j=1,6
write(unit_chpt) z(j,i)
enddo
enddo
write(unit_chpt) sigma_t
write(unit_chpt) mean_t
write(unit_chpt) N_ini
endif
!--- enter last turn in summary table
do i = 1,j_tot
tmp_d = i
call double_to_table_curr('tracksumm ', 'number ', tmp_d)
tmp_d = last_turn(i)
call double_to_table_curr('tracksumm ', 'turn ', tmp_d)
do j = 1, 6
tmp_d = last_orbit(j,i) - orbit0(j)
call double_to_table_curr('tracksumm ', vec_names(j), tmp_d)
enddo
!hbu
spos=last_pos(i)
!hbu
call double_to_table_curr('tracksumm ',vec_names(7),spos)
call augment_count('tracksumm ')
enddo
goto 101
100 call aafail('TRACK: Fatal ', 'checkpoint_restart file corrupted')
101 continue
end subroutine trrun
subroutine ttmap(switch,code,el,track,ktrack,dxt,dyt,sum,turn,part_id, &
last_turn,last_pos,last_orbit,aperflag,maxaper,al_errors, onepass)
use twtrrfi
use twiss0fi
use name_lenfi
implicit none
!----------------------------------------------------------------------*
! Purpose: *
! - Test apertures *
! - Track through thin lenses ONLY. *
! *
! Input/output: *
! switch (integer) 1: RUN, 2: DYNAP fastune *
! SUM (double) Accumulated length. *
! TRACK(6,*)(double) Track coordinates: (X, PX, Y, PY, T, PT). *
! NUMBER(*) (integer) Number of current track. *
! KTRACK (integer) number of surviving tracks. *
!----------------------------------------------------------------------*
logical aperflag,fmap,onepass
integer turn,code,ktrack,part_id(*),last_turn(*),nn,jtrk, &
get_option, optiondebug, switch, i
double precision ap1,ap2,ap3,ap4,el,sum,node_value,track(6,*), &
last_pos(*),last_orbit(6,*),parvec(26),get_value,ct,tmp, &
aperture(maxnaper),one,maxaper(6),al_errors(align_max),st, &
theta,ek(6),re(6,6),te(6,6,6),craporb(6),dxt(*),dyt(*),offset(2), &
offx,offy,min_double
character(name_len) aptype
parameter(one=1d0)
parameter(min_double = 1.e-36)
optiondebug = get_option('debug ')
fmap=.false.
call dzero(ek,6)
call dzero(craporb,6)
call dzero(re,36)
call dzero(te,216)
!---- Drift space
if(code.eq.1) then
call ttdrf(el,track,ktrack)
go to 502
endif
!---- Rotate trajectory before entry
theta = node_value('tilt ')
if (theta .ne. 0d0) then
st = sin(theta)
ct = cos(theta)
do jtrk = 1,ktrack
tmp = track(1,jtrk)
track(1,jtrk) = ct * tmp + st * track(3,jtrk)
track(3,jtrk) = ct * track(3,jtrk) - st * tmp
tmp = track(2,jtrk)
track(2,jtrk) = ct * tmp + st * track(4,jtrk)
track(4,jtrk) = ct * track(4,jtrk) - st * tmp
enddo
endif
!---- Translation of particles // AK 23.02.2006
!---- useful for combination of different transfer lines or rings
if(code.eq.36) then
if (onepass) then
call tttrans(track,ktrack)
go to 500
endif
endif
!---- Beam-beam, standard 4D, no aperture
if(code.eq.22) then
parvec(5)=get_value('probe ', 'arad ')
parvec(6)=node_value('charge ') * get_value('probe ', 'npart ')
parvec(7)=get_value('probe ','gamma ')
call ttbb(track, ktrack, parvec)
go to 500
endif
!---- 2015-Mar-19 09:07:37 ghislain:
! if(code.eq.20) then
! call aawarn('trrun: found deprecated ECOLLIMATOR element;',' should be replaced by COLLIMATOR')
! endif
! if(code.eq.21) then
! call aawarn('trrun: found deprecated RCOLLIMATOR element;',' should be replaced by COLLIMATOR')
! endif
!---- Test aperture. ALL ELEMENTS BUT DRIFTS
if(aperflag) then
nn=name_len
call node_string('apertype ',aptype,nn)
call dzero(aperture,maxnaper)
call get_node_vector('aperture ',nn,aperture)
call dzero(offset,2)
call get_node_vector('aper_offset ',nn,offset)
if (optiondebug .ne. 0) then
print *, " aperture type ", aptype
print *, " parameters ", (aperture(i),i=1,4)
print *, " offsets ", (offset(i),i=1,2)
print *, " alignment errors ", (al_errors(i),i=11,12)
print *, " maximum aperture ", (offset(i),i=1,2)
print *, " "
endif
call trcoll(aptype, aperture, offset, al_errors, maxaper, &
turn, sum, part_id, last_turn, last_pos, last_orbit, track, ktrack)
endif
!-- switch on element type BUT DRIFT, COLLIMATORS, BEAM_BEAM / 13.03.03
!-- 500 has been specified at the relevant places in go to below
!-- code =1 for drift, treated above, go to 500 directly
! print *," CODE ",code
go to ( 500, 20, 30, 40, 50, 60, 70, 80, 90, 100, &
110, 120, 130, 140, 150, 160, 170, 180, 190, 500, &
500, 500, 230, 240, 250, 260, 270, 280, 290, 300, &
310, 320, &
! 330, 500, 350, 360, 370,500,500,400,410, 500, 500, 500, 500), code
! Use this line to enable non-linear thin lens
! 330, 500, 350, 360, 370,500,500,400,410, 420, 500, 500, 500, 500), code
! Enable non-linear thin lens and RF-Multipole
330, 500, 350, 360, 370, 500, 500, 400, &
410, 420, 430, 500, 500, 500), code
!
!---- Make sure that nothing is executed if element is not known
go to 500
!
!---- Bending magnet.
!---- AL: use thick-dipole element
20 continue ! RBEND
30 continue ! SBEND
call tttdipole(track,ktrack)
go to 500
!---- Arbitrary matrix. OBSOLETE, to be kept for go to
40 continue
call tmarb(.false.,.false.,craporb,fmap,ek,re,te)
call tttrak(ek,re,track,ktrack)
go to 500
!---- Quadrupole. OBSOLETE, to be kept for go to
!---- AL: not so fast, here is the thick quadrupole
50 continue
call tttquad(track,ktrack)
go to 500
!---- Sextupole. OBSOLETE, to be kept for go to
60 continue
go to 500
!---- Octupole. OBSOLETE, to be kept for go to
70 continue
go to 500
!---- Monitors, beam instrument., MONITOR, HMONITOR, VMONITOR
170 continue
180 continue
190 continue
go to 500
!---- Multipole
80 continue
call ttmult(track,ktrack,dxt,dyt,turn)
go to 500
!---- Solenoid.
90 continue
call trsol(track, ktrack)
go to 500
!---- RF cavity.