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

#include <JetAnalyzer.h>

Inheritance diagram for JetAnalyzer:
JetAnalyzerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &, const reco::CaloJetCollection &caloJets, const int numPV)
 Get the analysis. More...
 
void beginJob (DQMStore *dbe)
 Inizialize parameters for histo binning. More...
 
void endJob ()
 Finish up a job. More...
 
int getLeadJetFlag ()
 
 JetAnalyzer (const edm::ParameterSet &)
 Constructor. More...
 
void setJetHiPass (int pass)
 
void setJetLoPass (int pass)
 
void setLeadJetFlag (int flag)
 
void setSource (std::string source)
 
virtual ~JetAnalyzer ()
 Destructor. More...
 
- Public Member Functions inherited from JetAnalyzerBase
void analyze (const edm::Event &, const edm::EventSetup &, reco::CaloJet &caloJet)
 Get the analysis of the muon properties. More...
 
 JetAnalyzerBase ()
 Constructor. More...
 
virtual ~JetAnalyzerBase ()
 Destructor. More...
 

Private Attributes

double _asymmetryThirdJetCut
 
double _balanceThirdJetCut
 
double _fHPDMax
 
double _fHPDMaxLoose
 
double _fHPDMaxTight
 
int _JetHiPass
 
int _JetLoPass
 
int _leadJetFlag
 
int _n90HitsMin
 
int _n90HitsMinLoose
 
int _n90HitsMinTight
 
double _ptThreshold
 
double _resEMFMin
 
double _resEMFMinLoose
 
double _resEMFMinTight
 
double _sigmaEtaMinTight
 
double _sigmaPhiMinTight
 
std::string _source
 
int _theend
 
int eBin
 
double eMax
 
double eMin
 
int etaBin
 
double etaMax
 
double etaMin
 
int fillJIDPassFrac
 
reco::helper::JetIDHelperjetID
 
MonitorElementjetME
 
std::string jetname
 
int makedijetselection
 
MonitorElementmConstituents
 
MonitorElementmConstituents_Barrel
 
MonitorElementmConstituents_Barrel_Hi
 
MonitorElementmConstituents_EndCap
 
MonitorElementmConstituents_EndCap_Hi
 
MonitorElementmConstituents_Forward
 
MonitorElementmConstituents_Forward_Hi
 
MonitorElementmConstituents_profile
 
MonitorElementmDijetAsymmetry
 
MonitorElementmDijetBalance
 
MonitorElementmDPhi
 
MonitorElementmEFrac
 
MonitorElementmEFrac_Barrel
 
MonitorElementmEFrac_EndCap
 
MonitorElementmEFrac_Forward
 
MonitorElementmEFrac_profile
 
MonitorElementmEmEnergyInEB
 
MonitorElementmEmEnergyInEE
 
MonitorElementmEmEnergyInHF
 
MonitorElementmEta
 
MonitorElementmEta_Hi
 
MonitorElementmEta_profile
 
MonitorElementmEtaFirst
 
MonitorElementmfHPD
 
MonitorElementmfRBX
 
MonitorElementmHadEnergyInHB
 
MonitorElementmHadEnergyInHE
 
MonitorElementmHadEnergyInHF
 
MonitorElementmHadEnergyInHO
 
MonitorElementmHFrac
 
MonitorElementmHFrac_Barrel
 
MonitorElementmHFrac_Barrel_Hi
 
MonitorElementmHFrac_EndCap
 
MonitorElementmHFrac_EndCap_Hi
 
MonitorElementmHFrac_Forward
 
MonitorElementmHFrac_Forward_Hi
 
MonitorElementmHFrac_profile
 
MonitorElementmLooseJIDPassFractionVSeta
 
MonitorElementmLooseJIDPassFractionVSpt
 
MonitorElementmMaxEInEmTowers
 
MonitorElementmMaxEInHadTowers
 
MonitorElementmN90Hits
 
MonitorElementmNJets
 
MonitorElementmNJets_profile
 
MonitorElementmPhi
 
MonitorElementmPhi_Barrel
 
MonitorElementmPhi_Barrel_Hi
 
MonitorElementmPhi_EndCap
 
MonitorElementmPhi_EndCap_Hi
 
MonitorElementmPhi_Forward
 
MonitorElementmPhi_Forward_Hi
 
MonitorElementmPhi_Hi
 
MonitorElementmPhi_Lo
 
MonitorElementmPhi_profile
 
MonitorElementmPhiFirst
 
MonitorElementmPhiVSEta
 
MonitorElementmPt
 
MonitorElementmPt_1
 
MonitorElementmPt_2
 
MonitorElementmPt_3
 
MonitorElementmPt_Barrel
 
MonitorElementmPt_Barrel_Hi
 
MonitorElementmPt_EndCap
 
MonitorElementmPt_EndCap_Hi
 
MonitorElementmPt_Forward
 
MonitorElementmPt_Forward_Hi
 
MonitorElementmPt_Hi
 
MonitorElementmPt_Lo
 
MonitorElementmPt_profile
 
MonitorElementmPtFirst
 
MonitorElementmresEMF
 
MonitorElementmTightJIDPassFractionVSeta
 
MonitorElementmTightJIDPassFractionVSpt
 
edm::ParameterSet parameters
 
int pBin
 
int phiBin
 
double phiMax
 
double phiMin
 
double pMax
 
double pMin
 
int ptBin
 
double ptMax
 
double ptMin
 
edm::InputTag theCaloJetCollectionLabel
 

Detailed Description

DQM monitoring source for Calo Jets

Author
F. Chlebana - Fermilab

Definition at line 39 of file JetAnalyzer.h.

Constructor & Destructor Documentation

JetAnalyzer::JetAnalyzer ( const edm::ParameterSet pSet)

Constructor.

Definition at line 20 of file JetAnalyzer.cc.

References Parameters::parameters.

20  {
21 
22  parameters = pSet;
23  _leadJetFlag = 0;
24  _JetLoPass = 0;
25  _JetHiPass = 0;
26  _ptThreshold = 20.;
28  _balanceThirdJetCut = 0.2;
29  _n90HitsMin =0;
30  _fHPDMax=1.;
31  _resEMFMin=0.;
33  _fHPDMaxLoose=1.;
34  _resEMFMinLoose=0.;
36  _fHPDMaxTight=1.;
37  _resEMFMinTight=0.;
38  _sigmaEtaMinTight=-999.;
39  _sigmaPhiMinTight=-999.;
40 
41 }
int _n90HitsMinTight
Definition: JetAnalyzer.h:110
int _n90HitsMinLoose
Definition: JetAnalyzer.h:107
double _fHPDMaxTight
Definition: JetAnalyzer.h:108
double _resEMFMin
Definition: JetAnalyzer.h:101
double _fHPDMax
Definition: JetAnalyzer.h:100
double _fHPDMaxLoose
Definition: JetAnalyzer.h:105
edm::ParameterSet parameters
Definition: JetAnalyzer.h:81
double _balanceThirdJetCut
Definition: JetAnalyzer.h:95
int _leadJetFlag
Definition: JetAnalyzer.h:90
double _resEMFMinTight
Definition: JetAnalyzer.h:109
double _sigmaEtaMinTight
Definition: JetAnalyzer.h:111
double _asymmetryThirdJetCut
Definition: JetAnalyzer.h:94
int _JetHiPass
Definition: JetAnalyzer.h:89
double _resEMFMinLoose
Definition: JetAnalyzer.h:106
int _JetLoPass
Definition: JetAnalyzer.h:88
double _sigmaPhiMinTight
Definition: JetAnalyzer.h:112
double _ptThreshold
Definition: JetAnalyzer.h:92
JetAnalyzer::~JetAnalyzer ( )
virtual

Destructor.

Definition at line 44 of file JetAnalyzer.cc.

44 { }

Member Function Documentation

void JetAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const reco::CaloJetCollection caloJets,
const int  numPV 
)

Get the analysis.

Definition at line 271 of file JetAnalyzer.cc.

References alignCSCRings::corr, diffTreeTool::diff, dPhi(), eta(), edm::EventID::event(), edm::EventBase::id(), metsig::jet, bTagSequences_cff::jetID, LogTrace, phi, RecoTauCleanerPlugins::pt, rand(), and mathSSE::sqrt().

