CMS 3D CMS Logo

PartonHadronDecayGenEvtSelector.cc
Go to the documentation of this file.
1 #include <iostream>
4 using namespace std;
5 
7 {
8 
9  hadronId_ = pset.getParameter<vector<int> >("hadrons");
10  hadronStatus_ = pset.getParameter<vector<int> >("hadronStatus");
11  hadronEtaMax_ = pset.getParameter<vector<double> >("hadronEtaMax");
12  hadronEtaMin_ = pset.getParameter<vector<double> >("hadronEtaMin");
13  hadronPMin_ = pset.getParameter<vector<double> >("hadronPMin");
14  hadronPtMax_ = pset.getParameter<vector<double> >("hadronPtMax");
15  hadronPtMin_ = pset.getParameter<vector<double> >("hadronPtMin");
16 
17  decayId_ = pset.getParameter<int>("decays");
18  decayStatus_ = pset.getParameter<int>("decayStatus");
19  decayEtaMax_ = pset.getParameter<double>("decayEtaMax");
20  decayEtaMin_ = pset.getParameter<double>("decayEtaMin");
21  decayPMin_ = pset.getParameter<double>("decayPMin");
22  decayPtMax_ = pset.getParameter<double>("decayPtMax");
23  decayPtMin_ = pset.getParameter<double>("decayPtMin");
24  decayNtrig_ = pset.getParameter<int>("decayNtrig");
25 
26  partonId_ = pset.getParameter<vector<int> >("partons");
27  partonStatus_ = pset.getParameter<vector<int> >("partonStatus");
28  partonEtaMax_ = pset.getParameter<vector<double> >("partonEtaMax");
29  partonPtMin_ = pset.getParameter<vector<double> >("partonPtMin");
30 
31  int id = hadronId_.size();
32  int st = hadronStatus_.size();
33  int etamax = hadronEtaMax_.size();
34  int etamin = hadronEtaMin_.size();
35  int pmin = hadronPMin_.size();
36  int ptmax = hadronPtMax_.size();
37  int ptmin = hadronPtMin_.size();
38 
39  if( id!=st || id!=etamax || id!=etamin || id!=ptmax || id!=ptmin || id!=pmin)
40  {
41  throw edm::Exception(edm::errors::LogicError)<<"Hadron selection parameters: "<<id<<st<<etamax<<etamin<<pmin<<ptmax<<ptmin<<endl;
42  }
43 
44 
45  id = partonId_.size();
46  st = partonStatus_.size();
47  etamax = partonEtaMax_.size();
48  ptmin = partonPtMin_.size();
49 
50  if( id!=st || id!=etamax || id!=ptmin )
51  {
52  throw edm::Exception(edm::errors::LogicError)<<"Parton selection parameters: "<<id<<st<<etamax<<ptmin<<endl;
53  }
54 
55 }
56 
57 
58 //____________________________________________________________________________________________
59 bool PartonHadronDecayGenEvtSelector::filter(HepMC::GenEvent *evt)
60 {
61  // loop over HepMC event, and search for products of interest
62 
63  HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
64  HepMC::GenEvent::particle_const_iterator end = evt->particles_end();
65 
66  bool foundHadron = false;
67  bool foundDecay = false;
68  bool foundParton = false;
69 
70  HepMC::GenEvent::particle_const_iterator it = begin;
71  while( !foundParton && it != end )
72  {
73  for(unsigned i = 0; i < partonId_.size(); ++i)
74  {
75  if( selectParticle(*it,
77  partonPtMin_[i],partonEtaMax_[i]) ) foundParton = true;
78  }
79  ++it;
80  }
81 
82  int foundtrig = 0;
83  HepMC::GenEvent::particle_const_iterator it2 = begin;
84 
85  if(foundParton)
86  {
87  while( (!foundHadron || !foundDecay) && it2 != end )
88  {
89 
90  for(unsigned i = 0; i < hadronId_.size(); ++i)
91  {
92  if( selectParticle(*it2,
95  hadronPMin_[i],
96  hadronPtMax_[i],hadronPtMin_[i]) ) foundHadron = true;
97  }
98 
99  if( selectParticle(*it2,
102  decayPMin_,
103  decayPtMax_,decayPtMin_) ) foundtrig++;
104  if(decayNtrig_ == foundtrig) foundDecay = true;
105 
106  ++it2;
107  }
108  }
109 
110  return (foundHadron && foundDecay && foundParton);
111 }
112 
113 
114 //____________________________________________________________________________________________
T getParameter(std::string const &) const
PartonHadronDecayGenEvtSelector(const edm::ParameterSet &pset)
#define end
Definition: vmac.h:37
double ptmin
Definition: HydjetWrapper.h:90
#define begin
Definition: vmac.h:30
bool selectParticle(HepMC::GenParticle *par, int status, int pdg, double ptMin, double etaMax)