SND@LHC Software
Loading...
Searching...
No Matches
Track.cc
Go to the documentation of this file.
1/* Copyright 2008-2010, Technische Universitaet Muenchen,
2 Authors: Christian Hoeppner & Sebastian Neubert & Johannes Rauch
3
4 This file is part of GENFIT.
5
6 GENFIT is free software: you can redistribute it and/or modify
7 it under the terms of the GNU Lesser General Public License as published
8 by the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 GENFIT is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU Lesser General Public License for more details.
15
16 You should have received a copy of the GNU Lesser General Public License
17 along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#include "Track.h"
21
22#include "Exception.h"
23#include "KalmanFitterInfo.h"
24#include "KalmanFitStatus.h"
25#include "PlanarMeasurement.h"
26
27#include <algorithm>
28#include <iostream>
29#include <map>
30
31#include <TDatabasePDG.h>
32#include <TMath.h>
33#include <TBuffer.h>
34
35//#include <glog/logging.h>
36
37//#define DEBUG
38
39
40namespace genfit {
41
43 cardinalRep_(0), fitStatuses_(), stateSeed_(6), covSeed_(6)
44{
45 ;
46}
47
48
50 cardinalRep_(0), fitStatuses_(), stateSeed_(6), covSeed_(6)
51{
52
53 if (rep != NULL)
54 addTrackRep(rep);
55
56 // create the measurements using the factory.
57 std::vector <genfit::AbsMeasurement*> factoryHits = factory.createMany(trackCand);
58
59 // create TrackPoints
60 PlanarMeasurement* lastPlanarMeas(NULL);
61 PlanarMeasurement* planarMeas(NULL);
62 for (unsigned int i=0; i<factoryHits.size(); ++i){
63 planarMeas = dynamic_cast<PlanarMeasurement*>(factoryHits[i]);
64
65 if (lastPlanarMeas != NULL && planarMeas != NULL &&
66 lastPlanarMeas->getDetId() == planarMeas->getDetId() &&
67 planarMeas->getPlaneId() != -1 && // -1 is default plane id
68 lastPlanarMeas->getPlaneId() == planarMeas->getPlaneId() ) {
69 trackPoints_.back()->addRawMeasurement(factoryHits[i]);
70 }
71 else
72 insertPoint(new genfit::TrackPoint(factoryHits[i], this));
73
74 lastPlanarMeas = dynamic_cast<PlanarMeasurement*>(factoryHits[i]);
75 }
76
77 // fill seed state
78 stateSeed_ = trackCand.getStateSeed();
79
80 // initial guess for cov
81 covSeed_ = trackCand.getCovSeed();
82
83 // fill cache
85
86 // self test
87 assert(checkConsistency());
88}
89
90
91Track::Track(AbsTrackRep* trackRep, const TVectorD& stateSeed) :
92 cardinalRep_(0), fitStatuses_(), stateSeed_(stateSeed),
93 covSeed_(TMatrixDSym::kUnit, TMatrixDSym(6))
94{
95 addTrackRep(trackRep);
96}
97
98
99Track::Track(AbsTrackRep* trackRep, const TVector3& posSeed, const TVector3& momSeed) :
100 cardinalRep_(0), fitStatuses_(), stateSeed_(6),
101 covSeed_(TMatrixDSym::kUnit, TMatrixDSym(6))
102{
103 addTrackRep(trackRep);
104 setStateSeed(posSeed, momSeed);
105}
106
107
108Track::Track(AbsTrackRep* trackRep, const TVectorD& stateSeed, const TMatrixDSym& covSeed) :
109 cardinalRep_(0), fitStatuses_(), stateSeed_(stateSeed), covSeed_(covSeed)
110{
111 addTrackRep(trackRep);
112}
113
114
115Track::Track(const Track& rhs) :
116 TObject(rhs),
117 cardinalRep_(rhs.cardinalRep_), stateSeed_(rhs.stateSeed_), covSeed_(rhs.covSeed_)
118{
119 assert(rhs.checkConsistency());
120
121 std::map<const AbsTrackRep*, AbsTrackRep*> oldRepNewRep;
122
123 for (std::vector<AbsTrackRep*>::const_iterator it=rhs.trackReps_.begin(); it!=rhs.trackReps_.end(); ++it) {
124 AbsTrackRep* newRep = (*it)->clone();
125 addTrackRep(newRep);
126 oldRepNewRep[(*it)] = newRep;
127 }
128
129 trackPoints_.reserve(rhs.trackPoints_.size());
130 for (std::vector<TrackPoint*>::const_iterator it=rhs.trackPoints_.begin(); it!=rhs.trackPoints_.end(); ++it) {
131 trackPoints_.push_back(new TrackPoint(**it, oldRepNewRep));
132 trackPoints_.back()->setTrack(this);
133 }
134
135 for (std::map< const AbsTrackRep*, FitStatus* >::const_iterator it=rhs.fitStatuses_.begin(); it!=rhs.fitStatuses_.end(); ++it) {
136 setFitStatus(it->second->clone(), oldRepNewRep[it->first]);
137 }
138
140
141 // self test
142 assert(checkConsistency());
143}
144
146 swap(other);
147
148 for (std::vector<TrackPoint*>::const_iterator it=trackPoints_.begin(); it!=trackPoints_.end(); ++it) {
149 trackPoints_.back()->setTrack(this);
150 }
151
153
154 // self test
155 assert(checkConsistency());
156
157 return *this;
158}
159
160void Track::swap(Track& other) {
161 std::swap(this->trackReps_, other.trackReps_);
162 std::swap(this->cardinalRep_, other.cardinalRep_);
163 std::swap(this->trackPoints_, other.trackPoints_);
165 std::swap(this->fitStatuses_, other.fitStatuses_);
166 std::swap(this->stateSeed_, other.stateSeed_);
167 std::swap(this->covSeed_, other.covSeed_);
168
169}
170
172 // causes problem for python, is it needed ? TR 2014 this->Clear();
173 this->Clear();
174}
175
176void Track::Clear(Option_t*)
177{
178 // This function is needed for TClonesArray embedding.
179 // FIXME: smarter containers or pointers needed ...
180 for (size_t i = 0; i < trackPoints_.size(); ++i)
181 delete trackPoints_[i];
182
183 trackPoints_.clear();
185
186 for (std::map< const AbsTrackRep*, FitStatus* >::iterator it = fitStatuses_.begin(); it!= fitStatuses_.end(); ++it)
187 delete it->second;
188 fitStatuses_.clear();
189
190 for (size_t i = 0; i < trackReps_.size(); ++i)
191 delete trackReps_[i];
192 trackReps_.clear();
193
194 cardinalRep_ = 0;
195
196 stateSeed_ *= 0;
197 covSeed_ *= 0;
198}
199
200
202 if (id < 0)
203 id += trackPoints_.size();
204
205 return trackPoints_.at(id);
206}
207
208
210 if (id < 0)
211 id += trackPointsWithMeasurement_.size();
212
213 return trackPointsWithMeasurement_.at(id);
214}
215
216
218 int i(0);
219 for (std::vector<TrackPoint*>::const_iterator it = trackPointsWithMeasurement_.begin(); it != trackPointsWithMeasurement_.end(); ++it) {
220 if ((*it)->hasFitterInfo(rep)) {
221 if (id == i)
222 return (*it);
223 ++i;
224 }
225 }
226
227 return NULL;
228}
229
230
231const MeasuredStateOnPlane& Track::getFittedState(int id, const AbsTrackRep* rep, bool biased) const {
232 if (rep == NULL)
233 rep = getCardinalRep();
234
236 if (point == NULL) {
237 Exception exc("Track::getFittedState ==> no trackPoint with fitterInfo for rep",__LINE__,__FILE__);
238 exc.setFatal();
239 throw exc;
240 }
241 return point->getFitterInfo(rep)->getFittedState(biased);
242}
243
244
245int Track::getIdForRep(const AbsTrackRep* rep) const
246{
247 for (size_t i = 0; i < trackReps_.size(); ++i)
248 if (trackReps_[i] == rep)
249 return i;
250
251 // std::cout << "track debug " << trackReps_.size() <<" "<< trackReps_[0] << " " << rep <<"\n";
252
253 //assert(0 == 1); // Cannot happen.
254 return 0; //RDM
255}
256
257
258bool Track::hasFitStatus(const AbsTrackRep* rep) const {
259 if (rep == NULL)
260 rep = getCardinalRep();
261
262 if (fitStatuses_.find(rep) == fitStatuses_.end())
263 return false;
264
265 return (fitStatuses_.at(rep) != NULL);
266}
267
268
270 if (rep == NULL)
271 rep = getCardinalRep();
272
273 if (fitStatuses_.find(rep) == fitStatuses_.end())
274 return false;
275
276 return (dynamic_cast<KalmanFitStatus*>(fitStatuses_.at(rep)) != NULL);
277}
278
279
281 return dynamic_cast<KalmanFitStatus*>(getFitStatus(rep));
282}
283
284
285void Track::setFitStatus(FitStatus* fitStatus, const AbsTrackRep* rep) {
286 if (fitStatuses_.find(rep) != fitStatuses_.end())
287 delete fitStatuses_.at(rep);
288
289 fitStatuses_[rep] = fitStatus;
290}
291
292
293void Track::setStateSeed(const TVector3& pos, const TVector3& mom) {
294 stateSeed_.ResizeTo(6);
295
296 stateSeed_(0) = pos.X();
297 stateSeed_(1) = pos.Y();
298 stateSeed_(2) = pos.Z();
299
300 stateSeed_(3) = mom.X();
301 stateSeed_(4) = mom.Y();
302 stateSeed_(5) = mom.Z();
303}
304
305
306
307void Track::insertPoint(TrackPoint* point, int id) {
308
309 point->setTrack(this);
310
311 #ifdef DEBUG
312 std::cout << "Track::insertPoint at position " << id << "\n";
313 #endif
314 assert(point!=NULL);
316
317 point->setTrack(this);
318
319 if (trackPoints_.size() == 0) {
320 trackPoints_.push_back(point);
321
322 if (point->hasRawMeasurements())
323 trackPointsWithMeasurement_.push_back(point);
324
325 return;
326 }
327
328 if (id == -1 || id == (int)trackPoints_.size()) {
329 trackPoints_.push_back(point);
330
331 if (point->hasRawMeasurements())
332 trackPointsWithMeasurement_.push_back(point);
333
334 deleteReferenceInfo(std::max(0, (int)trackPoints_.size()-2), (int)trackPoints_.size()-1);
335
336 // delete fitter infos if inserted point has a measurement
337 if (point->hasRawMeasurements()) {
338 deleteForwardInfo(-1, -1);
339 deleteBackwardInfo(0, -2);
340 }
341
342 return;
343 }
344
345 // [-size, size-1] is the allowed range
346 assert(id < (ssize_t)trackPoints_.size() || -id-1 <= (ssize_t)trackPoints_.size());
347
348 if (id < 0)
349 id += trackPoints_.size() + 1;
350
351 // insert
352 trackPoints_.insert(trackPoints_.begin() + id, point); // insert inserts BEFORE
353
354 // delete fitter infos if inserted point has a measurement
355 if (point->hasRawMeasurements()) {
356 deleteForwardInfo(id, -1);
357 deleteBackwardInfo(0, id);
358 }
359
360 // delete reference info of neighbouring points
361 deleteReferenceInfo(std::max(0, id-1), std::min((int)trackPoints_.size()-1, id+1));
362
364}
365
366
367void Track::insertPoints(std::vector<genfit::TrackPoint*> points, int id) {
368
369 int nBefore = getNumPoints();
370 int n = points.size();
371
372 if (n == 0)
373 return;
374 if (n == 1) {
375 insertPoint(points[0], id);
376 return;
377 }
378
379 for (std::vector<genfit::TrackPoint*>::iterator p = points.begin(); p != points.end(); ++p)
380 (*p)->setTrack(this);
381
382 if (id == -1 || id == (int)trackPoints_.size()) {
383 trackPoints_.insert(trackPoints_.end(), points.begin(), points.end());
384
385 deleteReferenceInfo(std::max(0, nBefore-1), nBefore);
386
387 deleteForwardInfo(nBefore, -1);
388 deleteBackwardInfo(0, std::max(0, nBefore-1));
389
391
392 return;
393 }
394
395
396 assert(id < (ssize_t)trackPoints_.size() || -id-1 <= (ssize_t)trackPoints_.size());
397
398 if (id < 0)
399 id += trackPoints_.size() + 1;
400
401
402 // insert
403 trackPoints_.insert(trackPoints_.begin() + id, points.begin(), points.end()); // insert inserts BEFORE
404
405 // delete fitter infos if inserted point has a measurement
406 deleteForwardInfo(id, -1);
407 deleteBackwardInfo(0, id+n);
408
409 // delete reference info of neighbouring points
410 deleteReferenceInfo(std::max(0, id-1), std::min((int)trackPoints_.size()-1, id));
411 deleteReferenceInfo(std::max(0, id+n-1), std::min((int)trackPoints_.size()-1, id+n));
412
414}
415
416
417void Track::deletePoint(int id) {
418
419 #ifdef DEBUG
420 std::cout << "Track::deletePoint at position " << id << "\n";
421 #endif
422
424
425 if (id < 0)
426 id += trackPoints_.size();
427 assert(id>0);
428
429
430 // delete forwardInfo after point (backwardInfo before point) if deleted point has a measurement
431 if (trackPoints_[id]->hasRawMeasurements()) {
432 deleteForwardInfo(id, -1);
433 deleteBackwardInfo(0, id-1);
434 }
435
436 // delete reference info of neighbouring points
437 deleteReferenceInfo(std::max(0, id-1), std::min((int)trackPoints_.size()-1, id+1));
438
439 // delete point
440 std::vector<TrackPoint*>::iterator it = std::find(trackPointsWithMeasurement_.begin(), trackPointsWithMeasurement_.end(), trackPoints_[id]);
441 if (it != trackPointsWithMeasurement_.end())
443
444 delete trackPoints_[id];
445 trackPoints_.erase(trackPoints_.begin()+id);
446
448
449}
450
451
452void Track::insertMeasurement(AbsMeasurement* measurement, int id) {
453 insertPoint(new TrackPoint(measurement, this), id);
454}
455
456
457void Track::mergeTrack(const Track* other, int id) {
458
459 #ifdef DEBUG
460 std::cout << "Track::mergeTrack\n";
461 #endif
462
463 if (other->getNumPoints() == 0)
464 return;
465
466 std::map<const AbsTrackRep*, AbsTrackRep*> otherRepThisRep;
467 std::vector<const AbsTrackRep*> otherRepsToRemove;
468
469 for (std::vector<AbsTrackRep*>::const_iterator otherRep=other->trackReps_.begin(); otherRep!=other->trackReps_.end(); ++otherRep) {
470 bool found(false);
471 for (std::vector<AbsTrackRep*>::const_iterator thisRep=trackReps_.begin(); thisRep!=trackReps_.end(); ++thisRep) {
472 if ((*thisRep)->isSame(*otherRep)) {
473 otherRepThisRep[*otherRep] = *thisRep;
474 #ifdef DEBUG
475 std::cout << " map other rep " << *otherRep << " to " << (*thisRep) << "\n";
476 #endif
477 if (found) {
478 Exception exc("Track::mergeTrack ==> more than one matching rep.",__LINE__,__FILE__);
479 exc.setFatal();
480 throw exc;
481 }
482 found = true;
483 break;
484 }
485 }
486 if (!found) {
487 otherRepsToRemove.push_back(*otherRep);
488 #ifdef DEBUG
489 std::cout << " remove other rep " << *otherRep << "\n";
490 #endif
491 }
492 }
493
494
495 std::vector<TrackPoint*> points;
496 points.reserve(other->getNumPoints());
497
498 for (std::vector<TrackPoint*>::const_iterator otherTp=other->trackPoints_.begin(); otherTp!=other->trackPoints_.end(); ++otherTp) {
499 points.push_back(new TrackPoint(**otherTp, otherRepThisRep, &otherRepsToRemove));
500 }
501
502 insertPoints(points, id);
503}
504
505
507 trackReps_.push_back(trackRep);
508 fitStatuses_[trackRep] = new FitStatus();
509}
510
511
513 if (id < 0)
514 id += trackReps_.size();
515
516 AbsTrackRep* rep = trackReps_.at(id);
517
518 // update cardinalRep_
519 if (int(cardinalRep_) == id)
520 cardinalRep_ = 0; // reset
521 else if (int(cardinalRep_) > id)
522 --cardinalRep_; // make cardinalRep_ point to the same TrackRep before and after deletion
523
524 // delete FitterInfos related to the deleted TrackRep
525 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin(); pointIt != trackPoints_.end(); ++pointIt) {
526 (*pointIt)->deleteFitterInfo(rep);
527 }
528
529 // delete fitStatus
530 delete fitStatuses_.at(rep);
531 fitStatuses_.erase(rep);
532
533 // delete rep
534 delete rep;
535 trackReps_.erase(trackReps_.begin()+id);
536}
537
538
540
541 if (id < 0)
542 id += trackReps_.size();
543
544 if (id >= 0 && (unsigned int)id < trackReps_.size())
545 cardinalRep_ = id;
546 else {
547 cardinalRep_ = 0;
548 std::cerr << "Track::setCardinalRep: Attempted to set cardinalRep_ to a value out of bounds. Resetting cardinalRep_ to 0." << std::endl;
549 }
550}
551
552
554
555 // Todo: test
556
557 if (trackReps_.size() <= 1)
558 return;
559
560 double minChi2(9.E99);
561 const AbsTrackRep* bestRep(NULL);
562
563 for (std::map< const AbsTrackRep*, FitStatus* >::const_iterator it = fitStatuses_.begin(); it != fitStatuses_.end(); ++it) {
564 if (it->second->isFitConverged()) {
565 if (it->second->getChi2() < minChi2) {
566 minChi2 = it->second->getChi2();
567 bestRep = it->first;
568 }
569 }
570 }
571
572 if (bestRep != NULL) {
573 setCardinalRep(getIdForRep(bestRep));
574 }
575}
576
577
579 #ifdef DEBUG
580 std::cout << "Track::sort \n";
581 #endif
582
583 int nPoints(trackPoints_.size());
584 // original order
585 std::vector<TrackPoint*> pointsBefore(trackPoints_);
586
587 // sort
588 std::stable_sort(trackPoints_.begin(), trackPoints_.end(), TrackPointComparator());
589
590 // see where order changed
591 int equalUntil(-1), equalFrom(nPoints);
592 for (int i = 0; i<nPoints; ++i) {
593 if (pointsBefore[i] == trackPoints_[i])
594 equalUntil = i;
595 else
596 break;
597 }
598
599 if (equalUntil == nPoints-1)
600 return false; // sorting did not change anything
601
602
604
605 for (int i = nPoints-1; i>equalUntil; --i) {
606 if (pointsBefore[i] == trackPoints_[i])
607 equalFrom = i;
608 else
609 break;
610 }
611
612 #ifdef DEBUG
613 std::cout << "Track::sort. Equal up to (including) hit " << equalUntil << " and from (including) hit " << equalFrom << " \n";
614 #endif
615
616 deleteForwardInfo(equalUntil+1, -1);
617 deleteBackwardInfo(0, equalFrom-1);
618
619 deleteReferenceInfo(std::max(0, equalUntil+1), std::min((int)trackPoints_.size()-1, equalFrom-1));
620
622
623 return true;
624}
625
626
627bool Track::udpateSeed(int id, AbsTrackRep* rep, bool biased) {
628 try {
629 const MeasuredStateOnPlane& fittedState = getFittedState(id, rep, biased);
630 setStateSeed(fittedState.get6DState());
631 setCovSeed(fittedState.get6DCov());
632
633 double fittedCharge = fittedState.getCharge();
634
635 for (unsigned int i = 0; i<trackReps_.size(); ++i) {
636 if (trackReps_[i]->getPDGCharge() * fittedCharge < 0) {
637 trackReps_[i]->switchPDGSign();
638 }
639 }
640 }
641 catch (Exception& e) {
642 // in this case the original track seed will be used
643 return false;
644 }
645 return true;
646}
647
648
650
651 std::reverse(trackPoints_.begin(),trackPoints_.end());
652
653 deleteForwardInfo(0, -1);
654 deleteBackwardInfo(0, -1);
655 deleteReferenceInfo(0, -1);
656
658}
659
660
662 if (rep != NULL) {
663 rep->switchPDGSign();
664 return;
665 }
666
667 for (unsigned int i = 0; i<trackReps_.size(); ++i) {
668 trackReps_[i]->switchPDGSign();
669 }
670}
671
672
674 udpateSeed(-1); // set fitted state of last hit as new seed
675 reverseMomSeed(); // flip momentum direction
677 reverseTrackPoints(); // also deletes all fitterInfos
678}
679
680
681void Track::deleteForwardInfo(int startId, int endId, const AbsTrackRep* rep) {
682 #ifdef DEBUG
683 std::cout << "Track::deleteForwardInfo from position " << startId << " to " << endId << "\n";
684 #endif
685
687
688 if (startId < 0)
689 startId += trackPoints_.size();
690 if (endId < 0)
691 endId += trackPoints_.size();
692 endId += 1;
693
694 assert (endId >= startId);
695
696 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
697 if (rep != NULL) {
698 if ((*pointIt)->hasFitterInfo(rep))
699 (*pointIt)->getFitterInfo(rep)->deleteForwardInfo();
700 }
701 else {
702 const std::vector<AbsFitterInfo*> fitterInfos = (*pointIt)->getFitterInfos();
703 for (std::vector<AbsFitterInfo*>::const_iterator fitterInfoIt = fitterInfos.begin(); fitterInfoIt != fitterInfos.end(); ++fitterInfoIt) {
704 (*fitterInfoIt)->deleteForwardInfo();
705 }
706 }
707 }
708}
709
710void Track::deleteBackwardInfo(int startId, int endId, const AbsTrackRep* rep) {
711
712 #ifdef DEBUG
713 std::cout << "Track::deleteBackwardInfo from position " << startId << " to " << endId << "\n";
714 #endif
715
717
718 if (startId < 0)
719 startId += trackPoints_.size();
720 if (endId < 0)
721 endId += trackPoints_.size();
722 endId += 1;
723
724 assert (endId >= startId);
725
726
727 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
728 if (rep != NULL) {
729 if ((*pointIt)->hasFitterInfo(rep))
730 (*pointIt)->getFitterInfo(rep)->deleteBackwardInfo();
731 }
732 else {
733 const std::vector<AbsFitterInfo*> fitterInfos = (*pointIt)->getFitterInfos();
734 for (std::vector<AbsFitterInfo*>::const_iterator fitterInfoIt = fitterInfos.begin(); fitterInfoIt != fitterInfos.end(); ++fitterInfoIt) {
735 (*fitterInfoIt)->deleteBackwardInfo();
736 }
737 }
738 }
739}
740
741void Track::deleteReferenceInfo(int startId, int endId, const AbsTrackRep* rep) {
742
743 #ifdef DEBUG
744 std::cout << "Track::deleteReferenceInfo from position " << startId << " to " << endId << "\n";
745 #endif
746
748
749 if (startId < 0)
750 startId += trackPoints_.size();
751 if (endId < 0)
752 endId += trackPoints_.size();
753 endId += 1;
754
755 assert (endId >= startId);
756
757 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
758 if (rep != NULL) {
759 if ((*pointIt)->hasFitterInfo(rep))
760 (*pointIt)->getFitterInfo(rep)->deleteReferenceInfo();
761 }
762 else {
763 std::vector<AbsFitterInfo*> fitterInfos = (*pointIt)->getFitterInfos();
764 for (std::vector<AbsFitterInfo*>::const_iterator fitterInfoIt = fitterInfos.begin(); fitterInfoIt != fitterInfos.end(); ++fitterInfoIt) {
765 (*fitterInfoIt)->deleteReferenceInfo();
766 }
767 }
768 }
769}
770
771void Track::deleteMeasurementInfo(int startId, int endId, const AbsTrackRep* rep) {
772
773 #ifdef DEBUG
774 std::cout << "Track::deleteMeasurementInfo from position " << startId << " to " << endId << "\n";
775 #endif
776
778
779 if (startId < 0)
780 startId += trackPoints_.size();
781 if (endId < 0)
782 endId += trackPoints_.size();
783 endId += 1;
784
785 assert (endId >= startId);
786
787 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
788 if (rep != NULL) {
789 if ((*pointIt)->hasFitterInfo(rep))
790 (*pointIt)->getFitterInfo(rep)->deleteMeasurementInfo();
791 }
792 else {
793 std::vector<AbsFitterInfo*> fitterInfos = (*pointIt)->getFitterInfos();
794 for (std::vector<AbsFitterInfo*>::const_iterator fitterInfoIt = fitterInfos.begin(); fitterInfoIt != fitterInfos.end(); ++fitterInfoIt) {
795 (*fitterInfoIt)->deleteMeasurementInfo();
796 }
797 }
798 }
799}
800
801void Track::deleteFitterInfo(int startId, int endId, const AbsTrackRep* rep) {
802
803 #ifdef DEBUG
804 std::cout << "Track::deleteFitterInfo from position " << startId << " to " << endId << "\n";
805 #endif
806
808
809 if (startId < 0)
810 startId += trackPoints_.size();
811 if (endId < 0)
812 endId += trackPoints_.size();
813 endId += 1;
814
815 assert (endId >= startId);
816
817 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
818 if (rep != NULL) {
819 if ((*pointIt)->hasFitterInfo(rep))
820 (*pointIt)->deleteFitterInfo(rep);
821 }
822 else {
823 for (std::vector<AbsTrackRep*>::const_iterator repIt = trackReps_.begin(); repIt != trackReps_.end(); ++repIt) {
824 if ((*pointIt)->hasFitterInfo(*repIt))
825 (*pointIt)->deleteFitterInfo(*repIt);
826 }
827 }
828 }
829}
830
831
832double Track::getTrackLen(AbsTrackRep* rep, int startId, int endId) const {
833
834 if (startId < 0)
835 startId += trackPoints_.size();
836 if (endId < 0)
837 endId += trackPoints_.size();
838
839 bool backwards(false);
840 if (startId > endId) {
841 double temp = startId;
842 startId = endId;
843 endId = temp;
844 backwards = true;
845 }
846
847 endId += 1;
848
849 if (rep == NULL)
850 rep = getCardinalRep();
851
852 double trackLen(0);
853 StateOnPlane state;
854
855 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
856 if (! (*pointIt)->hasFitterInfo(rep)) {
857 Exception e("Track::getTracklength: trackPoint has no fitterInfo", __LINE__,__FILE__);
858 throw e;
859 }
860
861 if (pointIt != trackPoints_.begin() + startId) {
862 trackLen += rep->extrapolateToPlane(state, (*pointIt)->getFitterInfo(rep)->getPlane());
863 }
864
865 state = (*pointIt)->getFitterInfo(rep)->getFittedState();
866 }
867
868 if (backwards)
869 trackLen *= -1.;
870
871 return trackLen;
872}
873
874
875double Track::getTOF(AbsTrackRep* rep, int startId, int endId) const {
876
877 if (startId < 0)
878 startId += trackPoints_.size();
879 if (endId < 0)
880 endId += trackPoints_.size();
881
882 bool backwards(false);
883 if (startId > endId) {
884 double temp = startId;
885 startId = endId;
886 endId = temp;
887 backwards = true;
888 }
889
890 endId += 1;
891
892 if (rep == NULL)
893 rep = getCardinalRep();
894
895 double tof(0);
896 StateOnPlane state;
897
898 for (std::vector<TrackPoint*>::const_iterator pointIt = trackPoints_.begin() + startId; pointIt != trackPoints_.begin() + endId; ++pointIt) {
899 if (! (*pointIt)->hasFitterInfo(rep)) {
900 Exception e("Track::getTOF: trackPoint has no fitterInfo", __LINE__,__FILE__);
901 throw e;
902 }
903
904 if (pointIt != trackPoints_.begin() + startId) {
905 rep->extrapolateToPlane(state, (*pointIt)->getFitterInfo(rep)->getPlane());
906 tof += rep->getTOF();
907 }
908
909 state = (*pointIt)->getFitterInfo(rep)->getFittedState();
910 }
911
912 if (backwards)
913 tof *= -1.;
914
915 return tof;
916}
917
918
919void Track::fixWeights(AbsTrackRep* rep, int startId, int endId) {
920
921 if (startId < 0)
922 startId += trackPoints_.size();
923 if (endId < 0)
924 endId += trackPoints_.size();
925
926 assert(startId >= 0);
927 assert(startId <= endId);
928 assert(endId <= (int)trackPoints_.size());
929
930 std::vector< genfit::AbsFitterInfo* > fis;
931
932 for (std::vector<TrackPoint*>::iterator tp = trackPoints_.begin() + startId; tp != trackPoints_.begin() + endId; ++tp) {
933 fis.clear();
934 if (rep == NULL) {
935 fis = (*tp)->getFitterInfos();
936 }
937 else if ((*tp)->hasFitterInfo(rep)) {
938 fis.push_back((*tp)->getFitterInfo(rep));
939 }
940
941 for (std::vector< genfit::AbsFitterInfo* >::iterator fi = fis.begin(); fi != fis.end(); ++fi) {
942 KalmanFitterInfo* kfi = dynamic_cast<KalmanFitterInfo*>(*fi);
943 if (kfi == NULL)
944 continue;
945
946 kfi->fixWeights();
947 }
948 }
949}
950
951
952void Track::prune(const Option_t* option) {
953
954 TString opt = option;
955 opt.ToUpper();
956
957 for (std::map< const AbsTrackRep*, FitStatus* >::const_iterator it=fitStatuses_.begin(); it!=fitStatuses_.end(); ++it) {
958 it->second->setIsTrackPruned();
959 }
960
961 // prune trackPoints
962 if (opt.Contains("F") || opt.Contains("L")) {
963 TrackPoint* firstPoint = getPointWithMeasurement(0);
964 TrackPoint* lastPoint = getPointWithMeasurement(-1);
965 for (unsigned int i = 0; i<trackPoints_.size(); ++i) {
966 if (trackPoints_[i] == firstPoint && opt.Contains("F"))
967 continue;
968
969 if (trackPoints_[i] == lastPoint && opt.Contains("L"))
970 continue;
971
972 delete trackPoints_[i];
973 trackPoints_.erase(trackPoints_.begin()+i);
974 --i;
975 }
976 }
977
978 // prune TrackReps
979 if (opt.Contains("C")) {
980 for (unsigned int i = 0; i < trackReps_.size(); ++i) {
981 if (i != cardinalRep_) {
983 --i;
984 }
985 }
986 }
987
988
989 // from remaining trackPoints: prune measurementsOnPlane, unneeded fitterInfoStuff
990 for (unsigned int i = 0; i<trackPoints_.size(); ++i) {
991 if (opt.Contains("W"))
992 trackPoints_[i]->deleteRawMeasurements();
993
994 std::vector< genfit::AbsFitterInfo* > fis = trackPoints_[i]->getFitterInfos();
995 for (unsigned int j = 0; j<fis.size(); ++j) {
996
997 if (i == 0 && opt.Contains("I") && opt.Contains("F") && opt.Contains("L"))
998 fis[j]->deleteForwardInfo();
999 else if (i == trackPoints_.size()-1 && opt.Contains("I") && opt.Contains("F") && opt.Contains("L"))
1000 fis[j]->deleteBackwardInfo();
1001 else if (opt.Contains("I") && opt.Contains("F"))
1002 fis[j]->deleteForwardInfo();
1003 else if (opt.Contains("I") && opt.Contains("L"))
1004 fis[j]->deleteBackwardInfo();
1005
1006 if (opt.Contains("U") && dynamic_cast<KalmanFitterInfo*>(fis[j]) != NULL) {
1007 static_cast<KalmanFitterInfo*>(fis[j])->deletePredictions();
1008 }
1009
1010 if (opt.Contains("R"))
1011 fis[j]->deleteReferenceInfo();
1012 if (opt.Contains("M"))
1013 fis[j]->deleteMeasurementInfo();
1014 }
1015 }
1016
1018
1019 #ifdef DEBUG
1020 std::cout << "pruned Track: "; Print();
1021 #endif
1022
1023}
1024
1025
1026void Track::Print(const Option_t* option) const {
1027 TString opt = option;
1028 opt.ToUpper();
1029 if (opt.Contains("C")) { // compact
1030
1031 std::cout << "\n ";
1032 for (unsigned int i=0; i<trackPoints_.size(); ++i) {
1033
1034 int color = 32*(size_t)(trackPoints_[i]) % 15;
1035 switch (color) {
1036 case 0:
1037 std::cout<<"\033[1;30m";
1038 break;
1039 case 1:
1040 std::cout<<"\033[0;34m";
1041 break;
1042 case 2:
1043 std::cout<<"\033[1;34m";
1044 break;
1045 case 3:
1046 std::cout<<"\033[0;32m";
1047 break;
1048 case 4:
1049 std::cout<<"\033[1;32m";
1050 break;
1051 case 5:
1052 std::cout<<"\033[0;36m";
1053 break;
1054 case 6:
1055 std::cout<<"\033[1;36m";
1056 break;
1057 case 7:
1058 std::cout<<"\033[0;31m";
1059 break;
1060 case 8:
1061 std::cout<<"\033[1;31m";
1062 break;
1063 case 9:
1064 std::cout<<"\033[0;35m";
1065 break;
1066 case 10:
1067 std::cout<<"\033[1;35m";
1068 break;
1069 case 11:
1070 std::cout<<"\033[0;33m";
1071 break;
1072 case 12:
1073 std::cout<<"\033[1;33m";
1074 break;
1075 case 13:
1076 std::cout<<"\033[0;37m";
1077 break;
1078 default:
1079 ;
1080 }
1081 std::cout << trackPoints_[i] << "\033[00m ";
1082 }
1083 std::cout << "\n";
1084
1085 std::cout << " ";
1086 for (unsigned int i=0; i<trackPoints_.size(); ++i) {
1087 printf("% -9.3g ", trackPoints_[i]->getSortingParameter());
1088 }
1089
1090 for (std::vector<AbsTrackRep*>::const_iterator rep = trackReps_.begin(); rep != trackReps_.end(); ++rep) {
1091 std::cout << "\n" << getIdForRep(*rep) << " ";
1092 for (unsigned int i=0; i<trackPoints_.size(); ++i) {
1093 if (! trackPoints_[i]->hasFitterInfo(*rep)) {
1094 std::cout << " ";
1095 continue;
1096 }
1097 AbsFitterInfo* fi = trackPoints_[i]->getFitterInfo(*rep);
1098 if (fi->hasMeasurements())
1099 std::cout << "M";
1100 else
1101 std::cout << " ";
1102
1103 if (fi->hasReferenceState())
1104 std::cout << "R";
1105 else
1106 std::cout << " ";
1107
1108 std::cout << " ";
1109 }
1110 std::cout << "\n";
1111
1112 std::cout << " -> ";
1113 for (unsigned int i=0; i<trackPoints_.size(); ++i) {
1114 if (! trackPoints_[i]->hasFitterInfo(*rep)) {
1115 std::cout << " ";
1116 continue;
1117 }
1118 AbsFitterInfo* fi = trackPoints_[i]->getFitterInfo(*rep);
1119 if (fi->hasForwardPrediction())
1120 std::cout << "P";
1121 else
1122 std::cout << " ";
1123
1124 if (fi->hasForwardUpdate())
1125 std::cout << "U";
1126 else
1127 std::cout << " ";
1128
1129 std::cout << " ";
1130 }
1131 std::cout << "\n";
1132
1133 std::cout << " <- ";
1134 for (unsigned int i=0; i<trackPoints_.size(); ++i) {
1135 if (! trackPoints_[i]->hasFitterInfo(*rep)) {
1136 std::cout << " ";
1137 continue;
1138 }
1139 AbsFitterInfo* fi = trackPoints_[i]->getFitterInfo(*rep);
1140 if (fi->hasBackwardPrediction())
1141 std::cout << "P";
1142 else
1143 std::cout << " ";
1144
1145 if (fi->hasBackwardUpdate())
1146 std::cout << "U";
1147 else
1148 std::cout << " ";
1149
1150 std::cout << " ";
1151 }
1152
1153 std::cout << "\n";
1154
1155 } //end loop over reps
1156
1157 std::cout << "\n";
1158 return;
1159 }
1160
1161
1162
1163 std::cout << "=======================================================================================\n";
1164 std::cout << "genfit::Track, containing " << trackPoints_.size() << " TrackPoints and " << trackReps_.size() << " TrackReps.\n";
1165 std::cout << " Seed state: "; stateSeed_.Print();
1166
1167 for (unsigned int i=0; i<trackReps_.size(); ++i) {
1168 std::cout << " TrackRep Nr. " << i;
1169 if (i == cardinalRep_)
1170 std::cout << " (This is the cardinal rep)";
1171 std::cout << "\n";
1172 trackReps_[i]->Print();
1173 }
1174
1175 std::cout << "---------------------------------------------------------------------------------------\n";
1176
1177 for (unsigned int i=0; i<trackPoints_.size(); ++i) {
1178 std::cout << "TrackPoint Nr. " << i << "\n";
1179 trackPoints_[i]->Print();
1180 std::cout << "..........................................................................\n";
1181 }
1182
1183 for (std::map< const AbsTrackRep*, FitStatus* >::const_iterator it=fitStatuses_.begin(); it!=fitStatuses_.end(); ++it) {
1184 it->second->Print();
1185 }
1186
1187 std::cout << "=======================================================================================\n";
1188
1189}
1190
1191
1193
1194 bool retVal(true);
1195
1196 std::map<const AbsTrackRep*, const KalmanFitterInfo*> prevFis;
1197
1198 // check if seed is 6D
1199 if (stateSeed_.GetNrows() != 6) {
1200 std::cerr << "Track::checkConsistency(): stateSeed_ dimension != 6" << std::endl;
1201 retVal = false;
1202 }
1203
1204 if (covSeed_.GetNrows() != 6) {
1205 std::cerr << "Track::checkConsistency(): covSeed_ dimension != 6" << std::endl;
1206 retVal = false;
1207 }
1208
1209 if (covSeed_.Max() == 0.) {
1210 std::cerr << "Track::checkConsistency(): Warning: covSeed_ is zero" << std::endl;
1211 //retVal = false;
1212 }
1213
1214 // check if cardinalRep_ is in range of trackReps_
1215 if (trackReps_.size() && cardinalRep_ >= trackReps_.size()) {
1216 std::cerr << "Track::checkConsistency(): cardinalRep id " << cardinalRep_ << " out of bounds" << std::endl;
1217 retVal = false;
1218 }
1219
1220 for (std::vector<AbsTrackRep*>::const_iterator rep = trackReps_.begin(); rep != trackReps_.end(); ++rep) {
1221 // check for NULL
1222 if ((*rep) == NULL) {
1223 std::cerr << "Track::checkConsistency(): TrackRep is NULL" << std::endl;
1224 retVal = false;
1225 }
1226
1227 // check for valid pdg code
1228 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle((*rep)->getPDG());
1229 if (particle == NULL) {
1230 std::cerr << "Track::checkConsistency(): TrackRep pdg ID " << (*rep)->getPDG() << " is not valid" << std::endl;
1231 retVal = false;
1232 }
1233
1234 }
1235
1236 // check TrackPoints
1237 for (std::vector<TrackPoint*>::const_iterator tp = trackPoints_.begin(); tp != trackPoints_.end(); ++tp) {
1238 // check for NULL
1239 if ((*tp) == NULL) {
1240 std::cerr << "Track::checkConsistency(): TrackPoint is NULL" << std::endl;
1241 retVal = false;
1242 }
1243 // check if trackPoint points back to this track
1244 if ((*tp)->getTrack() != this) {
1245 std::cerr << "Track::checkConsistency(): TrackPoint does not point back to this track" << std::endl;
1246 retVal = false;
1247 }
1248
1249 // check rawMeasurements
1250 const std::vector<AbsMeasurement*>& rawMeasurements = (*tp)->getRawMeasurements();
1251 for (std::vector<AbsMeasurement*>::const_iterator m = rawMeasurements.begin(); m != rawMeasurements.end(); ++m) {
1252 // check for NULL
1253 if ((*m) == NULL) {
1254 std::cerr << "Track::checkConsistency(): Measurement is NULL" << std::endl;
1255 retVal = false;
1256 }
1257 // check if measurement points back to TrackPoint
1258 if ((*m)->getTrackPoint() != *tp) {
1259 std::cerr << "Track::checkConsistency(): Measurement does not point back to correct TrackPoint" << std::endl;
1260 retVal = false;
1261 }
1262 }
1263
1264 // check fitterInfos
1265 std::vector<AbsFitterInfo*> fitterInfos = (*tp)->getFitterInfos();
1266 for (std::vector<AbsFitterInfo*>::const_iterator fi = fitterInfos.begin(); fi != fitterInfos.end(); ++fi) {
1267 // check for NULL
1268 if ((*fi) == NULL) {
1269 std::cerr << "Track::checkConsistency(): FitterInfo is NULL. TrackPoint: " << *tp << std::endl;
1270 retVal = false;
1271 }
1272
1273 if (!( (*fi)->checkConsistency() ) ) {
1274 std::cerr << "Track::checkConsistency(): FitterInfo not consistent. TrackPoint: " << *tp << std::endl;
1275 retVal = false;
1276 }
1277
1278 // check if fitterInfos point to valid TrackReps in trackReps_
1279 int mycount (0);
1280 for (std::vector<AbsTrackRep*>::const_iterator rep = trackReps_.begin(); rep != trackReps_.end(); ++rep) {
1281 if ( (*rep) == (*fi)->getRep() ) {
1282 ++mycount;
1283 }
1284 }
1285 if (mycount == 0) {
1286 std::cerr << "Track::checkConsistency(): fitterInfo points to TrackRep which is not in Track" << std::endl;
1287 retVal = false;
1288 }
1289
1290 if (dynamic_cast<KalmanFitterInfo*>(*fi) != NULL) {
1291 if (prevFis[(*fi)->getRep()] != NULL &&
1292 static_cast<KalmanFitterInfo*>(*fi)->hasReferenceState() &&
1293 prevFis[(*fi)->getRep()]->hasReferenceState() ) {
1294 double len = static_cast<KalmanFitterInfo*>(*fi)->getReferenceState()->getForwardSegmentLength();
1295 double prevLen = prevFis[(*fi)->getRep()]->getReferenceState()->getBackwardSegmentLength();
1296 if (fabs(prevLen + len) > 1E-10 ) {
1297 std::cerr << "Track::checkConsistency(): segment lengths of reference states for rep " << (*fi)->getRep() << " (id " << getIdForRep((*fi)->getRep()) << ") at TrackPoint " << (*tp) << " don't match" << std::endl;
1298 std::cerr << prevLen << " + " << len << " = " << prevLen + len << std::endl;
1299 std::cerr << "TrackPoint " << *tp << ", FitterInfo " << *fi << ", rep " << getIdForRep((*fi)->getRep()) << std::endl;
1300 retVal = false;
1301 }
1302 }
1303
1304 prevFis[(*fi)->getRep()] = static_cast<KalmanFitterInfo*>(*fi);
1305 }
1306 else
1307 prevFis[(*fi)->getRep()] = NULL;
1308
1309 } // end loop over FitterInfos
1310
1311 } // end loop over TrackPoints
1312
1313
1314 // check trackPointsWithMeasurement_
1315 std::vector<TrackPoint*> trackPointsWithMeasurement;
1316
1317 for (std::vector<TrackPoint*>::const_iterator it = trackPoints_.begin(); it != trackPoints_.end(); ++it) {
1318 if ((*it)->hasRawMeasurements()) {
1319 trackPointsWithMeasurement.push_back(*it);
1320 }
1321 }
1322
1323 if (trackPointsWithMeasurement.size() != trackPointsWithMeasurement_.size()) {
1324 std::cerr << "Track::checkConsistency(): trackPointsWithMeasurement_ has incorrect size" << std::endl;
1325 retVal = false;
1326 }
1327
1328 for (unsigned int i = 0; i < trackPointsWithMeasurement.size(); ++i) {
1329 if (trackPointsWithMeasurement[i] != trackPointsWithMeasurement_[i]) {
1330 std::cerr << "Track::checkConsistency(): trackPointsWithMeasurement_ is not correct" << std::endl;
1331 std::cerr << "has id " << i << ", adress " << trackPointsWithMeasurement_[i] << std::endl;
1332 std::cerr << "should have id " << i << ", adress " << trackPointsWithMeasurement[i] << std::endl;
1333 retVal = false;
1334 }
1335 }
1336
1337 return retVal;
1338}
1339
1340
1342
1343 #ifdef DEBUG
1344 std::cout << "Track::trackHasChanged \n";
1345 #endif
1346
1347 if (fitStatuses_.empty())
1348 return;
1349
1350 for (std::map< const AbsTrackRep*, FitStatus* >::const_iterator it=fitStatuses_.begin(); it!=fitStatuses_.end(); ++it) {
1351 it->second->setHasTrackChanged();
1352 }
1353}
1354
1355
1359
1360 for (std::vector<TrackPoint*>::const_iterator it = trackPoints_.begin(); it != trackPoints_.end(); ++it) {
1361 if ((*it)->hasRawMeasurements()) {
1362 trackPointsWithMeasurement_.push_back(*it);
1363 }
1364 }
1365}
1366
1367
1368void Track::Streamer(TBuffer &R__b)
1369{
1370 // Stream an object of class genfit::Track.
1371 const bool streamTrackPoints = true; // debugging option
1372 //This works around a msvc bug and should be harmless on other platforms
1373 typedef ::genfit::Track thisClass;
1374 UInt_t R__s, R__c;
1375 if (R__b.IsReading()) {
1376 this->Clear();
1377 Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
1378 TObject::Streamer(R__b);
1379 {
1380 std::vector<AbsTrackRep*> &R__stl = trackReps_;
1381 R__stl.clear();
1382 TClass *R__tcl1 = TBuffer::GetClass(typeid(genfit::AbsTrackRep));
1383 if (R__tcl1==0) {
1384 Error("trackReps_ streamer","Missing the TClass object for genfit::AbsTrackRep!");
1385 return;
1386 }
1387 int R__i, R__n;
1388 R__b >> R__n;
1389 R__stl.reserve(R__n);
1390 for (R__i = 0; R__i < R__n; R__i++) {
1391 genfit::AbsTrackRep* R__t;
1392 R__b >> R__t;
1393 R__stl.push_back(R__t);
1394 }
1395 }
1396 R__b >> cardinalRep_;
1397 if (streamTrackPoints)
1398 {
1399 std::vector<TrackPoint*> &R__stl = trackPoints_;
1400 R__stl.clear();
1401 TClass *R__tcl1 = TBuffer::GetClass(typeid(genfit::TrackPoint));
1402 if (R__tcl1==0) {
1403 Error("trackPoints_ streamer","Missing the TClass object for genfit::TrackPoint!");
1404 return;
1405 }
1406 int R__i, R__n;
1407 R__b >> R__n;
1408 R__stl.reserve(R__n);
1409 for (R__i = 0; R__i < R__n; R__i++) {
1410 genfit::TrackPoint* R__t;
1411 R__t = (genfit::TrackPoint*)R__b.ReadObjectAny(R__tcl1);
1412 R__t->setTrack(this);
1413 R__t->fixupRepsForReading();
1414 R__stl.push_back(R__t);
1415 }
1416 }
1417 {
1418 std::map<const AbsTrackRep*,FitStatus*> &R__stl = fitStatuses_;
1419 R__stl.clear();
1420 TClass *R__tcl2 = TBuffer::GetClass(typeid(genfit::FitStatus));
1421 if (R__tcl2==0) {
1422 Error("fitStatuses_ streamer","Missing the TClass object for genfit::FitStatus!");
1423 return;
1424 }
1425 int R__i, R__n;
1426 R__b >> R__n;
1427 for (R__i = 0; R__i < R__n; R__i++) {
1428 int id;
1429 R__b >> id;
1430 genfit::FitStatus* R__t2;
1431 R__t2 = (genfit::FitStatus*)R__b.ReadObjectAny(R__tcl2);
1432
1433 R__stl[this->getTrackRep(id)] = R__t2;
1434 }
1435 }
1436 stateSeed_.Streamer(R__b);
1437 covSeed_.Streamer(R__b);
1438 R__b.CheckByteCount(R__s, R__c, thisClass::IsA());
1439
1441 } else {
1442 R__c = R__b.WriteVersion(thisClass::IsA(), kTRUE);
1443 TObject::Streamer(R__b);
1444 {
1445 std::vector<AbsTrackRep*> &R__stl = trackReps_;
1446 int R__n=int(R__stl.size());
1447 R__b << R__n;
1448 if(R__n) {
1449 TClass *R__tcl1 = TBuffer::GetClass(typeid(genfit::AbsTrackRep));
1450 if (R__tcl1==0) {
1451 Error("trackReps_ streamer","Missing the TClass object for genfit::AbsTrackRep!");
1452 return;
1453 }
1454 std::vector<AbsTrackRep*>::iterator R__k;
1455 for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1456 R__b << *R__k;
1457 }
1458 }
1459 }
1460 R__b << cardinalRep_;
1461 if (streamTrackPoints)
1462 {
1463 std::vector<TrackPoint*> &R__stl = trackPoints_;
1464 int R__n=int(R__stl.size());
1465 R__b << R__n;
1466 if(R__n) {
1467 std::vector<TrackPoint*>::iterator R__k;
1468 for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1469 R__b << (*R__k);
1470 }
1471 }
1472 }
1473 {
1474 std::map<const AbsTrackRep*,FitStatus*> &R__stl = fitStatuses_;
1475 int R__n=int(R__stl.size());
1476 R__b << R__n;
1477 if(R__n) {
1478 TClass *R__tcl1 = TBuffer::GetClass(typeid(genfit::AbsTrackRep));
1479 if (R__tcl1==0) {
1480 Error("fitStatuses_ streamer","Missing the TClass object for genfit::AbsTrackRep!");
1481 return;
1482 }
1483 std::map<const AbsTrackRep*,FitStatus*>::iterator R__k;
1484/* this piece does not work for me TR July 2014
1485 for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1486 int id = this->getIdForRep((*R__k).first);
1487 R__b << id;
1488 R__b << ((*R__k).second);
1489 }
1490 use instead */
1491 int n=0;
1492 for (R__k = R__stl.begin(); n < 1; ++R__k){
1493 #ifdef DEBUG
1494 std::cout << "i am here " << (*R__k).first << " "<<n<< " "<<(*R__k).second << std::endl;
1495 #endif
1496 n+=1;
1497 if (n<100){
1498 int id = this->getIdForRep((*R__k).first);
1499 R__b << id;
1500 R__b << ((*R__k).second);
1501 }
1502 }
1503 }
1504 }
1505 stateSeed_.Streamer(R__b);
1506 covSeed_.Streamer(R__b);
1507 R__b.SetByteCount(R__c, kTRUE);
1508 }
1509}
1510
1511
1512
1513} /* End of namespace genfit */
Double_t m
This class collects all information needed and produced by a specific AbsFitter and is specific to on...
virtual bool hasForwardUpdate() const =0
virtual bool hasForwardPrediction() const =0
virtual bool hasBackwardUpdate() const =0
virtual const MeasuredStateOnPlane & getFittedState(bool biased=true) const =0
virtual bool hasReferenceState() const =0
virtual void Print(const Option_t *="") const
virtual bool hasBackwardPrediction() const =0
virtual bool hasMeasurements() const =0
Contains the measurement and covariance in raw detector coordinates.
Abstract base class for a track representation.
Definition AbsTrackRep.h:66
virtual double getTOF() const =0
Get the time of flight [ns] of the last extrapolation.
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,...
bool switchPDGSign()
try to multiply pdg code with -1. (Switch from particle to anti-particle and vice versa).
virtual AbsTrackRep * clone() const =0
Clone the trackRep.
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition Exception.h:48
void setFatal(bool b=true)
Set fatal flag.
Definition Exception.h:61
Class where important numbers and properties of a fit can be stored.
Definition FitStatus.h:36
FitStatus for use with AbsKalmanFitter implementations.
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
ReferenceStateOnPlane * getReferenceState() const
void fixWeights(bool arg=true)
StateOnPlane with additional covariance matrix.
Factory object to create AbsMeasurement objects from digitized and clustered data.
Measurement class implementing a planar hit geometry (1 or 2D).
A state with arbitrary dimension defined in a DetPlane.
double getCharge() const
TVectorD get6DState() const
const SharedPlanePtr & getPlane() const
Track candidate – seed values and indices.
Definition TrackCand.h:69
const TVectorD & getStateSeed() const
Returns the 6D seed state; should be in global coordinates.
Definition TrackCand.h:134
const TMatrixDSym & getCovSeed() const
get the covariance matrix seed (6D).
Definition TrackCand.h:131
Helper class for TrackPoint sorting, used in Track::sort().
Definition Track.h:45
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition TrackPoint.h:50
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=NULL) const
Get fitterInfo for rep. Per default, use cardinal rep.
bool hasRawMeasurements() const
Definition TrackPoint.h:96
void setTrack(Track *track)
Definition TrackPoint.h:91
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
Definition Track.h:71
TrackPoint * getPoint(int id) const
Definition Track.cc:201
void addTrackRep(AbsTrackRep *trackRep)
Definition Track.cc:506
std::map< const AbsTrackRep *, FitStatus * > fitStatuses_
helper
Definition Track.h:298
void setCardinalRep(int id)
Definition Track.cc:539
AbsTrackRep * getTrackRep(int id) const
Definition Track.h:127
void deleteFitterInfo(int startId=0, int endId=-1, const AbsTrackRep *rep=NULL)
Definition Track.cc:801
void deleteTrackRep(int id)
Delete a AbsTrackRep and all corresponding AbsFitterInfo objects in every TrackPoint.
Definition Track.cc:512
void deletePoint(int id)
Definition Track.cc:417
void reverseTrackPoints()
Flip the ordering of the TrackPoints.
Definition Track.cc:649
bool udpateSeed(int id=0, AbsTrackRep *rep=NULL, bool biased=true)
Definition Track.cc:627
void reverseMomSeed()
Flip direction of momentum seed.
Definition Track.h:225
double getTOF(AbsTrackRep *rep=NULL, int startId=0, int endId=-1) const
get time of flight in ns between to trackPoints (if NULL, for cardinal rep)
Definition Track.cc:875
std::vector< AbsTrackRep * > trackReps_
Definition Track.h:292
unsigned int getNumPoints() const
Definition Track.h:108
int getIdForRep(const AbsTrackRep *rep) const
This is used when streaming TrackPoints.
Definition Track.cc:245
virtual ~Track()
Definition Track.cc:171
void insertMeasurement(AbsMeasurement *measurement, int id=-1)
Creates a new TrackPoint containing the measurement, and adds it to the track.
Definition Track.cc:452
void deleteForwardInfo(int startId=0, int endId=-1, const AbsTrackRep *rep=NULL)
Definition Track.cc:681
TrackPoint * getPointWithMeasurement(int id) const
Definition Track.cc:209
void prune(const Option_t *="CFLWRMIU")
Delete unneeded information from the Track.
Definition Track.cc:952
KalmanFitStatus * getKalmanFitStatus(const AbsTrackRep *rep=NULL) const
If FitStatus is a KalmanFitStatus, return it. Otherwise return NULL.
Definition Track.cc:280
void setCovSeed(const TMatrixDSym &c)
Definition Track.h:163
const MeasuredStateOnPlane & getFittedState(int id=0, const AbsTrackRep *rep=NULL, bool biased=true) const
Shortcut to get FittedStates.
Definition Track.cc:231
TVectorD stateSeed_
Definition Track.h:301
virtual void Clear(Option_t *="")
Definition Track.cc:176
void deleteBackwardInfo(int startId=0, int endId=-1, const AbsTrackRep *rep=NULL)
Definition Track.cc:710
Track & operator=(Track)
Definition Track.cc:145
void insertPoint(TrackPoint *point, int id=-1)
Insert TrackPoint BEFORE TrackPoint with position id, if id >= 0.
Definition Track.cc:307
void insertPoints(std::vector< genfit::TrackPoint * > points, int id=-1)
Insert TrackPoints BEFORE TrackPoint with position id, if id >= 0.
Definition Track.cc:367
void fixWeights(AbsTrackRep *rep=NULL, int startId=0, int endId=-1)
Definition Track.cc:919
std::vector< TrackPoint * > trackPoints_
Definition Track.h:295
FitStatus * getFitStatus(const AbsTrackRep *rep=NULL) const
Get FitStatus for a AbsTrackRep. Per default, return FitStatus for cardinalRep.
Definition Track.h:149
void Print(const Option_t *="") const
Definition Track.cc:1026
AbsTrackRep * getCardinalRep() const
Get cardinal track representation.
Definition Track.h:140
void deleteReferenceInfo(int startId=0, int endId=-1, const AbsTrackRep *rep=NULL)
Definition Track.cc:741
void deleteMeasurementInfo(int startId=0, int endId=-1, const AbsTrackRep *rep=NULL)
Definition Track.cc:771
bool hasKalmanFitStatus(const AbsTrackRep *rep=NULL) const
Check if track has a KalmanFitStatus for given AbsTrackRep. Per default, check for cardinal rep.
Definition Track.cc:269
bool checkConsistency() const
Definition Track.cc:1192
void mergeTrack(const Track *other, int id=-1)
Merge two tracks.
Definition Track.cc:457
bool hasFitStatus(const AbsTrackRep *rep=NULL) const
Check if track has a FitStatus for given AbsTrackRep. Per default, check for cardinal rep.
Definition Track.cc:258
TrackPoint * getPointWithMeasurementAndFitterInfo(int id, const AbsTrackRep *rep) const
Definition Track.cc:217
void determineCardinalRep()
See with which AbsTrackRep the track was fitted best (converged fit w/ smallest chi2) and set the car...
Definition Track.cc:553
bool sort()
Sort TrackPoint and according to their sorting parameters.
Definition Track.cc:578
std::vector< TrackPoint * > trackPointsWithMeasurement_
Definition Track.h:296
void trackHasChanged()
Definition Track.cc:1341
void switchPDGSigns(AbsTrackRep *rep=NULL)
Switch the pdg signs of specified rep (of all reps if rep == NULL).
Definition Track.cc:661
unsigned int cardinalRep_
Definition Track.h:293
TMatrixDSym covSeed_
Definition Track.h:302
void swap(Track &other)
Definition Track.cc:160
void fillPointsWithMeasurement()
Definition Track.cc:1356
void reverseTrack()
Make track ready to be fitted in reverse direction.
Definition Track.cc:673
double getTrackLen(AbsTrackRep *rep=NULL, int startId=0, int endId=-1) const
get TrackLength between to trackPoints (if NULL, for cardinal rep)
Definition Track.cc:832
void setFitStatus(FitStatus *fitStatus, const AbsTrackRep *rep)
Definition Track.cc:285
void setStateSeed(const TVectorD &s)
Definition Track.h:159
Matrix inversion tools.
Definition AbsBField.h:29