CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PythiaFilterGammaJetWithBg Class Reference

#include <PythiaFilterGammaJetWithBg.h>

Inheritance diagram for PythiaFilterGammaJetWithBg:
edm::EDFilter edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 
 PythiaFilterGammaJetWithBg (const edm::ParameterSet &)
 
 ~PythiaFilterGammaJetWithBg ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

double cone
 
double deltaEB
 
double deltaEE
 
double detaMax
 
double dphiMin
 
double ebEtaMax
 
double etaMax
 
double etaPhotonCut2
 
int maxnumberofeventsinrun
 
double ptMax
 
double ptMin
 
double ptSeed
 
int theNumberOfSelected
 
edm::EDGetTokenT
< edm::HepMCProduct
token_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

PythiaFilterGammaJet filter implements generator-level preselections for photon+jet like events to be used in jet energy calibration. Ported from fortran code written by V.Konoplianikov.

Author
A.Ulyanov, ITEP

Definition at line 30 of file PythiaFilterGammaJetWithBg.h.

Constructor & Destructor Documentation

PythiaFilterGammaJetWithBg::PythiaFilterGammaJetWithBg ( const edm::ParameterSet iConfig)
explicit

Definition at line 40 of file PythiaFilterGammaJetWithBg.cc.

References deltaEB, deltaEE, and theNumberOfSelected.

40  :
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",10)){
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 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::HepMCProduct > token_
#define M_PI
PythiaFilterGammaJetWithBg::~PythiaFilterGammaJetWithBg ( )

Definition at line 58 of file PythiaFilterGammaJetWithBg.cc.

58 {}

Member Function Documentation

bool PythiaFilterGammaJetWithBg::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDFilter.

Definition at line 62 of file PythiaFilterGammaJetWithBg.cc.

References funct::abs(), cone, gather_cfg::cout, deltaEB, deltaEE, SiPixelRawToDigiRegional_cfi::deltaPhi, reco::deltaR2(), detaMax, ebEtaMax, HLT_25ns14e33_v1_cff::etaJet, etaMax, etaPhotonCut2, GenParticle::GenParticle, edm::Event::getByToken(), edm::EventSetup::getData(), i, AlCaHLTBitMon_ParallelJobs::p, sysUtil::pid, EnergyCorrector::pt, ptMax, ptSeed, theNumberOfSelected, and token_.