273  {
274  int numofjets=0;
275  double fstPhi=0.;
276  double sndPhi=0.;
277  double diff = 0.;
278  double corr = 0.;
279  double dphi = -999. ;
280  bool thiscleaned=false;
281  bool Loosecleaned=false;
282  bool Tightcleaned=false;
283  bool thisemfclean=true;
284  bool emfcleanLoose=true;
285  bool emfcleanTight=true;
286 
287  srand( iEvent.id().event() % 10000);
288 
289 
290  if (makedijetselection == 1) {
291  //Dijet selection - careful: the pT is uncorrected!
292  //if(makedijetselection==1 && caloJets.size()>=2){
293  if(caloJets.size()>=2){
294  double dphiDJ = -999. ;
295  bool emfcleanLooseFirstJet=true;
296  bool emfcleanLooseSecondJet=true;
297  bool emfcleanTightFirstJet=true;
298  bool emfcleanTightSecondJet=true;
299  bool LoosecleanedFirstJet = false;
300  bool LoosecleanedSecondJet = false;
301  bool TightcleanedFirstJet = false;
302  bool TightcleanedSecondJet = false;
303  //both jets pass pt threshold
304  if ((caloJets.at(0)).pt() > _ptThreshold && (caloJets.at(1)).pt() > _ptThreshold ) {
305  if(fabs((caloJets.at(0)).eta())<3. && fabs((caloJets.at(1)).eta())<3. ){
306  //calculate dphi
307  dphiDJ = fabs((caloJets.at(0)).phi()-(caloJets.at(1)).phi());
308  if (dphiDJ > 3.14) dphiDJ=fabs(dphiDJ -6.28 );
309  //fill DPhi histo (before cutting)
310  if (mDPhi) mDPhi->Fill (dphiDJ);
311  //dphi cut
312  if(fabs(dphiDJ)>2.1){
313  //JetID
314  emfcleanLooseFirstJet=true;
315  emfcleanTightFirstJet=true;
316  emfcleanLooseSecondJet=true;
317  emfcleanTightSecondJet=true;
318  //jetID for first jet
319  jetID->calculate(iEvent, (caloJets.at(0)));
320  if(jetID->restrictedEMF()<_resEMFMinLoose && fabs((caloJets.at(0)).eta())<2.6) emfcleanLooseFirstJet=false;
321  if(jetID->restrictedEMF()<_resEMFMinTight && fabs((caloJets.at(0)).eta())<2.6) emfcleanTightFirstJet=false;
322  if(jetID->n90Hits()>=_n90HitsMinLoose && jetID->fHPD()<_fHPDMaxLoose && emfcleanLooseFirstJet) LoosecleanedFirstJet=true;
323  if(jetID->n90Hits()>=_n90HitsMinTight && jetID->fHPD()<_fHPDMaxTight && sqrt((caloJets.at(0)).etaetaMoment())>_sigmaEtaMinTight && sqrt((caloJets.at(0)).phiphiMoment())>_sigmaPhiMinTight && emfcleanTightFirstJet) TightcleanedFirstJet=true;
324  //fill the JID variables histograms BEFORE you cut on them
325  if (mN90Hits) mN90Hits->Fill (jetID->n90Hits());
326  if (mfHPD) mfHPD->Fill (jetID->fHPD());
328  if (mfRBX) mfRBX->Fill (jetID->fRBX());
329 
330  //jetID for second jet
331  jetID->calculate(iEvent, (caloJets.at(1)));
332  if(jetID->restrictedEMF()<_resEMFMinLoose && fabs((caloJets.at(1)).eta())<2.6) emfcleanLooseSecondJet=false;
333  if(jetID->restrictedEMF()<_resEMFMinTight && fabs((caloJets.at(1)).eta())<2.6) emfcleanTightSecondJet=false;
334  if(jetID->n90Hits()>=_n90HitsMinLoose && jetID->fHPD()<_fHPDMaxLoose && emfcleanLooseSecondJet) LoosecleanedSecondJet=true;
335  if(jetID->n90Hits()>=_n90HitsMinTight && jetID->fHPD()<_fHPDMaxTight && sqrt((caloJets.at(1)).etaetaMoment())>_sigmaEtaMinTight && sqrt((caloJets.at(1)).phiphiMoment())>_sigmaPhiMinTight && emfcleanTightSecondJet) TightcleanedSecondJet=true;
336  //fill the JID variables histograms BEFORE you cut on them
337  if (mN90Hits) mN90Hits->Fill (jetID->n90Hits());
338  if (mfHPD) mfHPD->Fill (jetID->fHPD());
340  if (mfRBX) mfRBX->Fill (jetID->fRBX());
341 
342  if(LoosecleanedFirstJet && LoosecleanedSecondJet) { //only if both jets are (loose) cleaned
343  //fill histos for first jet
344  if (mPt) mPt->Fill ((caloJets.at(0)).pt());
345  if (mEta) mEta->Fill ((caloJets.at(0)).eta());
346  if (mPhi) mPhi->Fill ((caloJets.at(0)).phi());
347  if (mPhiVSEta) mPhiVSEta->Fill((caloJets.at(0)).eta(),(caloJets.at(0)).phi());
348  if (mConstituents) mConstituents->Fill ((caloJets.at(0)).nConstituents());
349  if (mHFrac) mHFrac->Fill ((caloJets.at(0)).energyFractionHadronic());
350  if (mEFrac) mEFrac->Fill ((caloJets.at(0)).emEnergyFraction());
351  //if (mE) mE->Fill ((caloJets.at(0)).energy());
352  //if (mP) mP->Fill ((caloJets.at(0)).p());
353  //if (mMass) mMass->Fill ((caloJets.at(0)).mass());
354  if (mMaxEInEmTowers) mMaxEInEmTowers->Fill ((caloJets.at(0)).maxEInEmTowers());
355  if (mMaxEInHadTowers) mMaxEInHadTowers->Fill ((caloJets.at(0)).maxEInHadTowers());
356  if (mN90Hits) mN90Hits->Fill (jetID->n90Hits());
357  if (mfHPD) mfHPD->Fill (jetID->fHPD());
359  if (mfRBX) mfRBX->Fill (jetID->fRBX());
360  //sigmaeta and sigmaphi only used in the tight selection.
361  //fill the histos for them AFTER the loose selection
362  // if (msigmaEta) msigmaEta->Fill(sqrt((caloJets.at(0)).etaetaMoment()));
363  // if (msigmaPhi) msigmaPhi->Fill(sqrt((caloJets.at(0)).phiphiMoment()));
364  //fill histos for second jet
365  if (mPt) mPt->Fill ((caloJets.at(1)).pt());
366  if (mEta) mEta->Fill ((caloJets.at(1)).eta());
367  if (mPhi) mPhi->Fill ((caloJets.at(1)).phi());
368  if (mPhiVSEta) mPhiVSEta->Fill((caloJets.at(1)).eta(),(caloJets.at(1)).phi());
369  if (mConstituents) mConstituents->Fill ((caloJets.at(1)).nConstituents());
370  if (mHFrac) mHFrac->Fill ((caloJets.at(1)).energyFractionHadronic());
371  if (mEFrac) mEFrac->Fill ((caloJets.at(1)).emEnergyFraction());
372  //if (mE) mE->Fill ((caloJets.at(1)).energy());
373  //if (mP) mP->Fill ((caloJets.at(1)).p());
374  //if (mMass) mMass->Fill ((caloJets.at(1)).mass());
375  if (mMaxEInEmTowers) mMaxEInEmTowers->Fill ((caloJets.at(1)).maxEInEmTowers());
376  if (mMaxEInHadTowers) mMaxEInHadTowers->Fill ((caloJets.at(1)).maxEInHadTowers());
377  //sigmaeta and sigmaphi only used in the tight selection.
378  //fill the histos for them AFTER the loose selection
379  // if (msigmaEta) msigmaEta->Fill(sqrt((caloJets.at(1)).etaetaMoment()));
380  // if (msigmaPhi) msigmaPhi->Fill(sqrt((caloJets.at(1)).phiphiMoment()));
381 
382 
383  // Fill NPV profiles
384  //----------------------------------------------------------------
385  for (int ijet=0; ijet<2; ijet++) {
386 
387  if (mPt_profile) mPt_profile ->Fill(numPV, (caloJets.at(ijet)).pt());
388  if (mEta_profile) mEta_profile ->Fill(numPV, (caloJets.at(ijet)).eta());
389  if (mPhi_profile) mPhi_profile ->Fill(numPV, (caloJets.at(ijet)).phi());
390  if (mConstituents_profile) mConstituents_profile->Fill(numPV, (caloJets.at(ijet)).nConstituents());
391  if (mHFrac_profile) mHFrac_profile ->Fill(numPV, (caloJets.at(ijet)).energyFractionHadronic());
392  if (mEFrac_profile) mEFrac_profile ->Fill(numPV, (caloJets.at(ijet)).emEnergyFraction());
393  }
394  }
395 
396 
397  //let's see how many of these jets passed the JetID cleaning
398  if(fillJIDPassFrac==1) {
399  if(LoosecleanedFirstJet) {
401  mLooseJIDPassFractionVSpt->Fill((caloJets.at(0)).pt(),1.);
402  } else {
404  mLooseJIDPassFractionVSpt->Fill((caloJets.at(0)).pt(),0.);
405  }
406  if(LoosecleanedSecondJet) {
408  mLooseJIDPassFractionVSpt->Fill((caloJets.at(1)).pt(),1.);
409  } else {
411  mLooseJIDPassFractionVSpt->Fill((caloJets.at(1)).pt(),0.);
412  }
413  //TIGHT JID
414  if(TightcleanedFirstJet) {
416  mTightJIDPassFractionVSpt->Fill((caloJets.at(0)).pt(),1.);
417  } else {
419  mTightJIDPassFractionVSpt->Fill((caloJets.at(0)).pt(),0.);
420  }
421  if(TightcleanedSecondJet) {
423  mTightJIDPassFractionVSpt->Fill((caloJets.at(1)).pt(),1.);
424  } else {
426  mTightJIDPassFractionVSpt->Fill((caloJets.at(1)).pt(),0.);
427  }
428 
429  }//if fillJIDPassFrac
430  }// FABS DPHI < 2.1
431  }// fabs eta < 3
432  }// pt jets > threshold
433  //now do the dijet balance and asymmetry calculations
434  if (fabs(caloJets.at(0).eta() < 1.4)) {
435  double pt_dijet = (caloJets.at(0).pt() + caloJets.at(1).pt())/2;
436 
437  double dPhi = fabs((caloJets.at(0)).phi()-(caloJets.at(1)).phi());
438  if (dPhi > 3.14) dPhi=fabs(dPhi -6.28 );
439 
440  if (dPhi > 2.7) {
441  double pt_probe;
442  double pt_barrel;
443  int jet1, jet2;
444 
445  int randJet = rand() % 2;
446 
447  if (fabs(caloJets.at(1).eta() < 1.4)) {
448  if (randJet) {
449  jet1 = 0;
450  jet2 = 1;
451  }
452  else {
453  jet1 = 1;
454  jet2 = 0;
455  }
456 
457  /***Di-Jet Asymmetry****
458  * leading jets eta < 1.4
459  * leading jets dphi > 2.7
460  * pt_third jet < threshold
461  * A = (pt_1 - pt_2)/(pt_1 + pt_2)
462  * jets 1 and two are randomly ordered
463  */
464  bool thirdJetCut = true;
465  for (unsigned int third = 2; third < caloJets.size(); ++third)
466  if (caloJets.at(third).pt() > _asymmetryThirdJetCut)
467  thirdJetCut = false;
468  if (thirdJetCut) {
469  double dijetAsymmetry = (caloJets.at(jet1).pt() - caloJets.at(jet2).pt()) / (caloJets.at(jet1).pt() + caloJets.at(jet2).pt());
470  mDijetAsymmetry->Fill(dijetAsymmetry);
471  }// end restriction on third jet pt in asymmetry calculation
472 
473  }
474  else {
475  jet1 = 0;
476  jet2 = 1;
477  }
478 
479  pt_barrel = caloJets.at(jet1).pt();
480  pt_probe = caloJets.at(jet2).pt();
481 
482  //dijet balance cuts
483  /***Di-Jet Balance****
484  * pt_dijet = (pt_probe+pt_barrel)/2
485  * leading jets dphi > 2.7
486  * reject evnets where pt_third/pt_dijet > 0.2
487  * pv selection
488  * B = (pt_probe - pt_barrel)/pt_dijet
489  * select probe randomly from 2 jets if both leading jets are in the barrel
490  */
491  bool thirdJetCut = true;
492  for (unsigned int third = 2; third < caloJets.size(); ++third)
493  if (caloJets.at(third).pt()/pt_dijet > _balanceThirdJetCut)
494  thirdJetCut = false;
495  if (thirdJetCut) {
496  double dijetBalance = (pt_probe - pt_barrel) / pt_dijet;
497  mDijetBalance->Fill(dijetBalance);
498  }// end restriction on third jet pt ratio in balance calculation
499  }// dPhi > 2.7
500  }// leading jet eta cut for asymmetry and balance calculations
501  }//jet size >= 2
502  }// do dijet selection
503  else {
504  for (reco::CaloJetCollection::const_iterator jet = caloJets.begin(); jet!=caloJets.end(); ++jet) {
505  LogTrace(jetname)<<"[JetAnalyzer] Analyze Calo Jet";
506  Loosecleaned=false;
507  Tightcleaned=false;
508  if (jet == caloJets.begin()) {
509  fstPhi = jet->phi();
510  _leadJetFlag = 1;
511  } else {
512  _leadJetFlag = 0;
513  }
514  if (jet == (caloJets.begin()+1)) sndPhi = jet->phi();
515  //jetID
516  jetID->calculate(iEvent, *jet);
517  //minimal (uncorrected!) pT cut
518  if (jet->pt() > _ptThreshold) {
519  // if (msigmaEta) msigmaEta->Fill(sqrt(jet->etaetaMoment()));
520  // if (msigmaPhi) msigmaPhi->Fill(sqrt(jet->phiphiMoment()));
521  //cleaning to use for filling histograms
522  thisemfclean=true;
523  if(jetID->restrictedEMF()<_resEMFMin && fabs(jet->eta())<2.6) thisemfclean=false;
524  if(jetID->n90Hits()>=_n90HitsMin && jetID->fHPD()<_fHPDMax && thisemfclean) thiscleaned=true;
525  //loose and tight cleaning, used to fill the JetIDPAssFraction histos
526  if(jetID->n90Hits()>=_n90HitsMinLoose && jetID->fHPD()<_fHPDMaxLoose && emfcleanLoose) Loosecleaned=true;
527  if(jetID->n90Hits()>=_n90HitsMinTight && jetID->fHPD()<_fHPDMaxTight && sqrt(jet->etaetaMoment())>_sigmaEtaMinTight && sqrt(jet->phiphiMoment())>_sigmaPhiMinTight && emfcleanTight) Tightcleaned=true;
528 
529  if(fillJIDPassFrac==1) {
530  if(Loosecleaned) {
531  mLooseJIDPassFractionVSeta->Fill(jet->eta(),1.);
533  } else {
534  mLooseJIDPassFractionVSeta->Fill(jet->eta(),0.);
536  }
537  //TIGHT
538  if(Tightcleaned) {
539  mTightJIDPassFractionVSeta->Fill(jet->eta(),1.);
541  } else {
542  mTightJIDPassFractionVSeta->Fill(jet->eta(),0.);
544  }
545  }
546  //eventually we could define the "cleaned" flag differently for e.g. HF
547  if(thiscleaned) {
548  numofjets++ ;
549  jetME->Fill(1);
550 
551  // Leading jet
552  // Histograms are filled once per event
553  if (_leadJetFlag == 1) {
554  if (mEtaFirst) mEtaFirst->Fill (jet->eta());
555  if (mPhiFirst) mPhiFirst->Fill (jet->phi());
556  //if (mEFirst) mEFirst->Fill (jet->energy());
557  if (mPtFirst) mPtFirst->Fill (jet->pt());
558  }
559  // --- Passed the low pt jet trigger
560  if (_JetLoPass == 1) {
561  /* if (fabs(jet->eta()) <= 1.3) {
562  if (mPt_Barrel_Lo) mPt_Barrel_Lo->Fill(jet->pt());
563  if (mEta_Lo) mEta_Lo->Fill(jet->eta());
564  if (mPhi_Barrel_Lo) mPhi_Barrel_Lo->Fill(jet->phi());
565  }
566  if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
567  if (mPt_EndCap_Lo) mPt_EndCap_Lo->Fill(jet->pt());
568  if (mEta_Lo) mEta_Lo->Fill(jet->eta());
569  if (mPhi_EndCap_Lo) mPhi_EndCap_Lo->Fill(jet->phi());
570  }
571  if (fabs(jet->eta()) > 3.0) {
572  if (mPt_Forward_Lo) mPt_Forward_Lo->Fill(jet->pt());
573  if (mEta_Lo) mEta_Lo->Fill(jet->eta());
574  if (mPhi_Forward_Lo) mPhi_Forward_Lo->Fill(jet->phi());
575  } */
576  //if (mEta_Lo) mEta_Lo->Fill (jet->eta());
577  if (mPhi_Lo) mPhi_Lo->Fill (jet->phi());
578  if (mPt_Lo) mPt_Lo->Fill (jet->pt());
579  }
580 
581  // --- Passed the high pt jet trigger
582  if (_JetHiPass == 1) {
583  if (fabs(jet->eta()) <= 1.3) {
584  if (mPt_Barrel_Hi && jet->pt()>100.) mPt_Barrel_Hi->Fill(jet->pt());
585  if (mEta_Hi && jet->pt()>100.) mEta_Hi->Fill(jet->eta());
586  if (mPhi_Barrel_Hi) mPhi_Barrel_Hi->Fill(jet->phi());
587  //if (mConstituents_Barrel_Hi) mConstituents_Barrel_Hi->Fill(jet->nConstituents());
588  //if (mHFrac_Barrel_Hi) mHFrac_Barrel_Hi->Fill(jet->energyFractionHadronic());
589  }
590  if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
591  if (mPt_EndCap_Hi && jet->pt()>100.) mPt_EndCap_Hi->Fill(jet->pt());
592  if (mEta_Hi && jet->pt()>100.) mEta_Hi->Fill(jet->eta());
593  if (mPhi_EndCap_Hi) mPhi_EndCap_Hi->Fill(jet->phi());
594  //if (mConstituents_EndCap_Hi) mConstituents_EndCap_Hi->Fill(jet->nConstituents());
595  //if (mHFrac_EndCap_Hi) mHFrac_EndCap_Hi->Fill(jet->energyFractionHadronic());
596  }
597  if (fabs(jet->eta()) > 3.0) {
598  if (mPt_Forward_Hi && jet->pt()>100.) mPt_Forward_Hi->Fill(jet->pt());
599  if (mEta_Hi && jet->pt()>100.) mEta_Hi->Fill(jet->eta());
600  if (mPhi_Forward_Hi) mPhi_Forward_Hi->Fill(jet->phi());
601  //if (mConstituents_Forward_Hi) mConstituents_Forward_Hi->Fill(jet->nConstituents());
602  //if (mHFrac_Forward_Hi) mHFrac_Forward_Hi->Fill(jet->energyFractionHadronic());
603  }
604 
605  if (mEta_Hi && jet->pt()>100.) mEta_Hi->Fill (jet->eta());
606  if (mPhi_Hi) mPhi_Hi->Fill (jet->phi());
607  if (mPt_Hi) mPt_Hi->Fill (jet->pt());
608  }
609 
610  if (mPt) mPt->Fill (jet->pt());
611  if (mPt_1) mPt_1->Fill (jet->pt());
612  if (mPt_2) mPt_2->Fill (jet->pt());
613  if (mPt_3) mPt_3->Fill (jet->pt());
614  if (mEta) mEta->Fill (jet->eta());
615  if (mPhi) mPhi->Fill (jet->phi());
616 
617  if (mPhiVSEta) mPhiVSEta->Fill(jet->eta(),jet->phi());
618 
619  if (mConstituents) mConstituents->Fill (jet->nConstituents());
620  if (mHFrac) mHFrac->Fill (jet->energyFractionHadronic());
621  if (mEFrac) mEFrac->Fill (jet->emEnergyFraction());
622 
623 
624  // Fill NPV profiles
625  //--------------------------------------------------------------------
626  if (mPt_profile) mPt_profile ->Fill(numPV, jet->pt());
627  if (mEta_profile) mEta_profile ->Fill(numPV, jet->eta());
628  if (mPhi_profile) mPhi_profile ->Fill(numPV, jet->phi());
629  if (mConstituents_profile) mConstituents_profile->Fill(numPV, jet->nConstituents());
630  if (mHFrac_profile) mHFrac_profile ->Fill(numPV, jet->energyFractionHadronic());
631  if (mEFrac_profile) mEFrac_profile ->Fill(numPV, jet->emEnergyFraction());
632 
633 
634  if (fabs(jet->eta()) <= 1.3) {
635  if (mPt_Barrel) mPt_Barrel->Fill (jet->pt());
636  if (mPhi_Barrel) mPhi_Barrel->Fill (jet->phi());
637  //if (mE_Barrel) mE_Barrel->Fill (jet->energy());
638  if (mConstituents_Barrel) mConstituents_Barrel->Fill(jet->nConstituents());
639  if (mHFrac_Barrel) mHFrac_Barrel->Fill(jet->energyFractionHadronic());
640  if (mEFrac_Barrel) mEFrac_Barrel->Fill(jet->emEnergyFraction());
641  }
642  if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
643  if (mPt_EndCap) mPt_EndCap->Fill (jet->pt());
644  if (mPhi_EndCap) mPhi_EndCap->Fill (jet->phi());
645  //if (mE_EndCap) mE_EndCap->Fill (jet->energy());
646  if (mConstituents_EndCap) mConstituents_EndCap->Fill(jet->nConstituents());
647  if (mHFrac_EndCap) mHFrac_EndCap->Fill(jet->energyFractionHadronic());
648  if (mEFrac_EndCap) mEFrac_EndCap->Fill(jet->emEnergyFraction());
649  }
650  if (fabs(jet->eta()) > 3.0) {
651  if (mPt_Forward) mPt_Forward->Fill (jet->pt());
652  if (mPhi_Forward) mPhi_Forward->Fill (jet->phi());
653  //if (mE_Forward) mE_Forward->Fill (jet->energy());
654  if (mConstituents_Forward) mConstituents_Forward->Fill(jet->nConstituents());
655  if (mHFrac_Forward) mHFrac_Forward->Fill(jet->energyFractionHadronic());
656  if (mEFrac_Forward) mEFrac_Forward->Fill(jet->emEnergyFraction());
657  }
658 
659  //if (mE) mE->Fill (jet->energy());
660  //if (mP) mP->Fill (jet->p());
661  // if (mMass) mMass->Fill (jet->mass());
662 
663  if (mMaxEInEmTowers) mMaxEInEmTowers->Fill (jet->maxEInEmTowers());
664  if (mMaxEInHadTowers) mMaxEInHadTowers->Fill (jet->maxEInHadTowers());
665 
666  if (mHadEnergyInHO) mHadEnergyInHO->Fill (jet->hadEnergyInHO());
667  if (mHadEnergyInHB) mHadEnergyInHB->Fill (jet->hadEnergyInHB());
668  if (mHadEnergyInHF) mHadEnergyInHF->Fill (jet->hadEnergyInHF());
669  if (mHadEnergyInHE) mHadEnergyInHE->Fill (jet->hadEnergyInHE());
670  if (mEmEnergyInEB) mEmEnergyInEB->Fill (jet->emEnergyInEB());
671  if (mEmEnergyInEE) mEmEnergyInEE->Fill (jet->emEnergyInEE());
672  if (mEmEnergyInHF) mEmEnergyInHF->Fill (jet->emEnergyInHF());
673 
674  if (mN90Hits) mN90Hits->Fill (jetID->n90Hits());
675  if (mfHPD) mfHPD->Fill (jetID->fHPD());
677  if (mfRBX) mfRBX->Fill (jetID->fRBX());
678 
679  //calculate correctly the dphi
680  if(numofjets>1) {
681  diff = fabs(fstPhi - sndPhi);
682  corr = 2*acos(-1.) - diff;
683  if(diff < acos(-1.)) {
684  dphi = diff;
685  } else {
686  dphi = corr;
687  }
688  }
689  }
690  }//pt cut
691  }
692 
693  if (mNJets) mNJets->Fill (numofjets);
694  if (mDPhi && dphi>-998.) mDPhi->Fill (dphi);
695 
696 
697  if (mNJets_profile) mNJets_profile->Fill(numPV, numofjets);
698 
699 
700  }//not dijet
701 }
MonitorElement * mEFrac_EndCap
Definition: JetAnalyzer.h:184
MonitorElement * mNJets
Definition: JetAnalyzer.h:213
EventNumber_t event() const
Definition: EventID.h:44
MonitorElement * mPt_Forward_Hi
Definition: JetAnalyzer.h:199
int _n90HitsMinTight
Definition: JetAnalyzer.h:110
MonitorElement * mHadEnergyInHB
Definition: JetAnalyzer.h:227
MonitorElement * mMaxEInEmTowers
Definition: JetAnalyzer.h:224
MonitorElement * mfRBX
Definition: JetAnalyzer.h:237
MonitorElement * mPhi_Forward
Definition: JetAnalyzer.h:173
int _n90HitsMinLoose
Definition: JetAnalyzer.h:107
MonitorElement * mPt_Lo
Definition: JetAnalyzer.h:249
MonitorElement * mEmEnergyInHF
Definition: JetAnalyzer.h:232
MonitorElement * mPhiFirst
Definition: JetAnalyzer.h:218
MonitorElement * mPt_Forward
Definition: JetAnalyzer.h:172
MonitorElement * mHFrac_Forward
Definition: JetAnalyzer.h:188
MonitorElement * mEmEnergyInEE
Definition: JetAnalyzer.h:231
int makedijetselection
Definition: JetAnalyzer.h:97
MonitorElement * jetME
Definition: JetAnalyzer.h:136
MonitorElement * mPhi_EndCap
Definition: JetAnalyzer.h:170
double fHPD() const
Definition: JetIDHelper.h:34
MonitorElement * mN90Hits
Definition: JetAnalyzer.h:235
MonitorElement * mEta_Hi
Definition: JetAnalyzer.h:251
MonitorElement * mEFrac_profile
Definition: JetAnalyzer.h:268
MonitorElement * mPt_3
Definition: JetAnalyzer.h:158
double restrictedEMF() const
Definition: JetIDHelper.h:51
T eta() const
MonitorElement * mEta
Definition: JetAnalyzer.h:159
double _fHPDMaxTight
Definition: JetAnalyzer.h:108
MonitorElement * mConstituents_Forward
Definition: JetAnalyzer.h:187
MonitorElement * mPt_Hi
Definition: JetAnalyzer.h:253
double _resEMFMin
Definition: JetAnalyzer.h:101
double _fHPDMax
Definition: JetAnalyzer.h:100
MonitorElement * mTightJIDPassFractionVSpt
Definition: JetAnalyzer.h:244
MonitorElement * mLooseJIDPassFractionVSpt
Definition: JetAnalyzer.h:242
double _fHPDMaxLoose
Definition: JetAnalyzer.h:105
MonitorElement * mDPhi
Definition: JetAnalyzer.h:214
void Fill(long long x)
MonitorElement * mresEMF
Definition: JetAnalyzer.h:238
MonitorElement * mConstituents_profile
Definition: JetAnalyzer.h:266
MonitorElement * mPt_Barrel_Hi
Definition: JetAnalyzer.h:191
MonitorElement * mEmEnergyInEB
Definition: JetAnalyzer.h:230
double _balanceThirdJetCut
Definition: JetAnalyzer.h:95
MonitorElement * mDijetBalance
Definition: JetAnalyzer.h:256
MonitorElement * mPhi_Barrel
Definition: JetAnalyzer.h:167
MonitorElement * mConstituents_EndCap
Definition: JetAnalyzer.h:182
MonitorElement * mPhi_EndCap_Hi
Definition: JetAnalyzer.h:196
int _leadJetFlag
Definition: JetAnalyzer.h:90
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
MonitorElement * mPt_1
Definition: JetAnalyzer.h:156
T sqrt(T t)
Definition: SSEVec.h:48
MonitorElement * mConstituents_Barrel
Definition: JetAnalyzer.h:177
std::string jetname
Definition: JetAnalyzer.h:83
MonitorElement * mPhi
Definition: JetAnalyzer.h:160
MonitorElement * mPt_profile
Definition: JetAnalyzer.h:263
MonitorElement * mEtaFirst
Definition: JetAnalyzer.h:217
MonitorElement * mEFrac_Forward
Definition: JetAnalyzer.h:189
double _resEMFMinTight
Definition: JetAnalyzer.h:109
MonitorElement * mPt_2
Definition: JetAnalyzer.h:157
MonitorElement * mHFrac_Barrel
Definition: JetAnalyzer.h:178
MonitorElement * mNJets_profile
Definition: JetAnalyzer.h:262
#define LogTrace(id)
MonitorElement * mEFrac
Definition: JetAnalyzer.h:163
double _sigmaEtaMinTight
Definition: JetAnalyzer.h:111
MonitorElement * mHadEnergyInHF
Definition: JetAnalyzer.h:228
MonitorElement * mfHPD
Definition: JetAnalyzer.h:236
MonitorElement * mDijetAsymmetry
Definition: JetAnalyzer.h:257
MonitorElement * mPt
Definition: JetAnalyzer.h:155
MonitorElement * mEFrac_Barrel
Definition: JetAnalyzer.h:179
double _asymmetryThirdJetCut
Definition: JetAnalyzer.h:94
MonitorElement * mPhi_Barrel_Hi
Definition: JetAnalyzer.h:192
MonitorElement * mHFrac_EndCap
Definition: JetAnalyzer.h:183
MonitorElement * mPhiVSEta
Definition: JetAnalyzer.h:164
int _JetHiPass
Definition: JetAnalyzer.h:89
MonitorElement * mTightJIDPassFractionVSeta
Definition: JetAnalyzer.h:243
MonitorElement * mLooseJIDPassFractionVSeta
Definition: JetAnalyzer.h:241
MonitorElement * mHFrac_profile
Definition: JetAnalyzer.h:267
reco::helper::JetIDHelper * jetID
Definition: JetAnalyzer.h:139
edm::EventID id() const
Definition: EventBase.h:56
MonitorElement * mPt_EndCap_Hi
Definition: JetAnalyzer.h:195
int fillJIDPassFrac
Definition: JetAnalyzer.h:104
double _resEMFMinLoose
Definition: JetAnalyzer.h:106
MonitorElement * mMaxEInHadTowers
Definition: JetAnalyzer.h:225
MonitorElement * mPhi_Hi
Definition: JetAnalyzer.h:252
Signal rand(Signal arg)
Definition: vlib.cc:442
MonitorElement * mPtFirst
Definition: JetAnalyzer.h:220
MonitorElement * mPhi_profile
Definition: JetAnalyzer.h:265
MonitorElement * mPt_Barrel
Definition: JetAnalyzer.h:166
MonitorElement * mHadEnergyInHE
Definition: JetAnalyzer.h:229
MonitorElement * mConstituents
Definition: JetAnalyzer.h:161
MonitorElement * mPhi_Forward_Hi
Definition: JetAnalyzer.h:200
MonitorElement * mEta_profile
Definition: JetAnalyzer.h:264
MonitorElement * mPt_EndCap
Definition: JetAnalyzer.h:169
int _JetLoPass
Definition: JetAnalyzer.h:88
double _sigmaPhiMinTight
Definition: JetAnalyzer.h:112
MonitorElement * mPhi_Lo
Definition: JetAnalyzer.h:248
MonitorElement * mHFrac
Definition: JetAnalyzer.h:162
void calculate(const edm::Event &event, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:86
double fRBX() const
Definition: JetIDHelper.h:35
double _ptThreshold
Definition: JetAnalyzer.h:92
Definition: DDAxes.h:10
MonitorElement * mHadEnergyInHO
Definition: JetAnalyzer.h:226
void JetAnalyzer::beginJob ( DQMStore dbe)
virtual

Inizialize parameters for histo binning.

Implements JetAnalyzerBase.

Definition at line 48 of file JetAnalyzer.cc.

References DQMStore::book1D(), DQMStore::book2D(), DQMStore::bookProfile(), jptDQMConfig_cff::eMax, jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, bTagSequences_cff::jetID, LogTrace, nbinsPV, Parameters::parameters, jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, jptDQMConfig_cff::pMax, jptDQMConfig_cff::ptMax, PtMinSelector_cfg::ptMin, PVlow, PVup, MonitorElement::setAxisTitle(), MonitorElement::setBinLabel(), and DQMStore::setCurrentFolder().

48  {
49 
50  jetname = "jetAnalyzer";
51 
52  LogTrace(jetname)<<"[JetAnalyzer] Parameters initialization";
53  dbe->setCurrentFolder("JetMET/Jet/"+_source);
54 
55  jetME = dbe->book1D("jetReco", "jetReco", 3, 1, 4);
56  jetME->setBinLabel(1,"CaloJets",1);
57 
58  //
60  //
61 
62  fillJIDPassFrac = parameters.getParameter<int>("fillJIDPassFrac");
63  makedijetselection = parameters.getParameter<int>("makedijetselection");
64 
65  // monitoring of eta parameter
66  etaBin = parameters.getParameter<int>("etaBin");
67  etaMin = parameters.getParameter<double>("etaMin");
68  etaMax = parameters.getParameter<double>("etaMax");
69 
70  // monitoring of phi paramater
71  phiBin = parameters.getParameter<int>("phiBin");
72  phiMin = parameters.getParameter<double>("phiMin");
73  phiMax = parameters.getParameter<double>("phiMax");
74 
75  // monitoring of the transverse momentum
76  ptBin = parameters.getParameter<int>("ptBin");
77  ptMin = parameters.getParameter<double>("ptMin");
78  ptMax = parameters.getParameter<double>("ptMax");
79 
80  //
81  eBin = parameters.getParameter<int>("eBin");
82  eMin = parameters.getParameter<double>("eMin");
83  eMax = parameters.getParameter<double>("eMax");
84 
85  //
86  pBin = parameters.getParameter<int>("pBin");
87  pMin = parameters.getParameter<double>("pMin");
88  pMax = parameters.getParameter<double>("pMax");
89 
90  //
91  _ptThreshold = parameters.getParameter<double>("ptThreshold");
92  _asymmetryThirdJetCut = parameters.getParameter<double>("asymmetryThirdJetCut");
93  _balanceThirdJetCut = parameters.getParameter<double>("balanceThirdJetCut");
94  _n90HitsMin = parameters.getParameter<int>("n90HitsMin");
95  _fHPDMax = parameters.getParameter<double>("fHPDMax");
96  _resEMFMin = parameters.getParameter<double>("resEMFMin");
97  _sigmaEtaMinTight = parameters.getParameter<double>("sigmaEtaMinTight");
98  _sigmaPhiMinTight = parameters.getParameter<double>("sigmaPhiMinTight");
99 
100  _n90HitsMinLoose = parameters.getParameter<int>("n90HitsMinLoose");
101  _fHPDMaxLoose = parameters.getParameter<double>("fHPDMaxLoose");
102  _resEMFMinLoose = parameters.getParameter<double>("resEMFMinLoose");
103  _n90HitsMinTight = parameters.getParameter<int>("n90HitsMinTight");
104  _fHPDMaxTight = parameters.getParameter<double>("fHPDMaxTight");
105  _resEMFMinTight = parameters.getParameter<double>("resEMFMinTight");
106 
107 
108  // Generic jet parameters
109  mPt = dbe->book1D("Pt", "pt", ptBin, ptMin, ptMax);
110  mEta = dbe->book1D("Eta", "eta", etaBin, etaMin, etaMax);
111  mPhi = dbe->book1D("Phi", "phi", phiBin, phiMin, phiMax);
112  mConstituents = dbe->book1D("Constituents", "# of constituents", 50, 0, 100);
113  mHFrac = dbe->book1D("HFrac", "HFrac", 120, -0.1, 1.1);
114  mEFrac = dbe->book1D("EFrac", "EFrac", 120, -0.1, 1.1);
115 
116 
117  // Book NPV profiles
118  //----------------------------------------------------------------------------
119  mPt_profile = dbe->bookProfile("Pt_profile", "pt", nbinsPV, PVlow, PVup, ptBin, ptMin, ptMax);
120  mEta_profile = dbe->bookProfile("Eta_profile", "eta", nbinsPV, PVlow, PVup, etaBin, etaMin, etaMax);
121  mPhi_profile = dbe->bookProfile("Phi_profile", "phi", nbinsPV, PVlow, PVup, phiBin, phiMin, phiMax);
122  mConstituents_profile = dbe->bookProfile("Constituents_profile", "# of constituents", nbinsPV, PVlow, PVup, 50, 0, 100);
123  mHFrac_profile = dbe->bookProfile("HFrac_profile", "HFrac", nbinsPV, PVlow, PVup, 120, -0.1, 1.1);
124  mEFrac_profile = dbe->bookProfile("EFrac_profile", "EFrac", nbinsPV, PVlow, PVup, 120, -0.1, 1.1);
125 
126  if (makedijetselection != 1)
127  mNJets_profile = dbe->bookProfile("NJets_profile", "number of jets", nbinsPV, PVlow, PVup, 100, 0, 100);
128 
129 
130  // Set NPV profiles x-axis title
131  //----------------------------------------------------------------------------
132  mPt_profile ->setAxisTitle("nvtx",1);
133  mEta_profile ->setAxisTitle("nvtx",1);
134  mPhi_profile ->setAxisTitle("nvtx",1);
136  mHFrac_profile ->setAxisTitle("nvtx",1);
137  mEFrac_profile ->setAxisTitle("nvtx",1);
138 
139  if (makedijetselection != 1) {
140  mNJets_profile->setAxisTitle("nvtx",1);
141  }
142 
143 
144  //mE = dbe->book1D("E", "E", eBin, eMin, eMax);
145  //mP = dbe->book1D("P", "P", pBin, pMin, pMax);
146  // mMass = dbe->book1D("Mass", "Mass", 100, 0, 25);
147  //
148  mPhiVSEta = dbe->book2D("PhiVSEta", "PhiVSEta", 50, etaMin, etaMax, 24, phiMin, phiMax);
149  if(makedijetselection!=1){
150  mPt_1 = dbe->book1D("Pt1", "Pt1", 20, 0, 100);
151  mPt_2 = dbe->book1D("Pt2", "Pt2", 60, 0, 300);
152  mPt_3 = dbe->book1D("Pt3", "Pt3", 100, 0, 5000);
153  // Low and high pt trigger paths
154  mPt_Lo = dbe->book1D("Pt_Lo", "Pt (Pass Low Pt Jet Trigger)", 20, 0, 100);
155  //mEta_Lo = dbe->book1D("Eta_Lo", "Eta (Pass Low Pt Jet Trigger)", etaBin, etaMin, etaMax);
156  mPhi_Lo = dbe->book1D("Phi_Lo", "Phi (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
157 
158  mPt_Hi = dbe->book1D("Pt_Hi", "Pt (Pass Hi Pt Jet Trigger)", 60, 0, 300);
159  mEta_Hi = dbe->book1D("Eta_Hi", "Eta (Pass Hi Pt Jet Trigger)", etaBin, etaMin, etaMax);
160  mPhi_Hi = dbe->book1D("Phi_Hi", "Phi (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
161  mNJets = dbe->book1D("NJets", "number of jets", 100, 0, 100);
162 
163  //mPt_Barrel_Lo = dbe->book1D("Pt_Barrel_Lo", "Pt Barrel (Pass Low Pt Jet Trigger)", 20, 0, 100);
164  //mPhi_Barrel_Lo = dbe->book1D("Phi_Barrel_Lo", "Phi Barrel (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
165  mConstituents_Barrel = dbe->book1D("Constituents_Barrel", "Constituents Barrel above", 50, 0, 100);
166  mHFrac_Barrel = dbe->book1D("HFrac_Barrel", "HFrac Barrel", 100, 0, 1);
167  mEFrac_Barrel = dbe->book1D("EFrac_Barrel", "EFrac Barrel", 110, -0.05, 1.05);
168 
169  //mPt_EndCap_Lo = dbe->book1D("Pt_EndCap_Lo", "Pt EndCap (Pass Low Pt Jet Trigger)", 20, 0, 100);
170  //mPhi_EndCap_Lo = dbe->book1D("Phi_EndCap_Lo", "Phi EndCap (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
171  mConstituents_EndCap = dbe->book1D("Constituents_EndCap", "Constituents EndCap", 50, 0, 100);
172  mHFrac_EndCap = dbe->book1D("HFrac_Endcap", "HFrac EndCap", 100, 0, 1);
173  mEFrac_EndCap = dbe->book1D("EFrac_Endcap", "EFrac EndCap", 110, -0.05, 1.05);
174 
175  //mPt_Forward_Lo = dbe->book1D("Pt_Forward_Lo", "Pt Forward (Pass Low Pt Jet Trigger)", 20, 0, 100);
176  //mPhi_Forward_Lo = dbe->book1D("Phi_Forward_Lo", "Phi Forward (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
177  mConstituents_Forward = dbe->book1D("Constituents_Forward", "Constituents Forward", 50, 0, 100);
178  mHFrac_Forward = dbe->book1D("HFrac_Forward", "HFrac Forward", 100, 0, 1);
179  mEFrac_Forward = dbe->book1D("EFrac_Forward", "EFrac Forward", 110, -0.05, 1.05);
180 
181  mPt_Barrel_Hi = dbe->book1D("Pt_Barrel_Hi", "Pt Barrel (Pass Hi Pt Jet Trigger)", 60, 0, 300);
182  mPhi_Barrel_Hi = dbe->book1D("Phi_Barrel_Hi", "Phi Barrel (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
183  //mConstituents_Barrel_Hi = dbe->book1D("Constituents_Barrel_Hi", "Constituents Barrel (Pass Hi Pt Jet Trigger)", 50, 0, 100);
184  //mHFrac_Barrel_Hi = dbe->book1D("HFrac_Barrel_Hi", "HFrac Barrel (Pass Hi Pt Jet Trigger)", 100, 0, 1);
185 
186  mPt_EndCap_Hi = dbe->book1D("Pt_EndCap_Hi", "Pt EndCap (Pass Hi Pt Jet Trigger)", 60, 0, 300);
187  mPhi_EndCap_Hi = dbe->book1D("Phi_EndCap_Hi", "Phi EndCap (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
188  //mConstituents_EndCap_Hi = dbe->book1D("Constituents_EndCap_Hi", "Constituents EndCap (Pass Hi Pt Jet Trigger)", 50, 0, 100);
189  //mHFrac_EndCap_Hi = dbe->book1D("HFrac_EndCap_Hi", "HFrac EndCap (Pass Hi Pt Jet Trigger)", 100, 0, 1);
190 
191  mPt_Forward_Hi = dbe->book1D("Pt_Forward_Hi", "Pt Forward (Pass Hi Pt Jet Trigger)", 60, 0, 300);
192  mPhi_Forward_Hi = dbe->book1D("Phi_Forward_Hi", "Phi Forward (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
193  //mConstituents_Forward_Hi = dbe->book1D("Constituents_Forward_Hi", "Constituents Forward (Pass Hi Pt Jet Trigger)", 50, 0, 100);
194  //mHFrac_Forward_Hi = dbe->book1D("HFrac_Forward_Hi", "HFrac Forward (Pass Hi Pt Jet Trigger)", 100, 0, 1);
195 
196  mPhi_Barrel = dbe->book1D("Phi_Barrel", "Phi_Barrel", phiBin, phiMin, phiMax);
197  //mE_Barrel = dbe->book1D("E_Barrel", "E_Barrel", eBin, eMin, eMax);
198  mPt_Barrel = dbe->book1D("Pt_Barrel", "Pt_Barrel", ptBin, ptMin, ptMax);
199 
200  mPhi_EndCap = dbe->book1D("Phi_EndCap", "Phi_EndCap", phiBin, phiMin, phiMax);
201  //mE_EndCap = dbe->book1D("E_EndCap", "E_EndCap", eBin, eMin, 2*eMax);
202  mPt_EndCap = dbe->book1D("Pt_EndCap", "Pt_EndCap", ptBin, ptMin, ptMax);
203 
204  mPhi_Forward = dbe->book1D("Phi_Forward", "Phi_Forward", phiBin, phiMin, phiMax);
205  //mE_Forward = dbe->book1D("E_Forward", "E_Forward", eBin, eMin, 4*eMax);
206  mPt_Forward = dbe->book1D("Pt_Forward", "Pt_Forward", ptBin, ptMin, ptMax);
207 
208  // Leading Jet Parameters
209  mEtaFirst = dbe->book1D("EtaFirst", "EtaFirst", 100, -5, 5);
210  mPhiFirst = dbe->book1D("PhiFirst", "PhiFirst", 70, -3.5, 3.5);
211  //mEFirst = dbe->book1D("EFirst", "EFirst", 100, 0, 1000);
212  mPtFirst = dbe->book1D("PtFirst", "PtFirst", 100, 0, 500);
213  if(fillJIDPassFrac==1){//fillJIDPassFrac defines a collection of cleaned jets, for which we will want to fill the cleaning passing fraction
214  mLooseJIDPassFractionVSeta = dbe->bookProfile("LooseJIDPassFractionVSeta","LooseJIDPassFractionVSeta",etaBin, etaMin, etaMax,0.,1.2);
215  mLooseJIDPassFractionVSpt = dbe->bookProfile("LooseJIDPassFractionVSpt","LooseJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
216  mTightJIDPassFractionVSeta = dbe->bookProfile("TightJIDPassFractionVSeta","TightJIDPassFractionVSeta",etaBin, etaMin, etaMax,0.,1.2);
217  mTightJIDPassFractionVSpt = dbe->bookProfile("TightJIDPassFractionVSpt","TightJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
218 
219 
220  }
221  }
222  // CaloJet specific
223  mMaxEInEmTowers = dbe->book1D("MaxEInEmTowers", "MaxEInEmTowers", 100, 0, 100);
224  mMaxEInHadTowers = dbe->book1D("MaxEInHadTowers", "MaxEInHadTowers", 100, 0, 100);
225  if(makedijetselection!=1) {
226  mHadEnergyInHO = dbe->book1D("HadEnergyInHO", "HadEnergyInHO", 100, 0, 10);
227  mHadEnergyInHB = dbe->book1D("HadEnergyInHB", "HadEnergyInHB", 100, 0, 50);
228  mHadEnergyInHF = dbe->book1D("HadEnergyInHF", "HadEnergyInHF", 100, 0, 50);
229  mHadEnergyInHE = dbe->book1D("HadEnergyInHE", "HadEnergyInHE", 100, 0, 100);
230  mEmEnergyInEB = dbe->book1D("EmEnergyInEB", "EmEnergyInEB", 100, 0, 50);
231  mEmEnergyInEE = dbe->book1D("EmEnergyInEE", "EmEnergyInEE", 100, 0, 50);
232  mEmEnergyInHF = dbe->book1D("EmEnergyInHF", "EmEnergyInHF", 120, -20, 100);
233  }
234  mDPhi = dbe->book1D("DPhi", "dPhi btw the two leading jets", 100, 0., acos(-1.));
235 
236  //JetID variables
237 
238  mresEMF = dbe->book1D("resEMF", "resEMF", 50, 0., 1.);
239  mN90Hits = dbe->book1D("N90Hits", "N90Hits", 50, 0., 50);
240  mfHPD = dbe->book1D("fHPD", "fHPD", 50, 0., 1.);
241  mfRBX = dbe->book1D("fRBX", "fRBX", 50, 0., 1.);
242 
243  // msigmaEta = dbe->book1D("sigmaEta", "sigmaEta", 50, 0., 1.);
244  // msigmaPhi = dbe->book1D("sigmaPhi", "sigmaPhi", 50, 0., 0.5);
245 
246  if(makedijetselection==1) {
247  mDijetAsymmetry = dbe->book1D("DijetAsymmetry", "DijetAsymmetry", 100, -1., 1.);
248  mDijetBalance = dbe->book1D("DijetBalance", "DijetBalance", 100, -2., 2.);
249  if (fillJIDPassFrac==1) {
250  mLooseJIDPassFractionVSeta = dbe->bookProfile("LooseJIDPassFractionVSeta","LooseJIDPassFractionVSeta",50, -3., 3.,0.,1.2);
251  mLooseJIDPassFractionVSpt = dbe->bookProfile("LooseJIDPassFractionVSpt","LooseJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
252  mTightJIDPassFractionVSeta = dbe->bookProfile("TightJIDPassFractionVSeta","TightJIDPassFractionVSeta",50, -3., 3.,0.,1.2);
253  mTightJIDPassFractionVSpt = dbe->bookProfile("TightJIDPassFractionVSpt","TightJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
254  }
255  }
256 }
MonitorElement * mEFrac_EndCap
Definition: JetAnalyzer.h:184
MonitorElement * mNJets
Definition: JetAnalyzer.h:213
T getParameter(std::string const &) const
MonitorElement * mPt_Forward_Hi
Definition: JetAnalyzer.h:199
int _n90HitsMinTight
Definition: JetAnalyzer.h:110
MonitorElement * mHadEnergyInHB
Definition: JetAnalyzer.h:227
MonitorElement * mMaxEInEmTowers
Definition: JetAnalyzer.h:224
MonitorElement * mfRBX
Definition: JetAnalyzer.h:237
MonitorElement * mPhi_Forward
Definition: JetAnalyzer.h:173
int _n90HitsMinLoose
Definition: JetAnalyzer.h:107
MonitorElement * mPt_Lo
Definition: JetAnalyzer.h:249
MonitorElement * mEmEnergyInHF
Definition: JetAnalyzer.h:232
MonitorElement * mPhiFirst
Definition: JetAnalyzer.h:218
MonitorElement * mPt_Forward
Definition: JetAnalyzer.h:172
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
double phiMax
Definition: JetAnalyzer.h:121
MonitorElement * mHFrac_Forward
Definition: JetAnalyzer.h:188
MonitorElement * mEmEnergyInEE
Definition: JetAnalyzer.h:231
double eMin
Definition: JetAnalyzer.h:128
int makedijetselection
Definition: JetAnalyzer.h:97
MonitorElement * jetME
Definition: JetAnalyzer.h:136
MonitorElement * mPhi_EndCap
Definition: JetAnalyzer.h:170
MonitorElement * mN90Hits
Definition: JetAnalyzer.h:235
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * mEta_Hi
Definition: JetAnalyzer.h:251
MonitorElement * mEFrac_profile
Definition: JetAnalyzer.h:268
MonitorElement * mPt_3
Definition: JetAnalyzer.h:158
MonitorElement * mEta
Definition: JetAnalyzer.h:159
double _fHPDMaxTight
Definition: JetAnalyzer.h:108
MonitorElement * mConstituents_Forward
Definition: JetAnalyzer.h:187
MonitorElement * mPt_Hi
Definition: JetAnalyzer.h:253
double _resEMFMin
Definition: JetAnalyzer.h:101
double ptMin
Definition: JetAnalyzer.h:124
const int nbinsPV
double _fHPDMax
Definition: JetAnalyzer.h:100
MonitorElement * mTightJIDPassFractionVSpt
Definition: JetAnalyzer.h:244
MonitorElement * mLooseJIDPassFractionVSpt
Definition: JetAnalyzer.h:242
double etaMax
Definition: JetAnalyzer.h:117
double _fHPDMaxLoose
Definition: JetAnalyzer.h:105
edm::ParameterSet parameters
Definition: JetAnalyzer.h:81
MonitorElement * mDPhi
Definition: JetAnalyzer.h:214
MonitorElement * mresEMF
Definition: JetAnalyzer.h:238
MonitorElement * mConstituents_profile
Definition: JetAnalyzer.h:266
MonitorElement * mPt_Barrel_Hi
Definition: JetAnalyzer.h:191
MonitorElement * mEmEnergyInEB
Definition: JetAnalyzer.h:230
double _balanceThirdJetCut
Definition: JetAnalyzer.h:95
MonitorElement * mDijetBalance
Definition: JetAnalyzer.h:256
const double PVup
MonitorElement * mPhi_Barrel
Definition: JetAnalyzer.h:167
MonitorElement * mConstituents_EndCap
Definition: JetAnalyzer.h:182
MonitorElement * mPhi_EndCap_Hi
Definition: JetAnalyzer.h:196
MonitorElement * mPt_1
Definition: JetAnalyzer.h:156
MonitorElement * mConstituents_Barrel
Definition: JetAnalyzer.h:177
std::string jetname
Definition: JetAnalyzer.h:83
MonitorElement * mPhi
Definition: JetAnalyzer.h:160
std::string _source
Definition: JetAnalyzer.h:84
MonitorElement * mPt_profile
Definition: JetAnalyzer.h:263
MonitorElement * bookProfile(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const char *option="s")
Definition: DQMStore.cc:1186
MonitorElement * mEtaFirst
Definition: JetAnalyzer.h:217
MonitorElement * mEFrac_Forward
Definition: JetAnalyzer.h:189
double _resEMFMinTight
Definition: JetAnalyzer.h:109
double pMax
Definition: JetAnalyzer.h:133
MonitorElement * mPt_2
Definition: JetAnalyzer.h:157
MonitorElement * mHFrac_Barrel
Definition: JetAnalyzer.h:178
double phiMin
Definition: JetAnalyzer.h:120
MonitorElement * mNJets_profile
Definition: JetAnalyzer.h:262
#define LogTrace(id)
const double PVlow
MonitorElement * mEFrac
Definition: JetAnalyzer.h:163
double _sigmaEtaMinTight
Definition: JetAnalyzer.h:111
double ptMax
Definition: JetAnalyzer.h:125
MonitorElement * mHadEnergyInHF
Definition: JetAnalyzer.h:228
MonitorElement * mfHPD
Definition: JetAnalyzer.h:236
MonitorElement * mDijetAsymmetry
Definition: JetAnalyzer.h:257
MonitorElement * mPt
Definition: JetAnalyzer.h:155
MonitorElement * mEFrac_Barrel
Definition: JetAnalyzer.h:179
double _asymmetryThirdJetCut
Definition: JetAnalyzer.h:94
MonitorElement * mPhi_Barrel_Hi
Definition: JetAnalyzer.h:192
MonitorElement * mHFrac_EndCap
Definition: JetAnalyzer.h:183
MonitorElement * mPhiVSEta
Definition: JetAnalyzer.h:164
MonitorElement * mTightJIDPassFractionVSeta
Definition: JetAnalyzer.h:243
double etaMin
Definition: JetAnalyzer.h:116
MonitorElement * mLooseJIDPassFractionVSeta
Definition: JetAnalyzer.h:241
MonitorElement * mHFrac_profile
Definition: JetAnalyzer.h:267
reco::helper::JetIDHelper * jetID
Definition: JetAnalyzer.h:139
MonitorElement * mPt_EndCap_Hi
Definition: JetAnalyzer.h:195
int fillJIDPassFrac
Definition: JetAnalyzer.h:104
double _resEMFMinLoose
Definition: JetAnalyzer.h:106
MonitorElement * mMaxEInHadTowers
Definition: JetAnalyzer.h:225
MonitorElement * mPhi_Hi
Definition: JetAnalyzer.h:252
MonitorElement * mPtFirst
Definition: JetAnalyzer.h:220
MonitorElement * mPhi_profile
Definition: JetAnalyzer.h:265
MonitorElement * mPt_Barrel
Definition: JetAnalyzer.h:166
MonitorElement * mHadEnergyInHE
Definition: JetAnalyzer.h:229
MonitorElement * mConstituents
Definition: JetAnalyzer.h:161
MonitorElement * mPhi_Forward_Hi
Definition: JetAnalyzer.h:200
MonitorElement * mEta_profile
Definition: JetAnalyzer.h:264
MonitorElement * mPt_EndCap
Definition: JetAnalyzer.h:169
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1000
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
double eMax
Definition: JetAnalyzer.h:129
double pMin
Definition: JetAnalyzer.h:132
double _sigmaPhiMinTight
Definition: JetAnalyzer.h:112
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
MonitorElement * mPhi_Lo
Definition: JetAnalyzer.h:248
MonitorElement * mHFrac
Definition: JetAnalyzer.h:162
double _ptThreshold
Definition: JetAnalyzer.h:92
MonitorElement * mHadEnergyInHO
Definition: JetAnalyzer.h:226
void JetAnalyzer::endJob ( void  )

Finish up a job.

Definition at line 265 of file JetAnalyzer.cc.

References bTagSequences_cff::jetID.

265  {
266  delete jetID;
267 }
reco::helper::JetIDHelper * jetID
Definition: JetAnalyzer.h:139
int JetAnalyzer::getLeadJetFlag ( )
inline

Definition at line 67 of file JetAnalyzer.h.

References _leadJetFlag.

67  {
68  return _leadJetFlag;
69  }
int _leadJetFlag
Definition: JetAnalyzer.h:90
void JetAnalyzer::setJetHiPass ( int  pass)
inline

Definition at line 74 of file JetAnalyzer.h.

References _JetHiPass.

74  {
75  _JetHiPass = pass;
76  }
int _JetHiPass
Definition: JetAnalyzer.h:89
void JetAnalyzer::setJetLoPass ( int  pass)
inline

Definition at line 70 of file JetAnalyzer.h.

References _JetLoPass.

70  {
71  _JetLoPass = pass;
72  }
int _JetLoPass
Definition: JetAnalyzer.h:88
void JetAnalyzer::setLeadJetFlag ( int  flag)
inline

Definition at line 64 of file JetAnalyzer.h.

References _leadJetFlag.

64  {
65  _leadJetFlag = flag;
66  }
int _leadJetFlag
Definition: JetAnalyzer.h:90
void JetAnalyzer::setSource ( std::string  source)
inline

Definition at line 60 of file JetAnalyzer.h.

References _source, and source.

60  {
61  _source = source;
62  }
std::string _source
Definition: JetAnalyzer.h:84
static std::string const source
Definition: EdmProvDump.cc:43

Member Data Documentation

double JetAnalyzer::_asymmetryThirdJetCut
private

Definition at line 94 of file JetAnalyzer.h.

double JetAnalyzer::_balanceThirdJetCut
private

Definition at line 95 of file JetAnalyzer.h.

double JetAnalyzer::_fHPDMax
private

Definition at line 100 of file JetAnalyzer.h.

double JetAnalyzer::_fHPDMaxLoose
private

Definition at line 105 of file JetAnalyzer.h.

double JetAnalyzer::_fHPDMaxTight
private

Definition at line 108 of file JetAnalyzer.h.

int JetAnalyzer::_JetHiPass
private

Definition at line 89 of file JetAnalyzer.h.

Referenced by setJetHiPass().

int JetAnalyzer::_JetLoPass
private

Definition at line 88 of file JetAnalyzer.h.

Referenced by setJetLoPass().

int JetAnalyzer::_leadJetFlag
private

Definition at line 90 of file JetAnalyzer.h.

Referenced by getLeadJetFlag(), and setLeadJetFlag().

int JetAnalyzer::_n90HitsMin
private

Definition at line 102 of file JetAnalyzer.h.

int JetAnalyzer::_n90HitsMinLoose
private

Definition at line 107 of file JetAnalyzer.h.

int JetAnalyzer::_n90HitsMinTight
private

Definition at line 110 of file JetAnalyzer.h.

double JetAnalyzer::_ptThreshold
private

Definition at line 92 of file JetAnalyzer.h.

double JetAnalyzer::_resEMFMin
private

Definition at line 101 of file JetAnalyzer.h.

double JetAnalyzer::_resEMFMinLoose
private

Definition at line 106 of file JetAnalyzer.h.

double JetAnalyzer::_resEMFMinTight
private

Definition at line 109 of file JetAnalyzer.h.

double JetAnalyzer::_sigmaEtaMinTight
private

Definition at line 111 of file JetAnalyzer.h.

double JetAnalyzer::_sigmaPhiMinTight
private

Definition at line 112 of file JetAnalyzer.h.

std::string JetAnalyzer::_source
private

Definition at line 84 of file JetAnalyzer.h.

Referenced by setSource().

int JetAnalyzer::_theend
private

Definition at line 91 of file JetAnalyzer.h.

int JetAnalyzer::eBin
private

Definition at line 127 of file JetAnalyzer.h.

double JetAnalyzer::eMax
private

Definition at line 129 of file JetAnalyzer.h.

double JetAnalyzer::eMin
private

Definition at line 128 of file JetAnalyzer.h.

int JetAnalyzer::etaBin
private

Definition at line 115 of file JetAnalyzer.h.

double JetAnalyzer::etaMax
private

Definition at line 117 of file JetAnalyzer.h.

double JetAnalyzer::etaMin
private

Definition at line 116 of file JetAnalyzer.h.

int JetAnalyzer::fillJIDPassFrac
private

Definition at line 104 of file JetAnalyzer.h.

reco::helper::JetIDHelper* JetAnalyzer::jetID
private

Definition at line 139 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::jetME
private

Definition at line 136 of file JetAnalyzer.h.

std::string JetAnalyzer::jetname
private

Definition at line 83 of file JetAnalyzer.h.

int JetAnalyzer::makedijetselection
private

Definition at line 97 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents
private

Definition at line 161 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_Barrel
private

Definition at line 177 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_Barrel_Hi
private

Definition at line 193 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_EndCap
private

Definition at line 182 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_EndCap_Hi
private

Definition at line 197 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_Forward
private

Definition at line 187 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_Forward_Hi
private

Definition at line 201 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mConstituents_profile
private

Definition at line 266 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mDijetAsymmetry
private

Definition at line 257 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mDijetBalance
private

Definition at line 256 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mDPhi
private

Definition at line 214 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEFrac
private

Definition at line 163 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEFrac_Barrel
private

Definition at line 179 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEFrac_EndCap
private

Definition at line 184 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEFrac_Forward
private

Definition at line 189 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEFrac_profile
private

Definition at line 268 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEmEnergyInEB
private

Definition at line 230 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEmEnergyInEE
private

Definition at line 231 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEmEnergyInHF
private

Definition at line 232 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEta
private

Definition at line 159 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEta_Hi
private

Definition at line 251 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEta_profile
private

Definition at line 264 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mEtaFirst
private

Definition at line 217 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mfHPD
private

Definition at line 236 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mfRBX
private

Definition at line 237 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHadEnergyInHB
private

Definition at line 227 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHadEnergyInHE
private

Definition at line 229 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHadEnergyInHF
private

Definition at line 228 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHadEnergyInHO
private

Definition at line 226 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac
private

Definition at line 162 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_Barrel
private

Definition at line 178 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_Barrel_Hi
private

Definition at line 194 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_EndCap
private

Definition at line 183 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_EndCap_Hi
private

Definition at line 198 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_Forward
private

Definition at line 188 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_Forward_Hi
private

Definition at line 202 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mHFrac_profile
private

Definition at line 267 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mLooseJIDPassFractionVSeta
private

Definition at line 241 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mLooseJIDPassFractionVSpt
private

Definition at line 242 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mMaxEInEmTowers
private

Definition at line 224 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mMaxEInHadTowers
private

Definition at line 225 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mN90Hits
private

Definition at line 235 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mNJets
private

Definition at line 213 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mNJets_profile
private

Definition at line 262 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi
private

Definition at line 160 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_Barrel
private

Definition at line 167 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_Barrel_Hi
private

Definition at line 192 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_EndCap
private

Definition at line 170 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_EndCap_Hi
private

Definition at line 196 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_Forward
private

Definition at line 173 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_Forward_Hi
private

Definition at line 200 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_Hi
private

Definition at line 252 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_Lo
private

Definition at line 248 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhi_profile
private

Definition at line 265 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhiFirst
private

Definition at line 218 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPhiVSEta
private

Definition at line 164 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt
private

Definition at line 155 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_1
private

Definition at line 156 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_2
private

Definition at line 157 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_3
private

Definition at line 158 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_Barrel
private

Definition at line 166 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_Barrel_Hi
private

Definition at line 191 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_EndCap
private

Definition at line 169 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_EndCap_Hi
private

Definition at line 195 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_Forward
private

Definition at line 172 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_Forward_Hi
private

Definition at line 199 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_Hi
private

Definition at line 253 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_Lo
private

Definition at line 249 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPt_profile
private

Definition at line 263 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mPtFirst
private

Definition at line 220 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mresEMF
private

Definition at line 238 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mTightJIDPassFractionVSeta
private

Definition at line 243 of file JetAnalyzer.h.

MonitorElement* JetAnalyzer::mTightJIDPassFractionVSpt
private

Definition at line 244 of file JetAnalyzer.h.

edm::ParameterSet JetAnalyzer::parameters
private
int JetAnalyzer::pBin
private

Definition at line 131 of file JetAnalyzer.h.

int JetAnalyzer::phiBin
private

Definition at line 119 of file JetAnalyzer.h.

double JetAnalyzer::phiMax
private

Definition at line 121 of file JetAnalyzer.h.

double JetAnalyzer::phiMin
private

Definition at line 120 of file JetAnalyzer.h.

double JetAnalyzer::pMax
private

Definition at line 133 of file JetAnalyzer.h.

double JetAnalyzer::pMin
private

Definition at line 132 of file JetAnalyzer.h.

int JetAnalyzer::ptBin
private

Definition at line 123 of file JetAnalyzer.h.

double JetAnalyzer::ptMax
private

Definition at line 125 of file JetAnalyzer.h.

double JetAnalyzer::ptMin
private

Definition at line 124 of file JetAnalyzer.h.

edm::InputTag JetAnalyzer::theCaloJetCollectionLabel
private

Definition at line 86 of file JetAnalyzer.h.