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_
 
edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtTopoToken_
 

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")),
46  tTopoToken_(iC.esConsumes()){};
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_

Member Function Documentation

◆ begin()

const_iterator HitPixelLayersTPSelector::begin ( void  ) const
inline

Definition at line 110 of file HitPixelLayersTPSelector.h.

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

110 { return selected_.begin(); }
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ end()

const_iterator HitPixelLayersTPSelector::end ( void  ) const
inline

Definition at line 113 of file HitPixelLayersTPSelector.h.

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

113 { return selected_.end(); }
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228

◆ goodHitPattern()

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

Definition at line 101 of file HitPixelLayersTPSelector.h.

Referenced by select().

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

◆ pixelHitPattern()

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

Definition at line 93 of file HitPixelLayersTPSelector.h.

Referenced by select().

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

◆ select()

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

Definition at line 50 of file HitPixelLayersTPSelector.h.

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

50  {
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  }
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
TrackingParticleCollection collection
T sqrt(T t)
Definition: SSEVec.h:23
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
void clear()
Clear the vector.
Definition: RefVector.h:142
bool goodHitPattern(const std::vector< bool > &hitpattern)
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67

◆ size()

size_t HitPixelLayersTPSelector::size ( void  ) const
inline

Definition at line 116 of file HitPixelLayersTPSelector.h.

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

Referenced by ntupleDataFormat._Collection::__iter__(), and ntupleDataFormat._Collection::__len__().

116 { return selected_.size(); }
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102

Member Data Documentation

◆ chargedOnly_

bool HitPixelLayersTPSelector::chargedOnly_

Definition at line 129 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ lip_

double HitPixelLayersTPSelector::lip_

Definition at line 126 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ maxRapidity_

double HitPixelLayersTPSelector::maxRapidity_

Definition at line 124 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ minHit_

int HitPixelLayersTPSelector::minHit_

Definition at line 127 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ minRapidity_

double HitPixelLayersTPSelector::minRapidity_

Definition at line 123 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ pdgId_

std::vector<int> HitPixelLayersTPSelector::pdgId_

Definition at line 132 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ primaryOnly_

bool HitPixelLayersTPSelector::primaryOnly_

Definition at line 130 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ ptMin_

double HitPixelLayersTPSelector::ptMin_

Definition at line 122 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ selected_

container HitPixelLayersTPSelector::selected_

Definition at line 120 of file HitPixelLayersTPSelector.h.

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

◆ signalOnly_

bool HitPixelLayersTPSelector::signalOnly_

Definition at line 128 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tip_

double HitPixelLayersTPSelector::tip_

Definition at line 125 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tpStatusBased_

bool HitPixelLayersTPSelector::tpStatusBased_

Definition at line 131 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tripletSeedOnly_

bool HitPixelLayersTPSelector::tripletSeedOnly_

Definition at line 121 of file HitPixelLayersTPSelector.h.

Referenced by select().

◆ tTopoToken_

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> HitPixelLayersTPSelector::tTopoToken_

Definition at line 133 of file HitPixelLayersTPSelector.h.

Referenced by select().