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 
9 
18 {
19 
20  public:
21  // input collection type
23 
24 
25  // output collection type
26  typedef std::vector<const TrackingParticle*> container;
27 
28  // iterator over result collection type.
29  typedef container::const_iterator const_iterator;
30 
31  // constructor from parameter set configurability
33  tripletSeedOnly_(iConfig.getParameter<bool>("tripletSeedOnly")),
34  ptMin_(iConfig.getParameter<double>("ptMin")),
35  minRapidity_(iConfig.getParameter<double>("minRapidity")),
36  maxRapidity_(iConfig.getParameter<double>("maxRapidity")),
37  tip_(iConfig.getParameter<double>("tip")),
38  lip_(iConfig.getParameter<double>("lip")),
39  minHit_(iConfig.getParameter<int>("minHit")),
40  signalOnly_(iConfig.getParameter<bool>("signalOnly")),
41  chargedOnly_(iConfig.getParameter<bool>("chargedOnly")),
42  primaryOnly_(iConfig.getParameter<bool>("primaryOnly")),
43  tpStatusBased_(iConfig.getParameter<bool>("tpStatusBased")),
44  pdgId_(iConfig.getParameter< std::vector<int> >("pdgId"))
45  {};
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  {
51  selected_.clear();
52 
53 
54  const collection & tpc = *(TPCH.product());
55 
56  for (TrackingParticleCollection::size_type i=0; i<tpc.size(); i++)
57  {
58  TrackingParticleRef tpr(TPCH, i);
59 
60  // quickly reject if it is from pile-up
61  if (signalOnly_ && !(tpr->eventId().bunchCrossing()==0 && tpr->eventId().event()==0) ) continue;
62  if (chargedOnly_ && tpr->charge()==0) continue; //select only if charge!=0
63  if (tpStatusBased_ && primaryOnly_ && tpr->status()!=1 ) continue; // TP status based sel primary
64  if ((!tpStatusBased_) && primaryOnly_ && tpr->parentVertex()->nSourceTracks()!=0 ) continue; // vertex based sel for primary
65 
66  // loop over specified PID values
67  bool testId = false;
68  unsigned int idSize = pdgId_.size();
69  if (idSize==0) testId = true;
70  else for (unsigned int it=0;it!=idSize;++it){
71  if (tpr->pdgId()==pdgId_[it]) testId = true;
72  }
73 
74  // selection criteria
75  if ( tpr->matchedHit() >= minHit_ &&
76  sqrt(tpr->momentum().perp2()) >= ptMin_ &&
77  tpr->momentum().eta() >= minRapidity_ && tpr->momentum().eta() <= maxRapidity_ &&
78  sqrt(tpr->vertex().perp2()) <= tip_ &&
79  fabs(tpr->vertex().z()) <= lip_ &&
80  testId)
81  {
82  if (tripletSeedOnly_ && !goodHitPattern(pixelHitPattern(tpr)) ) continue; //findable triplet seed
83  const TrackingParticle * trap = &(tpc[i]);
84  selected_.push_back(trap);
85  }
86 
87  }
88  }
89 
90  // return pixel layer hit pattern
91  std::vector<bool> pixelHitPattern( const TrackingParticleRef& simTrack )
92  {
93  std::vector<bool> hitpattern(5,false); // PXB 0,1,2 PXF 0,1
94 
95  for(std::vector<PSimHit>::const_iterator simHit = simTrack->pSimHit_begin();simHit!= simTrack->pSimHit_end();simHit++){
96 
97  DetId id = DetId(simHit->detUnitId());
98  uint32_t detid = id.det();
99  uint32_t subdet = id.subdetId();
100 
101  if (detid == DetId::Tracker) {
102  if (subdet == PixelSubdetector::PixelBarrel)
103  hitpattern[PXBDetId(id).layer()-1]=true;
104  else if (subdet == PixelSubdetector::PixelEndcap)
105  hitpattern[PXFDetId(id).disk()+2]=true;
106  }
107 
108  }// end simhit loop
109 
110  return hitpattern;
111  }
112 
113  // test whether hit pattern would give a pixel triplet seed
114  bool goodHitPattern( const std::vector<bool> hitpattern )
115  {
116  if( (hitpattern[0] && hitpattern[1] && hitpattern[2]) ||
117  (hitpattern[0] && hitpattern[1] && hitpattern[3]) ||
118  (hitpattern[0] && hitpattern[3] && hitpattern[4]) )
119  return true;
120  else
121  return false;
122  }
123 
124  // iterators over selected objects: collection begin
126  {
127  return selected_.begin();
128  }
129 
130  // iterators over selected objects: collection end
132  {
133  return selected_.end();
134  }
135 
136  // true if no object has been selected
137  size_t size() const
138  {
139  return selected_.size();
140  }
141 
142  //private:
143 
146  double ptMin_;
147  double minRapidity_;
148  double maxRapidity_;
149  double tip_;
150  double lip_;
151  int minHit_;
156  std::vector<int> pdgId_;
157 
158 
159 };
160 
161 
162 #endif
int i
Definition: DBlmapReader.cc:9
std::vector< const TrackingParticle * > container
uint16_t size_type
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
HitPixelLayersTPSelector(const edm::ParameterSet &iConfig)
int iEvent
Definition: GenABIO.cc:243
TrackingParticleCollection collection
bool goodHitPattern(const std::vector< bool > hitpattern)
void select(const edm::Handle< collection > &TPCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
T sqrt(T t)
Definition: SSEVec.h:28
std::vector< bool > pixelHitPattern(const TrackingParticleRef &simTrack)
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
container::const_iterator const_iterator
Definition: DetId.h:20
T const * product() const
Definition: Handle.h:74
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator begin() const