CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HitPixelLayersTPSelector.h
Go to the documentation of this file.
1 #ifndef HitPixelLayersTrackSelection_h
2 #define HitPixelLayersTrackSelection_h
3 
6 
11 
20 {
21 
22  public:
23  // input collection type
25 
26 
27  // output collection type
28  typedef std::vector<const TrackingParticle*> container;
29 
30  // iterator over result collection type.
31  typedef container::const_iterator const_iterator;
32 
33  // constructor from parameter set configurability
35  tripletSeedOnly_(iConfig.getParameter<bool>("tripletSeedOnly")),
36  ptMin_(iConfig.getParameter<double>("ptMin")),
37  minRapidity_(iConfig.getParameter<double>("minRapidity")),
38  maxRapidity_(iConfig.getParameter<double>("maxRapidity")),
39  tip_(iConfig.getParameter<double>("tip")),
40  lip_(iConfig.getParameter<double>("lip")),
41  minHit_(iConfig.getParameter<int>("minHit")),
42  signalOnly_(iConfig.getParameter<bool>("signalOnly")),
43  chargedOnly_(iConfig.getParameter<bool>("chargedOnly")),
44  primaryOnly_(iConfig.getParameter<bool>("primaryOnly")),
45  tpStatusBased_(iConfig.getParameter<bool>("tpStatusBased")),
46  pdgId_(iConfig.getParameter< std::vector<int> >("pdgId"))
47  {};
48 
49  // select object from a collection and
50  // possibly event content
51  void select( const edm::Handle<collection> & TPCH, const edm::Event & iEvent, const edm::EventSetup & iSetup)
52  {
53  selected_.clear();
54  //Retrieve tracker topology from geometry
56  iSetup.get<TrackerTopologyRcd>().get(tTopoHand);
57  const TrackerTopology *tTopo=tTopoHand.product();
58 
59 
60  const collection & tpc = *(TPCH.product());
61 
62  for (TrackingParticleCollection::size_type i=0; i<tpc.size(); i++)
63  {
64  TrackingParticleRef tpr(TPCH, i);
65 
66  // quickly reject if it is from pile-up
67  if (signalOnly_ && !(tpr->eventId().bunchCrossing()==0 && tpr->eventId().event()==0) ) continue;
68  if (chargedOnly_ && tpr->charge()==0) continue; //select only if charge!=0
69  if (tpStatusBased_ && primaryOnly_ && tpr->status()!=1 ) continue; // TP status based sel primary
70  if ((!tpStatusBased_) && primaryOnly_ && tpr->parentVertex()->nSourceTracks()!=0 ) continue; // vertex based sel for primary
71 
72  // loop over specified PID values
73  bool testId = false;
74  unsigned int idSize = pdgId_.size();
75  if (idSize==0) testId = true;
76  else for (unsigned int it=0;it!=idSize;++it){
77  if (tpr->pdgId()==pdgId_[it]) testId = true;
78  }
79 
80  // selection criteria
81  if ( tpr->numberOfTrackerLayers() >= minHit_ &&
82  sqrt(tpr->momentum().perp2()) >= ptMin_ &&
83  tpr->momentum().eta() >= minRapidity_ && tpr->momentum().eta() <= maxRapidity_ &&
84  sqrt(tpr->vertex().perp2()) <= tip_ &&
85  fabs(tpr->vertex().z()) <= lip_ &&
86  testId)
87  {
88  if (tripletSeedOnly_ && !goodHitPattern(pixelHitPattern(tpr,tTopo)) ) continue; //findable triplet seed
89  const TrackingParticle * trap = &(tpc[i]);
90  selected_.push_back(trap);
91  }
92 
93  }
94  }
95 
96  // return pixel layer hit pattern
97  std::vector<bool> pixelHitPattern( const TrackingParticleRef& simTrack, const TrackerTopology *tTopo )
98  {
99  std::vector<bool> hitpattern(5,false); // PXB 0,1,2 PXF 0,1
100  // 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
101 
102  return hitpattern;
103  }
104 
105  // test whether hit pattern would give a pixel triplet seed
106  bool goodHitPattern( const std::vector<bool>& hitpattern )
107  {
108  if( (hitpattern[0] && hitpattern[1] && hitpattern[2]) ||
109  (hitpattern[0] && hitpattern[1] && hitpattern[3]) ||
110  (hitpattern[0] && hitpattern[3] && hitpattern[4]) )
111  return true;
112  else
113  return false;
114  }
115 
116  // iterators over selected objects: collection begin
118  {
119  return selected_.begin();
120  }
121 
122  // iterators over selected objects: collection end
124  {
125  return selected_.end();
126  }
127 
128  // true if no object has been selected
129  size_t size() const
130  {
131  return selected_.size();
132  }
133 
134  //private:
135 
138  double ptMin_;
139  double minRapidity_;
140  double maxRapidity_;
141  double tip_;
142  double lip_;
143  int minHit_;
148  std::vector<int> pdgId_;
149 
150 
151 };
152 
153 
154 #endif
int i
Definition: DBlmapReader.cc:9
std::vector< TrackingParticle > TrackingParticleCollection
std::vector< bool > pixelHitPattern(const TrackingParticleRef &simTrack, const TrackerTopology *tTopo)
std::vector< const TrackingParticle * > container
uint16_t size_type
int iEvent
Definition: GenABIO.cc:230
TrackingParticleCollection collection
void select(const edm::Handle< collection > &TPCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
T sqrt(T t)
Definition: SSEVec.h:48
container::const_iterator const_iterator
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
bool goodHitPattern(const std::vector< bool > &hitpattern)
HitPixelLayersTPSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Monte Carlo truth information used for tracking validation.
const_iterator begin() const