CMS 3D CMS Logo

Public Member Functions | Private Attributes

PFJetAnalyzer Class Reference

#include <PFJetAnalyzer.h>

Inheritance diagram for PFJetAnalyzer:
PFJetAnalyzerBase

List of all members.

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &, const reco::PFJetCollection &pfJets)
 Get the analysis.
void beginJob (DQMStore *dbe)
 Inizialize parameters for histo binning.
int getLeadJetFlag ()
 PFJetAnalyzer (const edm::ParameterSet &)
 Constructor.
void setJetHiPass (int pass)
void setJetLoPass (int pass)
void setLeadJetFlag (int flag)
void setSource (std::string source)
virtual ~PFJetAnalyzer ()
 Destructor.

Private Attributes

int _JetHiPass
int _JetLoPass
int _leadJetFlag
double _LooseCEFMax
double _LooseCHFMin
double _LooseNEFMax
double _LooseNHFMax
double _ptThreshold
std::string _source
double _ThisCEFMax
double _ThisCHFMin
double _ThisNEFMax
double _ThisNHFMax
double _TightCEFMax
double _TightCHFMin
double _TightNEFMax
double _TightNHFMax
int eBin
double eMax
double eMin
int etaBin
double etaMax
double etaMin
int fillpfJIDPassFrac
MonitorElementjetME
MonitorElementmChargedEmEnergy
MonitorElementmChargedHadronEnergy
MonitorElementmChargedMuEnergy
MonitorElementmChargedMultiplicity
MonitorElementmConstituents
MonitorElementmConstituents_Barrel_Hi
MonitorElementmConstituents_Barrel_Lo
MonitorElementmConstituents_EndCap_Hi
MonitorElementmConstituents_EndCap_Lo
MonitorElementmConstituents_Forward_Hi
MonitorElementmConstituents_Forward_Lo
MonitorElementmDPhi
MonitorElementmE
MonitorElementmE_Barrel
MonitorElementmE_EndCap
MonitorElementmE_Forward
MonitorElementmEEffChargedFraction
MonitorElementmEEffNeutralFraction
MonitorElementmEFirst
MonitorElementmEFrac
MonitorElementmEResChargedFraction
MonitorElementmEResNeutralFraction
MonitorElementmEta
MonitorElementmEta_Hi
MonitorElementmEta_Lo
MonitorElementmEtaFirst
std::string metname
MonitorElementmHFrac
MonitorElementmHFrac_Barrel_Hi
MonitorElementmHFrac_Barrel_Lo
MonitorElementmHFrac_EndCap_Hi
MonitorElementmHFrac_EndCap_Lo
MonitorElementmHFrac_Forward_Hi
MonitorElementmHFrac_Forward_Lo
MonitorElementmLooseJIDPassFractionVSeta
MonitorElementmLooseJIDPassFractionVSpt
MonitorElementmMass
MonitorElementmMuonMultiplicity
MonitorElementmNeutralEmEnergy
MonitorElementmNeutralFraction
MonitorElementmNeutralFraction2
MonitorElementmNeutralHadronEnergy
MonitorElementmNeutralMultiplicity
MonitorElementmNJets
MonitorElementmP
MonitorElementmPhi
MonitorElementmPhi_Barrel
MonitorElementmPhi_Barrel_Hi
MonitorElementmPhi_Barrel_Lo
MonitorElementmPhi_EndCap
MonitorElementmPhi_EndCap_Hi
MonitorElementmPhi_EndCap_Lo
MonitorElementmPhi_Forward
MonitorElementmPhi_Forward_Hi
MonitorElementmPhi_Forward_Lo
MonitorElementmPhi_Hi
MonitorElementmPhi_Lo
MonitorElementmPhiFirst
MonitorElementmPhiVSEta
MonitorElementmPt
MonitorElementmPt_1
MonitorElementmPt_2
MonitorElementmPt_3
MonitorElementmPt_Barrel
MonitorElementmPt_Barrel_Hi
MonitorElementmPt_Barrel_Lo
MonitorElementmPt_EndCap
MonitorElementmPt_EndCap_Hi
MonitorElementmPt_EndCap_Lo
MonitorElementmPt_Forward
MonitorElementmPt_Forward_Hi
MonitorElementmPt_Forward_Lo
MonitorElementmPt_Hi
MonitorElementmPt_Lo
MonitorElementmPtFirst
MonitorElementmTightJIDPassFractionVSeta
MonitorElementmTightJIDPassFractionVSpt
MonitorElementnEEff
edm::ParameterSet parameters
int pBin
int phiBin
double phiMax
double phiMin
double pMax
double pMin
int ptBin
double ptMax
double ptMin
edm::InputTag thePFJetCollectionLabel

Detailed Description

DQM monitoring source for PFlow Jets

Date:
2010/10/15 13:49:54
Revision:
1.6
Author:
F. Chlebana - Fermilab

Definition at line 30 of file PFJetAnalyzer.h.


Constructor & Destructor Documentation

