CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BdecayFilter.cc
Go to the documentation of this file.
2 
3 using namespace edm;
4 using namespace std;
5 using namespace HepMC;
6 
7 
8 
10 {
11  token_ = consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"));
12  motherParticle = iConfig.getParameter< int >("motherParticle");
13 
14  firstDaughter.type = iConfig.getParameter< int >("firstDaughter");
15  firstDaughter.decayProduct = iConfig.getParameter< vector<int> >("firstDaughterDecay");
16  firstDaughter.etaMin = iConfig.getParameter<double>("firstDaughterDecayEtaMin");
17  firstDaughter.etaMax = iConfig.getParameter<double>("firstDaughterDecayEtaMax");
18  firstDaughter.ptMin = iConfig.getParameter<double>("firstDaughterDecayPtMin");
19 
20  secondDaughter.type = iConfig.getParameter< int >("secondDaughter");
21  secondDaughter.decayProduct = iConfig.getParameter< vector<int> >("secondDaughterDecay");
22  secondDaughter.etaMin = iConfig.getParameter<double>("secondDaughterDecayEtaMin");
23  secondDaughter.etaMax = iConfig.getParameter<double>("secondDaughterDecayEtaMax");
24  secondDaughter.ptMin = iConfig.getParameter<double>("secondDaughterDecayPtMin");
25 
26  noAccepted = 0;
27 }
28 
29 
31 {
32  std::cout << "Total number of accepted events = " << noAccepted << std::endl;
33 }
34 
35 /*
36 HepMC::GenParticle * BdecayFilter::findParticle(const GenPartVect genPartVect,
37  const int requested_id)
38 {
39  for (GenPartVectIt p = genPartVect.begin(); p != genPartVect.end(); p++)
40  {
41  int event_particle_id = abs( (*p)->pdg_id() );
42  cout << "isC "<<event_particle_id<<"\n";
43  if (requested_id == event_particle_id) return *p;
44  }
45  return 0;
46 }
47 */
48 
49 HepMC::GenParticle * BdecayFilter::findParticle(HepMC::GenVertex* vertex,
50  const int requested_id)
51 {
52  // for(std::set<GenParticle*>::const_iterator p = vertex->particles_out_const_begin();
53  for(GenVertex::particles_out_const_iterator p = vertex->particles_out_const_begin();
54  p != vertex->particles_out_const_end(); p++)
55  {
56  int event_particle_id = abs( (*p)->pdg_id() );
57  cout << "particle Id: "<<event_particle_id<<"\n";
58  if (requested_id == event_particle_id) return *p;
59  }
60  return 0;
61 }
62 
63 HepMC::GenEvent::particle_const_iterator BdecayFilter::getNextBs(const HepMC::GenEvent::particle_const_iterator start, const HepMC::GenEvent::particle_const_iterator end)
64 {
65  HepMC::GenEvent::particle_const_iterator p;
66  for (p = start; p != end; p++)
67  {
68  int event_particle_id = abs( (*p)->pdg_id() );
69 // cout << "search "<<event_particle_id<<"\n";
70  if (event_particle_id == motherParticle) return p;
71  }
72  return p;
73 }
74 
75 
77 {
79  iEvent.getByToken(token_, evt);
80 
81  const HepMC::GenEvent * generated_event = evt->GetEvent();
82  //cout << "Start\n";
83 
84  bool event_passed = false;
85  HepMC::GenEvent::particle_const_iterator bs = getNextBs(generated_event->particles_begin(), generated_event->particles_end());
86  while (bs!= generated_event->particles_end() ) {
87 
88  // vector< GenParticle * > bsChild = (*bs)->listChildren();
89 
90  //***
91  HepMC::GenVertex* outVertex = (*bs)->end_vertex();
92  //***
93 
94  GenParticle * jpsi = 0;
95  GenParticle * phi = 0;
96  // cout << "bs size "<<bsChild.size()<<endl;
97  //***
98  int numChildren = outVertex->particles_out_size();
99  cout<< "B children "<<numChildren<<endl;
100  //***
101 
102  /* if ((bsChild.size()==2) && ((jpsi = findParticle(bsChild, 443))!=0) &&
103  ((phi = findParticle(bsChild, 333))!=0)) {
104  cout << bsChild[0]->momentum()<<" "<<bsChild[0]->momentum().eta()
105  <<" "<<bsChild[1]->momentum()<<" "<<bsChild[1]->momentum().eta()<<endl;
106  */
107 
108  //***
109  if( (numChildren==2) && ((jpsi = findParticle(outVertex, firstDaughter.type))!=0) &&
110  ((phi = findParticle(outVertex, secondDaughter.type))!=0)) {
111 
112  cout << jpsi->momentum().rho() <<" "<<jpsi->momentum().eta() <<" "<<phi->momentum().rho()<<" "<<phi->momentum().eta()<<endl;
113  //cout <<"bs dec trouve"<<endl;
114  if (cuts(phi, secondDaughter) && cuts(jpsi, firstDaughter)) {
115  cout <<"decay found"<<endl;
116  event_passed = true;
117  break;
118  }
119  }
120  bs = getNextBs(++bs, generated_event->particles_end());
121  }
122 
123  if (event_passed) noAccepted++;
124  cout << "End filter\n";
125 
126  delete generated_event;
127 
128  return event_passed;
129 }
130 
131 
133 {
134  cout << "start cuts" << endl;
135  HepMC::GenVertex* myVertex = jpsi->end_vertex();
136  int numChildren = myVertex->particles_out_size();
137  int numDecProd = cut.decayProduct.size();
138  std::vector<HepMC::GenParticle*> psiChild;
139  // for(std::set<GenParticle*>::const_iterator p = myVertex->particles_out_const_begin();
140  for(GenVertex::particles_out_const_iterator p = myVertex->particles_out_const_begin();
141  p != myVertex->particles_out_const_end(); p++)
142  psiChild.push_back((*p));
143 
144  if (numChildren!=numDecProd) return false;
145  cout << psiChild[0]->pdg_id()<<" "<<psiChild[1]->pdg_id()<<endl;
146 
147  for (int i=0; i<numChildren; ++i) {
148  bool goodPart = false;
149  for (int j=0; j < numChildren; ++j){
150  if (abs(psiChild[i]->pdg_id()) == abs(cut.decayProduct[j])) goodPart = true;
151  }
152  cout << psiChild[i]->momentum().perp() << endl;
153  if ( !goodPart || (!etaInRange(psiChild[i]->momentum().eta(), cut.etaMin, cut.etaMax)) || (psiChild[i]->momentum().perp() < cut.ptMin) ) return false;
154  }
155  cout << "cuts true" << endl;
156  return true;
157 }
158 
159 bool BdecayFilter::etaInRange(float eta, float etamin, float etamax)
160 {
161  return ( (etamin < eta) && (eta < etamax) );
162 }
T getParameter(std::string const &) const
virtual bool filter(edm::Event &, const edm::EventSetup &)
Definition: BdecayFilter.cc:76
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
bool cuts(const HepMC::GenParticle *jpsi, const CutStruct &cut)
int iEvent
Definition: GenABIO.cc:230
HepMC::GenParticle * findParticle(HepMC::GenVertex *, const int requested_id)
Definition: BdecayFilter.cc:49
std::vector< int > decayProduct
Definition: BdecayFilter.h:41
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
HepMC::GenEvent::particle_const_iterator getNextBs(const HepMC::GenEvent::particle_const_iterator start, const HepMC::GenEvent::particle_const_iterator end)
Definition: BdecayFilter.cc:63
BdecayFilter(const edm::ParameterSet &)
Definition: BdecayFilter.cc:9
bool etaInRange(float eta, float etamin, float etamax)
Geom::Phi< T > phi() const
#define jpsi
tuple cout
Definition: gather_cfg.py:121