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 <AbsKalmanFitter.h>
16#include <KalmanFitter.h>
18#include <KalmanFitterInfo.h>
19#include <KalmanFitStatus.h>
20#include <DAF.h>
21#include <GFGbl.h>
23#include <MeasurementOnPlane.h>
24#include <FullMeasurement.h>
25#include <PlanarMeasurement.h>
29#include <SharedPlanePtr.h>
31#include <StateOnPlane.h>
32#include <Tools.h>
33#include <TrackCand.h>
34#include <TrackCandHit.h>
35#include <Track.h>
36#include <TrackPoint.h>
37#include <WireMeasurement.h>
39
40#include <MaterialEffects.h>
41#include <RKTools.h>
42#include <RKTrackRep.h>
43#include <StepLimits.h>
44#include <TGeoMaterialInterface.h>
45
46#include <EventDisplay.h>
47
48#include <HelixTrackModel.h>
49#include <MeasurementCreator.h>
50
51#include <TApplication.h>
52#include <TCanvas.h>
53#include <TDatabasePDG.h>
54#include <TEveManager.h>
55#include <TGeoManager.h>
56#include <TH1D.h>
57#include <TRandom.h>
58#include <TStyle.h>
59#include <TVector3.h>
60#include <vector>
61
62#include <TROOT.h>
63#include <TFile.h>
64#include <TTree.h>
65#include "TDatabasePDG.h"
66#include <TMath.h>
67#include <TString.h>
68
69#include <memory>
70
71
72void handler(int sig) {
73 void *array[10];
74 size_t size;
75
76 // get void*'s for all entries on the stack
77 size = backtrace(array, 10);
78
79 // print out all the frames to stderr
80 fprintf(stderr, "Error: signal %d:\n", sig);
81 backtrace_symbols_fd(array, size, 2);
82 exit(1);
83}
84
85int randomPdg() {
86 int pdg;
87
88 switch(int(gRandom->Uniform(8))) {
89 case 1:
90 pdg = -11; break;
91 case 2:
92 pdg = 11; break;
93 case 3:
94 pdg = 13; break;
95 case 4:
96 pdg = -13; break;
97 case 5:
98 pdg = 211; break;
99 case 6:
100 pdg = -211; break;
101 case 7:
102 pdg = 2212; break;
103 default:
104 pdg = 211;
105 }
106
107 return pdg;
108}
109
110
112 if (gRandom->Uniform(1) > 0.5)
113 return 1;
114 return -1;
115}
116
117//---------------------------------------------------------------------------------------------------------
118//---------------------------------------------------------------------------------------------------------
119//---------------------------------------------------------------------------------------------------------
120//---------------------------------------------------------------------------------------------------------
121//---------------------------------------------------------------------------------------------------------
122
123//#define VALGRIND
124
125#ifdef VALGRIND
126 #include <valgrind/callgrind.h>
127#else
128#define CALLGRIND_START_INSTRUMENTATION
129#define CALLGRIND_STOP_INSTRUMENTATION
130#define CALLGRIND_DUMP_STATS
131#endif
132
133int main() {
134 std::cout<<"main"<<std::endl;
135 gRandom->SetSeed(14);
136
137
138 const unsigned int nEvents = 100;
139 const unsigned int nMeasurements = 56;
140 const double BField = 15.; // kGauss
141 const double momentum = 0.4; // GeV
142 const double theta = 120; // degree
143 const double thetaDetPlane = 90; // degree
144 const double phiDetPlane = 0; // degree
145 const double pointDist = 5; // cm; approx. distance between measurements
146 const double resolution = 0.02; // cm; resolution of generated measurements
147
148 const double resolutionWire = 5*resolution; // cm; resolution of generated measurements
149 const TVector3 wireDir(0,0,1);
150 const double skewAngle(5);
151 const bool useSkew = true;
152 const int nSuperLayer = 10;
153 const double minDrift = 0.;
154 const double maxDrift = 2;
155 const bool idealLRResolution = false; // resolve the l/r ambiguities of the wire measurements
156
157 const double outlierProb = -0.1;
158 const double outlierRange = 2;
159
160 const double hitSwitchProb = -0.1; // probability to give hits to fit in wrong order (flip two hits)
161
162 const int splitTrack = -5; // for track merging testing.
163 const bool fullMeasurement = false; // put fit result of first tracklet as FullMeasurement into second tracklet, don't merge
164
166 //const genfit::eFitterType fitterId = genfit::RefKalman;
167 //const genfit::eFitterType fitterId = genfit::DafRef;
168 //const genfit::eFitterType fitterId = genfit::DafSimple;
169 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::weightedAverage;
170 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::unweightedClosestToReference;
171 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::unweightedClosestToPrediction;
172 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::weightedClosestToReference;
173 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::weightedClosestToPrediction;
174 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::unweightedClosestToReferenceWire;
176 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::weightedClosestToReferenceWire;
177 //const genfit::eMultipleMeasurementHandling mmHandling = genfit::weightedClosestToPredictionWire;
178 const int nIter = 20; // max number of iterations
179 const double dPVal = 1.E-3; // convergence criterion
180
181 const bool resort = true;
182 const bool prefit = false; // make a simple Kalman iteration before the actual fit
183 const bool refit = false; // if fit did not converge, try to fit again
184
185 const bool twoReps = false; // test if everything works with more than one rep in the tracks
186
187 const int pdg = 13; // particle pdg code
188
189 const bool smearPosMom = true; // init the Reps with smeared pos and mom
190 const double chargeSwitchProb = -0.1; // probability to seed with wrong charge sign
191 const double posSmear = 10*resolution; // cm
192 const double momSmear = 5. /180.*TMath::Pi(); // rad
193 const double momMagSmear = 0.2; // relative
194 const double zSmearFac = 2;
195
196
197 const bool matFX = false; // include material effects; can only be disabled for RKTrackRep!
198
199 const bool debug = false;
200 const bool onlyDisplayFailed = true; // only load non-converged tracks into the display
201
202 std::vector<genfit::eMeasurementType> measurementTypes;
203 measurementTypes.push_back(genfit::Pixel);
204 measurementTypes.push_back(genfit::Pixel);
205 measurementTypes.push_back(genfit::StripUV);
206 measurementTypes.push_back(genfit::StripUV);
207 measurementTypes.push_back(genfit::StripUV);
208 measurementTypes.push_back(genfit::StripUV);
209 for (unsigned int i = 0; i<50; ++i)
210 measurementTypes.push_back(genfit::Wire);
211
212
213 signal(SIGSEGV, handler); // install our handler
214
215 // init fitter
216 genfit::AbsKalmanFitter* fitter = 0;
217 switch (fitterId) {
219 fitter = new genfit::KalmanFitter(nIter, dPVal);
220 fitter->setMultipleMeasurementHandling(mmHandling);
221 break;
222
224 fitter = new genfit::KalmanFitterRefTrack(nIter, dPVal);
225 fitter->setMultipleMeasurementHandling(mmHandling);
226 break;
227
229 fitter = new genfit::DAF(false);
230 break;
231 case genfit::DafRef:
232 fitter = new genfit::DAF();
233 break;
234 }
235 fitter->setMaxIterations(nIter);
236 if (debug)
237 fitter->setDebugLvl(10);
238
239 /*if (dynamic_cast<genfit::DAF*>(fitter) != NULL) {
240 //static_cast<genfit::DAF*>(fitter)->setBetas(100, 50, 25, 12, 6, 3, 1, 0.5, 0.1);
241 //static_cast<genfit::DAF*>(fitter)->setBetas(81, 8, 4, 0.5, 0.1);
242 static_cast<genfit::DAF*>(fitter)->setAnnealingScheme(100, 0.1, 5);
243 //static_cast<genfit::DAF*>(fitter)->setConvergenceDeltaWeight(0.0001);
244 //fitter->setMaxIterations(nIter);
245 }*/
246
247
248 // init MeasurementCreator
249 genfit::MeasurementCreator measurementCreator;
250 measurementCreator.setResolution(resolution);
251 measurementCreator.setResolutionWire(resolutionWire);
252 measurementCreator.setOutlierProb(outlierProb);
253 measurementCreator.setOutlierRange(outlierRange);
254 measurementCreator.setThetaDetPlane(thetaDetPlane);
255 measurementCreator.setPhiDetPlane(phiDetPlane);
256 measurementCreator.setWireDir(wireDir);
257 measurementCreator.setMinDrift(minDrift);
258 measurementCreator.setMaxDrift(maxDrift);
259 measurementCreator.setIdealLRResolution(idealLRResolution);
260 measurementCreator.setUseSkew(useSkew);
261 measurementCreator.setSkewAngle(skewAngle);
262 measurementCreator.setNSuperLayer(nSuperLayer);
263 measurementCreator.setDebug(debug);
264
265
266 // init geometry and mag. field
267 new TGeoManager("Geometry", "Geane geometry");
268 TGeoManager::Import("genfitGeom.root");
272
273 const double charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/(3.);
274
275 // init event display
276#ifndef VALGRIND
278 display->reset();
279#endif
280
281
282#ifndef VALGRIND
283 // create histograms
284 gROOT->SetStyle("Plain");
285 gStyle->SetPalette(1);
286 gStyle->SetOptFit(1111);
287
288 TH1D *hmomRes = new TH1D("hmomRes","mom res",500,-10*resolution*momentum/nMeasurements,10*resolution*momentum/nMeasurements);
289 TH1D *hupRes = new TH1D("hupRes","u' res",500,-5*resolution/nMeasurements, 5*resolution/nMeasurements);
290 TH1D *hvpRes = new TH1D("hvpRes","v' res",500,-5*resolution/nMeasurements, 5*resolution/nMeasurements);
291 TH1D *huRes = new TH1D("huRes","u res",500,-5*resolution, 5*resolution);
292 TH1D *hvRes = new TH1D("hvRes","v res",500,-5*resolution, 5*resolution);
293
294 TH1D *hqopPu = new TH1D("hqopPu","q/p pull",200,-6.,6.);
295 TH1D *pVal = new TH1D("pVal","p-value",100,0.,1.00000001);
296 pVal->SetMinimum(0);
297 TH1D *hupPu = new TH1D("hupPu","u' pull",200,-6.,6.);
298 TH1D *hvpPu = new TH1D("hvpPu","v' pull",200,-6.,6.);
299 TH1D *huPu = new TH1D("huPu","u pull",200,-6.,6.);
300 TH1D *hvPu = new TH1D("hvPu","v pull",200,-6.,6.);
301
302 TH1D *weights = new TH1D("weights","Daf vs true weights",500,-1.01,1.01);
303
304 TH1D *trackLenRes = new TH1D("trackLenRes","(trueLen - FittedLen) / trueLen",500,-0.01,0.01);
305#endif
306
307 double maxWeight(0);
308 unsigned int nTotalIterConverged(0);
309 unsigned int nTotalIterNotConverged(0);
310 unsigned int nTotalIterSecondConverged(0);
311 unsigned int nTotalIterSecondNotConverged(0);
312 unsigned int nConvergedFits(0);
313 unsigned int nUNConvergedFits(0);
314 unsigned int nConvergedFitsSecond(0);
315 unsigned int nUNConvergedFitsSecond(0);
316
317
319
320 genfit::Track* fitTrack(NULL);
321 genfit::Track* secondTrack(NULL);
322
323 // main loop
324 for (unsigned int iEvent=0; iEvent<nEvents; ++iEvent){
325
326 /*measurementTypes.clear();
327 for (unsigned int i = 0; i < nMeasurements; ++i)
328 measurementTypes.push_back(genfit::eMeasurementType(gRandom->Uniform(8)));*/
329
330
331 if (debug || (iEvent)%10==0)
332 std::cout << "Event Nr. " << iEvent << " ";
333 else std::cout << ". ";
334 if (debug || (iEvent+1)%10==0)
335 std::cout << "\n";
336
337
338 // clean up
339 delete fitTrack;
340 fitTrack = NULL;
341 delete secondTrack;
342 secondTrack = NULL;
343
344
345 // true start values
346 TVector3 pos(0, 0, 0);
347 TVector3 mom(1.,0,0);
348 mom.SetPhi(gRandom->Uniform(0.,2*TMath::Pi()));
349 //mom.SetTheta(gRandom->Uniform(0.4*TMath::Pi(),0.6*TMath::Pi()));
350 mom.SetTheta(theta*TMath::Pi()/180);
351 mom.SetMag(momentum);
352 TMatrixDSym covM(6);
353 for (int i = 0; i < 3; ++i)
354 covM(i,i) = resolution*resolution;
355 for (int i = 3; i < 6; ++i)
356 covM(i,i) = pow(resolution / nMeasurements / sqrt(3), 2);
357
358 if (debug) {
359 std::cout << "start values \n";
360 pos.Print();
361 mom.Print();
362 }
363
364 // calc helix parameters
365 genfit::HelixTrackModel* helix = new genfit::HelixTrackModel(pos, mom, charge);
366 measurementCreator.setTrackModel(helix);
367
368 // smeared start values
369 TVector3 posM(pos);
370 TVector3 momM(mom);
371 if (smearPosMom) {
372 posM.SetX(gRandom->Gaus(posM.X(),posSmear));
373 posM.SetY(gRandom->Gaus(posM.Y(),posSmear));
374 posM.SetZ(gRandom->Gaus(posM.Z(),zSmearFac*posSmear));
375
376 momM.SetPhi(gRandom->Gaus(mom.Phi(),momSmear));
377 momM.SetTheta(gRandom->Gaus(mom.Theta(),momSmear));
378 momM.SetMag(gRandom->Gaus(mom.Mag(), momMagSmear*mom.Mag()));
379 }
380
381
382 // trackrep for creating measurements
383 double sign(1.);
384 if (chargeSwitchProb > gRandom->Uniform(1.))
385 sign = -1.;
386 genfit::AbsTrackRep* rep = new genfit::RKTrackRep(sign*pdg);
387 sign = 1.;
388 if (chargeSwitchProb > gRandom->Uniform(1.))
389 sign = -1.;
390 genfit::AbsTrackRep* secondRep = new genfit::RKTrackRep(sign*-211);
391 genfit::MeasuredStateOnPlane stateRef(rep);
392 rep->setPosMomCov(stateRef, pos, mom, covM);
393
394 // smeared start state
395 genfit::MeasuredStateOnPlane stateSmeared(rep);
396 rep->setPosMomCov(stateSmeared, posM, momM, covM);
397
398 //rep->setPropDir(1);
399
401
402 // remember original initial state
403 const genfit::StateOnPlane stateRefOrig(stateRef);
404
405 // create smeared measurements
406 std::vector< std::vector<genfit::AbsMeasurement*> > measurements;
407
408 std::vector<bool> outlierTrue;
409 bool outlier;
410 // true values for left right. 0 for non wire measurements
411 std::vector<int> leftRightTrue;
412 int lr;
413
414 double trueLen;
415
416 try{
417 for (unsigned int i=0; i<measurementTypes.size(); ++i){
418 trueLen = i*pointDist;
419
420 measurements.push_back(measurementCreator.create(measurementTypes[i], trueLen, outlier, lr));
421 outlierTrue.push_back(outlier);
422 leftRightTrue.push_back(lr);
423 }
424 assert(measurementTypes.size() == leftRightTrue.size());
425 assert(measurementTypes.size() == outlierTrue.size());
426 }
427 catch(genfit::Exception& e){
428 std::cerr<<"Exception, next track"<<std::endl;
429 std::cerr << e.what();
430 continue; // here is a memleak!
431 }
432
433 if (debug) std::cout << "... done creating measurements \n";
434
435
436
437 // create track
438 TVectorD seedState(6);
439 TMatrixDSym seedCov(6);
440 rep->get6DStateCov(stateSmeared, seedState, seedCov);
441 fitTrack = new genfit::Track(rep, seedState, seedCov); //initialized with smeared rep
442 secondTrack = new genfit::Track(rep->clone(), seedState, seedCov); //initialized with smeared rep
443 if (twoReps) {
444 fitTrack->addTrackRep(secondRep);
445 secondTrack->addTrackRep(secondRep->clone());
446 }
447 else
448 delete secondRep;
449 //if (debug) fitTrack->Print("C");
450
451 assert(fitTrack->checkConsistency());
452 //fitTrack->addTrackRep(rep->clone()); // check if everything works fine with more than one rep
453
454 // add measurements
455 for(unsigned int i=0; i<measurements.size(); ++i){
456 if (splitTrack > 0 && (int)i >= splitTrack)
457 break;
458 if (i>0 && hitSwitchProb > gRandom->Uniform(1.))
459 fitTrack->insertPoint(new genfit::TrackPoint(measurements[i], fitTrack), -2);
460 else
461 fitTrack->insertPoint(new genfit::TrackPoint(measurements[i], fitTrack));
462
463 assert(fitTrack->checkConsistency());
464 //if (debug) fitTrack->Print("C");
465 }
466
467 if (splitTrack > 0) {
468 for(unsigned int i=splitTrack; i<measurements.size(); ++i){
469 if (i>0 && hitSwitchProb > gRandom->Uniform(1.))
470 secondTrack->insertPoint(new genfit::TrackPoint(measurements[i], secondTrack), -2);
471 else
472 secondTrack->insertPoint(new genfit::TrackPoint(measurements[i], secondTrack));
473
474 //if (debug) fitTrack->Print("C");
475 }
476 }
477
478 assert(fitTrack->checkConsistency());
479 assert(secondTrack->checkConsistency());
480
481 //if (debug) fitTrack->Print();
482
483 // do the fit
484 try{
485 if (debug) std::cout<<"Starting the fitter"<<std::endl;
486
487 if (prefit) {
488 genfit::KalmanFitter prefitter(1, dPVal);
490 prefitter.processTrackWithRep(fitTrack, fitTrack->getCardinalRep());
491 }
492
493 fitter->processTrack(fitTrack, resort);
494 if (splitTrack > 0)
495 fitter->processTrack(secondTrack, resort);
496
497 if (debug) std::cout<<"fitter is finished!"<<std::endl;
498 }
499 catch(genfit::Exception& e){
500 std::cerr << e.what();
501 std::cerr << "Exception, next track" << std::endl;
502 continue;
503 }
504
505 if (splitTrack > 0) {
506 if (debug) fitTrack->Print("C");
507 if (debug) secondTrack->Print("C");
508
509 if (fullMeasurement) {
511 fitTrack->insertPoint(new genfit::TrackPoint(fullM, fitTrack));
512 }
513 else
514 fitTrack->mergeTrack(secondTrack);
515
516 if (debug) fitTrack->Print("C");
517
518 try{
519 if (debug) std::cout<<"Starting the fitter"<<std::endl;
520 fitter->processTrack(fitTrack, resort);
521 if (debug) std::cout<<"fitter is finished!"<<std::endl;
522 }
523 catch(genfit::Exception& e){
524 std::cerr << e.what();
525 std::cerr << "Exception, next track" << std::endl;
526 continue;
527 }
528 }
529
530
531 if (refit && !fitTrack->getFitStatus(rep)->isFitConverged()) {
532 std::cout<<"Trying to fit again "<<std::endl;
533 fitter->processTrack(fitTrack, resort);
534 }
535
536
537
538 if (debug) {
539 fitTrack->Print("C");
540 fitTrack->getFitStatus(rep)->Print();
541 }
542
543 assert(fitTrack->checkConsistency());
544 assert(secondTrack->checkConsistency());
545
546#ifndef VALGRIND
547 if (!onlyDisplayFailed && iEvent < 1000)
548 display->addEvent(fitTrack);
549 else if (onlyDisplayFailed &&
550 (!fitTrack->getFitStatus(rep)->isFitConverged() ||
551 fitTrack->getFitStatus(rep)->getPVal() < 0.01)) {
552 // add track to event display
553 display->addEvent(fitTrack);
554 }
555#endif
556
557
558 if (fitTrack->getFitStatus(rep)->isFitConverged()) {
559 nTotalIterConverged += static_cast<genfit::KalmanFitStatus*>(fitTrack->getFitStatus(rep))->getNumIterations();
560 nConvergedFits += 1;
561 }
562 else {
563 nTotalIterNotConverged += static_cast<genfit::KalmanFitStatus*>(fitTrack->getFitStatus(rep))->getNumIterations();
564 nUNConvergedFits += 1;
565 }
566
567 if (twoReps) {
568 if (fitTrack->getFitStatus(secondRep)->isFitConverged()) {
569 nTotalIterSecondConverged += static_cast<genfit::KalmanFitStatus*>(fitTrack->getFitStatus(secondRep))->getNumIterations();
570 nConvergedFitsSecond += 1;
571 }
572 else {
573 nTotalIterSecondNotConverged += static_cast<genfit::KalmanFitStatus*>(fitTrack->getFitStatus(secondRep))->getNumIterations();
574 nUNConvergedFitsSecond += 1;
575 }
576 }
577
578
579 // check if fit was successful
580 if (! fitTrack->getFitStatus(rep)->isFitConverged()) {
581 std::cout << "Track could not be fitted successfully! Fit is not converged! \n";
582 continue;
583 }
584
585
586 genfit::TrackPoint* tp = fitTrack->getPointWithMeasurementAndFitterInfo(0, rep);
587 if (tp == NULL) {
588 std::cout << "Track has no TrackPoint with fitterInfo! \n";
589 continue;
590 }
591 genfit::KalmanFittedStateOnPlane kfsop(*(static_cast<genfit::KalmanFitterInfo*>(tp->getFitterInfo(rep))->getBackwardUpdate()));
592 if (debug) {
593 std::cout << "state before extrapolating back to reference plane \n";
594 kfsop.Print();
595 }
596
597 // extrapolate back to reference plane.
598 try{
599 rep->extrapolateToPlane(kfsop, stateRefOrig.getPlane());;
600 }
601 catch(genfit::Exception& e){
602 std::cerr<<"Exception, next track"<<std::endl;
603 std::cerr << e.what();
604 continue;
605 }
606
607#ifndef VALGRIND
608 // calculate pulls
609 const TVectorD& referenceState = stateRefOrig.getState();
610
611 const TVectorD& state = kfsop.getState();
612 const TMatrixDSym& cov = kfsop.getCov();
613
614 double pval = fitter->getPVal(fitTrack, rep);
615 assert( fabs(pval - static_cast<genfit::KalmanFitStatus*>(fitTrack->getFitStatus(rep))->getBackwardPVal()) < 1E-10 );
616
617 hmomRes->Fill( (charge/state[0]-momentum));
618 hupRes->Fill( (state[1]-referenceState[1]));
619 hvpRes->Fill( (state[2]-referenceState[2]));
620 huRes->Fill( (state[3]-referenceState[3]));
621 hvRes->Fill( (state[4]-referenceState[4]));
622
623 hqopPu->Fill( (state[0]-referenceState[0]) / sqrt(cov[0][0]) );
624 pVal->Fill( pval);
625 hupPu->Fill( (state[1]-referenceState[1]) / sqrt(cov[1][1]) );
626 hvpPu->Fill( (state[2]-referenceState[2]) / sqrt(cov[2][2]) );
627 huPu->Fill( (state[3]-referenceState[3]) / sqrt(cov[3][3]) );
628 hvPu->Fill( (state[4]-referenceState[4]) / sqrt(cov[4][4]) );
629
630
631 try {
632 trackLenRes->Fill( (trueLen - fitTrack->getTrackLen(rep)) / trueLen );
633
634 if (debug) {
635 std::cout << "true track length = " << trueLen << "; fitted length = " << fitTrack->getTrackLen(rep) << "\n";
636 std::cout << "fitted tof = " << fitTrack->getTOF(rep) << " ns\n";
637 }
638 }
639 catch (genfit::Exception& e) {
640 std::cerr << e.what();
641 std::cout << "could not get TrackLen or TOF! \n";
642 }
643
644
645 // check l/r resolution and outlier rejection
646 if (dynamic_cast<genfit::DAF*>(fitter) != NULL) {
647 for (unsigned int i=0; i<fitTrack->getNumPointsWithMeasurement(); ++i){
648
649 if (! fitTrack->getPointWithMeasurement(i)->hasFitterInfo(rep))
650 continue;
651
652 if (debug) {
653 std::vector<double> dafWeights = dynamic_cast<genfit::KalmanFitterInfo*>(fitTrack->getPointWithMeasurement(i)->getFitterInfo(rep))->getWeights();
654 std::cout << "hit " << i;
655 if (outlierTrue[i]) std::cout << " is an OUTLIER";
656 std::cout << " weights: ";
657 for (unsigned int j=0; j<dafWeights.size(); ++j){
658 std::cout << dafWeights[j] << " ";
659 }
660 std::cout << " l/r: " << leftRightTrue[i];
661 std::cout << "\n";
662 }
663 int trueSide = leftRightTrue[i];
664 if (trueSide == 0) continue; // not a wire measurement
665 if (outlierTrue[i]) continue; // an outlier
666 std::vector<double> dafWeightLR = dynamic_cast<genfit::KalmanFitterInfo*>(fitTrack->getPointWithMeasurement(i)->getFitterInfo(rep))->getWeights();
667 if(dafWeightLR.size() != 2)
668 continue;
669
670 double weightCorrectSide, weightWrongSide;
671
672 if (trueSide < 0) {
673 weightCorrectSide = dafWeightLR[0];
674 weightWrongSide = dafWeightLR[1];
675 }
676 else {
677 weightCorrectSide = dafWeightLR[1];
678 weightWrongSide = dafWeightLR[0];
679 }
680 weightWrongSide -= 1.;
681
682 weights->Fill(weightCorrectSide);
683 weights->Fill(weightWrongSide);
684
685 if (weightCorrectSide>maxWeight) maxWeight = weightCorrectSide;
686 }
687
688 for (unsigned int i=0; i<fitTrack->getNumPointsWithMeasurement(); ++i){
689 if (! fitTrack->getPointWithMeasurement(i)->hasFitterInfo(rep))
690 continue;
691
692 std::vector<double> dafWeights = dynamic_cast<genfit::KalmanFitterInfo*>(fitTrack->getPointWithMeasurement(i)->getFitterInfo(rep))->getWeights();
693
694 if (outlierTrue[i]) { // an outlier
695 for (unsigned int j=0; j<dafWeights.size(); ++j){
696 weights->Fill(dafWeights[j]-1);
697 }
698 }
699 else if (leftRightTrue[i] == 0) { // only for non wire hits
700 for (unsigned int j=0; j<dafWeights.size(); ++j){
701 weights->Fill(dafWeights[j]);
702 }
703 }
704 }
705
706 }
707
708
709 fitTrack->prune();
710 if (debug) {
711 std::cout<<" pruned track: ";
712 fitTrack->Print();
713 }
714
715#endif
716
717
718 }// end loop over events
719
720 delete fitTrack;
721 delete secondTrack;
722 delete fitter;
723
726
727 std::cout<<"maxWeight = " << maxWeight << std::endl;
728 std::cout<<"avg nr iterations = " << (double)(nTotalIterConverged + nTotalIterNotConverged)/(double)(nConvergedFits + nUNConvergedFits) << std::endl;
729 std::cout<<"avg nr iterations of converged fits = " << (double)(nTotalIterConverged)/(double)(nConvergedFits) << std::endl;
730 std::cout<<"avg nr iterations of UNconverged fits = " << (double)(nTotalIterNotConverged)/(double)(nUNConvergedFits) << std::endl;
731 std::cout<<"fit efficiency = " << (double)nConvergedFits/nEvents << std::endl;
732
733 if (twoReps) {
734 std::cout<<"second rep: \navg nr iterations = " << (double)(nTotalIterSecondConverged + nTotalIterSecondNotConverged)/(double)(nConvergedFitsSecond + nUNConvergedFitsSecond) << std::endl;
735 std::cout<<"avg nr iterations of converged fits = " << (double)(nTotalIterSecondConverged)/(double)(nConvergedFitsSecond) << std::endl;
736 std::cout<<"avg nr iterations of UNconverged fits = " << (double)(nTotalIterSecondNotConverged)/(double)(nUNConvergedFitsSecond) << std::endl;
737 std::cout<<"fit efficiency = " << (double)nConvergedFitsSecond/nEvents << std::endl;
738 }
739
740
741 //std::cout<<"avg nr iterations (2nd rep) = " << (double)nTotalIterSecond/nSuccessfullFitsSecond << std::endl;
742 //std::cout<<"fit efficiency (2nd rep) = " << (double)nConvergedFitsSecond/nEvents << std::endl;
743
744
745#ifndef VALGRIND
746
747 if (debug) std::cout<<"Draw histograms ...";
748 // fit and draw histograms
749 TCanvas* c1 = new TCanvas();
750 c1->Divide(2,3);
751
752 c1->cd(1);
753 hmomRes->Fit("gaus");
754 hmomRes->Draw();
755
756 c1->cd(2);
757 weights->Draw();
758
759 c1->cd(3);
760 hupRes->Fit("gaus");
761 hupRes->Draw();
762
763 c1->cd(4);
764 hvpRes->Fit("gaus");
765 hvpRes->Draw();
766
767 c1->cd(5);
768 huRes->Fit("gaus");
769 huRes->Draw();
770
771 c1->cd(6);
772 hvRes->Fit("gaus");
773 hvRes->Draw();
774
775 c1->Write();
776
777 TCanvas* c2 = new TCanvas();
778 c2->Divide(2,3);
779
780 c2->cd(1);
781 hqopPu->Fit("gaus");
782 hqopPu->Draw();
783
784 c2->cd(2);
785 pVal->Fit("pol1");
786 pVal->Draw();
787 c2->cd(3);
788 hupPu->Fit("gaus");
789 hupPu->Draw();
790
791 c2->cd(4);
792 hvpPu->Fit("gaus");
793 hvpPu->Draw();
794
795 c2->cd(5);
796 huPu->Fit("gaus");
797 huPu->Draw();
798
799 c2->cd(6);
800 hvPu->Fit("gaus");
801 hvPu->Draw();
802
803 c2->Write();
804
805
806
807 TCanvas* c3 = new TCanvas();
808 //c3->Divide(2,3);
809
810 c3->cd(1);
811 trackLenRes->Fit("gaus");
812 trackLenRes->Draw();
813
814 c3->Write();
815
816 if (debug) std::cout<<"... done"<<std::endl;
817
818 // open event display
819 display->setOptions("ABDEFHMPT"); // G show geometry
820 if (matFX) display->setOptions("ABDEFGHMPT"); // G show geometry
821 display->open();
822
823
824#endif
825
826
827}
Abstract base class for Kalman fitter and derived fitting algorithms.
void setMultipleMeasurementHandling(eMultipleMeasurementHandling mmh)
How should multiple measurements be handled?
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 void setPosMomCov(MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const =0
Set position, momentum and covariance of state.
virtual void get6DStateCov(const MeasuredStateOnPlane &state, TVectorD &stateVec, TMatrixDSym &cov) const
Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance.
virtual AbsTrackRep * clone() const =0
Clone the trackRep.
Constant Magnetic field.
Definition ConstField.h:37
Determinstic Annealing Filter (DAF) implementation.
Definition DAF.h:48
Event display designed to run with Genfit.
void addEvent(std::vector< genfit::Track * > &tracks)
Add new event.
void reset()
Drop all events.
static EventDisplay * getInstance()
void setOptions(std::string opts)
Set the display options.
void open()
Open the event display.
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 useCache(bool opt=true, unsigned int nBuckets=8)
Cache last lookup positions, and use stored field values if a lookup at (almost) the same position is...
Measurement class implementing a measurement of all track parameters.
Helix track model for testing purposes.
FitStatus for use with AbsKalmanFitter implementations.
double getBackwardPVal() const
MeasuredStateOnPlane with additional info produced by a Kalman filter or DAF.
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
KalmanFittedStateOnPlane * getBackwardUpdate() const
Kalman filter implementation with linearization around a reference track.
Simple Kalman filter implementation.
void processTrackWithRep(Track *tr, const AbsTrackRep *rep, bool resortHits=false)
Hit resorting currently NOT supported.
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
Create different measurement types along a HelixTrackModel for testing purposes.
void setMaxDrift(double maxDrift)
void setSkewAngle(double skewAngle)
void setPhiDetPlane(double phiDetPlane)
void setResolutionWire(double resolutionWire)
void setOutlierRange(double outlierRange)
void setMinDrift(double minDrift)
void setWireDir(const TVector3 wireDir)
void setResolution(double resolution)
std::vector< genfit::AbsMeasurement * > create(eMeasurementType, double tracklength, bool &outlier, int &lr)
void setNSuperLayer(int nSuperLayer)
void setOutlierProb(double outlierProb)
void setTrackModel(const HelixTrackModel *model)
Takes ownership!
void setThetaDetPlane(double thetaDetPlane)
void setIdealLRResolution(bool idealLRResolution)
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.
const TVectorD & getState() const
const SharedPlanePtr & getPlane() const
AbsMaterialInterface implementation for use with ROOT's TGeoManager.
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition TrackPoint.h:50
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
Definition Track.h:71
void addTrackRep(AbsTrackRep *trackRep)
Definition Track.cc:506
const MeasuredStateOnPlane & getFittedState(int id=0, const AbsTrackRep *rep=NULL, bool biased=true) const
Shortcut to get FittedStates.
Definition Track.cc:231
void insertPoint(TrackPoint *point, int id=-1)
Insert TrackPoint BEFORE TrackPoint with position id, if id >= 0.
Definition Track.cc:307
void Print(const Option_t *="") const
Definition Track.cc:1026
bool checkConsistency() const
Definition Track.cc:1192
#define CALLGRIND_START_INSTRUMENTATION
Definition main.cc:128
void handler(int sig)
Definition main.cc:72
int randomSign()
Definition main.cc:111
int randomPdg()
Definition main.cc:85
#define CALLGRIND_DUMP_STATS
Definition main.cc:130
#define CALLGRIND_STOP_INSTRUMENTATION
Definition main.cc:129
int main()
Definition main.cc:133
eMultipleMeasurementHandling
@ weightedClosestToPrediction
@ unweightedClosestToPredictionWire
@ SimpleKalman