PFJetAnalyzer::PFJetAnalyzer ( const edm::ParameterSet pSet)

Constructor.

Definition at line 20 of file PFJetAnalyzer.cc.

References Parameters::parameters.

                                                        {

  parameters = pSet;
  _leadJetFlag = 0;
  _JetLoPass   = 0;
  _JetHiPass   = 0;
  _ptThreshold = 5.;
  _LooseCHFMin = -999.;
  _LooseNHFMax = -999.;
  _LooseCEFMax = -999.;
  _LooseNEFMax = -999.;
  _TightCHFMin = -999.;
  _TightNHFMax = -999.;
  _TightCEFMax = -999.;
  _TightNEFMax = -999.;
  _ThisCHFMin = -999.;
  _ThisNHFMax = -999.;
  _ThisCEFMax = -999.;
  _ThisNEFMax = -999.;

}
PFJetAnalyzer::~PFJetAnalyzer ( ) [virtual]

Destructor.

Definition at line 43 of file PFJetAnalyzer.cc.

{ }

Member Function Documentation

void PFJetAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const reco::PFJetCollection pfJets 
)

Get the analysis.

Definition at line 201 of file PFJetAnalyzer.cc.

References corr, diffTreeTool::diff, metsig::jet, LogTrace, and metname.

                                                                                                                  {

  int numofjets=0;
  double  fstPhi=0.;
  double  sndPhi=0.;
  double  diff = 0.;
  double  corr = 0.;
  double  dphi = -999. ;

  bool Thiscleaned=false; 
  bool Loosecleaned=false; 
  bool Tightcleaned=false; 
  bool ThisCHFcleaned=false;
  bool LooseCHFcleaned=false;
  bool TightCHFcleaned=false;

  for (reco::PFJetCollection::const_iterator jet = pfJets.begin(); jet!=pfJets.end(); ++jet){
    LogTrace(metname)<<"[JetAnalyzer] Analyze PFJet";
 
    Thiscleaned=false;
    Loosecleaned=false;
    Tightcleaned=false;

    if (jet == pfJets.begin()) {
      fstPhi = jet->phi();
      _leadJetFlag = 1;
    } else {
      _leadJetFlag = 0;
    }

    if (jet == (pfJets.begin()+1)) sndPhi = jet->phi();
    //  if (jet->pt() < _ptThreshold) return;
    if (jet->pt() > _ptThreshold) {
      numofjets++ ;
      jetME->Fill(2);
            
      //calculate the jetID
      ThisCHFcleaned=true;
      LooseCHFcleaned=true;
      TightCHFcleaned=true;
      if((jet->chargedHadronEnergy()/jet->energy())<=_ThisCHFMin && fabs(jet->eta())<2.4) ThisCHFcleaned=false; //apply CHF>0 only if |eta|<2.4
      if((jet->chargedHadronEnergy()/jet->energy())<=_LooseCHFMin && fabs(jet->eta())<2.4) LooseCHFcleaned=false; //apply CHF>0 only if |eta|<2.4
      if((jet->chargedHadronEnergy()/jet->energy())<=_TightCHFMin && fabs(jet->eta())<2.4) TightCHFcleaned=false; //apply CHF>0 only if |eta|<2.4
      if(ThisCHFcleaned && (jet->neutralHadronEnergy()/jet->energy())<_ThisNHFMax && (jet->chargedEmEnergy()/jet->energy())<_ThisCEFMax && (jet->neutralEmEnergy()/jet->energy())<_ThisNEFMax) Thiscleaned=true;
      if(LooseCHFcleaned && (jet->neutralHadronEnergy()/jet->energy())<_LooseNHFMax && (jet->chargedEmEnergy()/jet->energy())<_LooseCEFMax && (jet->neutralEmEnergy()/jet->energy())<_LooseNEFMax) Loosecleaned=true;
      if(TightCHFcleaned && (jet->neutralHadronEnergy()/jet->energy())<_TightNHFMax && (jet->chargedEmEnergy()/jet->energy())<_TightCEFMax && (jet->neutralEmEnergy()/jet->energy())<_TightNEFMax) Tightcleaned=true;
      
      if(fillpfJIDPassFrac==1) {
        //fill the profile for jid efficiency 
        if(Loosecleaned) {
          mLooseJIDPassFractionVSeta->Fill(jet->eta(),1.);
          mLooseJIDPassFractionVSpt->Fill(jet->pt(),1.);
        } else {
          mLooseJIDPassFractionVSeta->Fill(jet->eta(),0.);
          mLooseJIDPassFractionVSpt->Fill(jet->pt(),0.);
        }
        if(Tightcleaned) {
          mTightJIDPassFractionVSeta->Fill(jet->eta(),1.);
          mTightJIDPassFractionVSpt->Fill(jet->pt(),1.);
        } else {
          mTightJIDPassFractionVSeta->Fill(jet->eta(),0.);
          mTightJIDPassFractionVSpt->Fill(jet->pt(),0.);
        }
      }
      
      if(!Thiscleaned) continue;
      
      // Leading jet
      // Histograms are filled once per event
      if (_leadJetFlag == 1) { 
        
        if (mEtaFirst) mEtaFirst->Fill (jet->eta());
        if (mPhiFirst) mPhiFirst->Fill (jet->phi());
        if (mEFirst)   mEFirst->Fill (jet->energy());
        if (mPtFirst)  mPtFirst->Fill (jet->pt());
      }
      
      // --- Passed the low pt jet trigger
      if (_JetLoPass == 1) {
        if (fabs(jet->eta()) <= 1.3) {
          if (mPt_Barrel_Lo)           mPt_Barrel_Lo->Fill(jet->pt());
          if (mEta_Lo)          mEta_Lo->Fill(jet->eta());
          if (mPhi_Barrel_Lo)          mPhi_Barrel_Lo->Fill(jet->phi());
          if (mConstituents_Barrel_Lo) mConstituents_Barrel_Lo->Fill(jet->nConstituents());     
          if (mHFrac_Barrel_Lo)        mHFrac_Barrel_Lo->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction() );  
        }
        if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
          if (mPt_EndCap_Lo)           mPt_EndCap_Lo->Fill(jet->pt());
          if (mEta_Lo)          mEta_Lo->Fill(jet->eta());
          if (mPhi_EndCap_Lo)          mPhi_EndCap_Lo->Fill(jet->phi());
          if (mConstituents_EndCap_Lo) mConstituents_EndCap_Lo->Fill(jet->nConstituents());     
          if (mHFrac_EndCap_Lo)        mHFrac_EndCap_Lo->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());   
        }
        if (fabs(jet->eta()) > 3.0) {
          if (mPt_Forward_Lo)           mPt_Forward_Lo->Fill(jet->pt());
          if (mEta_Lo)          mEta_Lo->Fill(jet->eta());
          if (mPhi_Forward_Lo)          mPhi_Forward_Lo->Fill(jet->phi());
          if (mConstituents_Forward_Lo) mConstituents_Forward_Lo->Fill(jet->nConstituents());   
          if (mHFrac_Forward_Lo)        mHFrac_Forward_Lo->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction()); 
        }
        if (mEta_Lo) mEta_Lo->Fill (jet->eta());
        if (mPhi_Lo) mPhi_Lo->Fill (jet->phi());
        if (mPt_Lo)  mPt_Lo->Fill (jet->pt());
      }
      
      // --- Passed the high pt jet trigger
      if (_JetHiPass == 1) {
        if (fabs(jet->eta()) <= 1.3) {
          if (mPt_Barrel_Hi)           mPt_Barrel_Hi->Fill(jet->pt());
          if (mEta_Hi)          mEta_Hi->Fill(jet->eta());
          if (mPhi_Barrel_Hi)          mPhi_Barrel_Hi->Fill(jet->phi());
          if (mConstituents_Barrel_Hi) mConstituents_Barrel_Hi->Fill(jet->nConstituents());     
          if (mHFrac_Barrel_Hi)        mHFrac_Barrel_Hi->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());   
        }
        if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
          if (mPt_EndCap_Hi)           mPt_EndCap_Hi->Fill(jet->pt());
          if (mEta_Hi)          mEta_Hi->Fill(jet->eta());
          if (mPhi_EndCap_Hi)          mPhi_EndCap_Hi->Fill(jet->phi());
          if (mConstituents_EndCap_Hi) mConstituents_EndCap_Hi->Fill(jet->nConstituents());     
          if (mHFrac_EndCap_Hi)        mHFrac_EndCap_Hi->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());   
        }
        if (fabs(jet->eta()) > 3.0) {
          if (mPt_Forward_Hi)           mPt_Forward_Hi->Fill(jet->pt());
          if (mEta_Hi)          mEta_Hi->Fill(jet->eta());
          if (mPhi_Forward_Hi)          mPhi_Forward_Hi->Fill(jet->phi());
          if (mConstituents_Forward_Hi) mConstituents_Forward_Hi->Fill(jet->nConstituents());   
          if (mHFrac_Forward_Hi)        mHFrac_Forward_Hi->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction()); 
        }
        
        if (mEta_Hi) mEta_Hi->Fill (jet->eta());
        if (mPhi_Hi) mPhi_Hi->Fill (jet->phi());
        if (mPt_Hi)  mPt_Hi->Fill (jet->pt());
      }
      
      if (mPt)   mPt->Fill (jet->pt());
      if (mPt_1) mPt_1->Fill (jet->pt());
      if (mPt_2) mPt_2->Fill (jet->pt());
      if (mPt_3) mPt_3->Fill (jet->pt());
      if (mEta)  mEta->Fill (jet->eta());
      if (mPhi)  mPhi->Fill (jet->phi());
      if (mPhiVSEta) mPhiVSEta->Fill(jet->eta(),jet->phi());
      
      if (mConstituents) mConstituents->Fill (jet->nConstituents());
      if (mHFrac)        mHFrac->Fill (jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());
      if (mEFrac)        mEFrac->Fill (jet->chargedEmEnergyFraction() +jet->neutralEmEnergyFraction());
      
      if (fabs(jet->eta()) <= 1.3) {
        if (mPt_Barrel)   mPt_Barrel->Fill (jet->pt());
        if (mPhi_Barrel)  mPhi_Barrel->Fill (jet->phi());
        if (mE_Barrel)    mE_Barrel->Fill (jet->energy());
      }
      if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
        if (mPt_EndCap)   mPt_EndCap->Fill (jet->pt());
        if (mPhi_EndCap)  mPhi_EndCap->Fill (jet->phi());
        if (mE_EndCap)    mE_EndCap->Fill (jet->energy());
      }
      if (fabs(jet->eta()) > 3.0) {
        if (mPt_Forward)   mPt_Forward->Fill (jet->pt());
        if (mPhi_Forward)  mPhi_Forward->Fill (jet->phi());
        if (mE_Forward)    mE_Forward->Fill (jet->energy());
      }
      
      if (mE)    mE->Fill (jet->energy());
      if (mP)    mP->Fill (jet->p());
      if (mMass) mMass->Fill (jet->mass());
      
      
      if (mChargedHadronEnergy)  mChargedHadronEnergy->Fill (jet->chargedHadronEnergy());
      if (mNeutralHadronEnergy)  mNeutralHadronEnergy->Fill (jet->neutralHadronEnergy());
      if (mChargedEmEnergy) mChargedEmEnergy->Fill(jet->chargedEmEnergy());
      if (mChargedMuEnergy) mChargedMuEnergy->Fill (jet->chargedMuEnergy ());
      if (mNeutralEmEnergy) mNeutralEmEnergy->Fill(jet->neutralEmEnergy());
      if (mChargedMultiplicity ) mChargedMultiplicity->Fill(jet->chargedMultiplicity());
      if (mNeutralMultiplicity ) mNeutralMultiplicity->Fill(jet->neutralMultiplicity());
      if (mMuonMultiplicity )mMuonMultiplicity->Fill (jet-> muonMultiplicity());
      //_______________________________________________________
      if (mNeutralFraction) mNeutralFraction->Fill (jet->neutralMultiplicity()/jet->nConstituents());
      
      //calculate correctly the dphi
      if(numofjets>1) {
        diff = fabs(fstPhi - sndPhi);
        corr = 2*acos(-1.) - diff;
        if(diff < acos(-1.)) { 
          dphi = diff; 
        } else { 
          dphi = corr;
        }
      } // numofjets>1
    } // JetPt>_ptThreshold
  } // PF jet loop
  if (mNJets)    mNJets->Fill (numofjets);
  if (mDPhi)    mDPhi->Fill (dphi);
}
void PFJetAnalyzer::beginJob ( DQMStore dbe) [virtual]

