CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PythiaFilterGammaJet.cc
Go to the documentation of this file.
6 #include <iostream>
7 #include<list>
8 #include<vector>
9 #include<cmath>
10 
11 //using namespace edm;
12 //using namespace std;
13 
14 namespace{
15 
16  double deltaR2(double eta0, double phi0, double eta, double phi){
17  double dphi=phi-phi0;
18  if(dphi>M_PI) dphi-=2*M_PI;
19  else if(dphi<=-M_PI) dphi+=2*M_PI;
20  return dphi*dphi+(eta-eta0)*(eta-eta0);
21  }
22 
23  double deltaPhi(double phi0, double phi){
24  double dphi=phi-phi0;
25  if(dphi>M_PI) dphi-=2*M_PI;
26  else if(dphi<=-M_PI) dphi+=2*M_PI;
27  return dphi;
28  }
29 
30  class ParticlePtGreater{
31  public:
32  int operator()(const HepMC::GenParticle * p1,
33  const HepMC::GenParticle * p2) const{
34  return p1->momentum().perp() > p2->momentum().perp();
35  }
36  };
37 }
38 
39 
41 token_(consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared"))),
42 etaMax(iConfig.getUntrackedParameter<double>("MaxPhotonEta", 2.8)),
43 ptSeed(iConfig.getUntrackedParameter<double>("PhotonSeedPt", 5.)),
44 ptMin(iConfig.getUntrackedParameter<double>("MinPhotonPt")),
45 ptMax(iConfig.getUntrackedParameter<double>("MaxPhotonPt")),
46 dphiMin(iConfig.getUntrackedParameter<double>("MinDeltaPhi", -1)/180*M_PI),
47 detaMax(iConfig.getUntrackedParameter<double>("MaxDeltaEta", 10.)),
48 etaPhotonCut2(iConfig.getUntrackedParameter<double>("MinPhotonEtaForwardJet", 1.3)),
49 cone(0.5),ebEtaMax(1.479),
50 maxnumberofeventsinrun(iConfig.getUntrackedParameter<int>("MaxEvents",10000)){
51 
52  deltaEB=0.01745/2 *5; // delta_eta, delta_phi
53  deltaEE=2.93/317/2 *5; // delta_x/z, delta_y/z
55 }
56 
57 
59 
60 
61 // ------------ method called to produce the data ------------
63 
65  throw cms::Exception("endJob")<<"we have reached the maximum number of events ";
66  }
67 
68  bool accepted = false;
70  iEvent.getByToken(token_, evt);
71 
72  std::list<const HepMC::GenParticle *> seeds;
73  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
74 
75  if(myGenEvent->signal_process_id() == 14 || myGenEvent->signal_process_id() == 29) {
76 
77 
78  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
79 
80  if ( (*p)->pdg_id()==22 && (*p)->status()==1
81  && (*p)->momentum().perp() > ptSeed
82  && std::abs((*p)->momentum().eta()) < etaMax ) seeds.push_back(*p);
83 
84  }
85 
86  seeds.sort(ParticlePtGreater());
87  for(std::list<const HepMC::GenParticle *>::const_iterator is=
88  seeds.begin(); is!=seeds.end(); is++){
89 
90  double etaPhoton=(*is)->momentum().eta();
91  double phiPhoton=(*is)->momentum().phi();
92 
93  HepMC::GenEvent::particle_const_iterator ppp = myGenEvent->particles_begin();
94  for(int i=0;i<6;++i) ppp++;
95  HepMC::GenParticle* particle7 = (*ppp);
96  ppp++;
97  HepMC::GenParticle* particle8 = (*ppp);
98 
99  double dphi7=std::abs(deltaPhi(phiPhoton,
100  particle7->momentum().phi()));
101  double dphi=std::abs(deltaPhi(phiPhoton,
102  particle8->momentum().phi()));
103  int jetline=8;
104  if(dphi7>dphi) {
105  dphi=dphi7;
106  jetline=7;
107  }
108  if(dphi<dphiMin) continue;
109  //double etaJet= myGenEvent->particle(jetline)->momentum().eta();
110  double etaJet = 0.0;
111  if(jetline==8) etaJet = particle8->momentum().eta();
112  else etaJet = particle7->momentum().eta();
113 
114  double eta1=etaJet-detaMax;
115  double eta2=etaJet+detaMax;
116  if (eta1>etaPhotonCut2) eta1=etaPhotonCut2;
117  if (eta2<-etaPhotonCut2) eta2=-etaPhotonCut2;
118  if (etaPhoton<eta1 ||etaPhoton>eta2) continue;
119 
120  bool inEB(0);
121  double tgx(0);
122  double tgy(0);
123  if( std::abs(etaPhoton)<ebEtaMax) inEB=1;
124  else{
125  tgx=(*is)->momentum().px()/(*is)->momentum().pz();
126  tgy=(*is)->momentum().py()/(*is)->momentum().pz();
127  }
128 
129  double etPhoton=0;
130  double etPhotonCharged=0;
131  double etCone=0;
132  double etConeCharged=0;
133  double ptMaxHadron=0;
134 
135 
136  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
137 
138  if ( (*p)->status()!=1 ) continue;
139  int pid= (*p)->pdg_id();
140  int apid= std::abs(pid);
141  if (apid>11 && apid<20) continue; //get rid of muons and neutrinos
142  double eta=(*p)->momentum().eta();
143  double phi=(*p)->momentum().phi();
144  if (deltaR2(etaPhoton, phiPhoton, eta, phi)>cone*cone) continue;
145  double pt=(*p)->momentum().perp();
146 
147 
149  iSetup.getData( pdt );
150 
151 
152 // double charge=(*p)->particledata().charge();
153  //int charge3=(*p)->particleID().threeCharge();
154 
155  int charge3 = ((pdt->particle((*p)->pdg_id()))->ID().threeCharge());
156  etCone+=pt;
157  if(charge3 && pt<2) etConeCharged+=pt;
158 
159  //select particles matching a crystal array centered on photon
160  if(inEB) {
161  if( std::abs(eta-etaPhoton)> deltaEB ||
162  std::abs(deltaPhi(phi,phiPhoton)) > deltaEB) continue;
163  }
164  else if( std::abs((*p)->momentum().px()/(*p)->momentum().pz() - tgx)
165  > deltaEE ||
166  std::abs((*p)->momentum().py()/(*p)->momentum().pz() - tgy)
167  > deltaEE) continue;
168 
169  etPhoton+=pt;
170  if(charge3 && pt<2) etPhotonCharged+=pt;
171  if(apid>100 && apid!=310 && pt>ptMaxHadron) ptMaxHadron=pt;
172 
173  }
174 
175  if(etPhoton<ptMin ||etPhoton>ptMax) continue;
176 
177  //isolation cuts
178  if(etCone-etPhoton> 5+etPhoton/20-etPhoton*etPhoton/1e4) continue;
179  if(etCone-etPhoton-(etConeCharged-etPhotonCharged) >
180  3+etPhoton/20-etPhoton*etPhoton*etPhoton/1e6) continue;
181  if(ptMaxHadron > 4.5+etPhoton/40) continue;
182 
183  accepted=true;
184  break;
185 
186  } //loop over seeds
187 
188  } else {
189  // end of if(gammajetevent)
190  return true;
191  // accept all non-gammajet events
192  }
193 
194  if (accepted) {
196  return true;
197  }
198  else return false;
199 
200 }
201 
int i
Definition: DBlmapReader.cc:9
PythiaFilterGammaJet(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual bool filter(edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< edm::HepMCProduct > token_
void getData(T &iHolder) const
Definition: EventSetup.h:79
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
double p2[4]
Definition: TauolaWrapper.h:90
#define M_PI
bool operator()(const HepMC::GenParticle *a, const HepMC::GenParticle *b) const
tuple pid
Definition: sysUtil.py:22
Geom::Phi< T > phi() const
double p1[4]
Definition: TauolaWrapper.h:89