62  {
63 
64 // <<<<<<< PythiaFilterGammaJetWithBg.cc
65 // if(theNumberOfSelected>=maxnumberofeventsinrun) {
66 // throw cms::Exception("endJob")<<"we have reached the maximum number of events ";
67 // }
68 // =======
69 // >>>>>>> 1.4
70 
71  bool accepted = false;
73  iEvent.getByToken(token_, evt);
74 
75  std::list<const HepMC::GenParticle *> seeds;
76  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
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 
88 // std::cout<<" Number of seeds "<<seeds.size()<<" ProccessID "<<myGenEvent->signal_process_id()<<std::endl;
89 // std::cout<<" ParticleId 7= "<<myGenEvent->particle(7)->pdg_id()
90 // <<" pT "<<myGenEvent->particle(7)->momentum().perp()
91 // <<" Eta "<<myGenEvent->particle(7)->momentum().eta()
92 // <<" Phi "<<myGenEvent->particle(7)->momentum().phi()<<std::endl;
93 // std::cout<<" ParticleId 8= "<<myGenEvent->particle(8)->pdg_id()<<" pT "<<myGenEvent->particle(8)->momentum().perp()<<" Eta "<<myGenEvent->particle(8)->momentum().eta()<<" Phi "<<myGenEvent->particle(8)->momentum().phi()<<std::endl;
94 
95  for(std::list<const HepMC::GenParticle *>::const_iterator is=
96  seeds.begin(); is!=seeds.end(); is++){
97 
98  double etaPhoton=(*is)->momentum().eta();
99  double phiPhoton=(*is)->momentum().phi();
100 
101  /*
102  double dphi7=std::abs(deltaPhi(phiPhoton,
103  myGenEvent->particle(7)->momentum().phi()));
104  double dphi=std::abs(deltaPhi(phiPhoton,
105  myGenEvent->particle(8)->momentum().phi()));
106  */
107 
108  //***
109  HepMC::GenEvent::particle_const_iterator ppp = myGenEvent->particles_begin();
110  for(int i=0;i<6;++i) ppp++;
111  HepMC::GenParticle* particle7 = (*ppp);
112  ppp++;
113  HepMC::GenParticle* particle8 = (*ppp);
114 
115  double dphi7=std::abs(deltaPhi(phiPhoton,
116  particle7->momentum().phi()));
117  double dphi=std::abs(deltaPhi(phiPhoton,
118  particle8->momentum().phi()));
119  //***
120 
121  int jetline=8;
122  if(dphi7>dphi) {
123  dphi=dphi7;
124  jetline=7;
125  }
126 
127 // std::cout<<" Dphi "<<dphi<<" "<<dphiMin<<std::endl;
128 // if(dphi<dphiMin) {
129 // std::cout<<"Reject dphi"<<std::endl;
130 // continue;
131 // }
132 
133  //double etaJet= myGenEvent->particle(jetline)->momentum().eta();
134  //***
135  double etaJet = 0.0;
136  if(jetline==8) etaJet = particle8->momentum().eta();
137  else etaJet = particle7->momentum().eta();
138  //***
139 
140  double eta1=etaJet-detaMax;
141  double eta2=etaJet+detaMax;
142  if (eta1>etaPhotonCut2) eta1=etaPhotonCut2;
143  if (eta2<-etaPhotonCut2) eta2=-etaPhotonCut2;
144 // std::cout<<" Etaphoton "<<etaPhoton<<" "<<eta1<<" "<<eta2<<std::endl;
145  if (etaPhoton<eta1 ||etaPhoton>eta2) {
146 // std::cout<<"Reject eta"<<std::endl;
147  continue;
148  }
149  bool inEB(0);
150  double tgx(0);
151  double tgy(0);
152  if( std::abs(etaPhoton)<ebEtaMax) inEB=1;
153  else{
154  tgx=(*is)->momentum().px()/(*is)->momentum().pz();
155  tgy=(*is)->momentum().py()/(*is)->momentum().pz();
156  }
157 
158  double etPhoton=0;
159  double etPhotonCharged=0;
160  double etCone=0;
161  double etConeCharged=0;
162  double ptMaxHadron=0;
163 
164 
165  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
166 
167  if ( (*p)->status()!=1 ) continue;
168  int pid= (*p)->pdg_id();
169  int apid= std::abs(pid);
170  if (apid>11 && apid<21) continue; //get rid of muons and neutrinos
171  double eta=(*p)->momentum().eta();
172  double phi=(*p)->momentum().phi();
173  if (deltaR2(etaPhoton, phiPhoton, eta, phi)>cone*cone) continue;
174  double pt=(*p)->momentum().perp();
175 
176  //***
178  iSetup.getData( pdt );
179 
180 
181  // double charge=(*p)->particledata().charge();
182  // int charge3=(*p)->particleID().threeCharge();
183 
184  int charge3 = ((pdt->particle((*p)->pdg_id()))->ID().threeCharge());
185  //***
186 
187  etCone+=pt;
188  if(charge3 && pt<2) etConeCharged+=pt;
189 
190  //select particles matching a crystal array centered on photon
191  if(inEB) {
192  if( std::abs(eta-etaPhoton)> deltaEB ||
193  std::abs(deltaPhi(phi,phiPhoton)) > deltaEB) continue;
194  }
195  else if( std::abs((*p)->momentum().px()/(*p)->momentum().pz() - tgx)
196  > deltaEE ||
197  std::abs((*p)->momentum().py()/(*p)->momentum().pz() - tgy)
198  > deltaEE) continue;
199 
200  etPhoton+=pt;
201  if(charge3 && pt<2) etPhotonCharged+=pt;
202  if(apid>100 && apid!=310 && pt>ptMaxHadron) ptMaxHadron=pt;
203 
204  }
205 // std::cout<<" etPhoton "<<etPhoton<<" "<<ptMin<<" "<<ptMax<<std::endl;
206 
207  if(etPhoton<ptMin ||etPhoton>ptMax) {
208 // std::cout<<" Reject etPhoton "<<std::endl;
209  continue;
210  }
211  //isolation cuts
212 
213 // double isocut1 = 5+etPhoton/20-etPhoton*etPhoton/1e4;
214  double isocut2 = 3+etPhoton/20-etPhoton*etPhoton*etPhoton/1e6;
215  double isocut3 = 4.5+etPhoton/40;
216  if (etPhoton>165.)
217  {
218 // isocut1 = 5.+165./20.-165.*165./1e4;
219  isocut2 = 3.+165./20.-165.*165.*165./1e6;
220  isocut3 = 4.5+165./40.;
221  }
222 
223 // std::cout<<" etCone "<<etCone<<" "<<etPhoton<<" "<<etCone-etPhoton<<" "<<isocut1<<std::endl;
224 // std::cout<<"Second cut on iso "<<etCone-etPhoton-(etConeCharged-etPhotonCharged)<<" cut value "<<isocut2<<" etPhoton "<<etPhoton<<std::endl;
225 // std::cout<<" PtHadron "<<ptMaxHadron<<" "<<4.5+etPhoton/40<<std::endl;
226 
227  if(etCone-etPhoton> 5+etPhoton/20-etPhoton*etPhoton/1e4) continue;
228 // std::cout<<" Accept 1"<<std::endl;
229  if(etCone-etPhoton-(etConeCharged-etPhotonCharged) >
230  isocut2) continue;
231 // std::cout<<" Accept 2"<<std::endl;
232  if(ptMaxHadron > isocut3) continue;
233 
234 // std::cout<<"Accept event "<<std::endl;
235  accepted=true;
236  break;
237 
238  } //loop over seeds
239 
240 
241  if (accepted) {
243  std::cout<<" Event preselected "<<theNumberOfSelected<<" Proccess ID "<<myGenEvent->signal_process_id()<<std::endl;
244  return true;
245  }
246  else return false;
247 
248 }
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< edm::HepMCProduct > token_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
void getData(T &iHolder) const
Definition: EventSetup.h:79
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
tuple pid
Definition: sysUtil.py:22
tuple cout
Definition: gather_cfg.py:121

