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  tTopoToken_(iC.esConsumes()){};
47 
48  // select object from a collection and
49  // possibly event content
50  void select(const edm::Handle<collection>& TPCH, const edm::Event& iEvent, const edm::EventSetup& iSetup) {
51  selected_.clear();
52  //Retrieve tracker topology from geometry
53  const TrackerTopology* tTopo = &iSetup.getData(tTopoToken_);
54 
55  const collection& tpc = *(TPCH.product());
56 
57  for (TrackingParticleCollection::size_type i = 0; i < tpc.size(); i++) {
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))
62  continue;
63  if (chargedOnly_ && tpr->charge() == 0)
64  continue; //select only if charge!=0
65  if (tpStatusBased_ && primaryOnly_ && tpr->status() != 1)
66  continue; // TP status based sel primary
67  if ((!tpStatusBased_) && primaryOnly_ && tpr->parentVertex()->nSourceTracks() != 0)
68  continue; // vertex based sel for primary
69 
70  // loop over specified PID values
71  bool testId = false;
72  unsigned int idSize = pdgId_.size();
73  if (idSize == 0)
74  testId = true;
75  else
76  for (unsigned int it = 0; it != idSize; ++it) {
77  if (tpr->pdgId() == pdgId_[it])
78  testId = true;
79  }
80 
81  // selection criteria
82  if (tpr->numberOfTrackerLayers() >= minHit_ && sqrt(tpr->momentum().perp2()) >= ptMin_ &&
83  tpr->momentum().eta() >= minRapidity_ && tpr->momentum().eta() <= maxRapidity_ &&
84  sqrt(tpr->vertex().perp2()) <= tip_ && fabs(tpr->vertex().z()) <= lip_ && testId) {
85  if (tripletSeedOnly_ && !goodHitPattern(pixelHitPattern(tpr, tTopo)))
86  continue; //findable triplet seed
87  selected_.push_back(tpr);
88  }
89  }
90  }
91 
92  // return pixel layer hit pattern
93  std::vector<bool> pixelHitPattern(const TrackingParticleRef& simTrack, const TrackerTopology* tTopo) {
94  std::vector<bool> hitpattern(5, false); // PXB 0,1,2 PXF 0,1
95  // 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
96 
97  return hitpattern;
98  }
99 
100  // test whether hit pattern would give a pixel triplet seed
101  bool goodHitPattern(const std::vector<bool>& hitpattern) {
102  if ((hitpattern[0] && hitpattern[1] && hitpattern[2]) || (hitpattern[0] && hitpattern[1] && hitpattern[3]) ||
103  (hitpattern[0] && hitpattern[3] && hitpattern[4]))
104  return true;
105  else
106  return false;
107  }
108 
109  // iterators over selected objects: collection begin
110  const_iterator begin() const { return selected_.begin(); }
111 
112  // iterators over selected objects: collection end
113  const_iterator end() const { return selected_.end(); }
114 
115  // true if no object has been selected
116  size_t size() const { return selected_.size(); }
117 
118  //private:
119 
122  double ptMin_;
123  double minRapidity_;
124  double maxRapidity_;
125  double tip_;
126  double lip_;
127  int minHit_;
132  std::vector<int> pdgId_;
134 };
135 
136 #endif
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
T const * product() const
Definition: Handle.h:70
std::vector< bool > pixelHitPattern(const TrackingParticleRef &simTrack, const TrackerTopology *tTopo)
uint16_t size_type
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
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
void clear()
Clear the vector.
Definition: RefVector.h:142
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
bool goodHitPattern(const std::vector< bool > &hitpattern)
HitPixelLayersTPSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223