-
Notifications
You must be signed in to change notification settings - Fork 39
/
madx_ptc_twiss.f90
4968 lines (3803 loc) · 167 KB
/
madx_ptc_twiss.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
module madx_ptc_twiss_module
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! madx_ptc_distrib module
! Piotr K. Skowronski , Frank Schmidt (CERN)
!
! This module contains service for twiss distributions with PTC
use madx_ptc_module
use madx_ptc_intstate_module, only : getdebug
USE madx_ptc_setcavs_module
USE madx_ptc_knobs_module
USE madx_ptc_distrib_module
implicit none
save
!============================================================================================
! PUBLIC INTERFACE
public :: twiss,ptc_twiss
!============================================================================================
! PRIVATE
! datta structures
!PSk 2011.01.05 goes global to the modules so the slice tracking produces it for the summ table
type(probe_8) :: theTransferMap
type(probe_8) :: theRDTs
type(universal_taylor) :: unimap(6)
type twiss
logical(lp) nf
real(dp), dimension(3,3) :: beta,alfa,gama
real(dp), dimension(3,3) :: beta_p,alfa_p,gama_p ! derivatives of the above w.r.t delta_p
real(dp), dimension(3) :: mu
real(dp), dimension(6) :: disp
real(dp), dimension(6) :: disp_p ! derivatives of the dispersion w.r.t delta_p
real(dp), dimension(6) :: disp_p2 ! second derivatives of the dispersion w.r.t delta_p
real(dp), dimension(6) :: disp_p3 ! third order derivatives of dispersion w.r.t delta_p
real(dp), dimension(3) :: tune
real(dp), dimension(6,6) :: eigen
end type twiss
interface assignment (=)
module procedure equaltwiss
module procedure zerotwiss
module procedure normalform_normalform
end interface assignment (=)
interface alloc
module procedure alloctwiss
end interface alloc
interface kill
module procedure killtwiss
end interface kill
type(normalform) :: Normal
real(dp), private, dimension(2,ndim2) :: angp
real(dp), private, dimension(ndim2) :: dicu
integer, private, parameter :: ndd=ndim2
integer, private, dimension(4) :: iia,icoast
integer, private :: np, no
!new lattice function
real(dp), private, dimension(3) :: testold
real(dp), private, dimension(3) :: phase
real(dp), private, allocatable, dimension(:,:,:) :: savedTM
character(len=5), private, dimension(5), parameter :: str5 = (/'10000','01000','00100','00010','00001'/)
integer, private, dimension(6,6,3 ) :: Iaa ! for i=1,3 Ia(2*i-1,2*i-1,i) =1.d0; Ia(2*i,2*i,i) = 1.d0;
integer, private, allocatable :: J(:)
integer, private, dimension(6) :: j5 = (/0,0,0,0,1,0/)
integer, private, dimension(6) :: j6 = (/0,0,0,0,0,1/)
integer, private, dimension(6,6) :: fo = &
reshape( (/1,0,0,0,0,0,&
0,1,0,0,0,0,&
0,0,1,0,0,0,&
0,0,0,1,0,0,&
0,0,0,0,1,0,&
0,0,0,0,0,1 /), &
(/6,6/) )
logical :: slice_magnets, center_magnets, deltap_dependency, isRing
logical :: resetBetaExtrema;
real(dp) :: minBeta(3,3) ! jluc: to store extremas of Twiss functions (show-up in summary table
real(dp) :: maxBeta(3,3) ! jluc: to store extremas of Twiss functions (show-up in summary table)
real(dp) :: minBetX ! Edwards-Teng betas
real(dp) :: maxBetX ! Edwards-Teng betas
real(dp) :: minBetY ! Edwards-Teng betas
real(dp) :: maxBetY ! Edwards-Teng betas
real(dp) :: minDisp(4)
real(dp) :: maxDisp(4)
logical :: resetOrbitExtrema
real(dp) :: minOrbit(6)
real(dp) :: maxOrbit(6)
real(dp) :: sum2Orbit(6) ! sum of squares for rms calculatios
integer :: nobsOrbit ! counter of observation points for rms calculation
! slice tracking displays many points at the same position
! for rms only the last should be taken
real(dp) :: prevOrbit(6)
real(dp) :: prevS(6)
character(2000), private :: whymsg
character(48) :: nl_table_name='nonlin'
character(48) :: rdt_table_name='twissrdt'
!============================================================================================
! PRIVATE
! routines
private zerotwiss,equaltwiss,killtwiss
contains
!____________________________________________________________________________________________
subroutine initIaaMatrix()
implicit none
integer i
Iaa = 0;
do i=1,3
Iaa(2*i-1,2*i-1,i) =1; Iaa(2*i,2*i,i) = 1;
enddo
end subroutine initIaaMatrix
subroutine dispersion6D(A_script,disp)
implicit none
! type(probe_8), intent(in)::A_script_probe
type(real_8) ::A_script(6)
real(dp), intent(out) :: disp(4)
real(dp) :: amatrix(6,6)
real(dp) ::amatrix_inv(6,6),Ha(6,6,3)
integer i,j, inverr
! H based dispersion as in Chao-Sands paper
do i=1,6
do j=1,6
amatrix(i,j) = A_script(i).sub.fo(j,:)
enddo
enddo
inverr = 0
call matinv(amatrix, amatrix_inv,6,6,inverr) ! amap**(-1) ! and invert it and copy to 6:6 matrix
if (inverr .ne. 0) then
call fort_warn('ptc_twiss::dispersion6D ',' Can not invert A_script linear')
disp = zero;
endif
do i=1,c_%nd
Ha(1:6,1:6,i)=matmul(matmul(amatrix,Iaa(1:6,1:6,i)),amatrix_inv) ! (15b)
enddo
do i=1,4
disp(i) = Ha(i,5,3)/Ha(5,5,3)
enddo
end subroutine dispersion6D
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!version on polynomials to non-linear dispersion
subroutine dispersion6Dp(A_script,disp)
implicit none
type(real_8) ::A_script(6)
real(dp), intent(out) :: disp(4)
type(taylor) dispT(4)
type(damap) :: amap, dispMap
type(damap) :: Iaa_map
real(dp) :: amatrix(6,6)
real(dp) ::amatrix_inv(6,6),Ha(6,6,3), v
integer i,j,k
! H based dispersion as in Chao-Sands paper
call alloc(dispMap)
call alloc(dispT)
call alloc(amap)
amap = A_script ! move the transformation to type damap
print*,"++++++++++++++"
print*, "amap"
call print(amap,6)
amatrix =amap ! now we can copy to simple 6:6 matrix
amatrix_inv=amap**(-1) ! and invert it and copy to 6:6 matrixs
call alloc(Iaa_map)
do i=1,c_%nd
Ha(1:6,1:6,i)=matmul(matmul(amatrix,Iaa(1:6,1:6,i)),amatrix_inv) ! (15b)
enddo
do i=3,3 ! i=1,c_%nd
Iaa_map = zero;
do j=1,6
print*, Iaa(j,:,i)
do k=1,6
v = Iaa(j,k,i)
Iaa_map%v(j) = Iaa_map%v(j) + v*(1.0_dp.mono.fo(k,:))
enddo
enddo
print*,"++++++++++++++"
print*, "Iaa ", i
call print(Iaa_map,6)
dispMap = amap * Iaa_map * amap**(-1)
!do ii=1,6
! print*," ",Ha(1:6,ii,i)
!enddo
enddo
dispMap = dispMap * (1.0_dp / Ha(5,5,3));
! print*,"++++++++++++++"
! print*, "dispMap"
!call print(dispMap,6)
do i=1,4
disp(i) = Ha(i,5,3)/Ha(5,5,3)
dispT(i) = dispMap%v(i).par.fo(5,:)
if (getdebug() > 2) then
print*, 'Disp ', i, ' R', disp(i) , ' Map ', dispMap%v(i).sub.fo(5,:)
print*, "Taylor "
call print(dispT(i),6)
endif
enddo
call kill(amap)
end subroutine dispersion6Dp
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine equaltwiss(s1,A_script)
implicit none
type(twiss), intent(inout)::s1
type(real_8), intent(in) ::A_script(6)
real(dp) :: amatrix(6,6) ! first order A_script
integer jj,i,k, ndel
real(dp) :: lat(0:6,6,3)=0
real(dp) :: test, dph
real(dp) :: epsil=1e-12 !
integer :: J(lnv)
real(dp) :: beta(3)
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) ' check_stable ',check_stable,' c_%stable_da ',c_%stable_da,' PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('equaltwiss CHECK 0 : ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss CHECK 0 ",whymsg)
return
endif
lat = zero
ndel=0
if(c_%ndpt/=0 .and. isRing ) then
! in case of icase=56 the closed solution in long is not searched
!print*, "We are at the mode 6D + nocav"
ndel=1 !this is 6D without cavity (MADX icase=56) for a ring
endif
! calculation of alpha, beta, gamma
J=0;
!print*, "EqualTwiss nd=",c_%nd," nd2=",c_%nd2, " ndel=", ndel
beta(1) = (A_script(1).sub.'100000')**2 + (A_script(1).sub.'010000')**2
beta(2) = (A_script(3).sub.'001000')**2 + (A_script(3).sub.'000100')**2
if ( c_%nd > 2 ) beta(3) = (A_script(6).sub.'000010')**2 + (A_script(6).sub.'000001')**2
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) ' check_stable ',check_stable,' c_%stable_da ',c_%stable_da,' PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('equaltwiss CHECK 1 : ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss CHECK 1 ",whymsg)
return
endif
!print*, " Skowron 1 EqTwiss ", beta
do i=1,c_%nd2-2*ndel
do jj=i,c_%nd2-2*ndel
do k=1,c_%nd-ndel
J(2*k-1)=1
lat(i,jj,k)= (A_script(i).sub.J)*(A_script(jj).sub.J)
J(2*k-1)=0
J(2*k)=1
lat(i,jj,k)=lat(i,jj,k) + (A_script(i).sub.J)*(A_script(jj).sub.J)
lat(jj,i,k)=lat(i,jj,k)
J(2*k)=0
enddo
enddo
enddo
!print*,"BETZ=",(A_script(6).sub.'000010')**2 + (A_script(6).sub.'000001')**2
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) ' check_stable ',check_stable,' c_%stable_da ',c_%stable_da,' PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('equaltwiss CHECK 2 : ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss CHECK 2 ",whymsg)
return
endif
J=0
!here ND2=4 and delta is present nd2=6 and delta is a constant
! print*,"nv",c_%nv,"nd2",c_%nd2,"np",c_%np,"ndpt",c_%ndpt ,"=>",c_%nv-c_%nd2-c_%np
if( (c_%npara==5) .or. (c_%ndpt/=0) ) then
!when there is no cavity it gives us dispersions
do i=1,4
lat(0,i,1)=(A_script(i).sub.J5)
enddo
elseif (c_%nd2 == 6) then
call dispersion6D(A_script,lat(0,1:4,1))
else
do i=1,4
lat(0,i,1)=zero
enddo
endif
!when there is no cavity it gives us dispersions
do i=1,c_%nd2-2*ndel
s1%disp(i)=lat(0,i,1)
enddo
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) ' check_stable ',check_stable,' c_%stable_da ',c_%stable_da,' PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('equaltwiss CHECK 3 : ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss CHECK 3 ",whymsg)
return
endif
!!!!!!!!!!!!!!!!
! phase advance!
!!!!!!!!!!!!!!!!
k = 2
if(c_%nd2==6.and.c_%ndpt==0) k = 3
j=0
do i=1, k
jj=2*i -1
if (i == 3) then
TEST = ATAN2((A_script(6).SUB.fo(5,:)),(A_script(6).SUB.fo(6,:)))/TWOPI
else
TEST = ATAN2((A_script(2*i -1).SUB.fo(2*i,:)),(A_script(2*i-1).SUB.fo(2*i-1,:)))/TWOPI
endif
IF(TEST<zero.AND.abs(TEST)>EPSIL)TEST=TEST+one
DPH=TEST-TESTOLD(i)
IF(DPH<zero.AND.abs(DPH)>EPSIL) DPH=DPH+one
IF(DPH>half) DPH=DPH-one
PHASE(i)=PHASE(i)+DPH
TESTOLD(i)=TEST
enddo
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) ' check_stable ',check_stable,' c_%stable_da ',c_%stable_da,' PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('equaltwiss CHECK 4 : ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss CHECK 4 ",whymsg)
return
endif
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!print*," "
!print*, "beta 11 = ", lat(1,1,1)
!print*, "beta 22 = ", lat(3,3,2)
!print*, "beta 22 = ", lat(5,5,3)
!print*," manual 22 ", (A_script(3).sub.'0010')**2 + (A_script(3).sub.'0001')**2
do i=1,c_%nd
do k=1,c_%nd
s1%beta(k,i)= lat(2*k-1,2*k-1,i)
s1%alfa(k,i)=-lat(2*k-1,2*k ,i)
s1%gama(k,i)= lat(2*k ,2*k ,i)
enddo
enddo
!print*, " Skowron EqTwiss ", s1%beta
! --- derivatives of the Twiss parameters w.r.t delta_p
if (deltap_dependency) then
if( (c_%npara==5) .or. (c_%ndpt/=0) ) then ! condition to be checked
call computeDeltapDependency(A_script,s1)
endif
endif
! ---
!swap for longitudinal beta with gamma
if (c_%nd == 3) then
do i=1,c_%nd
test = s1%beta(3,i)
s1%beta(3,i) = s1%gama(3,i)
s1%gama(3,i) = test
enddo
endif
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) ' check_stable ',check_stable,' c_%stable_da ',c_%stable_da,' PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('equaltwiss CHECK 5 : ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss CHECK 5 ",whymsg)
return
endif
s1%mu=phase
do k=1,3
do i=1,6
s1%eigen(k*2-1,i) = A_script(k*2-1).sub.fo(i,:)
s1%eigen(k*2 ,i) = A_script(k*2 ).sub.fo(i,:)
enddo
enddo
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) 'DA got unstable: PTC msg: ',messagelost(:len_trim(messagelost))
call fort_warn('e qualtwiss: ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"equaltwiss: ",whymsg);
return
endif
end subroutine equaltwiss
subroutine alloctwiss(s1)
implicit none
type (twiss),intent(inout)::s1
s1=0
end subroutine alloctwiss
!_________________________________________________________________
subroutine killtwiss(s1)
implicit none
type (twiss),intent(inout)::s1
s1%nf=.false.
s1%beta(:,:)=zero
s1%beta_p(:,:)=zero
s1%alfa(:,:)=zero
s1%alfa_p(:,:)=zero
s1%gama(:,:)=zero
s1%gama_p(:,:)=zero
s1%mu(:)=zero
s1%disp(:)=zero
s1%disp_p(:)=zero
s1%disp_p2(:)=zero
s1%disp_p3(:)=zero
s1%tune(:)=zero
s1%eigen(:,:)=zero
end subroutine killtwiss
!_________________________________________________________________
!________________________________________________________________________________
subroutine zerotwiss(s1,i)
implicit none
type(twiss), intent(inout)::s1
integer, intent(in)::i
if(i==0) then
call liepeek(iia,icoast)
NP=iia(2)-c_%nd2
s1%nf=.false.
s1%beta(:,:)=zero
s1%beta_p(:,:)=zero
s1%alfa(:,:)=zero
s1%alfa_p(:,:)=zero
s1%gama(:,:)=zero
s1%gama_p(:,:)=zero
s1%mu(:)=zero
s1%disp(:)=zero
s1%disp_p(:)=zero
s1%disp_p2(:)=zero
s1%disp_p3(:)=zero
s1%tune(:)=zero
s1%eigen(:,:)=zero
dicu(:)=zero
angp(:,:)=zero
endif
end subroutine zerotwiss
!_________________________________________________________________
subroutine ptc_twiss(tab_name,summary_tab_name)
use twissafi
implicit none
logical(lp) :: closed_orbit,beta_flg, slice, goslice
integer :: k,i,ii
integer :: mynd2,npara,nda,icase,flag_index,why(9),my_nv,nv_min
integer :: ioptfun,iii,restart_sequ,advance_node,mf1,mf2
integer :: tab_name(*)
integer :: summary_tab_name(*)
real(dp) :: deltap0,deltap ,d_val
double precision :: get_value,suml,s
integer :: posstart, posnow
integer :: geterrorflag !C function that returns errorflag value
real(dp) :: orbit(6)=0.d0
type(probe) :: orbit_probe
type(probe_8) :: A_script_probe ! A_script == A**(-1); oneTurnMap = A_script o R o A
type(twiss) :: tw
type(fibre), POINTER :: current
type(integration_node), pointer :: nodePtr, stopNode
type(work) :: startfen !Fibre energy at the start
real(dp) :: relativisticBeta !current beta0, set by getdeltae()
real(dp) :: r,re(ndim2,ndim2),dt
logical(lp) :: initial_matrix_manual, initial_matrix_table, initial_map_manual
logical(lp) :: initial_distrib_manual, initial_ascript_manual, writetmap
logical(lp) :: maptable
logical(lp) :: ring_parameters !! forces isRing variable to true, i.e. calclulation of closed solution
logical(lp) :: doNormal !! do normal form analysis
logical(lp) :: doRDTtracking !!
logical(lp) :: isstochastic !! tempurary veriable used in switching off stochastic in closed orbit search
type(c_damap) :: AscriptInPhasor, dummyMap !! maps for RDTs calculations
type(c_vector_field) :: vectorField !! defined here to avoid every step alloc and kill
type(c_taylor) :: theRDTs !!
real(dp) :: emi(3)
logical(lp) :: isputdata ! in everystep mode (node by node) switch deciding if data are to be put in twiss table for a give node
logical(lp) :: rmatrix ! flag to mark that transfer matrix should be saved (otherwise we might not track theTransferMap)
logical(lp) :: isTMsave ! flag that TM was recorded during pre-run for closed solution search
logical(lp) :: doTMtrack ! true if rmatrix==true .and. isRing==true . do not track theTransferMap and save time
! .or. already tracked form closed solution search
logical(lp) :: usertableActive = .false. ! flag to mark that there was something requested with ptc_select
integer :: countSkipped
character(48) :: summary_table_name
character(12) :: tmfile='transfer.map'
character(48) :: charconv !routine
if(universe.le.0.or.EXCEPTION.ne.0) then
call fort_warn('return from ptc_twiss: ',' no universe created')
call seterrorflag(1,"ptc_twiss ","no universe created till now");
return
endif
if(index_mad.le.0.or.EXCEPTION.ne.0) then
call fort_warn('return from ptc_twiss: ',' no layout created')
call seterrorflag(2,"ptc_twiss ","no layout created till now");
return
endif
if (getdebug() > 1) then
if(.not.associated(MY_RING%t)) then
print*,"ptc_twiss: NODE LAYOUT ALREADY CREATED"
else
print*,"ptc_twiss: NODE LAYOUT NOT YET CREATED"
endif
endif
call resetBetaExtremas()
call initIaaMatrix()
!skipnormalform = my_false
countSkipped = 0
!all zeroing
testold = zero
phase = zero
do i=1,6
unimap(i) = zero
enddo
if (getdebug() > 1) then
print*,"ptc_twiss"
endif
call momfirstinit()
!------------------------------------------------------------------------------
table_name = charconv(tab_name)
summary_table_name = charconv(summary_tab_name)
if (getdebug() > 1) then
print*,"ptc_twiss: Table name is ",table_name
print*,"ptc_twiss: Summary table name is ", summary_table_name
endif
call cleartables() !defined in madx_ptc_knobs
call kill_para(my_ring) !removes all the previous parameters
nda = getnknobsall() !defined in madx_ptc_knobs
suml=zero
rmatrix = get_value('ptc_twiss ','rmatrix ') .ne. 0
if ( (getnknobsall() + getnpushes()) < 1) then
usertableActive = .false.
else
usertableActive = .true.
rmatrix = .true. ! for the time being force if ptc_select or ptc_knob was defined
endif
isTMsave = .false.
no = get_value('ptc_twiss ','no ')
if ( no .lt. 1 ) then
call fort_warn('madx_ptc_twiss.f90 <ptc_twiss>:','Order in twiss is smaller then 1')
print*, "Order is ", no
call tidy()
return
endif
icase = get_value('ptc_twiss ','icase ')
deltap0 = get_value('ptc_twiss ','deltap ')
deltap = zero
call my_state(icase,deltap,deltap0)
if (getdebug() > 2) then
print *, "ptc_twiss: internal state after my_state:"
call print(default,6)
endif
CALL UPDATE_STATES
! check that deltap-dependency selected if icase=5, which stands for
! 4 dimensions and deltap/p as an external parameter.
! indeed the simplified formulas we applied for derivation w.r.t deltap assume
! that deltap is an externally supplied constant parameter, which is no longer
! the case when icase = 6 as deltap then becomes a phase-space state-variable.
! and for icase=4, there is no dispersion.
deltap_dependency = get_value('ptc_twiss ','deltap_dependency ') .ne. 0
if (deltap_dependency) then
if (.not.((icase.eq.5) .or. (icase.eq.56))) then
call fort_warn('ptc_twiss: ','derivation w.r.t deltap assume deltap is fixed parameter')
call fort_warn('ptc_twiss: ','derivation w.r.t deltap assume icase=5 (formula differentiation) or icase=56 (from map)')
endif
if (no < 2) then
call fort_warn("ptc_twiss ","For dispersion's 1st order derivatives w.r.t delta-p there must be no>=2 ")
endif
if (no < 3) then
call fort_warn("ptc_twiss ","For dispersion's 2nd order derivatives w.r.t delta-p there must be no>=3 ")
endif
if (no < 4) then
call fort_warn("ptc_twiss ","For dispersion's 3rd order derivatives w.r.t delta-p there must be no>=4 ")
endif
endif
!############################################################################
!############################################################################
!############################################################################
slice_magnets = get_value('ptc_twiss ','slice_magnets ') .ne. 0
center_magnets = get_value('ptc_twiss ','center_magnets ') .ne. 0
slice = slice_magnets .or. center_magnets
if ( slice) then
call make_node_layout(my_ring)
call getBeamBeam()
endif
!############################################################################
!############################################################################
!############################################################################
orbit(:)=zero
! read the orbit
! if closed orbit is to be found pass it as the starting point for the searcher
orbit(1)=get_value('ptc_twiss ','x ')
orbit(2)=get_value('ptc_twiss ','px ')
orbit(3)=get_value('ptc_twiss ','y ')
orbit(4)=get_value('ptc_twiss ','py ')
orbit(6)=-get_value('ptc_twiss ','t ') ! swap of t sign
orbit(5)=orbit(5)+get_value('ptc_twiss ','pt ')
if(mytime) then
call Convert_dp_to_dt (deltap, dt)
else
dt=deltap
endif
if(icase.eq.5 .or. icase.eq.56) orbit(5)=dt
closed_orbit = get_value('ptc_twiss ','closed_orbit ') .ne. 0
if( closed_orbit .and. (icav .gt. 0) .and. (my_ring%closed .eqv. .false.)) then
call fort_warn('return from ptc_twiss: ',' Closed orbit requested on not closed layout.')
call seterrorflag(3,"ptc_twiss ","Closed orbit requested on not closed layout.")
call tidy()
return
endif
if(closed_orbit) then
if ( .not. c_%stable_da) then
call fort_warn('ptc_twiss: ','DA got unstable even before finding closed orbit')
call seterrorflag(10,"ptc_twiss ","DA got unstable even before finding closed orbit")
call aafail('ptc_twiss: ','DA got unstable even before finding closed orbit. program stops')
! return
endif
if (getdebug() > 2) then
print*, "Looking for orbit"
print*, "Init orbit ", orbit
call print(default,6)
endif
! disable stochastic for closed orbit seach
isstochastic = default%stochastic
default%stochastic = .false.
current=>my_ring%start
!global_verbose = .true.
call FIND_ORBIT_x(orbit,default,c_1d_8,fibre1=current)
!global_verbose = .false.
default%stochastic = isstochastic
if ( .not. check_stable) then
write(whymsg,*) 'DA got unstable during closed orbit search: PTC msg: ',messagelost(:len_trim(messagelost))
call fort_warn('ptc_twiss: ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"ptc_twiss ",whymsg);
call tidy()
return
endif
if (getdebug() > 1) then
CALL write_closed_orbit(icase,orbit)
endif
elseif((my_ring%closed .eqv. .true.) .and. (getdebug() > 1)) then
print*, "Closed orbit specified by the user!"
!CALL write_closed_orbit(icase,x) at this position it isn't read
endif
orbit_probe = orbit
mynd2 = 0
npara = 0
writetmap = get_value('ptc_twiss ','writetmap ') .ne. 0
!this must be before initialization of the Bertz
initial_distrib_manual = get_value('ptc_twiss ','initial_moments_manual ') .ne. 0
if (initial_distrib_manual) then
if (getdebug() > 1) then
print*,"Initializing map with initial_moments_manual=true"
endif
call readinitialdistrib()
endif
! old PTC
!call init(default,no,nda,BERZ,mynd2,npara)
!new complex PTC
call init_all(default,no,nda,BERZ,mynd2,npara,nclocks) ! need to add number of clocks
c_verbose=.false.
i_piotr(:) = 0
i_piotr(1)=1; i_piotr(2)=3; i_piotr(3)=5;
c_normal_auto=.true.;
! call init_all(default,no,nda)
! mynd2 and npara are outputs
if (getdebug() > 2) then
print *, "ptc_twiss: internal state after init:"
call print(default,6)
print*, "no,nda,BERZ,mynd2,npara"
print*, no,nda,BERZ,mynd2,npara
print*,"^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^"
print*,""
endif
if ( (cavsareset .eqv. .false.) .and. (my_ring%closed .eqv. .false.) ) then
call setcavities(my_ring,maxaccel)
if (geterrorflag() /= 0) then
call tidy()
return
endif
endif
call setknobs(my_ring)
call alloc(A_script_probe)
A_script_probe%u=my_false
A_script_probe%x=npara
A_script_probe%x=orbit
!This must be before init map
call alloc(theTransferMap)
theTransferMap%u = .false.
theTransferMap%x = npara
theTransferMap%x = orbit
!############################################################################
!############################################################################
!############################################################################
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! INIT A_script_probe that is tracked !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
call initmap(dt,slice)
if (geterrorflag() /= 0) then
!if arror occured then return
call tidy()
return
endif
if (( .not. check_stable ) .or. ( .not. c_%stable_da )) then
write(whymsg,*) 'DA got unstable during A_ initialization: The closed solution does not exist. PTC msg: ', &
messagelost(:len_trim(messagelost))
call fort_warn('ptc_twiss: ',whymsg(:len_trim(whymsg)))
call seterrorflag(10,"ptc_twiss INIT CHECK",whymsg)
call tidy()
return
endif
!############################################################################
!############################################################################
!############################################################################
doRDTtracking = get_value('ptc_twiss ','trackrdts ') .ne. 0
if (doRDTtracking) then
call alloc(theRDTs)
call alloc(vectorField)
call alloc(AscriptInPhasor)
call alloc(dummyMap)
endif
!############################################################################
!############################################################################
!############################################################################
! assume that we track the transfer map
doTMtrack = .true.
!and check now if we could skip it and save time
if ((usertableActive .eqv. .false.) ) then
!we do not need the full thing to fill usertable
! currently we save only the first order map (i.e. matrix)
! we can consider saving the map in universal taylor, if this still pays back
if ( isTMsave .and. ( (slice .eqv. .false.) .and. (rmatrix .eqv. .true.) ) ) then
!the matrix was tracked in the initmap, and slice is off (it is saved only for each element)
! in the future should do support for slicing
doTMtrack = .false.
endif
!another independent condition
if ( (isRing .eqv. .false.) .and. (rmatrix .eqv. .false.) ) then
!we do not have to do the normal form at the end (isRing false)
! and
!the user does not want transfer maps (rmatrix false)
!
doTMtrack = .false.
endif
endif
if (doTMtrack) then
!we will be tracking theTransferMap
! get clean initialization after initmap
call kill(theTransferMap)
call alloc(theTransferMap)
theTransferMap%u = .false.
theTransferMap%x = npara
theTransferMap%x = orbit
if (getdebug() > 1) then
print*, "doTMtrack=true, theTransferMap is new"
endif
elseif (getdebug() > 1) then
print*, "doTMtrack=false, theTransferMap stays as it was"
endif
!############################################################################
!############################################################################
!############################################################################
call alloc(tw)
!Y
!the initial twiss is needed to initialize propely calculation of some variables f.g. phase advance
tw = A_script_probe%x
if (geterrorflag() /= 0) then
call fort_warn('ptc_twiss: ','equaltwiss at the begining of the line returned with error')
call tidy()
return
endif
phase = zero !we have to do it after the very initial twiss params calculation above
current=>MY_RING%start
startfen = 0
! print*,'Check 5'
! call print(default,6)
! print*,'my_state'
startfen = current!setting up start energy for record
! print*,'Check 6'
! call print(default,6)
! print*,'my_state'
suml=zero
iii=restart_sequ()
print77=.false.
read77=.false.
if (getdebug() > 2) then
call kanalnummer(mf1)
open(unit=mf1,file='ptctwiss.txt')
print *, "ptc_twiss: internal state is:"
call print(default,6)
endif
call killsavedmaps() !delete all maps, if present
mapsorder = 0 !it is set at the end, so we are sure the twiss was successful
savemaps = get_value('ptc_twiss ','savemaps ') .ne. 0
if (savemaps) then
allocate(maps(MY_RING%n))
do i=1,MY_RING%n
do ii=1,6
call alloc(maps(i)%unimap(ii),0,0)
maps(i)%unimap(ii) = zero !this initializes and allocates the variables
enddo
enddo
else