Class Reference for E1039 Core & Analysis Software
UtilTrack.cc
Go to the documentation of this file.
1 #include <cassert>
2 #include <geom_svc/GeomSvc.h>
5 #include "UtilTrack.h"
6 using namespace std;
7 
9 
13 SQTrack* UtilTrack::FindTrackByID(const SQTrackVector* vec, const int id_trk, const bool do_assert)
14 {
15  for (SQTrackVector::ConstIter it = vec->begin(); it != vec->end(); it++) {
16  SQTrack* trk = *it;
17  if (trk->get_track_id() == id_trk) {
18  if (do_assert) assert(trk);
19  return trk;
20  }
21  }
22  return 0;
23 }
24 
26 
32 SQHitVector* UtilTrack::FindHitsOfTrack(const SQHitVector* vec_in, const int id_trk)
33 {
34  SQHitVector* vec = vec_in->Clone();
35  vec->clear();
36  for (SQHitVector::ConstIter it = vec_in->begin(); it != vec_in->end(); it++) {
37  SQHit* hit = *it;
38  if (hit->get_track_id() == id_trk) vec->push_back(hit);
39  }
40  return vec;
41 }
42 
44 
52 SQHitVector* UtilTrack::FindDetectorHitsOfTrack(const SQHitVector* vec_in, const int id_trk, const std::string det_name)
53 {
54  SQHitVector* vec = vec_in->Clone();
55  vec->clear();
56  for (SQHitVector::ConstIter it = vec_in->begin(); it != vec_in->end(); it++) {
57  SQHit* hit = *it;
58  if (hit->get_track_id() != id_trk) continue;
59  string name = GeomSvc::instance()->getDetectorName( hit->get_detector_id() );
60  if (name.substr(0, det_name.size()) != det_name) continue;
61  vec->push_back(hit);
62  }
63  return vec;
64 }
65 
67 SQHitVector* UtilTrack::FindDetectorHitsOfTrack(const SQHitVector* vec_in, const int id_trk, const char* det_name)
68 {
69  return FindDetectorHitsOfTrack(vec_in, id_trk, (const string)det_name);
70 }
71 
73 
79 SQHitVector* UtilTrack::FindHodoHitsOfTrack(const SQHitVector* vec_in, const int id_trk)
80 {
81  return FindDetectorHitsOfTrack(vec_in, id_trk, "H");
82 }
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
std::string getDetectorName(const int &detectorID) const
Definition: GeomSvc.h:223
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
virtual SQHitVector * Clone() const =0
std::vector< SQHit * >::const_iterator ConstIter
Definition: SQHitVector.h:37
virtual ConstIter end() const =0
virtual ConstIter begin() const =0
virtual void push_back(const SQHit *hit)=0
virtual void clear()=0
An SQ interface class to hold one detector hit.
Definition: SQHit.h:20
virtual short get_detector_id() const
Return the detector ID of this hit.
Definition: SQHit.h:42
virtual int get_track_id() const
Return the track ID associated with this hit. Probably the value is not properly set at present.
Definition: SQHit.h:66
An SQ interface class to hold a list of SQTrack objects.
Definition: SQTrackVector.h:19
virtual ConstIter begin() const =0
virtual ConstIter end() const =0
std::vector< SQTrack * >::const_iterator ConstIter
Definition: SQTrackVector.h:22
An SQ interface class to hold one true or reconstructed track.
Definition: SQTrack.h:8
virtual int get_track_id() const =0
Return the track ID, which is unique per event(?).
SQHitVector * FindDetectorHitsOfTrack(const SQHitVector *vec_in, const int id_trk, const std::string det_name)
Find track-associated hits whose detector name starts with 'det_name'.
Definition: UtilTrack.cc:52
SQHitVector * FindHodoHitsOfTrack(const SQHitVector *vec_in, const int id_trk)
Find all hodoscope hits hits associated with the given track.
Definition: UtilTrack.cc:79
SQHitVector * FindHitsOfTrack(const SQHitVector *vec_in, const int id_trk)
Find all hits associated with the given track.
Definition: UtilTrack.cc:32
SQTrack * FindTrackByID(const SQTrackVector *vec, const int id_trk, const bool do_assert=false)
Find a track by track ID in the given track list.
Definition: UtilTrack.cc:13