CMS 3D CMS Logo

PythiaFilterEMJet.cc
Go to the documentation of this file.
3 #include <iostream>
4 #include<list>
5 #include<map>
6 #include<vector>
7 #include<cmath>
8 
9 using namespace edm;
10 using namespace std;
11 
12 namespace{
13 
14  double deltaPhi(double phi0, double phi){
15  double dphi=phi-phi0;
16  if(dphi>M_PI) dphi-=2*M_PI;
17  else if(dphi<=-M_PI) dphi+=2*M_PI;
18  return dphi;
19  }
20 
21  class ParticlePtGreater{
22  public:
23  int operator()(const HepMC::GenParticle * p1,
24  const HepMC::GenParticle * p2) const{
25  return p1->momentum().perp() > p2->momentum().perp();
26  }
27  };
28 }
29 
30 
32 token_(consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"))),
33 etaMin(iConfig.getUntrackedParameter<double>("MinEMEta", 0)),
34 eTSumMin(iConfig.getUntrackedParameter<double>("ETSumMin", 50.)),
35 pTMin(iConfig.getUntrackedParameter<double>("MinEMpT", 5.)),
36 etaMax(iConfig.getUntrackedParameter<double>("MaxEMEta", 2.7)),
37 eTSumMax(iConfig.getUntrackedParameter<double>("ETSumMax", 100.)),
38 pTMax(iConfig.getUntrackedParameter<double>("MaxEMpT", 999999.)),
39 ebEtaMax(1.479),
40 maxnumberofeventsinrun(iConfig.getUntrackedParameter<int>("MaxEvents",3000000)){
41 
42  deltaEB=0.01745/2 *5; // delta_eta, delta_phi
43  deltaEE=2.93/317/2 *5; // delta_x/z, delta_y/z
46 
47  cout << " Max Events : " << maxnumberofeventsinrun << endl;
48  cout << " Cut Definition: " << endl;
49  cout << " MinEMEta = " << etaMin << endl;
50  cout << " ETSumMin = " << eTSumMin << endl;
51  cout << " MinEMpT = " << pTMin << endl;
52  cout << " MaxEMEta = " << etaMax << endl;
53  cout << " ETSumMax = " << eTSumMax << endl;
54  cout << " MaxEMpT = " << pTMax << endl;
55  cout << " 5x5 crystal cone around EM axis in ECAL Barrel = " << deltaEB << endl;
56  cout << " 5x5 crystal cone around EM axis in ECAL Endcap = " << deltaEE << endl;
57 
58 }
59 
61 {
62 std::cout << "Total number of tested events = " << theNumberOfTestedEvt << std::endl;
63 std::cout << "Total number of accepted events = " << theNumberOfSelected << std::endl;
64 }
65 
66 
67 // ------------ method called to produce the data ------------
69 
71  throw cms::Exception("endJob") << "we have reached the maximum number of events ";
72  }
73 
75  if(theNumberOfTestedEvt%1000 == 0) cout << "Number of tested events = " << theNumberOfTestedEvt << endl;
76 
77  bool accepted = false;
79  iEvent.getByToken(token_, evt);
80 
81  list<const HepMC::GenParticle *> EM_seeds;
82  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
83 
84  int particle_id = 1;
85 
86  //select e+/e-/gamma particles in the events
87  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
88  p != myGenEvent->particles_end();
89  ++p ) {
90 
91 
92  if ( (abs((*p)->pdg_id())==11 || (*p)->pdg_id()==22)
93  && (*p)->status()==1
94  && (*p)->momentum().perp() > pTMin
95  && (*p)->momentum().perp() < pTMax
96  && fabs((*p)->momentum().eta()) < etaMax
97  && fabs((*p)->momentum().eta()) > etaMin ) {
98  EM_seeds.push_back(*p);
99  }
100  particle_id++;
101  }
102 
103  EM_seeds.sort(ParticlePtGreater());
104 
105  double etaEMClus=0;
106  double phiEMClus=0;
107  double ptEMClus=0;
108  for(std::list<const HepMC::GenParticle *>::const_iterator is=EM_seeds.begin();
109  is!= EM_seeds.end();
110  ++is){
111  double etaEM=(*is)->momentum().eta();
112  double phiEM=(*is)->momentum().phi();
113  double ptEM=(*is)->momentum().perp();
114  //pass at the cluster the seed infos
115  etaEMClus = etaEM;
116  phiEMClus = phiEM;
117  ptEMClus = ptEM;
118 
119  //check if the EM particle is in the barrel
120  bool inEB(false);
121  double tgx(0);
122  double tgy(0);
123  if( std::abs(etaEM)<ebEtaMax ) inEB=true;
124  else{
125  tgx=(*is)->momentum().px()/(*is)->momentum().pz();
126  tgy=(*is)->momentum().py()/(*is)->momentum().pz();
127  }
128 
129 // std::vector<const HepMC::GenParticle*> takenEM ;
130 // std::vector<const HepMC::GenParticle*>::const_iterator itPart ;
131 
132  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
133  p != myGenEvent->particles_end();
134  ++p ) {
135 
136  if (((*p)->status()==1 && (*p)->pdg_id() == 22) || // gamma
137  ((*p)->status()==1 && abs((*p)->pdg_id()) == 11) ) // electron
138 // (*p)->pdg_id() == 111 || // pi0
139 // abs((*p)->pdg_id()) == 221 || // eta
140 // abs((*p)->pdg_id()) == 331 || // eta prime
141 // abs((*p)->pdg_id()) == 113 || // rho0
142 // abs((*p)->pdg_id()) == 223 ) // omega*/
143  {
144 // // check if found is daughter of one already taken
145 // bool isUsed = false ;
146 // const HepMC::GenParticle* mother = (*p)->production_vertex() ?
147 // *((*p)->production_vertex()->particles_in_const_begin()) : 0 ;
148 // const HepMC::GenParticle* motherMother = (mother != 0 && mother->production_vertex()) ?
149 // *(mother->production_vertex()->particles_in_const_begin()) : 0 ;
150 // const HepMC::GenParticle* motherMotherMother = (motherMother != 0 && motherMother->production_vertex()) ?
151 // *(motherMother->production_vertex()->particles_in_const_begin()) : 0 ;
152 // for(itPart = takenEM.begin(); itPart != takenEM.end(); ++itPart) {
153 // if ((*itPart) == mother ||
154 // (*itPart) == motherMother ||
155 // (*itPart) == motherMotherMother)
156 // {
157 // isUsed = true ;
158 // break ;
159 // }
160 // }
161 // if (!isUsed) takenEM.push_back(*p);
162 
163  double pt=(*p)->momentum().perp();
164  if (pt == ptEM) continue ; //discard the same particle of the seed
165  double eta=(*p)->momentum().eta();
166  double phi=(*p)->momentum().phi();
167 
168  if(inEB) {
169  if( std::abs(eta-etaEM)> deltaEB ||
170  std::abs(deltaPhi(phi,phiEM)) > deltaEB) continue;
171  }
172  else if( std::abs((*p)->momentum().px()/(*p)->momentum().pz() - tgx)
173  > deltaEE ||
174  std::abs((*p)->momentum().py()/(*p)->momentum().pz() - tgy)
175  > deltaEE) continue;
176  ptEMClus += pt ;
177  if(inEB)
178  {
179  etaEMClus += (eta-etaEMClus)*pt/ptEMClus ;
180  phiEMClus += deltaPhi(phi,phiEM)*pt/ptEMClus;
181  }
182  else
183  {
184  etaEMClus += ((*p)->momentum().px()/(*p)->momentum().pz() - tgx)*pt/ptEMClus ;
185  phiEMClus += ((*p)->momentum().py()/(*p)->momentum().pz() - tgy)*pt/ptEMClus;
186  }
187  }
188  }
189  if( ptEMClus > eTSumMin)
190  accepted = true ;
191  }
192 
193  if (accepted) {
195  cout << "========> Event preselected " << theNumberOfSelected
196  << " Proccess ID " << myGenEvent->signal_process_id() << endl;
197  return true;
198  }
199  else return false;
200 }
201 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< edm::HepMCProduct > token_
~PythiaFilterEMJet() override
int iEvent
Definition: GenABIO.cc:224
PythiaFilterEMJet(const edm::ParameterSet &)
bool filter(edm::Event &, const edm::EventSetup &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
#define M_PI
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
bool accepted(std::vector< std::string_view > const &, std::string_view)
HLT enums.
double p1[4]
Definition: TauolaWrapper.h:89