CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BsJpsiPhiFilter.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  label_ = iConfig.getUntrackedParameter("moduleLabel",std::string("generator"));
12  hadronCuts.type = iConfig.getParameter< int >("hadronType");
13  hadronCuts.etaMin = iConfig.getParameter<double>("hadronEtaMin");
14  hadronCuts.etaMax = iConfig.getParameter<double>("hadronEtaMax");
15  hadronCuts.ptMin = iConfig.getParameter<double>("hadronPtMin");
16  leptonCuts.type = iConfig.getParameter< int >("leptonType");
17  leptonCuts.etaMin = iConfig.getParameter<double>("leptonEtaMin");
18  leptonCuts.etaMax = iConfig.getParameter<double>("leptonEtaMax");
19  leptonCuts.ptMin = iConfig.getParameter<double>("leptonPtMin");
20 
21  noAccepted = 0;
22 }
23 
24 
26 {
27  std::cout << "Total number of accepted events = " << noAccepted << std::endl;
28 }
29 
30 /*
31 HepMC::GenParticle * BsJpsiPhiFilter::findParticle(const GenPartVect genPartVect,
32  const int requested_id)
33 {
34  for (GenPartVectIt p = genPartVect.begin(); p != genPartVect.end(); p++)
35  {
36  int event_particle_id = abs( (*p)->pdg_id() );
37  cout << "isC "<<event_particle_id<<"\n";
38  if (requested_id == event_particle_id) return *p;
39  }
40  return 0;
41 }
42 */
43 
45  const int requested_id)
46 {
47  for(std::vector<GenParticle*>::const_iterator p = vertex->particles_out_const_begin();
48  p != vertex->particles_out_const_end(); p++)
49  {
50  int event_particle_id = abs( (*p)->pdg_id() );
51  cout << "isC "<<event_particle_id<<"\n";
52  if (requested_id == event_particle_id) return *p;
53  }
54  return 0;
55 }
56 
57 HepMC::GenEvent::particle_const_iterator
58 BsJpsiPhiFilter::getNextBs(const HepMC::GenEvent::particle_const_iterator start,
59  const HepMC::GenEvent::particle_const_iterator end)
60 {
61  HepMC::GenEvent::particle_const_iterator p;
62  for (p = start; p != end; p++)
63  {
64  int event_particle_id = abs( (*p)->pdg_id() );
65 // cout << "search "<<event_particle_id<<"\n";
66  if (event_particle_id == 531) return p;
67  }
68  return p;
69 }
70 
71 
73 {
75  iEvent.getByLabel(label_, evt);
76 
77  const HepMC::GenEvent * generated_event = evt->GetEvent();
78  //cout << "Start\n";
79 
80  bool event_passed = false;
81  HepMC::GenEvent::particle_const_iterator bs = getNextBs(generated_event->particles_begin(),
82  generated_event->particles_end());
83  while (bs!= generated_event->particles_end() ) {
84 
85  // vector< GenParticle * > bsChild = (*bs)->listChildren();
86 
87  //***
88  HepMC::GenVertex* outVertex = (*bs)->end_vertex();
89  //***
90 
91  GenParticle * jpsi = 0;
92  GenParticle * phi = 0;
93  // cout << "bs size "<<bsChild.size()<<endl;
94  //***
95  int numChildren = outVertex->particles_out_size();
96  cout<< "bs size "<<numChildren<<endl;
97  //***
98 
99  /* if ((bsChild.size()==2) && ((jpsi = findParticle(bsChild, 443))!=0) &&
100  ((phi = findParticle(bsChild, 333))!=0)) {
101  cout << bsChild[0]->momentum()<<" "<<bsChild[0]->momentum().eta()
102  <<" "<<bsChild[1]->momentum()<<" "<<bsChild[1]->momentum().eta()<<endl;
103  */
104 
105  //***
106  if( (numChildren==2) && ((jpsi = findParticle(outVertex, 443))!=0) &&
107  ((phi = findParticle(outVertex, 333))!=0)) {
108 
109  cout << jpsi->momentum().rho()<<" "<<jpsi->momentum().eta()
110  <<" "<<phi->momentum().rho() <<" "<<phi->momentum().eta()<<endl;
111  cout <<"bs dec trouve"<<endl;
112  if (cuts(phi, hadronCuts) && cuts(jpsi, leptonCuts)) {
113  cout <<"OK trouve"<<endl;
114  event_passed = true;
115  break;
116  }
117  }
118  bs = getNextBs(++bs, generated_event->particles_end());
119  }
120 
121  if (event_passed) noAccepted++;
122  cout << "End filter\n";
123 
124  delete generated_event;
125 
126  return event_passed;
127 }
128 
129 /*
130 bool BsJpsiPhiFilter::cuts(const GenParticle * jpsi, const CutStruct cut)
131 {
132  GenPartVect psiChild = jpsi->listChildren();
133  cout << psiChild[0]->pdg_id()<<" "<<psiChild[1]->pdg_id()<<endl;
134  if (psiChild.size()==2 && (abs(psiChild[0]->pdg_id()) == cut.type) &&
135  (abs(psiChild[1]->pdg_id()) == cut.type))
136  {
137  cout << psiChild[0]->momentum()<<" "<<psiChild[0]->momentum().eta()
138  <<" "<<psiChild[1]->momentum()<<" "<<psiChild[1]->momentum().eta()<<endl;
139  return ( (etaInRange(psiChild[0]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
140  (etaInRange(psiChild[1]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
141  (psiChild[0]->momentum().perp()> cut.ptMin) &&
142  (psiChild[1]->momentum().perp()> cut.ptMin));
143  }
144  return false;
145 }
146 */
147 
149 {
150  HepMC::GenVertex* myVertex = jpsi->end_vertex();
151  int numChildren = myVertex->particles_out_size();
152  std::vector<HepMC::GenParticle*> psiChild;
153  for(std::vector<GenParticle*>::const_iterator p = myVertex->particles_out_const_begin();
154  p != myVertex->particles_out_const_end(); p++)
155  psiChild.push_back((*p));
156 
157  if(numChildren>1) {
158  cout << psiChild[0]->pdg_id()<<" "<<psiChild[1]->pdg_id()<<endl;
159 
160  if (psiChild.size()==2 && (abs(psiChild[0]->pdg_id()) == cut.type) &&
161  (abs(psiChild[1]->pdg_id()) == cut.type))
162  {
163  cout << psiChild[0]->momentum().rho()<<" "<<psiChild[0]->momentum().eta()
164  <<" "<<psiChild[1]->momentum().rho()<<" "<<psiChild[1]->momentum().eta()<<endl;
165  return ( (etaInRange(psiChild[0]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
166  (etaInRange(psiChild[1]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
167  (psiChild[0]->momentum().perp()> cut.ptMin) &&
168  (psiChild[1]->momentum().perp()> cut.ptMin));
169  }
170  return false;
171  }
172  return false;
173 }
174 
175 bool BsJpsiPhiFilter::etaInRange(float eta, float etamin, float etamax)
176 {
177  return ( (etamin < eta) && (eta < etamax) );
178 }
virtual bool filter(edm::Event &, const edm::EventSetup &)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
#define abs(x)
Definition: mlp_lapack.h:159
HepMC::GenEvent::particle_const_iterator getNextBs(const HepMC::GenEvent::particle_const_iterator start, const HepMC::GenEvent::particle_const_iterator end)
T eta() const
HepMC::GenParticle * findParticle(HepMC::GenVertex *, const int requested_id)
BsJpsiPhiFilter(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:243
#define end
Definition: vmac.h:38
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
bool etaInRange(float eta, float etamin, float etamax)
#define jpsi
tuple cout
Definition: gather_cfg.py:121
bool cuts(const HepMC::GenParticle *jpsi, const CutStruct cut)
Definition: DDAxes.h:10