SND@LHC Software
Loading...
Searching...
No Matches
main.cc
Go to the documentation of this file.
1#include <iostream>
2#include <execinfo.h>
3#include <signal.h>
4#include <stdlib.h>
5
6#include <AbsFinitePlane.h>
7#include <AbsFitterInfo.h>
8#include <AbsMeasurement.h>
9#include <AbsTrackRep.h>
10#include <ConstField.h>
11#include <DetPlane.h>
12#include <Exception.h>
13#include <FieldManager.h>
15#include <KalmanFitterInfo.h>
17#include <MeasurementOnPlane.h>
18#include <PlanarMeasurement.h>
22#include <SharedPlanePtr.h>
24#include <StateOnPlane.h>
25#include <Tools.h>
26#include <TrackCand.h>
27#include <TrackCandHit.h>
28#include <Track.h>
29#include <TrackPoint.h>
30#include <WireMeasurement.h>
32
33#include <MaterialEffects.h>
34#include <RKTools.h>
35#include <RKTrackRep.h>
36#include <StepLimits.h>
37#include <TGeoMaterialInterface.h>
38
39#include <TApplication.h>
40#include <TCanvas.h>
41#include <TDatabasePDG.h>
42#include <TEveManager.h>
43#include <TGeoManager.h>
44#include <TH1D.h>
45#include <TRandom.h>
46#include <TStyle.h>
47#include <TVector3.h>
48#include <vector>
49
50#include <TROOT.h>
51#include <TFile.h>
52#include <TTree.h>
53#include "TDatabasePDG.h"
54#include <TMath.h>
55#include <TString.h>
56
57
58void handler(int sig) {
59 void *array[10];
60 size_t size;
61
62 // get void*'s for all entries on the stack
63 size = backtrace(array, 10);
64
65 // print out all the frames to stderr
66 fprintf(stderr, "Error: signal %d:\n", sig);
67 backtrace_symbols_fd(array, size, 2);
68 exit(1);
69}
70
71int randomPdg() {
72 int pdg;
73
74 switch(int(gRandom->Uniform(8))) {
75 case 1:
76 pdg = -11; break;
77 case 2:
78 pdg = 11; break;
79 case 3:
80 pdg = 13; break;
81 case 4:
82 pdg = -13; break;
83 case 5:
84 pdg = 211; break;
85 case 6:
86 pdg = -211; break;
87 case 7:
88 pdg = 2212; break;
89 default:
90 pdg = 211;
91 }
92
93 return pdg;
94}
95
96
98 if (gRandom->Uniform(1) > 0.5)
99 return 1;
100 return -1;
101}
102
103
104bool compareMatrices(const TMatrixTBase<double>& A, const TMatrixTBase<double>& B, double maxRelErr) {
105 bool retVal = true;
106
107 // search max abs value
108 double max(0);
109 for (int i=0; i<A.GetNrows(); ++i) {
110 for (int j=0; j<A.GetNcols(); ++j) {
111 if (fabs(A(i,j)) > max)
112 max = fabs(A(i,j));
113 }
114 }
115
116 double maxAbsErr = maxRelErr*max;
117
118 for (int i=0; i<A.GetNrows(); ++i) {
119 for (int j=0; j<A.GetNcols(); ++j) {
120 double absErr = A(i,j) - B(i,j);
121 if ( fabs(absErr) > maxAbsErr ) {
122 double relErr = A(i,j)/B(i,j) - 1;
123 if ( fabs(relErr) > maxRelErr ) {
124 std::cout << "compareMatrices: A("<<i<<","<<j<<") = " << A(i,j) << " B("<<i<<","<<j<<") = " << B(i,j) << " absErr = " << absErr << " relErr = " << relErr << "\n";
125 retVal = false;
126 }
127 }
128 }
129 }
130 return retVal;
131}
132
133bool isCovMatrix(TMatrixTBase<double>& cov) {
134
135 if (!(cov.IsSymmetric())) {
136 std::cout << "isCovMatrix: not symmetric\n";
137 return false;
138 }
139
140 for (int i=0; i<cov.GetNrows(); ++i) {
141 for (int j=0; j<cov.GetNcols(); ++j) {
142 if (isnan(cov(i,j))) {
143 std::cout << "isCovMatrix: element isnan\n";
144 return false;
145 }
146 if (i==j && cov(i,j) < 0) {
147 std::cout << "isCovMatrix: negative diagonal element\n";
148 return false;
149 }
150 }
151 }
152
153 return true;
154}
155
156
157
159
160 double epsilonLen = 1.E-10;
161 double epsilonMom = 1.E-10;
162
163 int pdg = randomPdg();
165 rep = new genfit::RKTrackRep(pdg);
166
167 //TVector3 pos(0,0,0);
168 TVector3 pos(gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1));
169 TVector3 mom(0,0.5,gRandom->Gaus(0,0.3));
170 mom.SetMag(0.5);
171 mom *= randomSign();
172
173
174 genfit::StateOnPlane state(rep);
175 rep->setPosMom(state, pos, mom);
176
177 // check if we can set another position in the same plane
178 if (randomSign() == 1) {
179 genfit::SharedPlanePtr plane = state.getPlane();
180 const TVector3& u = plane->getU();
181 const TVector3& v = plane->getV();
182
183 // random position on plane
184 pos += gRandom->Gaus() * u;
185 pos += gRandom->Gaus() * v;
186
187 // new random momentum
188 mom.SetXYZ(0,0.5,gRandom->Gaus(0,0.3));
189 mom.SetMag(0.5);
190 mom *= randomSign();
191
192 rep->setPosMom(state, pos, mom);
193
194 // check if plane has changed
195 if (state.getPlane() != plane) {
196 std::cout << "plane has changed unexpectedly! \n";
197 delete rep;
198 return false;
199 }
200 }
201
202
203 // compare
204 if ((pos - rep->getPos(state)).Mag() > epsilonLen ||
205 (mom - rep->getMom(state)).Mag() > epsilonMom) {
206
207 state.Print();
208
209 std::cout << "pos difference = " << (pos - rep->getPos(state)).Mag() << "\n";
210 std::cout << "mom difference = " << (mom - rep->getMom(state)).Mag() << "\n";
211
212 std::cout << std::endl;
213
214 delete rep;
215 return false;
216 }
217
218 delete rep;
219 return true;
220
221}
222
223
225
226 double epsilonLen = 5.E-5; // 0.5 mu
227 double epsilonMom = 1.E-4; // 100 keV
228
229 int pdg = randomPdg();
231 rep = new genfit::RKTrackRep(pdg);
232
233 //TVector3 pos(0,0,0);
234 TVector3 pos(gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1));
235 TVector3 mom(0,0.5,gRandom->Gaus(0,0.3));
236 mom.SetMag(0.5);
237 mom *= randomSign();
238
239
240 genfit::StateOnPlane state(rep);
241 rep->setPosMom(state, pos, mom);
242
243 genfit::SharedPlanePtr origPlane = state.getPlane();
244 genfit::SharedPlanePtr plane(new genfit::DetPlane(TVector3(0,randomSign()*10,0), TVector3(0,randomSign()*1,0)));
245
246 genfit::StateOnPlane origState(state);
247
248 // forth
249 double extrapLen(0);
250 try {
251 extrapLen = rep->extrapolateToPlane(state, plane);
252 }
253 catch (genfit::Exception& e) {
254 std::cerr << e.what();
255
256 delete rep;
257 return false;
258 }
259
260 // back
261 double backExtrapLen(0);
262 try {
263 backExtrapLen = rep->extrapolateToPlane(state, origPlane);
264 }
265 catch (genfit::Exception& e) {
266 std::cerr << e.what();
267
268 delete rep;
269 return false;
270 }
271
272 // compare
273 if ((rep->getPos(origState) - rep->getPos(state)).Mag() > epsilonLen ||
274 (rep->getMom(origState) - rep->getMom(state)).Mag() > epsilonMom ||
275 fabs(extrapLen + backExtrapLen) > epsilonLen) {
276
277 origState.Print();
278 state.Print();
279
280 std::cout << "pos difference = " << (rep->getPos(origState) - rep->getPos(state)).Mag() << "\n";
281 std::cout << "mom difference = " << (rep->getMom(origState) - rep->getMom(state)).Mag() << "\n";
282 std::cout << "len difference = " << extrapLen + backExtrapLen << "\n";
283
284 std::cout << std::endl;
285
286 delete rep;
287 return false;
288 }
289
290 delete rep;
291 return true;
292
293}
294
295
297
298 bool retVal(true);
299
300 bool fx( randomSign() > 0);
302
303 double deltaJac = 0.005; // relative
304 double deltaNoise = 0.00001;
305 double deltaState = 3.E-6; // absolute
306
307 if (fx) {
308 deltaJac = 0.1; // relative
309 deltaNoise = 0.1;
310 deltaState = 5.E-4; // absolute
311 }
312
313 int pdg = randomPdg();
315 rep = new genfit::RKTrackRep(pdg);
316
317 //TVector3 pos(0,0,0);
318 //TVector3 mom(0,1,2);
319 TVector3 pos(gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1));
320 TVector3 mom(0, 0.5, gRandom->Gaus(0, 1));
321 mom *= randomSign();
322 mom.SetMag(gRandom->Uniform(2)+0.3);
323 //mom.SetMag(3);
324
325 TMatrixD jac_f, jac_fi, jac_b, jac_bi;
326 TMatrixDSym noise_f, noise_fi, noise_b, noise_bi;
327 TVectorD c_f, c_fi, c_b, c_bi;
328 TVectorD state_man, stateOrig_man;
329
330 // original state and plane
332 rep->setPosMom(state, pos, mom);
333
334 static const double smear = 0.2;
335 TVector3 normal(mom);
336 normal.SetMag(1);
337 normal.SetXYZ(gRandom->Gaus(normal.X(), smear),
338 gRandom->Gaus(normal.Y(), smear),
339 gRandom->Gaus(normal.Z(), smear));
340 genfit::DetPlane* origPlanePtr = new genfit::DetPlane (pos, normal);
341 //genfit::DetPlane* origPlanePtr = new genfit::DetPlane (pos, TVector3(1,0,0), TVector3(0,0,1));
342 double rotAngleOrig = gRandom->Uniform(2.*TMath::Pi());
343 origPlanePtr->rotate(rotAngleOrig);
344 genfit::SharedPlanePtr origPlane(origPlanePtr);
345 //genfit::SharedPlanePtr origPlane = state.getPlane();
346 rep->extrapolateToPlane(state, origPlane);
347
348 const genfit::StateOnPlane origState(state);
349
350
351 // dest plane
352 normal = mom;
353 normal.SetMag(1);
354 normal.SetXYZ(gRandom->Gaus(normal.X(), smear),
355 gRandom->Gaus(normal.Y(), smear),
356 gRandom->Gaus(normal.Z(), smear));
357 TVector3 dest(mom);
358 dest.SetMag(10);
359 genfit::DetPlane* planePtr = new genfit::DetPlane (dest, normal);
360 //genfit::DetPlane* planePtr = new genfit::DetPlane (dest, TVector3(1,0,0), TVector3(0,0,1));
361 double rotAngle = gRandom->Uniform(2.*TMath::Pi());
362 planePtr->rotate(rotAngle);
363 genfit::SharedPlanePtr plane(planePtr);
364
365 /* genfit::DetPlane* planePtr = new genfit::DetPlane (*origPlane);
366 planePtr->setO(TVector3(0,randomSign()*10,0));
367 //planePtr->rotate(rotAngle);
368 genfit::SharedPlanePtr plane(planePtr);
369*/
370
371 // numerical calculation
372 TMatrixD jac_f_num;
373 rep->calcJacobianNumerically(origState, plane, jac_f_num);
374
375
376 // forth
377 genfit::StateOnPlane extrapolatedState;
378 try {
379 //std::cout << "DO FORTH EXTRAPOLATION \n";
380 rep->extrapolateToPlane(state, plane);
381 //std::cout << "GET INFO FOR FORTH EXTRAPOLATION \n";
382 extrapolatedState = state;
383 rep->getForwardJacobianAndNoise(jac_f, noise_f, c_f);
384 rep->getBackwardJacobianAndNoise(jac_fi, noise_fi, c_fi);
385 }
386 catch (genfit::Exception& e) {
387 std::cerr << e.what();
388
389 delete rep;
391 return false;
392 }
393
394 // back
395 try {
396 //std::cout << "DO BACK EXTRAPOLATION \n";
397 rep->extrapolateToPlane(state, origPlane);
398 //std::cout << "GET INFO FOR BACK EXTRAPOLATION \n";
399 rep->getForwardJacobianAndNoise(jac_b, noise_b, c_b);
400 rep->getBackwardJacobianAndNoise(jac_bi, noise_bi, c_bi);
401 }
402 catch (genfit::Exception& e) {
403 std::cerr << e.what();
404
405 delete rep;
407 return false;
408 }
409
410
411 // compare
412 if (!((origState.getState() - state.getState()).Abs() < deltaState) ) {
413 std::cout << "(origState.getState() - state.getState()) ";
414 (origState.getState() - state.getState()).Print();
415
416 retVal = false;
417 }
418
419 // check c
420 if (!((jac_f * origState.getState() + c_f - extrapolatedState.getState()).Abs() < deltaState) ||
421 !((jac_bi * origState.getState() + c_bi - extrapolatedState.getState()).Abs() < deltaState) ||
422 !((jac_b * extrapolatedState.getState() + c_b - origState.getState()).Abs() < deltaState) ||
423 !((jac_fi * extrapolatedState.getState() + c_fi - origState.getState()).Abs() < deltaState) ) {
424
425 std::cout << "(jac_f * origState.getState() + c_f - extrapolatedState.getState()) ";
426 (jac_f * origState.getState() + c_f - extrapolatedState.getState()).Print();
427 std::cout << "(jac_bi * origState.getState() + c_bi - extrapolatedState.getState()) ";
428 (jac_bi * origState.getState() + c_bi - extrapolatedState.getState()).Print();
429 std::cout << "(jac_b * extrapolatedState.getState() + c_b - origState.getState()) ";
430 (jac_b * extrapolatedState.getState() + c_b - origState.getState()).Print();
431 std::cout << "(jac_fi * extrapolatedState.getState() + c_fi - origState.getState()) ";
432 (jac_fi * extrapolatedState.getState() + c_fi - origState.getState()).Print();
433
434 retVal = false;
435 }
436
437 if (!isCovMatrix(state.getCov())) {
438 retVal = false;
439 }
440
441
442 // compare
443 if (!compareMatrices(jac_f, jac_bi, deltaJac)) {
444 std::cout << "jac_f = "; jac_f.Print();
445 std::cout << "jac_bi = "; jac_bi.Print();
446 std::cout << std::endl;
447
448 retVal = false;
449 }
450
451 // compare
452 if (!compareMatrices(jac_b, jac_fi, deltaJac)) {
453 std::cout << "jac_b = "; jac_b.Print();
454 std::cout << "jac_fi = "; jac_fi.Print();
455 std::cout << std::endl;
456
457 retVal = false;
458 }
459
460 // compare
461 if (!compareMatrices(noise_f, noise_bi, deltaNoise)) {
462 std::cout << "noise_f = "; noise_f.Print();
463 std::cout << "noise_bi = "; noise_bi.Print();
464 std::cout << std::endl;
465
466 retVal = false;
467 }
468
469 // compare
470 if (!compareMatrices(noise_b, noise_fi, deltaNoise)) {
471 std::cout << "noise_b = "; noise_b.Print();
472 std::cout << "noise_fi = "; noise_fi.Print();
473 std::cout << std::endl;
474
475 retVal = false;
476 }
477
478
479 if (!fx) {
480 // compare
481 if (!compareMatrices(jac_f, jac_f_num, deltaJac)) {
482 std::cout << "jac_f = "; jac_f.Print();
483 std::cout << "jac_f_num = "; jac_f_num.Print();
484 std::cout << std::endl;
485
486 retVal = false;
487 }
488 }
489
490 delete rep;
492
493 return retVal;
494}
495
496
498
499 double epsilonLen = 1.E-4; // 1 mu
500 //double epsilonMom = 1.E-5; // 10 keV
501
502 double matRadius(1.);
503
504 int pdg = randomPdg();
506 rep = new genfit::RKTrackRep(pdg);
507
508 //TVector3 pos(0,0,0);
509 TVector3 pos(gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1));
510 TVector3 mom(0,0.5,gRandom->Gaus(0,0.1));
511 mom *= randomSign();
512
513
514 genfit::StateOnPlane state(rep);
515 rep->setPosMom(state, pos, mom);
516
517 genfit::SharedPlanePtr origPlane = state.getPlane();
518 genfit::SharedPlanePtr plane(new genfit::DetPlane(TVector3(0,randomSign()*10,0), TVector3(0,randomSign()*1,0)));
519
520 genfit::StateOnPlane origState(state);
521
522 // forth
523 try {
524 rep->extrapolateToPlane(state, plane, true);
525 }
526 catch (genfit::Exception& e) {
527 std::cerr << e.what();
528
529 delete rep;
530 return false;
531 }
532
533
534 // compare
535 if (fabs(rep->getPos(state).Perp() - matRadius) > epsilonLen) {
536
537 origState.Print();
538 state.Print();
539
540 std::cerr << "radius difference = " << rep->getPos(state).Perp() - matRadius << "\n";
541
542 std::cerr << std::endl;
543
544 delete rep;
545 return false;
546 }
547
548 delete rep;
549 return true;
550
551}
552
553
555
556 //double epsilonLen = 1.E-4; // 1 mu
557 //double epsilonMom = 1.E-5; // 10 keV
558
559 int pdg = randomPdg();
561 rep = new genfit::RKTrackRep(pdg);
562
563 //TVector3 pos(0,0,0);
564 TVector3 pos(gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1));
565 TVector3 mom(0,0.5,gRandom->Gaus(0,0.1));
566 mom *= randomSign();
567
568
570 rep->setPosMom(state, pos, mom);
571
572 genfit::SharedPlanePtr origPlane = state.getPlane();
573 genfit::SharedPlanePtr plane(new genfit::DetPlane(TVector3(0,randomSign()*50,0), TVector3(0,randomSign()*1,0)));
574
575 genfit::MeasuredStateOnPlane origState(state);
576
577 // forth
578 try {
579 rep->extrapolateToPlane(state, plane);
580 }
581 catch (genfit::Exception& e) {
582 std::cerr << e.what();
583
584 delete rep;
585 return false;
586 }
587
588
589 // check
590 if (!isCovMatrix(state.getCov())) {
591
592 origState.Print();
593 state.Print();
594
595 delete rep;
596 return false;
597 }
598
599 delete rep;
600 return true;
601
602}
603
604
606
607 double epsilonLen = 1.E-4; // 1 mu
608 double epsilonMom = 1.E-5; // 10 keV
609
610 int pdg = randomPdg();
612 rep = new genfit::RKTrackRep(pdg);
613
614 //TVector3 pos(0,0,0);
615 TVector3 pos(0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1));
616 TVector3 mom(0,0.5,0.);
617 mom *= randomSign();
618
619
620 genfit::StateOnPlane state(rep);
621 rep->setPosMom(state, pos, mom);
622
623 genfit::SharedPlanePtr origPlane = state.getPlane();
624 genfit::StateOnPlane origState(state);
625
626 TVector3 linePoint(gRandom->Gaus(),randomSign()*10+gRandom->Gaus(),gRandom->Gaus());
627 TVector3 lineDirection(gRandom->Gaus(),gRandom->Gaus(),randomSign()*10+gRandom->Gaus());
628
629 // forth
630 try {
631 rep->extrapolateToLine(state, linePoint, lineDirection, false);
632 }
633 catch (genfit::Exception& e) {
634 std::cerr << e.what();
635
636 delete rep;
637 return false;
638 }
639
640
641 // compare
642 if (fabs(state.getPlane()->distance(linePoint)) > epsilonLen ||
643 fabs(state.getPlane()->distance(linePoint+lineDirection)) > epsilonLen ||
644 (rep->getMom(state).Unit() * state.getPlane()->getNormal()) > epsilonMom) {
645
646 origState.Print();
647 state.Print();
648
649 std::cout << "distance of linePoint to plane = " << state.getPlane()->distance(linePoint) << "\n";
650 std::cout << "distance of linePoint+lineDirection to plane = " << state.getPlane()->distance(linePoint+lineDirection) << "\n";
651 std::cout << "direction * plane normal = " << rep->getMom(state).Unit() * state.getPlane()->getNormal() << "\n";
652
653 delete rep;
654 return false;
655 }
656
657 delete rep;
658 return true;
659
660}
661
662
664
665 double epsilonLen = 1.E-4; // 1 mu
666 double epsilonMom = 1.E-5; // 10 keV
667
668 int pdg = randomPdg();
670 rep = new genfit::RKTrackRep(pdg);
671
672 //TVector3 pos(0,0,0);
673 TVector3 pos(gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1),gRandom->Gaus(0,0.1));
674 TVector3 mom(0,0.5,gRandom->Gaus(0,0.1));
675 mom *= randomSign();
676
677
678 genfit::StateOnPlane state(rep);
679 rep->setPosMom(state, pos, mom);
680
681 genfit::SharedPlanePtr origPlane = state.getPlane();
682 genfit::StateOnPlane origState(state);
683
684 TVector3 point(gRandom->Gaus(),randomSign()*10+gRandom->Gaus(),gRandom->Gaus());
685
686 // forth
687 try {
688 rep->extrapolateToPoint(state, point, false);
689 }
690 catch (genfit::Exception& e) {
691 std::cerr << e.what();
692
693 delete rep;
694 return false;
695 }
696
697
698 // compare
699 if (fabs(state.getPlane()->distance(point)) > epsilonLen ||
700 fabs((rep->getMom(state).Unit() * state.getPlane()->getNormal())) - 1 > epsilonMom) {
701
702 origState.Print();
703 state.Print();
704
705 std::cout << "distance of point to plane = " << state.getPlane()->distance(point) << "\n";
706 std::cout << "direction * plane normal = " << rep->getMom(state).Unit() * state.getPlane()->getNormal() << "\n";
707
708 delete rep;
709 return false;
710 }
711
712 delete rep;
713 return true;
714
715}
716
717
719
720 double epsilonLen = 1.E-4; // 1 mu
721 //double epsilonMom = 1.E-5; // 10 keV
722
723 int pdg = randomPdg();
725 rep = new genfit::RKTrackRep(pdg);
726
727 //TVector3 pos(0,0,0);
728 TVector3 pos(0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1));
729 TVector3 mom(0,0.5,0.);
730 mom *= randomSign();
731
732
733 genfit::StateOnPlane state(rep);
734 rep->setPosMom(state, pos, mom);
735
736 genfit::SharedPlanePtr origPlane = state.getPlane();
737 genfit::StateOnPlane origState(state);
738
739 const TVector3 linePoint(gRandom->Gaus(0,5), gRandom->Gaus(0,5), gRandom->Gaus(0,5));
740 const TVector3 lineDirection(gRandom->Gaus(),gRandom->Gaus(),2+gRandom->Gaus());
741 const double radius = gRandom->Uniform(10);
742
743 // forth
744 try {
745 rep->extrapolateToCylinder(state, radius, linePoint, lineDirection, false);
746 }
747 catch (genfit::Exception& e) {
748 std::cerr << e.what();
749
750 delete rep;
751
752 static const char* bla = "cannot solve";
753 const char* what = e.what();
754 if (strstr(what, bla))
755 return true;
756 return false;
757 }
758
759 TVector3 pocaOnLine(lineDirection);
760 double t = 1./(pocaOnLine.Mag2()) * ((rep->getPos(state)*pocaOnLine) - (linePoint*pocaOnLine));
761 pocaOnLine *= t;
762 pocaOnLine += linePoint;
763
764 TVector3 radiusVec = rep->getPos(state) - pocaOnLine;
765
766 // compare
767 if (fabs(state.getPlane()->getNormal()*radiusVec.Unit())-1 > epsilonLen ||
768 fabs(lineDirection*radiusVec) > epsilonLen ||
769 fabs(radiusVec.Mag()-radius) > epsilonLen) {
770
771 origState.Print();
772 state.Print();
773
774 std::cout << "lineDirection*radiusVec = " << lineDirection*radiusVec << "\n";
775 std::cout << "radiusVec.Mag()-radius = " << radiusVec.Mag()-radius << "\n";
776
777 delete rep;
778 return false;
779 }
780
781 delete rep;
782 return true;
783
784}
785
786
788
789 double epsilonLen = 1.E-4; // 1 mu
790 //double epsilonMom = 1.E-5; // 10 keV
791
792 int pdg = randomPdg();
794 rep = new genfit::RKTrackRep(pdg);
795
796 //TVector3 pos(0,0,0);
797 TVector3 pos(0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1));
798 TVector3 mom(0,0.5,0.);
799 mom *= randomSign();
800
801
802 genfit::StateOnPlane state(rep);
803 rep->setPosMom(state, pos, mom);
804
805 genfit::SharedPlanePtr origPlane = state.getPlane();
806 genfit::StateOnPlane origState(state);
807
808 const TVector3 centerPoint(gRandom->Gaus(0,10), gRandom->Gaus(0,10), gRandom->Gaus(0,10));
809 const double radius = gRandom->Uniform(10);
810
811 // forth
812 try {
813 rep->extrapolateToSphere(state, radius, centerPoint, false);
814 }
815 catch (genfit::Exception& e) {
816 std::cerr << e.what();
817
818 delete rep;
819
820 static const char* bla = "cannot solve";
821 const char* what = e.what();
822 if (strstr(what, bla))
823 return true;
824 return false;
825 }
826
827
828 TVector3 radiusVec = rep->getPos(state) - centerPoint;
829
830 // compare
831 if (fabs(state.getPlane()->getNormal()*radiusVec.Unit())-1 > epsilonLen ||
832 fabs(radiusVec.Mag()-radius) > epsilonLen) {
833
834 origState.Print();
835 state.Print();
836
837 std::cout << "state.getPlane()->getNormal()*radiusVec = " << state.getPlane()->getNormal()*radiusVec << "\n";
838 std::cout << "radiusVec.Mag()-radius = " << radiusVec.Mag()-radius << "\n";
839
840 delete rep;
841 return false;
842 }
843
844 delete rep;
845 return true;
846
847}
848
849
851
852 double epsilonLen = 1.E-3; // 10 mu
853
854 int pdg = randomPdg();
856 rep = new genfit::RKTrackRep(pdg);
857
858 //TVector3 pos(0,0,0);
859 TVector3 pos(0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1),0+gRandom->Gaus(0,0.1));
860 TVector3 mom(0,0.5,0.);
861 mom *= randomSign();
862
863
864 genfit::StateOnPlane state(rep);
865 rep->setPosMom(state, pos, mom);
866
867 TVector3 posOrig(state.getPos());
868
869 genfit::SharedPlanePtr origPlane = state.getPlane();
870 genfit::StateOnPlane origState(state);
871
872 double step = gRandom->Uniform(-15.,15.);
873 double extrapolatedLen(0);
874
875 // forth
876 try {
877 extrapolatedLen = rep->extrapolateBy(state, step, false);
878 }
879 catch (genfit::Exception& e) {
880 return false;
881 }
882
883 TVector3 posExt(state.getPos());
884
885
886
887
888 // compare
889 if (fabs(extrapolatedLen-step) > epsilonLen ||
890 (posOrig - posExt).Mag() > fabs(step)) {
891
892 origState.Print();
893 state.Print();
894
895 std::cout << "extrapolatedLen-step = " << extrapolatedLen-step << "\n";
896 std::cout << "started extrapolation from: "; posOrig.Print();
897 std::cout << "extrapolated to "; posExt.Print();
898 std::cout << "difference = " << (posOrig - posExt).Mag() << "; step = " << step << "\n";
899
900 delete rep;
901 return false;
902 }
903
904 delete rep;
905 return true;
906
907}
908//=====================================================================================================================
909//=====================================================================================================================
910//=====================================================================================================================
911//=====================================================================================================================
912
913
914int main() {
915
916 const double BField = 15.; // kGauss
917 //const bool debug = true;
918
919 gRandom->SetSeed(10);
920 signal(SIGSEGV, handler); // install our handler
921
922 // init geometry and mag. field
923 new TGeoManager("Geometry", "Geane geometry");
924 TGeoManager::Import("genfitGeom.root");
927
928 TDatabasePDG::Instance()->GetParticle(211);
929
930
931 unsigned int nFailed(0);
932 const unsigned int nTests(100);
933
934 for (unsigned int i=0; i<nTests; ++i) {
935
936 if (!checkSetGetPosMom()) {
937 std::cout << "failed checkSetGetPosMom nr" << i << "\n";
938 ++nFailed;
939 }
940
942 std::cout << "failed compareForthBackExtrapolation nr" << i << "\n";
943 ++nFailed;
944 }
945
946 if (!checkStopAtBoundary()) {
947 std::cout << "failed checkStopAtBoundary nr" << i << "\n";
948 ++nFailed;
949 }
950
951 if (!checkErrorPropagation()) {
952 std::cout << "failed checkErrorPropagation nr" << i << "\n";
953 ++nFailed;
954 }
955
956 if (!checkExtrapolateToLine()) {
957 std::cout << "failed checkExtrapolateToLine nr" << i << "\n";
958 ++nFailed;
959 }
960
962 std::cout << "failed checkExtrapolateToPoint nr" << i << "\n";
963 ++nFailed;
964 }
965
967 std::cout << "failed checkExtrapolateToCylinder nr" << i << "\n";
968 ++nFailed;
969 }
970
972 std::cout << "failed checkExtrapolateToSphere nr" << i << "\n";
973 ++nFailed;
974 }
975
976 if (!checkExtrapolateBy()) {
977 std::cout << "failed checkExtrapolateBy nr" << i << "\n";
978 ++nFailed;
979 }
980
982 std::cout << "failed compareForthBackJacNoise nr" << i << "\n";
983 ++nFailed;
984 }
985
986 }
987
988 std::cout << "failed " << nFailed << " of " << nTests << " Tests." << std::endl;
989 if (nFailed == 0) {
990 std::cout << "passed all tests!" << std::endl;
991 }
992
993
994
995
996 return 0;
997}
998
999
Abstract base class for a track representation.
Definition AbsTrackRep.h:66
virtual double extrapolateToPlane(StateOnPlane &state, const genfit::SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to plane, and returns the extrapolation length and, via reference,...
virtual double extrapolateToCylinder(StateOnPlane &state, double radius, const TVector3 &linePoint=TVector3(0., 0., 0.), const TVector3 &lineDirection=TVector3(0., 0., 1.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to the cylinder surface, and returns the extrapolation length and,...
virtual double extrapolateToLine(StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to the POCA to a line, and returns the extrapolation length and,...
void calcJacobianNumerically(const genfit::StateOnPlane &origState, const genfit::SharedPlanePtr destPlane, TMatrixD &jacobian) const
Calculate Jacobian of transportation numerically. Slow but accurate. Can be used to validate (semi)an...
virtual TVector3 getMom(const StateOnPlane &state) const =0
Get the cartesian momentum vector of a state.
virtual void setPosMom(StateOnPlane &state, const TVector3 &pos, const TVector3 &mom) const =0
Set position and momentum of state.
virtual double extrapolateToSphere(StateOnPlane &state, double radius, const TVector3 &point=TVector3(0., 0., 0.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to the sphere surface, and returns the extrapolation length and,...
virtual TVector3 getPos(const StateOnPlane &state) const =0
Get the cartesian position of a state.
virtual void getForwardJacobianAndNoise(TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const =0
Get the jacobian and noise matrix of the last extrapolation.
virtual void getBackwardJacobianAndNoise(TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const =0
Get the jacobian and noise matrix of the last extrapolation if it would have been done in opposite di...
virtual double extrapolateToPoint(StateOnPlane &state, const TVector3 &point, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to the POCA to a point, and returns the extrapolation length and,...
virtual double extrapolateBy(StateOnPlane &state, double step, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state by step (cm) and returns the extrapolation length and, via reference,...
Constant Magnetic field.
Definition ConstField.h:37
Detector plane.
Definition DetPlane.h:61
void rotate(double angle)
rotate u and v around normal. Angle is in rad. More for debugging than for actual use.
Definition DetPlane.cc:320
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition Exception.h:48
virtual const char * what() const
Standard error message handling for exceptions. use like "std::cerr << e.what();".
Definition Exception.cc:51
void init(AbsBField *b)
set the magnetic field here. Magnetic field classes must be derived from AbsBField.
static FieldManager * getInstance()
Get singleton instance.
void setNoEffects(bool opt=true)
static MaterialEffects * getInstance()
void init(AbsMaterialInterface *matIfc)
set the material interface here. Material interface classes must be derived from AbsMaterialInterface...
StateOnPlane with additional covariance matrix.
const TMatrixDSym & getCov() const
virtual void Print(Option_t *option="") const
AbsTrackRep with 5D track parameterization in plane coordinates: (q/p, u', v', u, v)
Definition RKTrackRep.h:70
A state with arbitrary dimension defined in a DetPlane.
TVector3 getPos() const
const TVectorD & getState() const
virtual void Print(Option_t *option="") const
const SharedPlanePtr & getPlane() const
AbsMaterialInterface implementation for use with ROOT's TGeoManager.
void handler(int sig)
Definition main.cc:72
int randomSign()
Definition main.cc:111
int randomPdg()
Definition main.cc:85
int main()
Definition main.cc:133
boost::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.
bool isCovMatrix(TMatrixTBase< double > &cov)
Definition main.cc:133
bool compareMatrices(const TMatrixTBase< double > &A, const TMatrixTBase< double > &B, double maxRelErr)
Definition main.cc:104
bool checkStopAtBoundary()
Definition main.cc:497
bool checkExtrapolateToLine()
Definition main.cc:605
bool checkExtrapolateToPoint()
Definition main.cc:663
bool compareForthBackJacNoise()
Definition main.cc:296
bool checkSetGetPosMom()
Definition main.cc:158
bool checkExtrapolateToCylinder()
Definition main.cc:718
bool checkErrorPropagation()
Definition main.cc:554
bool checkExtrapolateToSphere()
Definition main.cc:787
bool checkExtrapolateBy()
Definition main.cc:850
bool compareForthBackExtrapolation()
Definition main.cc:224