CMS 3D CMS Logo

TrackingParticleSelector.h
Go to the documentation of this file.
1 #ifndef SimTracker_Common_TrackingParticleSelector_h
2 #define SimTracker_Common_TrackingParticleSelector_h
3 /* \class TrackingParticleSelector
4  *
5  * \author Giuseppe Cerati, INFN
6  *
7  * $Date: 2013/05/14 15:46:46 $
8  * $Revision: 1.5.4.2 $
9  *
10  */
15 
17 public:
20  double ptMax,
21  double minRapidity,
22  double maxRapidity,
23  double tip,
24  double lip,
25  int minHit,
26  bool signalOnly,
27  bool intimeOnly,
28  bool chargedOnly,
29  bool stableOnly,
30  const std::vector<int> &pdgId = std::vector<int>(),
31  bool invertRapidityCut = false,
32  double minPhi = -3.2,
33  double maxPhi = 3.2)
34  : ptMin2_(ptMin * ptMin),
35  ptMax2_(ptMax * ptMax),
38  meanPhi_((minPhi + maxPhi) / 2.),
39  rangePhi_((maxPhi - minPhi) / 2.),
40  tip2_(tip * tip),
41  lip_(lip),
42  minHit_(minHit),
47  pdgId_(pdgId),
49  if (minPhi >= maxPhi) {
50  throw cms::Exception("Configuration")
51  << "TrackingParticleSelector: minPhi (" << minPhi << ") must be smaller than maxPhi (" << maxPhi
52  << "). The range is constructed from minPhi to maxPhi around their "
53  "average.";
54  }
55  if (minPhi >= M_PI) {
56  throw cms::Exception("Configuration") << "TrackingParticleSelector: minPhi (" << minPhi
57  << ") must be smaller than PI. The range is constructed from minPhi "
58  "to maxPhi around their average.";
59  }
60  if (maxPhi <= -M_PI) {
61  throw cms::Exception("Configuration") << "TrackingParticleSelector: maxPhi (" << maxPhi
62  << ") must be larger than -PI. The range is constructed from minPhi "
63  "to maxPhi around their average.";
64  }
65  }
66 
68  bool operator()(const TrackingParticle &tp) const {
69  // signal only means no PU particles
70  if (signalOnly_ && !(tp.eventId().bunchCrossing() == 0 && tp.eventId().event() == 0))
71  return false;
72  // intime only means no OOT PU particles
73  if (intimeOnly_ && !(tp.eventId().bunchCrossing() == 0))
74  return false;
75 
76  auto pdgid = tp.pdgId();
77  if (!pdgId_.empty()) {
78  bool testId = false;
79  for (auto id : pdgId_) {
80  if (id == pdgid) {
81  testId = true;
82  break;
83  }
84  }
85  if (!testId)
86  return false;
87  }
88 
89  if (chargedOnly_ && tp.charge() == 0)
90  return false; // select only if charge!=0
91 
92  // select only stable particles
93  if (stableOnly_) {
94  for (TrackingParticle::genp_iterator j = tp.genParticle_begin(); j != tp.genParticle_end(); ++j) {
95  if (j->get() == nullptr || j->get()->status() != 1) {
96  return false;
97  }
98  }
99  // test for remaining unstabled due to lack of genparticle pointer
100  if (tp.status() == -99 && (std::abs(pdgid) != 11 && std::abs(pdgid) != 13 && std::abs(pdgid) != 211 &&
101  std::abs(pdgid) != 321 && std::abs(pdgid) != 2212 && std::abs(pdgid) != 3112 &&
102  std::abs(pdgid) != 3222 && std::abs(pdgid) != 3312 && std::abs(pdgid) != 3334))
103  return false;
104  }
105 
106  auto etaOk = [&](const TrackingParticle &p) -> bool {
107  float eta = etaFromXYZ(p.px(), p.py(), p.pz());
108  if (!invertRapidityCut_)
109  return (eta >= minRapidity_) && (eta <= maxRapidity_);
110  else
111  return (eta < minRapidity_ || eta > maxRapidity_);
112  };
113  auto phiOk = [&](const TrackingParticle &p) {
114  float dphi = deltaPhi(atan2f(p.py(), p.px()), meanPhi_);
115  return dphi >= -rangePhi_ && dphi <= rangePhi_;
116  };
117  auto ptOk = [&](const TrackingParticle &p) {
118  double pt2 = tp.p4().perp2();
119  return pt2 >= ptMin2_ && pt2 <= ptMax2_;
120  };
121  return (tp.numberOfTrackerLayers() >= minHit_ && ptOk(tp) && etaOk(tp) && phiOk(tp) &&
122  std::abs(tp.vertex().z()) <= lip_ && // vertex last to avoid to load it if not striclty
123  // necessary...
124  tp.vertex().perp2() <= tip2_);
125  }
126 
127 private:
128  double ptMin2_;
129  double ptMax2_;
132  float meanPhi_;
133  float rangePhi_;
134  double tip2_;
135  double lip_;
136  int minHit_;
141  std::vector<int> pdgId_;
143 };
144 
147 
148 namespace reco {
149  namespace modules {
150 
151  template <>
154  return make(cfg);
155  }
156 
158  return TrackingParticleSelector(cfg.getParameter<double>("ptMin"),
159  cfg.getParameter<double>("ptMax"),
160  cfg.getParameter<double>("minRapidity"),
161  cfg.getParameter<double>("maxRapidity"),
162  cfg.getParameter<double>("tip"),
163  cfg.getParameter<double>("lip"),
164  cfg.getParameter<int>("minHit"),
165  cfg.getParameter<bool>("signalOnly"),
166  cfg.getParameter<bool>("intimeOnly"),
167  cfg.getParameter<bool>("chargedOnly"),
168  cfg.getParameter<bool>("stableOnly"),
169  cfg.getParameter<std::vector<int>>("pdgId"),
170  cfg.getParameter<bool>("invertRapidityCut"),
171  cfg.getParameter<double>("minPhi"),
172  cfg.getParameter<double>("maxPhi"));
173  }
174  };
175 
176  } // namespace modules
177 } // namespace reco
178 
179 #endif
trackingParticleSelector_cfi.signalOnly
signalOnly
Definition: trackingParticleSelector_cfi.py:9
TrackingParticleSelector::TrackingParticleSelector
TrackingParticleSelector()
Definition: TrackingParticleSelector.h:18
HLT_FULL_cff.invertRapidityCut
invertRapidityCut
Definition: HLT_FULL_cff.py:52982
TrackingParticleSelector::chargedOnly_
bool chargedOnly_
Definition: TrackingParticleSelector.h:139
qcdUeDQM_cfi.minRapidity
minRapidity
Definition: qcdUeDQM_cfi.py:24
reco::modules::ParameterAdapter< TrackingParticleSelector >::make
static TrackingParticleSelector make(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
Definition: TrackingParticleSelector.h:153
reco::modules::ParameterAdapter< TrackingParticleSelector >::make
static TrackingParticleSelector make(const edm::ParameterSet &cfg)
Definition: TrackingParticleSelector.h:157
deltaPhi.h
modules
Definition: MuonCleanerBySegments.cc:35
qcdUeDQM_cfi.maxRapidity
maxRapidity
Definition: qcdUeDQM_cfi.py:27
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
ptMin
constexpr float ptMin
Definition: PhotonIDValueMapProducer.cc:155
TrackingParticleSelector
Definition: TrackingParticleSelector.h:16
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
TrackingParticleSelector::signalOnly_
bool signalOnly_
Definition: TrackingParticleSelector.h:137
qcdUeDQM_cfi.lip
lip
Definition: qcdUeDQM_cfi.py:25
GenParticle.h
HLT_FULL_cff.maxPhi
maxPhi
Definition: HLT_FULL_cff.py:52983
TrackingParticleSelector::pdgId_
std::vector< int > pdgId_
Definition: TrackingParticleSelector.h:141
AlignmentTrackSelector_cfi.ptMax
ptMax
Definition: AlignmentTrackSelector_cfi.py:12
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
PVValHelper::eta
Definition: PVValidationHelpers.h:70
TrackingParticleSelector::maxRapidity_
float maxRapidity_
Definition: TrackingParticleSelector.h:131
TrackingParticle
Monte Carlo truth information used for tracking validation.
Definition: TrackingParticle.h:29
reco::modules::TrackingParticleSelector
SingleObjectSelector< TrackingParticleCollection, ::TrackingParticleSelector > TrackingParticleSelector
Definition: TrackingParticleSelector.cc:17
TrackingParticleSelector::ptMax2_
double ptMax2_
Definition: TrackingParticleSelector.h:129
cmsswSequenceInfo.tp
tp
Definition: cmsswSequenceInfo.py:17
qcdUeDQM_cfi.tip
tip
Definition: qcdUeDQM_cfi.py:23
TrackingParticleSelector::meanPhi_
float meanPhi_
Definition: TrackingParticleSelector.h:132
TrackingParticleSelector::invertRapidityCut_
bool invertRapidityCut_
Definition: TrackingParticleSelector.h:142
trackingParticleSelector_cfi.intimeOnly
intimeOnly
Definition: trackingParticleSelector_cfi.py:10
edm::ParameterSet
Definition: ParameterSet.h:47
genCandidates_cfi.stableOnly
stableOnly
Definition: genCandidates_cfi.py:6
TrackingParticleSelector::minHit_
int minHit_
Definition: TrackingParticleSelector.h:136
TrackingParticleSelector::rangePhi_
float rangePhi_
Definition: TrackingParticleSelector.h:133
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
reco::modules::ParameterAdapter::make
static S make(const edm::ParameterSet &cfg)
Definition: ParameterAdapter.h:13
TrackingParticleSelector::TrackingParticleSelector
TrackingParticleSelector(double ptMin, double ptMax, double minRapidity, double maxRapidity, double tip, double lip, int minHit, bool signalOnly, bool intimeOnly, bool chargedOnly, bool stableOnly, const std::vector< int > &pdgId=std::vector< int >(), bool invertRapidityCut=false, double minPhi=-3.2, double maxPhi=3.2)
Definition: TrackingParticleSelector.h:19
HLT_FULL_cff.pt2
pt2
Definition: HLT_FULL_cff.py:9872
qcdUeDQM_cfi.minHit
minHit
Definition: qcdUeDQM_cfi.py:33
SingleObjectSelectorBase
Definition: SingleObjectSelector.h:26
TrackingParticleSelector::lip_
double lip_
Definition: TrackingParticleSelector.h:135
looper.cfg
cfg
Definition: looper.py:297
TrackingParticleSelector::ptMin2_
double ptMin2_
Definition: TrackingParticleSelector.h:128
HLT_FULL_cff.minPhi
minPhi
Definition: HLT_FULL_cff.py:52971
TrackingParticle.h
PtEtaPhiMass.h
Exception
Definition: hltDiff.cc:245
edm::RefVectorIterator
Definition: EDProductfwd.h:33
TrackingParticleSelector::tip2_
double tip2_
Definition: TrackingParticleSelector.h:134
TrackingParticleSelector::stableOnly_
bool stableOnly_
Definition: TrackingParticleSelector.h:140
cosmictrackingParticleSelector_cfi.chargedOnly
chargedOnly
Definition: cosmictrackingParticleSelector_cfi.py:5
ParameterAdapter.h
ConsumesCollector.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TrackingParticleSelector::minRapidity_
float minRapidity_
Definition: TrackingParticleSelector.h:130
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
EgammaValidation_cff.pdgid
pdgid
Definition: EgammaValidation_cff.py:30
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
TrackingParticleSelector::operator()
bool operator()(const TrackingParticle &tp) const
Operator() performs the selection: e.g. if (tPSelector(tp)) {...}.
Definition: TrackingParticleSelector.h:68
reco::modules::ParameterAdapter
Definition: ParameterAdapter.h:12
TrackingParticleSelector::intimeOnly_
bool intimeOnly_
Definition: TrackingParticleSelector.h:138