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 
18 public:
20  TrackingParticleSelector ( double ptMin, double ptMax, double minRapidity,double maxRapidity,
21  double tip,double lip,int minHit, bool signalOnly, bool intimeOnly, bool chargedOnly, bool stableOnly,
22  const std::vector<int>& pdgId = std::vector<int>(),
23  double minPhi=-3.2, double maxPhi=3.2) :
24  ptMin2_( ptMin*ptMin ), ptMax2_( ptMax*ptMax ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
26  tip2_( tip*tip ), lip_( lip ), minHit_( minHit ), signalOnly_(signalOnly), intimeOnly_(intimeOnly), chargedOnly_(chargedOnly), stableOnly_(stableOnly), pdgId_( pdgId ) {
27  if(minPhi >= maxPhi) {
28  throw cms::Exception("Configuration") << "TrackingParticleSelector: minPhi (" << minPhi << ") must be smaller than maxPhi (" << maxPhi << "). The range is constructed from minPhi to maxPhi around their average.";
29  }
30  if(minPhi >= M_PI) {
31  throw cms::Exception("Configuration") << "TrackingParticleSelector: minPhi (" << minPhi << ") must be smaller than PI. The range is constructed from minPhi to maxPhi around their average.";
32  }
33  if(maxPhi <= -M_PI) {
34  throw cms::Exception("Configuration") << "TrackingParticleSelector: maxPhi (" << maxPhi << ") must be larger than -PI. The range is constructed from minPhi to maxPhi around their average.";
35  }
36  }
37 
39  bool operator()( const TrackingParticle & tp ) const {
40  // signal only means no PU particles
41  if (signalOnly_ && !(tp.eventId().bunchCrossing()== 0 && tp.eventId().event() == 0)) return false;
42  // intime only means no OOT PU particles
43  if (intimeOnly_ && !(tp.eventId().bunchCrossing()==0)) return false;
44 
45  auto pdgid = tp.pdgId();
46  if(!pdgId_.empty()) {
47  bool testId = false;
48  for(auto id: pdgId_) {
49  if(id == pdgid) { testId = true; break;}
50  }
51  if(!testId) return false;
52  }
53 
54  if (chargedOnly_ && tp.charge()==0) return false;//select only if charge!=0
55 
56  // select only stable particles
57  if (stableOnly_) {
58  for( TrackingParticle::genp_iterator j = tp.genParticle_begin(); j != tp.genParticle_end(); ++ j ) {
59  if (j->get()==0 || j->get()->status() != 1) {
60  return false;
61  }
62  }
63  // test for remaining unstabled due to lack of genparticle pointer
64  if( tp.status() == -99 &&
65  (std::abs(pdgid) != 11 && std::abs(pdgid) != 13 && std::abs(pdgid) != 211 &&
66  std::abs(pdgid) != 321 && std::abs(pdgid) != 2212 && std::abs(pdgid) != 3112 &&
67  std::abs(pdgid) != 3222 && std::abs(pdgid) != 3312 && std::abs(pdgid) != 3334))
68  return false;
69  }
70 
71  auto etaOk = [&](const TrackingParticle& p)->bool{ float eta= etaFromXYZ(p.px(),p.py(),p.pz()); return (eta>= minRapidity_) & (eta<=maxRapidity_);};
72  auto phiOk = [&](const TrackingParticle& p) { float dphi = deltaPhi(atan2f(p.py(),p.px()), meanPhi_); return dphi >= -rangePhi_ && dphi <= rangePhi_; };
73  auto ptOk = [&](const TrackingParticle& p) { double pt2 = tp.p4().perp2(); return pt2 >= ptMin2_ && pt2 <= ptMax2_; };
74  return (
76  ptOk(tp) &&
77  etaOk(tp) &&
78  phiOk(tp) &&
79  std::abs(tp.vertex().z()) <= lip_ && // vertex last to avoid to load it if not striclty necessary...
80  tp.vertex().perp2() <= tip2_
81  );
82  }
83 
84 private:
85  double ptMin2_;
86  double ptMax2_;
87  float minRapidity_;
88  float maxRapidity_;
89  float meanPhi_;
90  float rangePhi_;
91  double tip2_;
92  double lip_;
93  int minHit_;
98  std::vector<int> pdgId_;
99 
100 };
101 
104 
105 namespace reco {
106  namespace modules {
107 
108  template<>
112  cfg.getParameter<double>( "ptMin" ),
113  cfg.getParameter<double>( "ptMax" ),
114  cfg.getParameter<double>( "minRapidity" ),
115  cfg.getParameter<double>( "maxRapidity" ),
116  cfg.getParameter<double>( "tip" ),
117  cfg.getParameter<double>( "lip" ),
118  cfg.getParameter<int>( "minHit" ),
119  cfg.getParameter<bool>( "signalOnly" ),
120  cfg.getParameter<bool>( "intimeOnly" ),
121  cfg.getParameter<bool>( "chargedOnly" ),
122  cfg.getParameter<bool>( "stableOnly" ),
123  cfg.getParameter<std::vector<int> >( "pdgId" ),
124  cfg.getParameter<double>( "minPhi" ),
125  cfg.getParameter<double>( "maxPhi" ));
126  }
127  };
128 
129  }
130 }
131 
132 #endif
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
T getParameter(std::string const &) const
genp_iterator genParticle_begin() const
iterators
int event() const
get the contents of the subdetector field (should be protected?)
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 >(), double minPhi=-3.2, double maxPhi=3.2)
int pdgId() const
PDG ID.
int status() const
Status word.
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
int bunchCrossing() const
get the detector field from this detid
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool operator()(const TrackingParticle &tp) const
Operator() performs the selection: e.g. if (tPSelector(tp)) {...}.
#define M_PI
genp_iterator genParticle_end() const
Point vertex() const
Parent vertex position.
EncodedEventId eventId() const
Signal source, crossing number.
fixed size matrix
Monte Carlo truth information used for tracking validation.
static TrackingParticleSelector make(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)