CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Public Attributes
HitPixelLayersTPSelector Class Reference

#include <HitPixelLayersTPSelector.h>

Public Types

typedef TrackingParticleCollection collection
 
typedef container::const_iterator const_iterator
 
typedef TrackingParticleRefVector container
 

Public Member Functions

const_iterator begin () const
 
const_iterator end () const
 
bool goodHitPattern (const std::vector< bool > &hitpattern)
 
 HitPixelLayersTPSelector (const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
 
std::vector< bool > pixelHitPattern (const TrackingParticleRef &simTrack, const TrackerTopology *tTopo)
 
void select (const edm::Handle< collection > &TPCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
size_t size () const
 

Public Attributes

bool chargedOnly_
 
double lip_
 
double maxRapidity_
 
int minHit_
 
double minRapidity_
 
std::vector< int > pdgId_
 
bool primaryOnly_
 
double ptMin_
 
container selected_
 
bool signalOnly_
 
double tip_
 
bool tpStatusBased_
 
bool tripletSeedOnly_
 

Detailed Description

Selector to select only tracking particles that leave hits in three pixel layers Additional selection done on pt, rapidity, impact parameter, min hits, pdg id, etc.

Inspired by CommonTools.RecoAlgos.TrackingParticleSelector.h

Definition at line 21 of file HitPixelLayersTPSelector.h.

Member Typedef Documentation

◆ collection

Definition at line 24 of file HitPixelLayersTPSelector.h.

◆ const_iterator

Definition at line 30 of file HitPixelLayersTPSelector.h.

◆ container

Definition at line 27 of file HitPixelLayersTPSelector.h.

Constructor & Destructor Documentation

◆ HitPixelLayersTPSelector()

HitPixelLayersTPSelector::HitPixelLayersTPSelector ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 33 of file HitPixelLayersTPSelector.h.

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")){};

Member Function Documentation

◆ begin()

const_iterator HitPixelLayersTPSelector::begin ( void  ) const
inline

Definition at line 111 of file HitPixelLayersTPSelector.h.

111 { return selected_.begin(); }

References edm::RefVector< C, T, F >::begin(), and selected_.

◆ end()

const_iterator HitPixelLayersTPSelector::end ( void  ) const
inline

Definition at line 114 of file HitPixelLayersTPSelector.h.

114 { return selected_.end(); }

References edm::RefVector< C, T, F >::end(), and selected_.

◆ goodHitPattern()

bool HitPixelLayersTPSelector::goodHitPattern ( const std::vector< bool > &  hitpattern)
inline

Definition at line 102 of file HitPixelLayersTPSelector.h.

102  {
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  }

Referenced by select().

◆ pixelHitPattern()

std::vector<bool> HitPixelLayersTPSelector::pixelHitPattern ( const TrackingParticleRef simTrack,
const TrackerTopology tTopo 
)
inline

Definition at line 94 of file HitPixelLayersTPSelector.h.

94  {
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  }

Referenced by select().

◆ select()

void HitPixelLayersTPSelector::select ( const edm::Handle< collection > &  TPCH,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)
inline

Definition at line 49 of file HitPixelLayersTPSelector.h.

49  {
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  }

References chargedOnly_, edm::RefVector< C, T, F >::clear(), edm::EventSetup::get(), get, goodHitPattern(), mps_fire::i, lip_, maxRapidity_, minHit_, minRapidity_, pdgId_, pixelHitPattern(), primaryOnly_, edm::Handle< T >::product(), edm::ESHandle< T >::product(), ptMin_, edm::RefVector< C, T, F >::push_back(), selected_, signalOnly_, mathSSE::sqrt(), tip_, tpStatusBased_, and tripletSeedOnly_.

◆ size()

size_t HitPixelLayersTPSelector::size ( void  ) const
inline

Member Data Documentation

◆ chargedOnly_

bool HitPixelLayersTPSelector::chargedOnly_

Definition at line 130 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ lip_

double HitPixelLayersTPSelector::lip_

Definition at line 127 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ maxRapidity_

double HitPixelLayersTPSelector::maxRapidity_

Definition at line 125 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ minHit_

int HitPixelLayersTPSelector::minHit_

Definition at line 128 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ minRapidity_

double HitPixelLayersTPSelector::minRapidity_

Definition at line 124 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ pdgId_

std::vector<int> HitPixelLayersTPSelector::pdgId_

Definition at line 133 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ primaryOnly_

bool HitPixelLayersTPSelector::primaryOnly_

Definition at line 131 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ ptMin_

double HitPixelLayersTPSelector::ptMin_

Definition at line 123 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ selected_

container HitPixelLayersTPSelector::selected_

Definition at line 121 of file HitPixelLayersTPSelector.h.

Referenced by begin(), end(), select(), and size().

◆ signalOnly_

bool HitPixelLayersTPSelector::signalOnly_

Definition at line 129 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tip_

double HitPixelLayersTPSelector::tip_

Definition at line 126 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tpStatusBased_

bool HitPixelLayersTPSelector::tpStatusBased_

Definition at line 132 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tripletSeedOnly_

bool HitPixelLayersTPSelector::tripletSeedOnly_

Definition at line 122 of file HitPixelLayersTPSelector.h.

Referenced by select().

HitPixelLayersTPSelector::lip_
double lip_
Definition: HitPixelLayersTPSelector.h:127
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
edm::RefVector::clear
void clear()
Clear the vector.
Definition: RefVector.h:142
mps_fire.i
i
Definition: mps_fire.py:428
edm::Handle::product
T const * product() const
Definition: Handle.h:70
HitPixelLayersTPSelector::minRapidity_
double minRapidity_
Definition: HitPixelLayersTPSelector.h:124
HitPixelLayersTPSelector::ptMin_
double ptMin_
Definition: HitPixelLayersTPSelector.h:123
TrackerTopology
Definition: TrackerTopology.h:16
edm::RefVector::begin
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
HitPixelLayersTPSelector::primaryOnly_
bool primaryOnly_
Definition: HitPixelLayersTPSelector.h:131
HitPixelLayersTPSelector::tpStatusBased_
bool tpStatusBased_
Definition: HitPixelLayersTPSelector.h:132
HitPixelLayersTPSelector::chargedOnly_
bool chargedOnly_
Definition: HitPixelLayersTPSelector.h:130
edm::Ref< TrackingParticleCollection >
edm::RefVector::end
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
HitPixelLayersTPSelector::pixelHitPattern
std::vector< bool > pixelHitPattern(const TrackingParticleRef &simTrack, const TrackerTopology *tTopo)
Definition: HitPixelLayersTPSelector.h:94
trigger::size_type
uint16_t size_type
Definition: TriggerTypeDefs.h:18
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
HitPixelLayersTPSelector::maxRapidity_
double maxRapidity_
Definition: HitPixelLayersTPSelector.h:125
edm::ESHandle< TrackerTopology >
HitPixelLayersTPSelector::tripletSeedOnly_
bool tripletSeedOnly_
Definition: HitPixelLayersTPSelector.h:122
HitPixelLayersTPSelector::collection
TrackingParticleCollection collection
Definition: HitPixelLayersTPSelector.h:24
HitPixelLayersTPSelector::goodHitPattern
bool goodHitPattern(const std::vector< bool > &hitpattern)
Definition: HitPixelLayersTPSelector.h:102
get
#define get
edm::RefVector::push_back
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
HitPixelLayersTPSelector::tip_
double tip_
Definition: HitPixelLayersTPSelector.h:126
HitPixelLayersTPSelector::selected_
container selected_
Definition: HitPixelLayersTPSelector.h:121
HitPixelLayersTPSelector::pdgId_
std::vector< int > pdgId_
Definition: HitPixelLayersTPSelector.h:133
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
TrackerTopologyRcd
Definition: TrackerTopologyRcd.h:10
edm::RefVector::size
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
HitPixelLayersTPSelector::signalOnly_
bool signalOnly_
Definition: HitPixelLayersTPSelector.h:129
HitPixelLayersTPSelector::minHit_
int minHit_
Definition: HitPixelLayersTPSelector.h:128