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
PythiaFilterGammaJet Class Reference

#include <PythiaFilterGammaJet.h>

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

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 
 PythiaFilterGammaJet (const edm::ParameterSet &)
 
 ~PythiaFilterGammaJet ()
 
- 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 PythiaFilterGammaJet.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file PythiaFilterGammaJet.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",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 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::HepMCProduct > token_
#define M_PI
PythiaFilterGammaJet::~PythiaFilterGammaJet ( )

Definition at line 58 of file PythiaFilterGammaJet.cc.

58 {}

Member Function Documentation

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

Implements edm::EDFilter.

Definition at line 62 of file PythiaFilterGammaJet.cc.

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

62  {
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 }
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< edm::HepMCProduct > token_
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

Member Data Documentation

double PythiaFilterGammaJet::cone
private

Definition at line 48 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::deltaEB
private

Definition at line 50 of file PythiaFilterGammaJet.h.

Referenced by filter(), and PythiaFilterGammaJet().

double PythiaFilterGammaJet::deltaEE
private

Definition at line 51 of file PythiaFilterGammaJet.h.

Referenced by filter(), and PythiaFilterGammaJet().

double PythiaFilterGammaJet::detaMax
private

Definition at line 45 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::dphiMin
private

Definition at line 44 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::ebEtaMax
private

Definition at line 49 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::etaMax
private

Definition at line 40 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::etaPhotonCut2
private

Definition at line 46 of file PythiaFilterGammaJet.h.

Referenced by filter().

int PythiaFilterGammaJet::maxnumberofeventsinrun
private

Definition at line 54 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::ptMax
private

Definition at line 43 of file PythiaFilterGammaJet.h.

Referenced by filter().

double PythiaFilterGammaJet::ptMin
private

Definition at line 42 of file PythiaFilterGammaJet.h.

double PythiaFilterGammaJet::ptSeed
private

Definition at line 41 of file PythiaFilterGammaJet.h.

Referenced by filter().

int PythiaFilterGammaJet::theNumberOfSelected
private

Definition at line 53 of file PythiaFilterGammaJet.h.

Referenced by filter(), and PythiaFilterGammaJet().

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

Definition at line 39 of file PythiaFilterGammaJet.h.

Referenced by filter().