-
Notifications
You must be signed in to change notification settings - Fork 0
/
FitHelix.C
1454 lines (1385 loc) · 53.3 KB
/
FitHelix.C
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
//
// Class definitions and tests for a kinematic 6-parameter track fit
// Original author: David Brown (LBNL), Dec. 2018
//
#include "TH1F.h"
#include "TSystem.h"
#include "THelix.h"
#include "TPolyLine3D.h"
#include "TArrow.h"
#include "TAxis3D.h"
#include "TCanvas.h"
#include "TStyle.h"
#include "TVector3.h"
#include "TPolyLine3D.h"
#include "TPolyMarker3D.h"
#include "TLegend.h"
#include "TGraph.h"
#include "TRandom3.h"
#include "TH2F.h"
#include "TF1.h"
#include "TDirectory.h"
#include "TProfile.h"
#include "TProfile2D.h"
#include "TMath.h"
#include "Math/SVector.h"
#include "Math/SMatrix.h"
//#include "Math/PxPyPzE4D.h"
#include <math.h>
#include <iostream>
#include <string>
#include <vector>
using namespace std;
using namespace ROOT::Math;
// global values: these should be obtained from services
enum effect {momfrac=0,theta1,theta2};
enum pars {_rad=0, _lambda=1, _cx=2, _cy=3, _phi0=4, _t0=5, _npars=6};
vector<string> parnames = {"Radius", "Lambda", "Cx", "Cy", "Phi0", "T0"};
vector<string> parunits = {"(mm)", "(mm)", "(mm)", "(mm)", "", "(nsec)"};
// _rad; // transverse radius
// _lambda; // longiduinal wavelength
// _cx, _cy; // cylinder center transverse coordinates
// _phi0; // phi angle at z=0 plane
// _t0; // time helix passes z=0 plane
double c=299.792;
double twopi = 6.283185;
double pi = twopi/2.0;
// set B to 1.0
double B=1.0;
// drift velocity
double vdrift(0.06); // drift velocity
double timez = 1600.0; // location of time hit
double minDOCA = 0.0; // don't use drift for DOCA below this
double rmin(350.0), rmax(700.0); // radial limits
TRandom* gRandom = new TRandom3(324913);
typedef SVector<double,6> PVec;
typedef SMatrix<double,6,6,MatRepSym<double,6> > PMat;
typedef SMatrix<double,6,1> PDer;
struct FourV {
double _x, _y, _z, _t;
FourV(): _x(0), _y(0), _z(0), _t(0) {}
FourV(double x, double y, double z, double t): _x(x), _y(y), _z(z), _t(t) {}
double dot3(FourV& other) {
return _x*other._x + _y*other._y + _z*other._z; }
double dot4(FourV& other) {
return dot3(other) - _t*other._t; }
double mag3() { return sqrt(dot3(*this)); }
double mag4() { return sqrt(dot4(*this)); }
FourV& normalize() { double norm = 1.0/mag3(); _x *= norm; _y *= norm; _z *= norm; return *this; }
FourV& operator += (FourV const& other) { _x += other._x; _y += other._y; _z += other._z; _t += other._t; return *this; }
FourV& operator -= (FourV const& other) { _x -= other._x; _y -= other._y; _z -= other._z; _t -= other._t; return *this; }
friend FourV operator + (FourV left, FourV const& right) { return left += right; }
friend FourV operator - (FourV left, FourV const& right) { return left -= right; }
// set energy given mass in MeV/c^2
void SetEnergy(double mass) { _t = sqrt(dot3(*this) + mass*mass); }
};
ostream& operator << (ostream& os, FourV const& vec) {
os << vec._x << " "
<< vec._y << " "
<< vec._z << " "
<< vec._t << endl;
return os;
}
struct HelixPars {
PVec _pars; // parameter vector
// the following is not a true parameter, but it defines the connection between time as a parametric variable
// and the geometric helix description, so must be included
double _mbar; // reduced mass in units of mm
double _mass; // in units of MeV
double pbar() const { return sqrt(_pars[_rad]*_pars[_rad] + _pars[_lambda]*_pars[_lambda] ); } // momentum in mm
double ebar() const { return sqrt(_pars[_rad]*_pars[_rad] + _pars[_lambda]*_pars[_lambda] + _mbar*_mbar); } // energy in mm
// angular rotation frequency
double omega() const { return copysign(c,_mbar)/ebar(); } // rotational velocity, sign set by magnetic force
double beta() const { return pbar()/ebar(); }
double phi(double t) const { return omega()*(t - _pars[_t0]) + _pars[_phi0]; }
double time(double zpos) const { return _pars[_t0] + zpos/(omega()*_pars[_lambda]); }
void invert() {
_mbar *= -1.0;
_pars[_t0] *= -1.0;
}
};
// drift hit along a helix; generated from a particular helix
struct DHit {
FourV _hpos; // helix position and time
FourV _wpos; // wire position and time
double _weta; // hit wire axis azimumth angle
double _rdrift; // hit true drift radius, signed by ambiguity
double _t; // hit measurement time; depends only on generation parameters
double _sigt; // hit time smearing sigma
};
// timing measurement hit
struct THit {
FourV _hpos; // helix position and time
double _t; // hit measurement time
double _sigt; // hit time smearing sigma
PVec _dtoca; // time derivatives
};
// DOCA information: this of for a particular helix and hit pair
struct HDOCA {
int _ambig; // reco ambig
double _doca; // DOCA
double _toca; // time of closest approach = propagation time + drift time
PVec _ddoca; // DOCA derivatives
PVec _dtoca; // time derivatives
};
void CrossProduct(FourV const& a, FourV const& b, FourV& c) {
c._t = a._t;
c._x = a._y*b._z - a._z*b._y;
c._y = a._z*b._x - a._x*b._z;
c._z = a._x*b._y - a._y*b._x;
}
// set the space components of the vector given the time
void HelixPos(HelixPars const& pars, FourV& pos) {
// compute rotational frequency
double omega = pars.omega();
// relative time
double dt = pos._t - pars._pars[_t0];
// compute azimuthal angle
double phi = omega*dt + pars._pars[_phi0];
// now compute position
pos._x = pars._pars[_cx] + pars._pars[_rad]*sin(phi);
pos._y = pars._pars[_cy] - pars._pars[_rad]*cos(phi);
pos._z = omega*pars._pars[_lambda]*dt;
}
void HelixVel(HelixPars const& pars, FourV& vel){
// compute rotational frequency
double omega = pars.omega();
// relative time
double dt = vel._t - pars._pars[_t0];
// compute azimuthal angle
double phi = omega*dt + pars._pars[_phi0];
// now compute position
vel._x = omega*pars._pars[_rad]*cos(phi);
vel._y = omega*pars._pars[_rad]*sin(phi);
vel._z = omega*pars._pars[_lambda];
}
void HelixMom(HelixPars const& pars, double time, FourV& mom){
// compute velocity vector
double phi = pars.phi(time);
double factor = pars._mass/pars._mbar;
mom._x = factor * pars._pars[_rad] * cos(phi);
mom._y = factor * pars._pars[_rad] * sin(phi);
mom._z = factor * pars._pars[_lambda];
mom.SetEnergy(pars._mass);
}
// convert 4-momentum positions and momentum, plus particle charge, into a helix.
// Z is assumed to be along Z
void HelixFromMom(FourV pos, FourV mom, double charge, double Bz, HelixPars& pars){
// speed of light in mm/msec
// compute some simple useful parameters
double pt = sqrt(mom._x*mom._x+mom._y*mom._y);
double phibar = atan2(mom._y,mom._x);
// translation factor from MeV/c to curvature radius in mm
double momToRad = 1000.0/(charge*Bz*c);
// mass in units of mev/c^2
pars._mass = sqrt(mom._t*mom._t - mom._x*mom._x - mom._y*mom._y - mom._z*mom._z);
// reduced mass
pars._mbar = -pars._mass*momToRad;
// transverse radius of the helix
pars._pars[_rad] = -pt*momToRad;
//tan dip
pars._pars[_lambda] = -mom._z*momToRad;
// time at z=0
pars._pars[_t0] = pos._t - pos._z/(pars.omega()*pars._pars[_lambda]);
// compute winding that miminimizes z1
double nwind = rint((pos._z/(pars._pars[_lambda]) - phibar)/twopi);
// cout << "winding number = " << nwind << endl;
// azimuth at z=0
pars._pars[_phi0] = phibar - pars.omega()*(pos._t-pars._pars[_t0]) + twopi*nwind;
// circle center
pars._pars[_cx] = pos._x + mom._y*momToRad;
pars._pars[_cy] = pos._y - mom._x*momToRad;
}
void HelixChangeDir(HelixPars const& pars, effect eff,FourV& unit) {
double phi = pars.omega()*(unit._t - pars._pars[_t0]) + pars._pars[_phi0];
double pbar = pars.pbar();
switch ( eff ) {
case theta1:
unit._x = pars._pars[_lambda]*cos(phi)/pbar;
unit._y = pars._pars[_lambda]*sin(phi)/pbar;
unit._z = pars._pars[_rad]/pbar;
break;
case theta2:
unit._x = sin(phi);
unit._y = -cos(phi);
unit._z = 0.0;
break;
case momfrac:
unit._x = pars._pars[_rad]*cos(phi)/pbar;
unit._y = pars._pars[_rad]*sin(phi)/pbar;
unit._z = -pars._pars[_lambda]/pbar;
break;
default:
cout << "Error: unknown effect!" << eff << endl;
}
}
void HelixChange(HelixPars const& pars,HelixPars& newpars, double time, double delta, effect eff) {
// initialize to no effect
newpars = pars;
// compute some useful quantities
double beta = pars.beta();
double omega = pars.omega();
double pbar = pars.pbar();
double phi = omega*(time - pars._pars[_t0]) + pars._pars[_phi0];
// cases
switch ( eff ) {
case theta1:
// polar bending: only momentum and position are unchanged
newpars._pars[_rad] += delta*pars._pars[_lambda];
newpars._pars[_lambda] += -delta*pars._pars[_rad];
newpars._pars[_t0] += delta*(time-pars._pars[_t0])/pars._pars[_lambda];
newpars._pars[_phi0] += delta*omega*(time-pars._pars[_t0])/pars._pars[_lambda];
newpars._pars[_cx] += -delta*pars._pars[_lambda]*sin(phi);
newpars._pars[_cy] += delta*pars._pars[_lambda]*cos(phi);
break;
case theta2:
// Azimuthal bending: R, Lambda, t0 are unchanged
newpars._pars[_phi0] += copysign(delta,delta*omega)*pbar/pars._pars[_rad];
newpars._pars[_cx] += -copysign(delta,delta*omega)*pbar*cos(phi);
newpars._pars[_cy] += -copysign(delta,delta*omega)*pbar*sin(phi);
break;
case momfrac:
// fractional momentum change: position and direction are unchanged
newpars._pars[_lambda] += delta*pars._pars[_lambda];
newpars._pars[_rad] += delta*pars._pars[_rad];
newpars._pars[_t0] += delta*(time-pars._pars[_t0])*(1.0-beta*beta);
newpars._pars[_phi0] += delta*omega*(time-pars._pars[_t0]);
newpars._pars[_cx] += -delta*pars._pars[_rad]*sin(phi);
newpars._pars[_cy] += +delta*pars._pars[_rad]*cos(phi);
break;
default:
cout << "Error: unknown effect!" << eff << endl;
}
}
void TestHelix( double charge, double x, double y, double z, double cost, double momphi, double momval=100.0, double mass=0.5,double tmin=-8, double tmax=8) {
// position and momentum vectors
FourV pos;
pos._x = x;
pos._y = y;
pos._z = z;
pos._t = 0.5*tmin+0.5*tmax;
FourV mom;
double sint = sqrt(1-cost*cost);
mom._x = momval*sint*cos(momphi);
mom._y = momval*sint*sin(momphi);
mom._z = momval*cost;
mom.SetEnergy(mass);
// construct the helix from these
HelixPars pars;
HelixFromMom(pos,mom,charge,B,pars);
cout << "Helix parameter:" << endl
<< "Transverse Radius = " << pars._pars[_rad] << endl
<< "Lambda = " << pars._pars[_lambda] << endl
<< "Time-0 = " << pars._pars[_t0] << endl
<< "phi-0 = " << pars._pars[_phi0] << endl
<< "Center x = " << pars._pars[_cx] << endl
<< "Center y = " << pars._pars[_cy] << endl
<< "mass = " << pars._mass << endl
<< "Reduced mass = " << pars._mbar << endl;
FourV testmom;
HelixMom(pars,pos._t,testmom);
cout << "Original mom " << mom << endl << "Helix mom " << testmom << endl;
HelixPars invpars = pars;
invpars.invert();
// test velocity
FourV vel;
vel._t = pos._t;
HelixVel(pars,vel);
double dot = testmom.dot3(vel)/c;
cout << "velocity dot mom = " << dot << endl;
double velmag = sqrt(vel._x*vel._x + vel._y*vel._y + vel._z*vel._z);
double beta = velmag/c;
double gamma = 1.0/sqrt(1.0 - beta*beta);
double energy = sqrt(momval*momval + mass*mass);
cout << "beta = " << beta << " gamma = " << gamma << " E/m = " << energy/mass << " P/E = " << momval/energy << endl;
// now graph this as a polyline over the specified time range.
double tstep = 0.1; // nanoseconds
double trange = tmax-tmin;
int nsteps = (int)rint(trange/tstep);
// create Canvase
TCanvas* hcan = new TCanvas("hcan","Helix",1000,1000);
//TPolyLine to graph the result
TPolyLine3D* hel = new TPolyLine3D(nsteps+1);
TPolyLine3D* invhel = new TPolyLine3D(nsteps+1);
FourV hpos;
for(int istep=0;istep<nsteps+1;++istep){
// compute the position from the time
hpos._t = tmin + tstep*istep;
HelixPos(pars,hpos);
// add these positions to the TPolyLine3D
hel->SetPoint(istep, hpos._x, hpos._y, hpos._z);
HelixPos(invpars,hpos);
invhel->SetPoint(istep, hpos._x, hpos._y, hpos._z);
}
// draw the helix
if(charge > 0)
hel->SetLineColor(kBlue);
else
hel->SetLineColor(kRed);
hel->Draw();
// invhel->SetLineColor(kCyan);
// invhel->Draw();
// draw the origin and axes
TAxis3D* rulers = new TAxis3D();
rulers->GetXaxis()->SetAxisColor(kBlue);
rulers->GetXaxis()->SetLabelColor(kBlue);
rulers->GetYaxis()->SetAxisColor(kCyan);
rulers->GetYaxis()->SetLabelColor(kCyan);
rulers->GetZaxis()->SetAxisColor(kOrange);
rulers->GetZaxis()->SetLabelColor(kOrange);
// TAxis3D::ToggleRulers();
// TAxis3D::ToggleRulers();
// cout <<"axis = " << rulers << endl;
rulers->Draw();
// TEveArrow* start = new TEveArrow(pos._x,pos._y,pos._z,0.0,0.0,0.0);
// start->Draw();
// TVector3* start =new TVector3(pos._x,pos._y,pos._z);
// start->Draw();
TPolyLine3D* refmom = new TPolyLine3D(2);
refmom->SetPoint(0,pos._x,pos._y,pos._z);
double mommag = sqrt(mom._x*mom._x + mom._y*mom._y + mom._z*mom._z);
double momscale = fabs(pars._pars[_rad])/mommag;
refmom->SetPoint(1,pos._x + mom._x*momscale,pos._y + mom._y*momscale ,pos._z + mom._z*momscale);
int momcolor;
if(charge>0.0)
momcolor = kRed;
else
momcolor = kBlack;
refmom->SetLineColor(momcolor);
refmom->Draw();
TPolyMarker3D* refp = new TPolyMarker3D(1,24);
refp->SetMarkerColor(momcolor);
refp->SetPoint(0,pos._x,pos._y,pos._z);
refp->Draw();
TPolyMarker3D* refmomp = new TPolyMarker3D(1,22);
refmomp->SetPoint(0,pos._x + mom._x*momscale,pos._y + mom._y*momscale ,pos._z + mom._z*momscale);
refmomp->SetMarkerColor(momcolor);
refmomp->Draw();
TPolyMarker3D* helixp = new TPolyMarker3D(1,3);
helixp->SetMarkerColor(kGreen);
hpos._t = pos._t;
HelixPos(pars,hpos);
helixp->SetPoint(0,hpos._x,hpos._y,hpos._z);
helixp->Draw();
TPolyMarker3D* testp = new TPolyMarker3D(1,25);
testp->SetMarkerColor(kOrange);
testp->SetPoint(0,hpos._x,hpos._y,hpos._z+2*M_PI*pars._pars[_lambda]);
testp->Draw();
TPolyMarker3D* startp = new TPolyMarker3D(1,21);
startp->SetMarkerColor(kBlue);
hpos._t = tmin;
HelixPos(pars,hpos);
startp->SetPoint(0,hpos._x,hpos._y,hpos._z);
startp->Draw();
TPolyMarker3D* endp = new TPolyMarker3D(1,22);
endp->SetMarkerColor(kBlue);
hpos._t = tmax;
HelixPos(pars,hpos);
endp->SetPoint(0,hpos._x,hpos._y,hpos._z);
endp->Draw();
TLegend* leg = new TLegend(0.8,0.8,1.0,1.0);
char title[80];
snprintf(title,80,"Helix, mass=%3.1g MeV/c^{2}, q=%1.1f",mass,charge);
leg->AddEntry(hel,title,"L");
snprintf(title,80,"Initial Momentum =%3.1g MeV/c",momval);
leg->AddEntry(refmom,title,"L");
leg->AddEntry(refp,"Initial Position","P");
snprintf(title,80,"Helix, t=%4.2g ns",pos._t);
leg->AddEntry(helixp,title,"P");
snprintf(title,80,"Helix, t=%4.2g ns",pos._t+tmin);
leg->AddEntry(startp,title,"P");
snprintf(title,80,"Helix, t=%4.2g ns",pos._t+tmax);
leg->AddEntry(endp,title,"P");
leg->Draw();
}
void TestHelixDerivs( double dmin, double dmax, double time, effect ieff, double charge=1.0, double momval=100.0, double momphi=1.5, double cost=0.7, double mass=0.5,double t0=0.0){
// position and momentum vectors
FourV pos;
pos._x = 0.0;
pos._y = 0.0;
pos._z = 0.0;
pos._t = t0;
FourV mom;
double sint = sqrt(1-cost*cost);
mom._x = momval*sint*cos(momphi);
mom._y = momval*sint*sin(momphi);
mom._z = momval*cost;
mom.SetEnergy(mass);
// construct original helix from these
HelixPars pars;
HelixFromMom(pos,mom,charge,B,pars);
// position at time of change
FourV testpos;
testpos._t = time;
HelixPos(pars,testpos);
FourV testmom;
HelixMom(pars,time,testmom);
cout << "initial momentum = " << mom._x << " " << mom._y << " " << mom._z << endl;
cout << "derived momentum = " << testmom._x << " " << testmom._y << " " << testmom._z << endl;
// check
HelixPars testpars;
HelixFromMom(testpos,testmom,charge,B,testpars);
cout << "Delta radius = " << testpars._pars[_rad] << " " << pars._pars[_rad] <<
" Delta lambda = " << testpars._pars[_lambda] << " " << pars._pars[_lambda] <<
" Delta t0 = " << testpars._pars[_t0] << " " << pars._pars[_t0] <<
" Delta phi0 = " << testpars._pars[_phi0] << " " << pars._pars[_phi0] <<
" Delta Cx = " << testpars._pars[_cx] << " " << pars._pars[_cx] <<
" Delta Cy = " << testpars._pars[_cy] << " " << pars._pars[_cy] << endl;
// momentum direction at time of change
double testphi = atan2(testmom._y,testmom._x);
double testcost = testmom._z/sqrt(testmom._x*testmom._x + testmom._y*testmom._y + testmom._z*testmom._z);
cout << "old phi = " << momphi << " new phi = " << testphi << endl;
cout << "old cost = " << cost << " new cost = " << testcost << endl;
// make graphs to compare each parameter
int ndel(100);
TGraph* radgraph = new TGraph(ndel+1);
radgraph->SetTitle("Radius;exact;1st derivative");
TGraph* lambdagraph = new TGraph(ndel+1);
lambdagraph->SetTitle("Lambda;exact;1st derivative");
TGraph* t0graph = new TGraph(ndel+1);
t0graph->SetTitle("t_{0};exact;1st derivative");
TGraph* phi0graph = new TGraph(ndel+1);
phi0graph->SetTitle("phi_{0};exact;1st derivative");
TGraph* cxgraph = new TGraph(ndel+1);
cxgraph->SetTitle("C_{x};exact;1st derivative");
TGraph* cygraph = new TGraph(ndel+1);
cygraph->SetTitle("C_{y};exact;1st derivative");
// graphs to compare momentum change
TGraph* mom0graph = new TGraph(ndel+1);
mom0graph->SetTitle("Momentum Direction;exact;1st derivative");
TGraph* mom1graph = new TGraph(ndel+1);
mom1graph->SetTitle("Theta Direction;exact;1st derivative");
TGraph* mom2graph = new TGraph(ndel+1);
mom2graph->SetTitle("Phi Direction;exact;1st derivative");
TGraph* gapgraph = new TGraph(ndel+1);
gapgraph->SetTitle("Gap;change;gap value (mm)");
// scan range of change
double del = (dmax-dmin)/ndel;
for(int id=0;id<=ndel+1;++id){
double delta = dmin + del*id;
// compute exact altered params
HelixPars xpars;
FourV newmom;
if(ieff == theta1){
double newcost = cos(acos(cost) + delta);
double newsint = sqrt(1.0-newcost*newcost);
newmom._x = momval*newsint*cos(testphi);
newmom._y = momval*newsint*sin(testphi);
newmom._z = momval*newcost;
newmom.SetEnergy(mass);
HelixFromMom(testpos,newmom,charge,B,xpars);
} else if(ieff == theta2){
double dphi = delta/sint;
double newphi = testphi + dphi;
newmom._x = momval*sint*cos(newphi);
newmom._y = momval*sint*sin(newphi);
newmom._z = momval*cost;
newmom.SetEnergy(mass);
HelixFromMom(testpos,newmom,charge,B,xpars);
} else if(ieff == momfrac) {
double newmomval = momval*(1.0+delta);
newmom._x = newmomval*sint*cos(testphi);
newmom._y = newmomval*sint*sin(testphi);
newmom._z = newmomval*cost;
newmom.SetEnergy(mass);
HelixFromMom(testpos,newmom,charge,B,xpars);
}
// now, also compute 1st order change in parameters
HelixPars d1pars;
HelixChange(testpars,d1pars, testpos._t, delta, ieff);
FourV d1pos;
d1pos._t = testpos._t;
HelixPos(d1pars,d1pos);
FourV gap = d1pos - testpos;
gapgraph->SetPoint(id,delta,gap.mag3());
FourV d1mom;
HelixMom(d1pars,testpos._t,d1mom);
radgraph->SetPoint(id,xpars._pars[_rad],d1pars._pars[_rad]);
lambdagraph->SetPoint(id,xpars._pars[_lambda],d1pars._pars[_lambda]);
t0graph->SetPoint(id,xpars._pars[_t0],d1pars._pars[_t0]);
phi0graph->SetPoint(id,xpars._pars[_phi0],d1pars._pars[_phi0]);
cxgraph->SetPoint(id,xpars._pars[_cx],d1pars._pars[_cx]);
cygraph->SetPoint(id,xpars._pars[_cy],d1pars._pars[_cy]);
// compare momenta after change
FourV dxmom = newmom - testmom;
FourV dd1mom = d1mom - testmom;
FourV changedir;
changedir._t = time;
HelixChangeDir(testpars,momfrac,changedir);
mom0graph->SetPoint(id,dxmom.dot3(changedir),dd1mom.dot3(changedir));
HelixChangeDir(pars, theta1, changedir );
mom1graph->SetPoint(id,dxmom.dot3(changedir),dd1mom.dot3(changedir));
HelixChangeDir(pars, theta2, changedir );
mom2graph->SetPoint(id,dxmom.dot3(changedir),dd1mom.dot3(changedir));
}
// draw comparisons
TCanvas* dhcan = new TCanvas("dhcan","Helix Change",1200,800);
dhcan->Divide(3,2);
dhcan->cd(1);
radgraph->Draw("AC*");
dhcan->cd(2);
lambdagraph->Draw("AC*");
dhcan->cd(3);
t0graph->Draw("AC*");
dhcan->cd(4);
phi0graph->Draw("AC*");
dhcan->cd(5);
cxgraph->Draw("AC*");
dhcan->cd(6);
cygraph->Draw("AC*");
TCanvas* dmomcan = new TCanvas("dmomcan","Mom Change",800,800);
dmomcan->Divide(2,2);
dmomcan->cd(1);
mom0graph->Draw("AC*");
dmomcan->cd(2);
mom1graph->Draw("AC*");
dmomcan->cd(3);
mom2graph->Draw("AC*");
dmomcan->cd(4);
gapgraph->Draw("AC*");
}
void TestHelixAccuracy(double mommin, double mommax, unsigned nbins, double mass, double xfrac,double costmin, double costmax, double xrange, unsigned nsamples, double rmax, double zmax) {
// set B to 1.0
double charge = 1.0;
// random generator
// diagnostics
double maxgap(100.0);
TProfile2D* mfgapP = new TProfile2D("mfgapP","Energy straggling gap;Momentum (MeV/c);cos(#theta);gap (mm)",nbins,mommin,mommax,nbins,costmin,costmax,0.0,maxgap);
TProfile2D* tsgapP = new TProfile2D("tsgapP","Theta scatter gap;Momentum (MeV/c);cos(#theta);gap (mm)",nbins,mommin,mommax,nbins,costmin,costmax,0.0,maxgap);
TProfile2D* psgapP = new TProfile2D("psgapP","Phi scatter gap;Momentum (MeV/c);cos(#theta);gap (mm)",nbins,mommin,mommax,nbins,costmin,costmax,0.0,maxgap);
TProfile2D* allgapP = new TProfile2D("allgapP","Total gap;Momentum (MeV/c);cos(#theta);gap (mm)",nbins,mommin,mommax,nbins,costmin,costmax,0.0,maxgap);
mfgapP->SetStats(0);
tsgapP->SetStats(0);
psgapP->SetStats(0);
allgapP->SetStats(0);
mfgapP->SetMaximum(0.1);
tsgapP->SetMaximum(2.0);
psgapP->SetMaximum(2.0);
allgapP->SetMaximum(2.0);
mfgapP->SetMinimum(0.0);
tsgapP->SetMinimum(0.0);
psgapP->SetMinimum(0.0);
allgapP->SetMinimum(0.0);
TProfile2D* effP[4] = {mfgapP,tsgapP,psgapP,allgapP};
// loop over samples
for (unsigned isamp=0;isamp<nsamples;++isamp){
double momval = gRandom->Uniform(mommin,mommax);
// initial position is random inside range (+-)
FourV pos;
pos._x = gRandom->Uniform(-xrange,xrange);
pos._y = gRandom->Uniform(-xrange,xrange);
pos._z = gRandom->Uniform(-xrange,xrange);
// initial time is arbitrary: set to 0
pos._t = 0.0;
FourV mom;
// random momentum directions
double cost = gRandom->Uniform(costmin,costmax);
double sint = sqrt(1-cost*cost);
double momphi = gRandom->Uniform(-pi,pi);
mom._x = momval*sint*cos(momphi);
mom._y = momval*sint*sin(momphi);
mom._z = momval*cost;
mom.SetEnergy(mass);
// construct original helix from these
HelixPars pars;
HelixFromMom(pos,mom,charge,B,pars);
// relativistic beta
double beta = momval/mom._t;
// time range is assumed to be symmetric accross the measurement range. Take
// the smaller of transverse or Z
double tmax(10);
double radius = pars._pars[_rad];
double drange;
if(2*radius < rmax) // spiral
drange = zmax/cost;
else
drange = 2.0*radius*asin(rmax/(2.0*radius));
double trange = drange/beta*c;
trange = trange < tmax ? trange : tmax;
// choose a random time for the effect
FourV testpos;
testpos._t = gRandom->Uniform(pars._pars[_t0]-trange,pars._pars[_t0]+trange);
HelixPos(pars,testpos);
// scattering sigma
double scatsig = 14*sqrt(xfrac)/(momval*beta);
// energy (momentum) loss fraction RMS
double momfsig = 3.0*xfrac/(momval*beta*beta);
// compute 1st order change in parameters
HelixPars allpars = pars;
// loop over scattering directions
for(int ieff=0;ieff<3;++ieff){
// choose correct sigma
double sigma = ieff == 0 ? momfsig : scatsig;
// compute random sample of this as a Gaussian (should be Moliere/landau!!, FIXME!!)
double delta = gRandom->Gaus(0.0,sigma);
HelixPars d1pars;
HelixChange(pars,d1pars, testpos._t, delta, (effect)ieff);
FourV d1pos;
d1pos._t = testpos._t;
HelixPos(d1pars,d1pos);
// compute the gap between positions
FourV gap = d1pos - testpos;
effP[ieff]->Fill(momval,cost,gap.mag3());
// accumulate the differences
HelixPars apars = allpars;
HelixChange(apars,allpars, testpos._t, delta, (effect)ieff);
}
// total
FourV d1pos;
d1pos._t = testpos._t;
HelixPos(allpars,d1pos);
// compute the gap between positions
FourV allgap = d1pos - testpos;
effP[3]->Fill(momval,cost,allgap.mag3());
}
// draw comparisons
TCanvas* hacc = new TCanvas("hacc","Helix Accuracy",1000,600);
hacc->Divide(2,2);
hacc->cd(1);
mfgapP->Draw("colorz");
hacc->cd(2);
tsgapP->Draw("colorz");
hacc->cd(3);
psgapP->Draw("colorz");
hacc->cd(4);
allgapP->Draw("colorz");
}
// change in residual for a small change in parameter.
void ParamDerivs(HelixPars const& pars, double time,double dvec[6]){
// small increments for each parameter
double dpar[6] = {0.1,0.1,0.1,0.1,0.001,0.1};
// position at the test point
FourV opos;
opos._t = time;
HelixPos(pars,opos);
// directions at test point
FourV momdir;
momdir._t = time;
FourV thetadir;
thetadir._t = time;
HelixChangeDir(pars,theta1,thetadir);
FourV phidir;
phidir._t = time;
HelixChangeDir(pars,theta2,phidir);
// now increment parameter
HelixPars cpars = pars;
for(unsigned ipar=0;ipar<5;++ipar){
switch(ipar) {
case 0:
cpars._pars[_rad] += dpar[ipar];
break;
case 1:
cpars._pars[_lambda] += dpar[ipar];
break;
case 2:
cpars._pars[_cx] += dpar[ipar];
break;
case 3:
cpars._pars[_cy] += dpar[ipar];
break;
case 4:
cpars._pars[_phi0] += dpar[ipar];
break;
case 5:
cpars._pars[_t0] += dpar[ipar];
break;
}
FourV cpos;
cpos._t = time;
HelixPos(cpars,cpos);
FourV delta = opos-cpos;
// projection perpendicular to momentum
double dphi = delta.dot3(phidir);
double dtheta = delta.dot3(thetadir);
// arbitrary sign
double del = copysign(sqrt(dphi*dphi+dtheta*dtheta),dphi);
// derivative
dvec[ipar] = del/dpar[ipar];
}
// t0 is a special case; derivative is constant, sign is arbitrary and unrelated to geometry (LR ambiguity).
dvec[5] = copysign(vdrift,cos(time/1.0e-5));
}
void ParamCovar(double mommag, double mass, double cost, unsigned nmeasure, double trange, double xfrac, double msig) {
// create helix
FourV pos(0.0,0.0,0.0,0.0);
double sint = sqrt(1.0 - cost*cost);
FourV mom(mommag*sint,0.0,mommag*cost,0.0);
mom.SetEnergy(mass);
HelixPars pars;
HelixFromMom(pos, mom, -1.0, 1.0, pars);
// scattering per measrurement
double xfmeas = xfrac/nmeasure;
// scattering sigma per measurement
double beta = pars.beta();
double scatsig = 14.0*sqrt(xfmeas)/(mommag*beta);
// measurement weight
double wm = 1.0/msig;
// time step per measurement; time assumed symmetric
double tstep = 2*trange/nmeasure;
// initialize covariance
// loop over measurements
for(unsigned imeas = 0;imeas<nmeasure;++imeas){
double time = imeas*tstep-trange;
// compute weight for this measurement
// scattering derivatives
HelixPars dt1pars;
HelixChange(pars,dt1pars,time,scatsig,theta1);
HelixPars dt2pars;
HelixChange(pars,dt2pars,time,scatsig,theta2);
// update covariance
}
}
void ParamCovarTest(double minmom, double maxmom, double mass, double mincost, double maxcost, unsigned nbins, unsigned nmeasure,double trange, double xfrac, double msig) {
double momstep = (maxmom-minmom)/(nbins-1);
double coststep = (maxcost-mincost)/(nbins-1);
// loop over momenta
for(unsigned mombin=0;mombin<nbins;++mombin){
double momval = minmom + mombin*momstep;
for(unsigned costbin=0;costbin<nbins;++costbin){
double cost = mincost + costbin*coststep;
// compute covariance for these values
ParamCovar(momval,mass,cost,nmeasure,trange,xfrac,msig);
}
}
}
void ComputeDOCA(DHit const& hit, HelixPars const& hpars, HDOCA& hdoca) {
PVec const& pars = hpars._pars;
// compute helix azimuth angle from the hit z position
double phi = hit._wpos._z/pars[_lambda] + pars[_phi0];
// denomiator
// angle difference
double phibar = phi-hit._weta;
double sinphibar = sin(phibar);
double cosphibar = cos(phibar);
double sineta = sin(hit._weta);
double coseta = cos(hit._weta);
double l2 = pars[_lambda]*pars[_lambda];
double r2 = pars[_rad]*pars[_rad];
double s2 = sinphibar*sinphibar;
double m2 = hpars._mbar*hpars._mbar;
double denom = l2 + s2*r2;
double Factor = pars[_lambda]/sqrt(denom);
double dx = pars[_cx] - hit._wpos._x;
double dy = pars[_cy] - hit._wpos._y;
double ddot = -sineta*dx + coseta*dy;
// doca; this approximation ignores 2nd order terms. It becomes
// innacurate when tandip is small and z large
double doca = -Factor*(pars[_rad]*cosphibar - ddot );
// cout << " doca " << " Factor " << Factor << " ddot " << ddot << " r*cosphi " << pars[_rad]*cosphibar << endl;
hdoca._doca = doca;
// propagation time
double tprop = hpars.time(hit._wpos._z);
// time estimate for this hit. This is equivalent to DOCA in the time domain
hdoca._toca = tprop + fabs(doca)/vdrift;
// cout << "wire t = " << hit._wpos._t << " hit t " << doca._toca << " tprop " << tprop << " htime " << hit._hpos._t << " hit z " << hit._wpos._z << endl;
double ambig = doca > 0.0 ? 1.0 : -1.0;
if(fabs(doca) > minDOCA)
hdoca._ambig = ambig;
else
hdoca._ambig = 0.0;
// cout << "tresid = " << dt << " tpull = " << dt/hit._sigt << endl;
// derivatives of DOCA
PVec ddoca;
ddoca[_t0] = 0.0; // no t0 dependence, purely geometric
ddoca[_cx] = -Factor*sineta;
ddoca[_cy] = Factor*coseta;
// components that depend on phi; this ignores terms of magnitude DOCA/Radius
ddoca[_phi0] = Factor*pars[_rad]*sinphibar;
ddoca[_rad] = -Factor*(l2*cosphibar + pars[_rad]*s2*ddot)/denom;
// this becomes inaccurate when z=0, but should have no effect on the fit
ddoca[_lambda] = -Factor*pars[_rad]*sinphibar*hit._wpos._z/l2;
// the extra terms are negligible, relative order DOCA/R and can be ignored
//ddoca[_lambda] = Factor*sinphibar*hit._wpos._z/l2 + sinphibar*(sinphibar + cosphibar*hit._wpos._z/(pars[_rad]*pars[_lambda]))*doca/(pars[_lambda]*denom);
hdoca._ddoca = ddoca;
// drift time residual. This is signed by the ambiguity. Use the true ambiguity
// here, but in a real experiment it would need to be determined from data
PVec dtdrift = ambig*(1.0/vdrift)*ddoca;
// propagation time residual
PVec dtprop;
// add explicit dependence on parameters due to propagation
dtprop[_t0] = 1.0; // direct dependence on t0
// double f2 = -hit._wpos._z/(c*pars[_lambda]*ebar);
// the following terms are negligible and contribute and can be ignored
// dtprop[_lambda] = f2*(pars[_rad]*pars[_rad] + mbar*mbar)/(c*pars[_lambda]);
// dtprop[_rad] = f2*mbar*mbar/pars[_rad];
hdoca._dtoca = dtprop + dtdrift;
// cout << " dtprop[tlambda] " << dtprop[_lambda] << " dtdrift[tlambda] " << dtdrift[_lambda]
// << " dtprop[rad] " << dtprop[_rad] << " dtdrift[rad] " << dtdrift[_rad] << endl;
}
void GenerateTimeHit(HelixPars const& pars, THit& thit, double tz) {
// compute the time for this hit location
thit._hpos._t = pars.time(tz);
// find helix position and momentum at this z
HelixPos(pars,thit._hpos);
// smear the time
thit._t = gRandom->Gaus(thit._hpos._t,thit._sigt);
// derivatives are trivial
thit._dtoca[_t0] = 1.0;
}
void GenerateDriftHits(HelixPars const& pars, std::vector<DHit>& hits, unsigned nz, double deltaz, double sigdrift=0.2, double rmax=2.5) {
hits.clear();
double pbar = pars.pbar();
double ebar = pars.ebar();
double beta = pbar/ebar;
double omega = pars.omega();
// cout << "omega " << omega << "omega * lambda " << omega * pars._pars[_lambda] << endl;
// divide z range around 0
for(int iz=-nz; iz<=(int)nz; ++iz){
FourV hpos;
double z = iz*deltaz;
// compute the time from this position
hpos._t = pars.time(z);
// find helix position and momentum at this z
HelixPos(pars,hpos);
double rho = sqrt(hpos._x*hpos._x+hpos._y*hpos._y);
// if(rho > rmin && rho < rmax) {
FourV mom;
HelixMom(pars, hpos._t, mom);
// generate a random azimuth angle for the wire
double eta = gRandom->Uniform(-pi,pi);
// define wire direction
FourV wdir(cos(eta),sin(eta),0.0,0.0);
// define perp direction to wire and track
FourV udir;
CrossProduct(wdir,mom,udir);
udir.normalize();
// generate random drift; this is along the U direction
double rdrift = gRandom->Uniform(-rmax,rmax); // sign defines ambiguity
// define wire position
FourV wpos = hpos;
wpos._x -= rdrift*udir._x;
wpos._y -= rdrift*udir._y;
wpos._z -= rdrift*udir._z;
double tdrift = fabs(rdrift/vdrift);
wpos._t = hpos._t + tdrift;
// create hit
DHit hit;
hit._hpos = hpos;
hit._wpos = wpos;
hit._weta = eta;
hit._rdrift = rdrift;
hit._sigt = sigdrift/vdrift;
// generate random time smearing
hit._t = gRandom->Gaus(wpos._t,hit._sigt);
hits.push_back(hit);
// }
}
}
void PrintHit(DHit const& hit ) {
cout << "Wire x " << hit._wpos._x << " Helix x " << hit._hpos._x
<< " Wire y " << hit._wpos._y << " Helix y " << hit._hpos._y
<< " Wire z " << hit._wpos._z << " Helix z " << hit._hpos._z
<< " hit t " << hit._t << " helix t " << hit._hpos._t
<< " eta " << hit._weta << " rdrift " << hit._rdrift << endl;
}
void PrintDOCA(HDOCA const& doca ) {
cout << "DOCA " << doca._doca << " TOCA " << doca._toca << endl;
}
void PrintHits(std::vector<DHit>const& hits) {
for(auto hit: hits) {
PrintHit(hit);
}
}
void TestHelixFitDerivs(double cost=0.7, double momval=105.0,double charge=1.0,double mass=0.5,double deltaz=75.0) {
double zrange = 3000.0;
double momphi = pi/4.0;
// parameter change for derivatives
vector<double> delpars { 0.1, 0.025, 0.1, 0.1, 0.001, 0.1}; // small parameter changes for derivative calcs
// position and momentum vectors
FourV pos;
FourV mom;
double sint = sqrt(1-cost*cost);
mom._x = momval*sint*cos(momphi);
mom._y = momval*sint*sin(momphi);
mom._z = momval*cost;
mom.SetEnergy(mass);
// construct the helix from these
HelixPars pars;
HelixFromMom(pos,mom,charge,B,pars);
pars._pars[_t0] = 2.8173; // random
cout << "helix pars " << pars._pars << " mass " << pars._mass << " mbar " << pars._mbar << endl;
// generate hits
std::vector<DHit> hits;
double pbar = pars.pbar();
double ebar = pars.ebar();
double beta = pbar/ebar;
unsigned nz = ceil(zrange/(2*deltaz));
GenerateDriftHits(pars,hits,nz,deltaz);
// print them
// PrintHits(hits);
// test derivatives; first setup parameter changes
// positive change in helix
vector<HelixPars> dpars(6,pars);
for(size_t ipar=0;ipar < _npars;ipar++)
dpars[ipar]._pars[ipar] += delpars[ipar];
// histogram differences
vector<TH1F*> dderivh(_npars);
vector<TGraph*> dderivg(_npars);
vector<TH1F*> dtderivh(_npars);
vector<TGraph*> dtderivg(_npars);
string htitle, gtitle, hname, gname;
for(int ipar=0;ipar<_npars;++ipar){
hname = string("dderiv_")+string(parnames[ipar]) + string("_H");
htitle = string("DOCA Change - Derivative*#Delta ")+string(parnames[ipar])+string(";#Delta DOCA (mm)");
gname = string("dderiv_")+string(parnames[ipar]) + string("_G");
gtitle = string("Derivative Change vs DOCA change: ") + string(parnames[ipar]) + string(";#Delta DOCA (mm);#Delta Dervi*#Delta")+string(parnames[ipar]) +string("(mm)");
dderivh[ipar] = new TH1F(hname.c_str(),htitle.c_str(),100,-0.01,0.01);
dderivg[ipar] = new TGraph(hits.size());
dderivg[ipar]->SetTitle(gtitle.c_str());
// time deriv
hname = string("dtderiv_")+string(parnames[ipar]) + string("_H");
htitle = string("Time Change - Time Derivative*#Delta ")+string(parnames[ipar])+string(";#Delta Time (ns)");
gname = string("dtderiv_")+string(parnames[ipar]) + string("_G");
gtitle = string("Time Derivative Change vs Time change: ") + string(parnames[ipar]) + string(";#Delta Time (ns);#Delta TDervi*#Delta")+string(parnames[ipar]) +string("(ns)");
dtderivh[ipar] = new TH1F(hname.c_str(),htitle.c_str(),100,-0.1,0.1);
dtderivg[ipar] = new TGraph(hits.size());
dtderivg[ipar]->SetTitle(gtitle.c_str());
}
// time residuals, etc
TH1F* tresid = new TH1F("tresid","Time Residuals;#Delta t(ns)",100,-25.0,25.0);
TH1F* tpull = new TH1F("tpull","Time Pulls",100,-10.0,10.0);
TH1F* tdrift = new TH1F("tdrift","Drift Time;ns",100,-1.0,50.0);
TH1F* tres = new TH1F("tres","Time Resolution;ns",100,-0.2,0.2);
unsigned ihit(0);
for(auto hit : hits) {
HDOCA hdoca; // nominal doca