Inizialize parameters for histo binning.

Implements PFJetAnalyzerBase.

Definition at line 46 of file PFJetAnalyzer.cc.

References DQMStore::book1D(), DQMStore::book2D(), DQMStore::bookProfile(), jptDQMConfig_cff::eMax, jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, LogTrace, metname, Parameters::parameters, jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, jptDQMConfig_cff::pMax, jptDQMConfig_cff::ptMax, PtMinSelector_cfg::ptMin, MonitorElement::setBinLabel(), and DQMStore::setCurrentFolder().

                                           {

  metname = "pFJetAnalyzer";

  LogTrace(metname)<<"[PFJetAnalyzer] Parameters initialization";
  //dbe->setCurrentFolder("JetMET/Jet/PFJets");//old version, now name set to source, which 
  //can be set for each instance of PFJetAnalyzer called inside JetMETAnalyzer. Useful, e.g., to 
  //name differently the dir for all jets and cleaned jets 
  dbe->setCurrentFolder("JetMET/Jet/"+_source);
  // dbe->setCurrentFolder("JetMET/Jet/PFJets");
                        
  jetME = dbe->book1D("jetReco", "jetReco", 3, 1, 4);
  jetME->setBinLabel(2,"PFJets",1);

  // monitoring of eta parameter
  etaBin = parameters.getParameter<int>("etaBin");
  etaMin = parameters.getParameter<double>("etaMin");
  etaMax = parameters.getParameter<double>("etaMax");

  // monitoring of phi paramater
  phiBin = parameters.getParameter<int>("phiBin");
  phiMin = parameters.getParameter<double>("phiMin");
  phiMax = parameters.getParameter<double>("phiMax");

  // monitoring of the transverse momentum
  ptBin = parameters.getParameter<int>("ptBin");
  ptMin = parameters.getParameter<double>("ptMin");
  ptMax = parameters.getParameter<double>("ptMax");

  // 
  eBin = parameters.getParameter<int>("eBin");
  eMin = parameters.getParameter<double>("eMin");
  eMax = parameters.getParameter<double>("eMax");

  // 
  pBin = parameters.getParameter<int>("pBin");
  pMin = parameters.getParameter<double>("pMin");
  pMax = parameters.getParameter<double>("pMax");

  _ptThreshold = parameters.getParameter<double>("ptThreshold");

  _TightCHFMin = parameters.getParameter<double>("TightCHFMin");
  _TightNHFMax = parameters.getParameter<double>("TightNHFMax");
  _TightCEFMax = parameters.getParameter<double>("TightCEFMax");
  _TightNEFMax = parameters.getParameter<double>("TightNEFMax");
  _LooseCHFMin = parameters.getParameter<double>("LooseCHFMin");
  _LooseNHFMax = parameters.getParameter<double>("LooseNHFMax");
  _LooseCEFMax = parameters.getParameter<double>("LooseCEFMax");
  _LooseNEFMax = parameters.getParameter<double>("LooseNEFMax");

  fillpfJIDPassFrac = parameters.getParameter<int>("fillpfJIDPassFrac");

  _ThisCHFMin = parameters.getParameter<double>("ThisCHFMin");
  _ThisNHFMax = parameters.getParameter<double>("ThisNHFMax");
  _ThisCEFMax = parameters.getParameter<double>("ThisCEFMax");
  _ThisNEFMax = parameters.getParameter<double>("ThisNEFMax");

  // Generic Jet Parameters
  mPt                      = dbe->book1D("Pt",  "Pt", ptBin, ptMin, ptMax);
  mPt_1                    = dbe->book1D("Pt1", "Pt1", 100, 0, 100);
  mPt_2                    = dbe->book1D("Pt2", "Pt2", 100, 0, 300);
  mPt_3                    = dbe->book1D("Pt3", "Pt3", 100, 0, 5000);
  mEta                     = dbe->book1D("Eta", "Eta", etaBin, etaMin, etaMax);
  mPhi                     = dbe->book1D("Phi", "Phi", phiBin, phiMin, phiMax);
  mConstituents            = dbe->book1D("Constituents", "# of Constituents", 100, 0, 100);
  mHFrac                   = dbe->book1D("HFrac", "HFrac", 120, -0.1, 1.1);
  mEFrac                   = dbe->book1D("EFrac", "EFrac", 120, -0.1, 1.1);
 //
  mPhiVSEta                     = dbe->book2D("PhiVSEta", "PhiVSEta", 50, etaMin, etaMax, 24, phiMin, phiMax);

  // Low and high pt trigger paths
  mPt_Lo                  = dbe->book1D("Pt_Lo", "Pt (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mEta_Lo                 = dbe->book1D("Eta_Lo", "Eta (Pass Low Pt Jet Trigger)", etaBin, etaMin, etaMax);
  mPhi_Lo                 = dbe->book1D("Phi_Lo", "Phi (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);

  mPt_Hi                  = dbe->book1D("Pt_Hi", "Pt (Pass Hi Pt Jet Trigger)", 100, 0, 300);
  mEta_Hi                 = dbe->book1D("Eta_Hi", "Eta (Pass Hi Pt Jet Trigger)", etaBin, etaMin, etaMax);
  mPhi_Hi                 = dbe->book1D("Phi_Hi", "Phi (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);

  mE                       = dbe->book1D("E", "E", eBin, eMin, eMax);
  mP                       = dbe->book1D("P", "P", pBin, pMin, pMax);
  mMass                    = dbe->book1D("Mass", "Mass", 100, 0, 25);
  mNJets                   = dbe->book1D("NJets", "Number of Jets", 100, 0, 100);

  mPt_Barrel_Lo            = dbe->book1D("Pt_Barrel_Lo", "Pt Barrel (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mPhi_Barrel_Lo           = dbe->book1D("Phi_Barrel_Lo", "Phi Barrel (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
  mConstituents_Barrel_Lo  = dbe->book1D("Constituents_Barrel_Lo", "Constituents Barrel (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mHFrac_Barrel_Lo         = dbe->book1D("HFrac_Barrel_Lo", "HFrac Barrel (Pass Low Pt Jet Trigger)", 100, 0, 1);

  mPt_EndCap_Lo            = dbe->book1D("Pt_EndCap_Lo", "Pt EndCap (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mPhi_EndCap_Lo           = dbe->book1D("Phi_EndCap_Lo", "Phi EndCap (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
  mConstituents_EndCap_Lo  = dbe->book1D("Constituents_EndCap_Lo", "Constituents EndCap (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mHFrac_EndCap_Lo         = dbe->book1D("HFrac_Endcap_Lo", "HFrac EndCap (Pass Low Pt Jet Trigger)", 100, 0, 1);

  mPt_Forward_Lo           = dbe->book1D("Pt_Forward_Lo", "Pt Forward (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mPhi_Forward_Lo          = dbe->book1D("Phi_Forward_Lo", "Phi Forward (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
  mConstituents_Forward_Lo = dbe->book1D("Constituents_Forward_Lo", "Constituents Forward (Pass Low Pt Jet Trigger)", 100, 0, 100);
  mHFrac_Forward_Lo        = dbe->book1D("HFrac_Forward_Lo", "HFrac Forward (Pass Low Pt Jet Trigger)", 100, 0, 1);

  mPt_Barrel_Hi            = dbe->book1D("Pt_Barrel_Hi", "Pt Barrel (Pass Hi Pt Jet Trigger)", 100, 0, 300);
  mPhi_Barrel_Hi           = dbe->book1D("Phi_Barrel_Hi", "Phi Barrel (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
  mConstituents_Barrel_Hi  = dbe->book1D("Constituents_Barrel_Hi", "Constituents Barrel (Pass Hi Pt Jet Trigger)", 100, 0, 100);
  mHFrac_Barrel_Hi         = dbe->book1D("HFrac_Barrel_Hi", "HFrac Barrel (Pass Hi Pt Jet Trigger)", 100, 0, 1);

  mPt_EndCap_Hi            = dbe->book1D("Pt_EndCap_Hi", "Pt EndCap (Pass Hi Pt Jet Trigger)", 100, 0, 300);
  mPhi_EndCap_Hi           = dbe->book1D("Phi_EndCap_Hi", "Phi EndCap (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
  mConstituents_EndCap_Hi  = dbe->book1D("Constituents_EndCap_Hi", "Constituents EndCap (Pass Hi Pt Jet Trigger)", 100, 0, 100);
  mHFrac_EndCap_Hi         = dbe->book1D("HFrac_EndCap_Hi", "HFrac EndCap (Pass Hi Pt Jet Trigger)", 100, 0, 1);

  mPt_Forward_Hi           = dbe->book1D("Pt_Forward_Hi", "Pt Forward (Pass Hi Pt Jet Trigger)", 100, 0, 300);
  mPhi_Forward_Hi          = dbe->book1D("Phi_Forward_Hi", "Phi Forward (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
  mConstituents_Forward_Hi = dbe->book1D("Constituents_Forward_Hi", "Constituents Forward (Pass Hi Pt Jet Trigger)", 100, 0, 100);
  mHFrac_Forward_Hi        = dbe->book1D("HFrac_Forward_Hi", "HFrac Forward (Pass Hi Pt Jet Trigger)", 100, 0, 1);

  mPhi_Barrel              = dbe->book1D("Phi_Barrel", "Phi_Barrel", phiBin, phiMin, phiMax);
  mE_Barrel                = dbe->book1D("E_Barrel", "E_Barrel", eBin, eMin, eMax);
  mPt_Barrel               = dbe->book1D("Pt_Barrel", "Pt_Barrel", ptBin, ptMin, ptMax);

  mPhi_EndCap              = dbe->book1D("Phi_EndCap", "Phi_EndCap", phiBin, phiMin, phiMax);
  mE_EndCap                = dbe->book1D("E_EndCap", "E_EndCap", eBin, eMin, eMax);
  mPt_EndCap               = dbe->book1D("Pt_EndCap", "Pt_EndCap", ptBin, ptMin, ptMax);

  mPhi_Forward             = dbe->book1D("Phi_Forward", "Phi_Forward", phiBin, phiMin, phiMax);
  mE_Forward               = dbe->book1D("E_Forward", "E_Forward", eBin, eMin, eMax);
  mPt_Forward              = dbe->book1D("Pt_Forward", "Pt_Forward", ptBin, ptMin, ptMax);

  // Leading Jet Parameters
  mEtaFirst                = dbe->book1D("EtaFirst", "EtaFirst", 100, -5, 5);
  mPhiFirst                = dbe->book1D("PhiFirst", "PhiFirst", 70, -3.5, 3.5);
  mEFirst                  = dbe->book1D("EFirst", "EFirst", 100, 0, 1000);
  mPtFirst                 = dbe->book1D("PtFirst", "PtFirst", 100, 0, 500);

  mDPhi                   = dbe->book1D("DPhi", "dPhi btw the two leading jets", 100, 0., acos(-1.));

  //
  mChargedHadronEnergy = dbe->book1D("mChargedHadronEnergy", "mChargedHadronEnergy", 100, 0, 100);
  mNeutralHadronEnergy = dbe->book1D("mNeutralHadronEnergy", "mNeutralHadronEnergy", 100, 0, 100);
  mChargedEmEnergy= dbe->book1D("mChargedEmEnergy ", "mChargedEmEnergy ", 100, 0, 100);
  mChargedMuEnergy = dbe->book1D("mChargedMuEnergy", "mChargedMuEnergy", 100, 0, 100);
  mNeutralEmEnergy= dbe->book1D("mNeutralEmEnergy", "mNeutralEmEnergy", 100, 0, 100);
  mChargedMultiplicity= dbe->book1D("mChargedMultiplicity ", "mChargedMultiplicity ", 100, 0, 100);
  mNeutralMultiplicity = dbe->book1D(" mNeutralMultiplicity", "mNeutralMultiplicity", 100, 0, 100);
  mMuonMultiplicity= dbe->book1D("mMuonMultiplicity", "mMuonMultiplicity", 100, 0, 100);
  //__________________________________________________
  mNeutralFraction = dbe->book1D("NeutralFraction","Neutral Fraction",100,0,1);
  if(fillpfJIDPassFrac==1) {
    mLooseJIDPassFractionVSeta= dbe->bookProfile("LooseJIDPassFractionVSeta","LooseJIDPassFractionVSeta",etaBin, etaMin, etaMax,0.,1.2);
    mLooseJIDPassFractionVSpt= dbe->bookProfile("LooseJIDPassFractionVSpt","LooseJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
    mTightJIDPassFractionVSeta= dbe->bookProfile("TightJIDPassFractionVSeta","TightJIDPassFractionVSeta",etaBin, etaMin, etaMax,0.,1.2);
    mTightJIDPassFractionVSpt= dbe->bookProfile("TightJIDPassFractionVSpt","TightJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
  }

  
}
int PFJetAnalyzer::getLeadJetFlag ( ) [inline]

Definition at line 52 of file PFJetAnalyzer.h.

References _leadJetFlag.

                       {
    return  _leadJetFlag;
  }
void PFJetAnalyzer::setJetHiPass ( int  pass) [inline]

Definition at line 59 of file PFJetAnalyzer.h.

References _JetHiPass.

                              {
    _JetHiPass = pass;
  }
void PFJetAnalyzer::setJetLoPass ( int  pass) [inline]

Definition at line 55 of file PFJetAnalyzer.h.

References _JetLoPass.

                              {
    _JetLoPass = pass;
  }
void PFJetAnalyzer::setLeadJetFlag ( int  flag) [inline]

Definition at line 49 of file PFJetAnalyzer.h.

References _leadJetFlag.

void PFJetAnalyzer::setSource ( std::string  source) [inline]

Definition at line 45 of file PFJetAnalyzer.h.

References _source, and LaserTracksInput_cfi::source.

                                   {
    _source = source;
  }

Member Data Documentation

Definition at line 68 of file PFJetAnalyzer.h.

Referenced by setJetHiPass().

Definition at line 67 of file PFJetAnalyzer.h.

Referenced by setJetLoPass().

Definition at line 69 of file PFJetAnalyzer.h.

Referenced by getLeadJetFlag(), and setLeadJetFlag().

double PFJetAnalyzer::_LooseCEFMax [private]

Definition at line 100 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_LooseCHFMin [private]

Definition at line 98 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_LooseNEFMax [private]

Definition at line 101 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_LooseNHFMax [private]

Definition at line 99 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_ptThreshold [private]

Definition at line 70 of file PFJetAnalyzer.h.

std::string PFJetAnalyzer::_source [private]

Definition at line 115 of file PFJetAnalyzer.h.

Referenced by setSource().

double PFJetAnalyzer::_ThisCEFMax [private]

Definition at line 96 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_ThisCHFMin [private]

Definition at line 94 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_ThisNEFMax [private]

Definition at line 97 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_ThisNHFMax [private]

Definition at line 95 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_TightCEFMax [private]

Definition at line 104 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_TightCHFMin [private]

Definition at line 102 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_TightNEFMax [private]

Definition at line 105 of file PFJetAnalyzer.h.

double PFJetAnalyzer::_TightNHFMax [private]

Definition at line 103 of file PFJetAnalyzer.h.

int PFJetAnalyzer::eBin [private]

Definition at line 84 of file PFJetAnalyzer.h.

double PFJetAnalyzer::eMax [private]

Definition at line 86 of file PFJetAnalyzer.h.

double PFJetAnalyzer::eMin [private]

Definition at line 85 of file PFJetAnalyzer.h.

int PFJetAnalyzer::etaBin [private]

Definition at line 72 of file PFJetAnalyzer.h.

double PFJetAnalyzer::etaMax [private]

Definition at line 74 of file PFJetAnalyzer.h.

double PFJetAnalyzer::etaMin [private]

Definition at line 73 of file PFJetAnalyzer.h.

Definition at line 92 of file PFJetAnalyzer.h.

Definition at line 108 of file PFJetAnalyzer.h.

Definition at line 207 of file PFJetAnalyzer.h.

Definition at line 205 of file PFJetAnalyzer.h.

Definition at line 208 of file PFJetAnalyzer.h.

Definition at line 210 of file PFJetAnalyzer.h.

Definition at line 137 of file PFJetAnalyzer.h.

Definition at line 166 of file PFJetAnalyzer.h.

Definition at line 153 of file PFJetAnalyzer.h.

Definition at line 170 of file PFJetAnalyzer.h.

Definition at line 157 of file PFJetAnalyzer.h.

Definition at line 174 of file PFJetAnalyzer.h.

Definition at line 161 of file PFJetAnalyzer.h.

Definition at line 187 of file PFJetAnalyzer.h.

Definition at line 183 of file PFJetAnalyzer.h.

Definition at line 179 of file PFJetAnalyzer.h.

Definition at line 180 of file PFJetAnalyzer.h.

Definition at line 181 of file PFJetAnalyzer.h.

Definition at line 220 of file PFJetAnalyzer.h.

Definition at line 219 of file PFJetAnalyzer.h.

Definition at line 192 of file PFJetAnalyzer.h.

Definition at line 139 of file PFJetAnalyzer.h.

Definition at line 222 of file PFJetAnalyzer.h.

Definition at line 221 of file PFJetAnalyzer.h.

Definition at line 135 of file PFJetAnalyzer.h.

Definition at line 200 of file PFJetAnalyzer.h.

Definition at line 196 of file PFJetAnalyzer.h.

Definition at line 190 of file PFJetAnalyzer.h.

std::string PFJetAnalyzer::metname [private]

Definition at line 112 of file PFJetAnalyzer.h.

Definition at line 138 of file PFJetAnalyzer.h.

Definition at line 167 of file PFJetAnalyzer.h.

Definition at line 154 of file PFJetAnalyzer.h.

Definition at line 171 of file PFJetAnalyzer.h.

Definition at line 158 of file PFJetAnalyzer.h.

Definition at line 175 of file PFJetAnalyzer.h.

Definition at line 162 of file PFJetAnalyzer.h.

Definition at line 225 of file PFJetAnalyzer.h.

Definition at line 226 of file PFJetAnalyzer.h.

Definition at line 185 of file PFJetAnalyzer.h.

Definition at line 212 of file PFJetAnalyzer.h.

Definition at line 209 of file PFJetAnalyzer.h.

Definition at line 216 of file PFJetAnalyzer.h.

Definition at line 217 of file PFJetAnalyzer.h.

Definition at line 206 of file PFJetAnalyzer.h.

Definition at line 211 of file PFJetAnalyzer.h.

Definition at line 186 of file PFJetAnalyzer.h.

Definition at line 184 of file PFJetAnalyzer.h.

Definition at line 136 of file PFJetAnalyzer.h.

Definition at line 143 of file PFJetAnalyzer.h.

Definition at line 165 of file PFJetAnalyzer.h.

Definition at line 152 of file PFJetAnalyzer.h.

Definition at line 146 of file PFJetAnalyzer.h.

Definition at line 169 of file PFJetAnalyzer.h.

Definition at line 156 of file PFJetAnalyzer.h.

Definition at line 149 of file PFJetAnalyzer.h.

Definition at line 173 of file PFJetAnalyzer.h.

Definition at line 160 of file PFJetAnalyzer.h.

Definition at line 201 of file PFJetAnalyzer.h.

Definition at line 197 of file PFJetAnalyzer.h.

Definition at line 191 of file PFJetAnalyzer.h.

Definition at line 140 of file PFJetAnalyzer.h.

Definition at line 131 of file PFJetAnalyzer.h.

Definition at line 132 of file PFJetAnalyzer.h.

Definition at line 133 of file PFJetAnalyzer.h.

Definition at line 134 of file PFJetAnalyzer.h.

Definition at line 142 of file PFJetAnalyzer.h.

Definition at line 164 of file PFJetAnalyzer.h.

Definition at line 151 of file PFJetAnalyzer.h.

Definition at line 145 of file PFJetAnalyzer.h.

Definition at line 168 of file PFJetAnalyzer.h.

Definition at line 155 of file PFJetAnalyzer.h.

Definition at line 148 of file PFJetAnalyzer.h.

Definition at line 172 of file PFJetAnalyzer.h.

Definition at line 159 of file PFJetAnalyzer.h.

Definition at line 202 of file PFJetAnalyzer.h.

Definition at line 198 of file PFJetAnalyzer.h.

Definition at line 193 of file PFJetAnalyzer.h.

Definition at line 227 of file PFJetAnalyzer.h.

Definition at line 228 of file PFJetAnalyzer.h.

Definition at line 223 of file PFJetAnalyzer.h.

Definition at line 110 of file PFJetAnalyzer.h.

int PFJetAnalyzer::pBin [private]

Definition at line 88 of file PFJetAnalyzer.h.

int PFJetAnalyzer::phiBin [private]

Definition at line 76 of file PFJetAnalyzer.h.

double PFJetAnalyzer::phiMax [private]

Definition at line 78 of file PFJetAnalyzer.h.

double PFJetAnalyzer::phiMin [private]

Definition at line 77 of file PFJetAnalyzer.h.

double PFJetAnalyzer::pMax [private]

Definition at line 90 of file PFJetAnalyzer.h.

double PFJetAnalyzer::pMin [private]

Definition at line 89 of file PFJetAnalyzer.h.

int PFJetAnalyzer::ptBin [private]

Definition at line 80 of file PFJetAnalyzer.h.

double PFJetAnalyzer::ptMax [private]

Definition at line 82 of file PFJetAnalyzer.h.

double PFJetAnalyzer::ptMin [private]

Definition at line 81 of file PFJetAnalyzer.h.

Definition at line 114 of file PFJetAnalyzer.h.