Member Data Documentation

double PythiaFilterGammaJetWithBg::cone
private

Definition at line 48 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

double PythiaFilterGammaJetWithBg::deltaEB
private

Definition at line 50 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter(), and PythiaFilterGammaJetWithBg().

double PythiaFilterGammaJetWithBg::deltaEE
private

Definition at line 51 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter(), and PythiaFilterGammaJetWithBg().

double PythiaFilterGammaJetWithBg::detaMax
private

Definition at line 45 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

double PythiaFilterGammaJetWithBg::dphiMin
private

Definition at line 44 of file PythiaFilterGammaJetWithBg.h.

double PythiaFilterGammaJetWithBg::ebEtaMax
private

Definition at line 49 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

double PythiaFilterGammaJetWithBg::etaMax
private

Definition at line 40 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

double PythiaFilterGammaJetWithBg::etaPhotonCut2
private

Definition at line 46 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

int PythiaFilterGammaJetWithBg::maxnumberofeventsinrun
private

Definition at line 54 of file PythiaFilterGammaJetWithBg.h.

double PythiaFilterGammaJetWithBg::ptMax
private

Definition at line 43 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

double PythiaFilterGammaJetWithBg::ptMin
private

Definition at line 42 of file PythiaFilterGammaJetWithBg.h.

double PythiaFilterGammaJetWithBg::ptSeed
private

Definition at line 41 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().

int PythiaFilterGammaJetWithBg::theNumberOfSelected
private

Definition at line 53 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter(), and PythiaFilterGammaJetWithBg().

edm::EDGetTokenT<edm::HepMCProduct> PythiaFilterGammaJetWithBg::token_
private

Definition at line 39 of file PythiaFilterGammaJetWithBg.h.

Referenced by filter().