CMS 3D CMS Logo

HadronDecayGenEvtSelector.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 
27  int id = hadronId_.size();
28  int st = hadronStatus_.size();
29  int etamax = hadronEtaMax_.size();
30  int etamin = hadronEtaMin_.size();
31  int pmin = hadronPMin_.size();
32  int ptmax = hadronPtMax_.size();
33  int ptmin = hadronPtMin_.size();
34 
35  if( id!=st || id!=etamax || id!=etamin || id!=ptmax || id!=ptmin || id!=pmin)
36  {
37  throw edm::Exception(edm::errors::LogicError)<<"Hadron selection parameters: "<<id<<st<<etamax<<etamin<<pmin<<ptmax<<ptmin<<endl;
38  }
39 
40 
41 }
42 
43 
44 //____________________________________________________________________________________________
45 bool HadronDecayGenEvtSelector::filter(HepMC::GenEvent *evt)
46 {
47  // loop over HepMC event, and search for products of interest
48 
49  HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
50  HepMC::GenEvent::particle_const_iterator end = evt->particles_end();
51 
52  bool foundHadron = false;
53  bool foundDecay = false;
54 
55  int foundtrig = 0;
56  HepMC::GenEvent::particle_const_iterator it = begin;
57  while( (!foundHadron || !foundDecay) && it != end )
58  {
59  for(unsigned i = 0; i < hadronId_.size(); ++i)
60  {
61  if( selectParticle(*it,
64  hadronPMin_[i],
65  hadronPtMax_[i],hadronPtMin_[i]) ) foundHadron = true;
66  }
67 
68  if( selectParticle(*it,
71  decayPMin_,
72  decayPtMax_,decayPtMin_) ) foundtrig++;
73  if(decayNtrig_ == foundtrig) foundDecay = true;
74 
75  ++it;
76  }
77 
78  return (foundHadron && foundDecay);
79 }
80 
81 
82 //____________________________________________________________________________________________
T getParameter(std::string const &) const
bool selectParticle(HepMC::GenParticle *par, int status, int pdg, double etaMax, double etaMin, double pMin, double ptMax, double ptMin)
std::vector< double > hadronEtaMax_
HadronDecayGenEvtSelector(const edm::ParameterSet &pset)
std::vector< double > hadronEtaMin_
#define end
Definition: vmac.h:39
std::vector< double > hadronPtMax_
std::vector< double > hadronPtMin_
double ptmin
Definition: HydjetWrapper.h:90
#define begin
Definition: vmac.h:32
bool filter(HepMC::GenEvent *) override