CMS 3D CMS Logo

HitPixelLayersTPSelector.h
Go to the documentation of this file.
1 #ifndef HitPixelLayersTrackSelection_h
2 #define HitPixelLayersTrackSelection_h
3 
6 
13 
22 public:
23  // input collection type
25 
26  // output collection type
28 
29  // iterator over result collection type.
31 
32  // constructor from parameter set configurability
34  : tripletSeedOnly_(iConfig.getParameter<bool>("tripletSeedOnly")),
35  ptMin_(iConfig.getParameter<double>("ptMin")),
36  minRapidity_(iConfig.getParameter<double>("minRapidity")),
37  maxRapidity_(iConfig.getParameter<double>("maxRapidity")),
38  tip_(iConfig.getParameter<double>("tip")),
39  lip_(iConfig.getParameter<double>("lip")),
40  minHit_(iConfig.getParameter<int>("minHit")),
41  signalOnly_(iConfig.getParameter<bool>("signalOnly")),
42  chargedOnly_(iConfig.getParameter<bool>("chargedOnly")),
43  primaryOnly_(iConfig.getParameter<bool>("primaryOnly")),
44  tpStatusBased_(iConfig.getParameter<bool>("tpStatusBased")),
45  pdgId_(iConfig.getParameter<std::vector<int> >("pdgId")){};
46 
47  // select object from a collection and
48  // possibly event content
49  void select(const edm::Handle<collection>& TPCH, const edm::Event& iEvent, const edm::EventSetup& iSetup) {
50  selected_.clear();
51  //Retrieve tracker topology from geometry
53  iSetup.get<TrackerTopologyRcd>().get(tTopoHand);
54  const TrackerTopology* tTopo = tTopoHand.product();
55 
56  const collection& tpc = *(TPCH.product());
57 
58  for (TrackingParticleCollection::size_type i = 0; i < tpc.size(); i++) {
59  TrackingParticleRef tpr(TPCH, i);
60 
61  // quickly reject if it is from pile-up
62  if (signalOnly_ && !(tpr->eventId().bunchCrossing() == 0 && tpr->eventId().event() == 0))
63  continue;
64  if (chargedOnly_ && tpr->charge() == 0)
65  continue; //select only if charge!=0
66  if (tpStatusBased_ && primaryOnly_ && tpr->status() != 1)
67  continue; // TP status based sel primary
68  if ((!tpStatusBased_) && primaryOnly_ && tpr->parentVertex()->nSourceTracks() != 0)
69  continue; // vertex based sel for primary
70 
71  // loop over specified PID values
72  bool testId = false;
73  unsigned int idSize = pdgId_.size();
74  if (idSize == 0)
75  testId = true;
76  else
77  for (unsigned int it = 0; it != idSize; ++it) {
78  if (tpr->pdgId() == pdgId_[it])
79  testId = true;
80  }
81 
82  // selection criteria
83  if (tpr->numberOfTrackerLayers() >= minHit_ && sqrt(tpr->momentum().perp2()) >= ptMin_ &&
84  tpr->momentum().eta() >= minRapidity_ && tpr->momentum().eta() <= maxRapidity_ &&
85  sqrt(tpr->vertex().perp2()) <= tip_ && fabs(tpr->vertex().z()) <= lip_ && testId) {
86  if (tripletSeedOnly_ && !goodHitPattern(pixelHitPattern(tpr, tTopo)))
87  continue; //findable triplet seed
88  selected_.push_back(tpr);
89  }
90  }
91  }
92 
93  // return pixel layer hit pattern
94  std::vector<bool> pixelHitPattern(const TrackingParticleRef& simTrack, const TrackerTopology* tTopo) {
95  std::vector<bool> hitpattern(5, false); // PXB 0,1,2 PXF 0,1
96  // This currently will always return false, since we can no loger use the sim hits to check for triplets. This would need to be fixed if we want to enable this feature, but it's not being used at the moment, since tripletSeedOnly is always set to False - Matt Nguyen, 24/7/2013
97 
98  return hitpattern;
99  }
100 
101  // test whether hit pattern would give a pixel triplet seed
102  bool goodHitPattern(const std::vector<bool>& hitpattern) {
103  if ((hitpattern[0] && hitpattern[1] && hitpattern[2]) || (hitpattern[0] && hitpattern[1] && hitpattern[3]) ||
104  (hitpattern[0] && hitpattern[3] && hitpattern[4]))
105  return true;
106  else
107  return false;
108  }
109 
110  // iterators over selected objects: collection begin
111  const_iterator begin() const { return selected_.begin(); }
112 
113  // iterators over selected objects: collection end
114  const_iterator end() const { return selected_.end(); }
115 
116  // true if no object has been selected
117  size_t size() const { return selected_.size(); }
118 
119  //private:
120 
121  container selected_;
123  double ptMin_;
124  double minRapidity_;
125  double maxRapidity_;
126  double tip_;
127  double lip_;
128  int minHit_;
133  std::vector<int> pdgId_;
134 };
135 
136 #endif
std::vector< TrackingParticle > TrackingParticleCollection
std::vector< bool > pixelHitPattern(const TrackingParticleRef &simTrack, const TrackerTopology *tTopo)
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
uint16_t size_type
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
int iEvent
Definition: GenABIO.cc:224
TrackingParticleCollection collection
void select(const edm::Handle< collection > &TPCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
TrackingParticleRefVector container
T sqrt(T t)
Definition: SSEVec.h:19
container::const_iterator const_iterator
T const * product() const
Definition: Handle.h:69
void clear()
Clear the vector.
Definition: RefVector.h:142
bool goodHitPattern(const std::vector< bool > &hitpattern)
HitPixelLayersTPSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
T get() const
Definition: EventSetup.h:73
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
T const * product() const
Definition: ESHandle.h:86
const_iterator begin() const