CMS 3D CMS Logo

GenParticleCustomSelector.h
Go to the documentation of this file.
1 #ifndef RecoSelectors_GenParticleCustomSelector_h
2 #define RecoSelectors_GenParticleCustomSelector_h
3 /* \class GenParticleCustomSelector
4  *
5  * \author Giuseppe Cerati, UCSD
6  *
7  */
8 
10 
12 
13 public:
16  double tip,double lip, bool chargedOnly, int status,
17  const std::vector<int>& pdgId = std::vector<int>()) :
18  ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
19  tip_( tip ), lip_( lip ), chargedOnly_(chargedOnly), status_(status), pdgId_( pdgId ) { }
20 
22  bool operator()( const reco::GenParticle & tp ) const {
23 
24  if (chargedOnly_ && tp.charge()==0) return false;//select only if charge!=0
25  bool testId = false;
26  unsigned int idSize = pdgId_.size();
27  if (idSize==0) testId = true;
28  else for (unsigned int it=0;it!=idSize;++it){
29  if (tp.pdgId()==pdgId_[it]) testId = true;
30  }
31 
32  return (
33  tp.pt() >= ptMin_ &&
34  tp.eta() >= minRapidity_ && tp.eta() <= maxRapidity_ &&
35  sqrt(tp.vertex().perp2()) <= tip_ &&
36  fabs(tp.vertex().z()) <= lip_ &&
37  tp.status() == status_ &&
38  testId
39  );
40  }
41 
42 private:
43  double ptMin_;
44  double minRapidity_;
45  double maxRapidity_;
46  double tip_;
47  double lip_;
49  int status_;
50  std::vector<int> pdgId_;
51 
52 };
53 
56 
57 namespace reco {
58  namespace modules {
59 
60  template<>
64  cfg.getParameter<double>( "ptMin" ),
65  cfg.getParameter<double>( "minRapidity" ),
66  cfg.getParameter<double>( "maxRapidity" ),
67  cfg.getParameter<double>( "tip" ),
68  cfg.getParameter<double>( "lip" ),
69  cfg.getParameter<bool>( "chargedOnly" ),
70  cfg.getParameter<int>( "status" ),
71  cfg.getParameter<std::vector<int> >( "pdgId" ));
72  }
73  };
74 
75  }
76 }
77 
78 #endif
T getParameter(std::string const &) const
int pdgId() const final
PDG identifier.
double eta() const final
momentum pseudorapidity
double pt() const final
transverse momentum
int charge() const final
electric charge
Definition: LeafCandidate.h:91
GenParticleCustomSelector(double ptMin, double minRapidity, double maxRapidity, double tip, double lip, bool chargedOnly, int status, const std::vector< int > &pdgId=std::vector< int >())
T sqrt(T t)
Definition: SSEVec.h:18
const Point & vertex() const override
vertex position (overwritten by PF...)
static GenParticleCustomSelector make(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
bool operator()(const reco::GenParticle &tp) const
Operator() performs the selection: e.g. if (tPSelector(tp)) {...}.
fixed size matrix
int status() const final
status word