CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/DQMOffline/JetMET/src/PFJetAnalyzer.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2011/08/12 15:29:40 $
00005  *  $Revision: 1.15.6.1 $
00006  *  \author F. Chlebana - Fermilab
00007  */
00008 
00009 #include "DQMOffline/JetMET/interface/PFJetAnalyzer.h"
00010 #include "DataFormats/Common/interface/Handle.h"
00011 
00012 #include "DataFormats/JetReco/interface/PFJet.h"
00013 
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 
00016 #include <string>
00017 using namespace edm;
00018 
00019 
00020 PFJetAnalyzer::PFJetAnalyzer(const edm::ParameterSet& pSet) {
00021 
00022   parameters = pSet;
00023   _leadJetFlag = 0;
00024   _JetLoPass   = 0;
00025   _JetHiPass   = 0;
00026   _ptThreshold = 5.;
00027   _LooseCHFMin = -999.;
00028   _LooseNHFMax = -999.;
00029   _LooseCEFMax = -999.;
00030   _LooseNEFMax = -999.;
00031   _TightCHFMin = -999.;
00032   _TightNHFMax = -999.;
00033   _TightCEFMax = -999.;
00034   _TightNEFMax = -999.;
00035   _ThisCHFMin = -999.;
00036   _ThisNHFMax = -999.;
00037   _ThisCEFMax = -999.;
00038   _ThisNEFMax = -999.;
00039 
00040 }
00041 
00042 
00043 PFJetAnalyzer::~PFJetAnalyzer() { }
00044 
00045 
00046 void PFJetAnalyzer::beginJob(DQMStore * dbe) {
00047 
00048   metname = "pFJetAnalyzer";
00049 
00050   LogTrace(metname)<<"[PFJetAnalyzer] Parameters initialization";
00051   //dbe->setCurrentFolder("JetMET/Jet/PFJets");//old version, now name set to source, which 
00052   //can be set for each instance of PFJetAnalyzer called inside JetMETAnalyzer. Useful, e.g., to 
00053   //name differently the dir for all jets and cleaned jets 
00054   dbe->setCurrentFolder("JetMET/Jet/"+_source);
00055   // dbe->setCurrentFolder("JetMET/Jet/PFJets");
00056                         
00057   jetME = dbe->book1D("jetReco", "jetReco", 3, 1, 4);
00058   jetME->setBinLabel(2,"PFJets",1);
00059 
00060   // monitoring of eta parameter
00061   etaBin = parameters.getParameter<int>("etaBin");
00062   etaMin = parameters.getParameter<double>("etaMin");
00063   etaMax = parameters.getParameter<double>("etaMax");
00064 
00065   // monitoring of phi paramater
00066   phiBin = parameters.getParameter<int>("phiBin");
00067   phiMin = parameters.getParameter<double>("phiMin");
00068   phiMax = parameters.getParameter<double>("phiMax");
00069 
00070   // monitoring of the transverse momentum
00071   ptBin = parameters.getParameter<int>("ptBin");
00072   ptMin = parameters.getParameter<double>("ptMin");
00073   ptMax = parameters.getParameter<double>("ptMax");
00074 
00075   // 
00076   eBin = parameters.getParameter<int>("eBin");
00077   eMin = parameters.getParameter<double>("eMin");
00078   eMax = parameters.getParameter<double>("eMax");
00079 
00080   // 
00081   pBin = parameters.getParameter<int>("pBin");
00082   pMin = parameters.getParameter<double>("pMin");
00083   pMax = parameters.getParameter<double>("pMax");
00084 
00085   _ptThreshold = parameters.getParameter<double>("ptThreshold");
00086 
00087   _TightCHFMin = parameters.getParameter<double>("TightCHFMin");
00088   _TightNHFMax = parameters.getParameter<double>("TightNHFMax");
00089   _TightCEFMax = parameters.getParameter<double>("TightCEFMax");
00090   _TightNEFMax = parameters.getParameter<double>("TightNEFMax");
00091   _LooseCHFMin = parameters.getParameter<double>("LooseCHFMin");
00092   _LooseNHFMax = parameters.getParameter<double>("LooseNHFMax");
00093   _LooseCEFMax = parameters.getParameter<double>("LooseCEFMax");
00094   _LooseNEFMax = parameters.getParameter<double>("LooseNEFMax");
00095 
00096   fillpfJIDPassFrac = parameters.getParameter<int>("fillpfJIDPassFrac");
00097 
00098   _ThisCHFMin = parameters.getParameter<double>("ThisCHFMin");
00099   _ThisNHFMax = parameters.getParameter<double>("ThisNHFMax");
00100   _ThisCEFMax = parameters.getParameter<double>("ThisCEFMax");
00101   _ThisNEFMax = parameters.getParameter<double>("ThisNEFMax");
00102 
00103   // Generic Jet Parameters
00104   mPt                      = dbe->book1D("Pt",  "Pt", ptBin, ptMin, ptMax);
00105   mPt_1                    = dbe->book1D("Pt1", "Pt1", 50, 0, 100);
00106   mPt_2                    = dbe->book1D("Pt2", "Pt2", 100, 0, 500);
00107   mPt_3                    = dbe->book1D("Pt3", "Pt3", 100, 0, 5000);
00108   mEta                     = dbe->book1D("Eta", "Eta", etaBin, etaMin, etaMax);
00109   mPhi                     = dbe->book1D("Phi", "Phi", phiBin, phiMin, phiMax);
00110   mConstituents            = dbe->book1D("Constituents", "# of Constituents", 50, 0, 100);
00111   mHFrac                   = dbe->book1D("HFrac", "HFrac", 120, -0.1, 1.1);
00112   mEFrac                   = dbe->book1D("EFrac", "EFrac", 120, -0.1, 1.1);
00113  //
00114   mPhiVSEta                     = dbe->book2D("PhiVSEta", "PhiVSEta", 50, etaMin, etaMax, 24, phiMin, phiMax);
00115 
00116   // Low and high pt trigger paths
00117   mPt_Lo                  = dbe->book1D("Pt_Lo", "Pt (Pass Low Pt Jet Trigger)", 20, 0, 100);
00118   mEta_Lo                 = dbe->book1D("Eta_Lo", "Eta (Pass Low Pt Jet Trigger)", etaBin, etaMin, etaMax);
00119   mPhi_Lo                 = dbe->book1D("Phi_Lo", "Phi (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
00120 
00121   mPt_Hi                  = dbe->book1D("Pt_Hi", "Pt (Pass Hi Pt Jet Trigger)", 100, 0, 500);
00122   mEta_Hi                 = dbe->book1D("Eta_Hi", "Eta (Pass Hi Pt Jet Trigger)", etaBin, etaMin, etaMax);
00123   mPhi_Hi                 = dbe->book1D("Phi_Hi", "Phi (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
00124 
00125   //removed for optimizations//mE                       = dbe->book1D("E", "E", eBin, eMin, eMax);
00126   //removed for optimizations//mP                       = dbe->book1D("P", "P", pBin, pMin, pMax);
00127   //removed for optimizations//mMass                    = dbe->book1D("Mass", "Mass", 100, 0, 25);
00128   mNJets                   = dbe->book1D("NJets", "Number of Jets", 100, 0, 100);
00129 
00130   mConstituents_Barrel     = dbe->book1D("Constituents_Barrel", "Constituents Barrel", 50, 0, 100);
00131   mHFrac_Barrel            = dbe->book1D("HFrac_Barrel", "HFrac Barrel", 100, 0, 1);
00132   mEFrac_Barrel            = dbe->book1D("EFrac_Barrel", "EFrac Barrel", 110, -0.05, 1.05);
00133   //removed for optimization//mPt_Barrel_Lo            = dbe->book1D("Pt_Barrel_Lo", "Pt Barrel (Pass Low Pt Jet Trigger)", 20, 0, 100);
00134   //removed for optimization//mPhi_Barrel_Lo           = dbe->book1D("Phi_Barrel_Lo", "Phi Barrel (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
00135   //removed for optimization//mConstituents_Barrel_Lo  = dbe->book1D("Constituents_Barrel_Lo", "Constituents Barrel (Pass Low Pt Jet Trigger)", 50, 0, 100);
00136   //removed for optimization//mHFrac_Barrel_Lo         = dbe->book1D("HFrac_Barrel_Lo", "HFrac Barrel (Pass Low Pt Jet Trigger)", 100, 0, 1);
00137 
00138   mConstituents_EndCap     = dbe->book1D("Constituents_EndCap", "Constituents EndCap", 50, 0, 100);
00139   mHFrac_EndCap            = dbe->book1D("HFrac_EndCap", "HFrac EndCap", 100, 0, 1);
00140   mEFrac_EndCap            = dbe->book1D("EFrac_EndCap", "EFrac EndCap", 110, -0.05, 1.05);
00141   //removed for optimization//mPt_EndCap_Lo            = dbe->book1D("Pt_EndCap_Lo", "Pt EndCap (Pass Low Pt Jet Trigger)", 20, 0, 100);
00142   //removed for optimization//mPhi_EndCap_Lo           = dbe->book1D("Phi_EndCap_Lo", "Phi EndCap (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
00143   //removed for optimization//mConstituents_EndCap_Lo  = dbe->book1D("Constituents_EndCap_Lo", "Constituents EndCap (Pass Low Pt Jet Trigger)", 50, 0, 100);
00144   //removed for optimization//mHFrac_EndCap_Lo         = dbe->book1D("HFrac_Endcap_Lo", "HFrac EndCap (Pass Low Pt Jet Trigger)", 100, 0, 1);
00145 
00146   mConstituents_Forward     = dbe->book1D("Constituents_Forward", "Constituents Forward", 50, 0, 100);
00147   mHFrac_Forward            = dbe->book1D("HFrac_Forward", "HFrac Forward", 100, 0, 1);
00148   mEFrac_Forward            = dbe->book1D("EFrac_Forward", "EFrac Forward", 110, -0.05, 1.05);
00149   //removed for optimization//mPt_Forward_Lo           = dbe->book1D("Pt_Forward_Lo", "Pt Forward (Pass Low Pt Jet Trigger)", 20, 0, 100);
00150   //removed for optimization//mPhi_Forward_Lo          = dbe->book1D("Phi_Forward_Lo", "Phi Forward (Pass Low Pt Jet Trigger)", phiBin, phiMin, phiMax);
00151   //removed for optimization//mConstituents_Forward_Lo = dbe->book1D("Constituents_Forward_Lo", "Constituents Forward (Pass Low Pt Jet Trigger)", 50, 0, 100);
00152   //removed for optimization//mHFrac_Forward_Lo        = dbe->book1D("HFrac_Forward_Lo", "HFrac Forward (Pass Low Pt Jet Trigger)", 100, 0, 1);
00153 
00154   mPt_Barrel_Hi            = dbe->book1D("Pt_Barrel_Hi", "Pt Barrel (Pass Hi Pt Jet Trigger)", 100, 0, 500);
00155   mPhi_Barrel_Hi           = dbe->book1D("Phi_Barrel_Hi", "Phi Barrel (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
00156   //removed for optimization//mConstituents_Barrel_Hi  = dbe->book1D("Constituents_Barrel_Hi", "Constituents Barrel (Pass Hi Pt Jet Trigger)", 50, 0, 100);
00157   //removed for optimization//mHFrac_Barrel_Hi         = dbe->book1D("HFrac_Barrel_Hi", "HFrac Barrel (Pass Hi Pt Jet Trigger)", 100, 0, 1);
00158 
00159   mPt_EndCap_Hi            = dbe->book1D("Pt_EndCap_Hi", "Pt EndCap (Pass Hi Pt Jet Trigger)", 100, 0, 500);
00160   mPhi_EndCap_Hi           = dbe->book1D("Phi_EndCap_Hi", "Phi EndCap (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
00161   //removed for optimization//mConstituents_EndCap_Hi  = dbe->book1D("Constituents_EndCap_Hi", "Constituents EndCap (Pass Hi Pt Jet Trigger)", 50, 0, 100);
00162   //removed for optimization//mHFrac_EndCap_Hi         = dbe->book1D("HFrac_EndCap_Hi", "HFrac EndCap (Pass Hi Pt Jet Trigger)", 100, 0, 1);
00163 
00164   mPt_Forward_Hi           = dbe->book1D("Pt_Forward_Hi", "Pt Forward (Pass Hi Pt Jet Trigger)", 100, 0, 500);
00165   mPhi_Forward_Hi          = dbe->book1D("Phi_Forward_Hi", "Phi Forward (Pass Hi Pt Jet Trigger)", phiBin, phiMin, phiMax);
00166   //removed for optimization//mConstituents_Forward_Hi = dbe->book1D("Constituents_Forward_Hi", "Constituents Forward (Pass Hi Pt Jet Trigger)", 50, 0, 100);
00167   //removed for optimization//mHFrac_Forward_Hi        = dbe->book1D("HFrac_Forward_Hi", "HFrac Forward (Pass Hi Pt Jet Trigger)", 100, 0, 1);
00168 
00169   mPhi_Barrel              = dbe->book1D("Phi_Barrel", "Phi_Barrel", phiBin, phiMin, phiMax);
00170   //removed for optimizations//mE_Barrel                = dbe->book1D("E_Barrel", "E_Barrel", eBin, eMin, eMax);
00171   mPt_Barrel               = dbe->book1D("Pt_Barrel", "Pt_Barrel", ptBin, ptMin, ptMax);
00172 
00173   mPhi_EndCap              = dbe->book1D("Phi_EndCap", "Phi_EndCap", phiBin, phiMin, phiMax);
00174   //removed for optimizations//mE_EndCap                = dbe->book1D("E_EndCap", "E_EndCap", eBin, eMin, eMax);
00175   mPt_EndCap               = dbe->book1D("Pt_EndCap", "Pt_EndCap", ptBin, ptMin, ptMax);
00176 
00177   mPhi_Forward             = dbe->book1D("Phi_Forward", "Phi_Forward", phiBin, phiMin, phiMax);
00178   //removed for optimizations//mE_Forward               = dbe->book1D("E_Forward", "E_Forward", eBin, eMin, eMax);
00179   mPt_Forward              = dbe->book1D("Pt_Forward", "Pt_Forward", ptBin, ptMin, ptMax);
00180 
00181   // Leading Jet Parameters
00182   mEtaFirst                = dbe->book1D("EtaFirst", "EtaFirst", 100, -5, 5);
00183   mPhiFirst                = dbe->book1D("PhiFirst", "PhiFirst", 70, -3.5, 3.5);
00184   //removed for optimizations//mEFirst                  = dbe->book1D("EFirst", "EFirst", 100, 0, 1000);
00185   mPtFirst                 = dbe->book1D("PtFirst", "PtFirst", 100, 0, 500);
00186 
00187   mDPhi                   = dbe->book1D("DPhi", "dPhi btw the two leading jets", 100, 0., acos(-1.));
00188 
00189   //
00190   mChargedHadronEnergy = dbe->book1D("mChargedHadronEnergy", "mChargedHadronEnergy", 100, 0, 100);
00191   mNeutralHadronEnergy = dbe->book1D("mNeutralHadronEnergy", "mNeutralHadronEnergy", 100, 0, 100);
00192   mChargedEmEnergy= dbe->book1D("mChargedEmEnergy ", "mChargedEmEnergy ", 100, 0, 100);
00193   mChargedMuEnergy = dbe->book1D("mChargedMuEnergy", "mChargedMuEnergy", 100, 0, 100);
00194   mNeutralEmEnergy= dbe->book1D("mNeutralEmEnergy", "mNeutralEmEnergy", 100, 0, 100);
00195   mChargedMultiplicity= dbe->book1D("mChargedMultiplicity ", "mChargedMultiplicity ", 100, 0, 100);
00196   mNeutralMultiplicity = dbe->book1D(" mNeutralMultiplicity", "mNeutralMultiplicity", 100, 0, 100);
00197   mMuonMultiplicity= dbe->book1D("mMuonMultiplicity", "mMuonMultiplicity", 100, 0, 100);
00198   //__________________________________________________
00199   mNeutralFraction = dbe->book1D("NeutralFraction","Neutral Fraction",100,0,1);
00200   if(fillpfJIDPassFrac==1) {
00201     mLooseJIDPassFractionVSeta= dbe->bookProfile("LooseJIDPassFractionVSeta","LooseJIDPassFractionVSeta",etaBin, etaMin, etaMax,0.,1.2);
00202     mLooseJIDPassFractionVSpt= dbe->bookProfile("LooseJIDPassFractionVSpt","LooseJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
00203     mTightJIDPassFractionVSeta= dbe->bookProfile("TightJIDPassFractionVSeta","TightJIDPassFractionVSeta",etaBin, etaMin, etaMax,0.,1.2);
00204     mTightJIDPassFractionVSpt= dbe->bookProfile("TightJIDPassFractionVSpt","TightJIDPassFractionVSpt",ptBin, ptMin, ptMax,0.,1.2);
00205   }
00206 
00207   
00208 }
00209 
00210 void PFJetAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::PFJetCollection& pfJets) {
00211 
00212   int numofjets=0;
00213   double  fstPhi=0.;
00214   double  sndPhi=0.;
00215   double  diff = 0.;
00216   double  corr = 0.;
00217   double  dphi = -999. ;
00218 
00219   bool Thiscleaned=false; 
00220   bool Loosecleaned=false; 
00221   bool Tightcleaned=false; 
00222   bool ThisCHFcleaned=false;
00223   bool LooseCHFcleaned=false;
00224   bool TightCHFcleaned=false;
00225 
00226   for (reco::PFJetCollection::const_iterator jet = pfJets.begin(); jet!=pfJets.end(); ++jet){
00227     LogTrace(metname)<<"[JetAnalyzer] Analyze PFJet";
00228  
00229     Thiscleaned=false;
00230     Loosecleaned=false;
00231     Tightcleaned=false;
00232 
00233     if (jet == pfJets.begin()) {
00234       fstPhi = jet->phi();
00235       _leadJetFlag = 1;
00236     } else {
00237       _leadJetFlag = 0;
00238     }
00239 
00240     if (jet == (pfJets.begin()+1)) sndPhi = jet->phi();
00241     //  if (jet->pt() < _ptThreshold) return;
00242     if (jet->pt() > _ptThreshold) {
00243       numofjets++ ;
00244       jetME->Fill(2);
00245             
00246       //calculate the jetID
00247       ThisCHFcleaned=true;
00248       LooseCHFcleaned=true;
00249       TightCHFcleaned=true;
00250       if((jet->chargedHadronEnergy()/jet->energy())<=_ThisCHFMin && fabs(jet->eta())<2.4) ThisCHFcleaned=false; //apply CHF>0 only if |eta|<2.4
00251       if((jet->chargedHadronEnergy()/jet->energy())<=_LooseCHFMin && fabs(jet->eta())<2.4) LooseCHFcleaned=false; //apply CHF>0 only if |eta|<2.4
00252       if((jet->chargedHadronEnergy()/jet->energy())<=_TightCHFMin && fabs(jet->eta())<2.4) TightCHFcleaned=false; //apply CHF>0 only if |eta|<2.4
00253       if(ThisCHFcleaned && (jet->neutralHadronEnergy()/jet->energy())<_ThisNHFMax && (jet->chargedEmEnergy()/jet->energy())<_ThisCEFMax && (jet->neutralEmEnergy()/jet->energy())<_ThisNEFMax) Thiscleaned=true;
00254       if(LooseCHFcleaned && (jet->neutralHadronEnergy()/jet->energy())<_LooseNHFMax && (jet->chargedEmEnergy()/jet->energy())<_LooseCEFMax && (jet->neutralEmEnergy()/jet->energy())<_LooseNEFMax) Loosecleaned=true;
00255       if(TightCHFcleaned && (jet->neutralHadronEnergy()/jet->energy())<_TightNHFMax && (jet->chargedEmEnergy()/jet->energy())<_TightCEFMax && (jet->neutralEmEnergy()/jet->energy())<_TightNEFMax) Tightcleaned=true;
00256       
00257       if(fillpfJIDPassFrac==1) {
00258         //fill the profile for jid efficiency 
00259         if(Loosecleaned) {
00260           mLooseJIDPassFractionVSeta->Fill(jet->eta(),1.);
00261           mLooseJIDPassFractionVSpt->Fill(jet->pt(),1.);
00262         } else {
00263           mLooseJIDPassFractionVSeta->Fill(jet->eta(),0.);
00264           mLooseJIDPassFractionVSpt->Fill(jet->pt(),0.);
00265         }
00266         if(Tightcleaned) {
00267           mTightJIDPassFractionVSeta->Fill(jet->eta(),1.);
00268           mTightJIDPassFractionVSpt->Fill(jet->pt(),1.);
00269         } else {
00270           mTightJIDPassFractionVSeta->Fill(jet->eta(),0.);
00271           mTightJIDPassFractionVSpt->Fill(jet->pt(),0.);
00272         }
00273       }
00274       
00275       if(!Thiscleaned) continue;
00276       
00277       // Leading jet
00278       // Histograms are filled once per event
00279       if (_leadJetFlag == 1) { 
00280         
00281         if (mEtaFirst) mEtaFirst->Fill (jet->eta());
00282         if (mPhiFirst) mPhiFirst->Fill (jet->phi());
00283         //removed for optimizations//if (mEFirst)   mEFirst->Fill (jet->energy());
00284         if (mPtFirst)  mPtFirst->Fill (jet->pt());
00285       }
00286       
00287       // --- Passed the low pt jet trigger
00288       if (_JetLoPass == 1) {
00289         /*if (fabs(jet->eta()) <= 1.3) {
00290           if (mPt_Barrel_Lo)           mPt_Barrel_Lo->Fill(jet->pt());
00291           //if (mEta_Lo)                 mEta_Lo->Fill(jet->eta());
00292           if (mPhi_Barrel_Lo)          mPhi_Barrel_Lo->Fill(jet->phi());
00293           //removed for optimization//if (mConstituents_Barrel_Lo) mConstituents_Barrel_Lo->Fill(jet->nConstituents()); 
00294           //removed for optimization//if (mHFrac_Barrel_Lo)        mHFrac_Barrel_Lo->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction() );      
00295         }
00296         if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
00297           if (mPt_EndCap_Lo)           mPt_EndCap_Lo->Fill(jet->pt());
00298           //if (mEta_Lo)               mEta_Lo->Fill(jet->eta());
00299           if (mPhi_EndCap_Lo)          mPhi_EndCap_Lo->Fill(jet->phi());
00300           //removed for optimization//if (mConstituents_EndCap_Lo) mConstituents_EndCap_Lo->Fill(jet->nConstituents()); 
00301           //removed for optimization//if (mHFrac_EndCap_Lo)        mHFrac_EndCap_Lo->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());       
00302         }
00303         if (fabs(jet->eta()) > 3.0) {
00304           if (mPt_Forward_Lo)           mPt_Forward_Lo->Fill(jet->pt());
00305           //if (mEta_Lo)                mEta_Lo->Fill(jet->eta());
00306           if (mPhi_Forward_Lo)          mPhi_Forward_Lo->Fill(jet->phi());
00307           //removed for optimization//if (mConstituents_Forward_Lo) mConstituents_Forward_Lo->Fill(jet->nConstituents());       
00308           //removed for optimization//if (mHFrac_Forward_Lo)        mHFrac_Forward_Lo->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());     
00309           }*/
00310         if (mEta_Lo) mEta_Lo->Fill (jet->eta());
00311         if (mPhi_Lo) mPhi_Lo->Fill (jet->phi());
00312         if (mPt_Lo)  mPt_Lo->Fill (jet->pt());
00313       }
00314       
00315       // --- Passed the high pt jet trigger
00316       if (_JetHiPass == 1) {
00317         if (fabs(jet->eta()) <= 1.3) {
00318           if (mPt_Barrel_Hi)           mPt_Barrel_Hi->Fill(jet->pt());
00319           //if (mEta_Hi)                 mEta_Hi->Fill(jet->eta());
00320           if (mPhi_Barrel_Hi)          mPhi_Barrel_Hi->Fill(jet->phi());
00321           //removed for optimization//if (mConstituents_Barrel_Hi) mConstituents_Barrel_Hi->Fill(jet->nConstituents()); 
00322           //removed for optimization//if (mHFrac_Barrel_Hi)        mHFrac_Barrel_Hi->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());       
00323         }
00324         if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
00325           if (mPt_EndCap_Hi)           mPt_EndCap_Hi->Fill(jet->pt());
00326           //if (mEta_Hi)                 mEta_Hi->Fill(jet->eta());
00327           if (mPhi_EndCap_Hi)          mPhi_EndCap_Hi->Fill(jet->phi());
00328           //removed for optimization//if (mConstituents_EndCap_Hi) mConstituents_EndCap_Hi->Fill(jet->nConstituents()); 
00329           //removed for optimization//if (mHFrac_EndCap_Hi)        mHFrac_EndCap_Hi->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());       
00330         }
00331         if (fabs(jet->eta()) > 3.0) {
00332           if (mPt_Forward_Hi)           mPt_Forward_Hi->Fill(jet->pt());
00333           //if (mEta_Hi)                  mEta_Hi->Fill(jet->eta());
00334           if (mPhi_Forward_Hi)          mPhi_Forward_Hi->Fill(jet->phi());
00335           //removed for optimization//if (mConstituents_Forward_Hi) mConstituents_Forward_Hi->Fill(jet->nConstituents());       
00336           //removed for optimization//if (mHFrac_Forward_Hi)        mHFrac_Forward_Hi->Fill(jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());     
00337         }
00338         
00339         if (mEta_Hi) mEta_Hi->Fill (jet->eta());
00340         if (mPhi_Hi) mPhi_Hi->Fill (jet->phi());
00341         if (mPt_Hi)  mPt_Hi->Fill (jet->pt());
00342       }
00343       
00344       if (mPt)   mPt->Fill (jet->pt());
00345       if (mPt_1) mPt_1->Fill (jet->pt());
00346       if (mPt_2) mPt_2->Fill (jet->pt());
00347       if (mPt_3) mPt_3->Fill (jet->pt());
00348       if (mEta)  mEta->Fill (jet->eta());
00349       if (mPhi)  mPhi->Fill (jet->phi());
00350       if (mPhiVSEta) mPhiVSEta->Fill(jet->eta(),jet->phi());
00351       
00352       if (mConstituents) mConstituents->Fill (jet->nConstituents());
00353       if (mHFrac)        mHFrac->Fill (jet->chargedHadronEnergyFraction()+jet->neutralHadronEnergyFraction());
00354       if (mEFrac)        mEFrac->Fill (jet->chargedEmEnergyFraction() +jet->neutralEmEnergyFraction());
00355       
00356       if (fabs(jet->eta()) <= 1.3) {
00357         if (mPt_Barrel)   mPt_Barrel->Fill (jet->pt());
00358         if (mPhi_Barrel)  mPhi_Barrel->Fill (jet->phi());
00359         //removed for optimizations//if (mE_Barrel)    mE_Barrel->Fill (jet->energy());
00360         if (mConstituents_Barrel)    mConstituents_Barrel->Fill(jet->nConstituents());  
00361         if (mHFrac_Barrel)           mHFrac_Barrel->Fill(jet->chargedHadronEnergyFraction() + jet->neutralHadronEnergyFraction() );
00362         if (mEFrac_Barrel)           mEFrac->Fill (jet->chargedEmEnergyFraction() + jet->neutralEmEnergyFraction());    
00363       }
00364       if ( (fabs(jet->eta()) > 1.3) && (fabs(jet->eta()) <= 3) ) {
00365         if (mPt_EndCap)   mPt_EndCap->Fill (jet->pt());
00366         if (mPhi_EndCap)  mPhi_EndCap->Fill (jet->phi());
00367         //removed for optimizations//if (mE_EndCap)    mE_EndCap->Fill (jet->energy());
00368         if (mConstituents_EndCap)    mConstituents_EndCap->Fill(jet->nConstituents());  
00369         if (mHFrac_EndCap)           mHFrac_EndCap->Fill(jet->chargedHadronEnergyFraction() + jet->neutralHadronEnergyFraction() );
00370         if (mEFrac_EndCap)           mEFrac->Fill (jet->chargedEmEnergyFraction() + jet->neutralEmEnergyFraction());    
00371       }
00372       if (fabs(jet->eta()) > 3.0) {
00373         if (mPt_Forward)   mPt_Forward->Fill (jet->pt());
00374         if (mPhi_Forward)  mPhi_Forward->Fill (jet->phi());
00375         //removed for optimizations//if (mE_Forward)    mE_Forward->Fill (jet->energy());
00376         if (mConstituents_Forward)    mConstituents_Forward->Fill(jet->nConstituents());        
00377         if (mHFrac_Forward)           mHFrac_Forward->Fill(jet->chargedHadronEnergyFraction() + jet->neutralHadronEnergyFraction() );
00378         if (mEFrac_Forward)           mEFrac->Fill (jet->chargedEmEnergyFraction() + jet->neutralEmEnergyFraction());   
00379       }
00380       
00381       //removed for optimizations//if (mE)    mE->Fill (jet->energy());
00382       //removed for optimizations//if (mP)    mP->Fill (jet->p());
00383       //removed for optimizations//if (mMass) mMass->Fill (jet->mass());
00384       
00385       
00386       if (mChargedHadronEnergy)  mChargedHadronEnergy->Fill (jet->chargedHadronEnergy());
00387       if (mNeutralHadronEnergy)  mNeutralHadronEnergy->Fill (jet->neutralHadronEnergy());
00388       if (mChargedEmEnergy) mChargedEmEnergy->Fill(jet->chargedEmEnergy());
00389       if (mChargedMuEnergy) mChargedMuEnergy->Fill (jet->chargedMuEnergy ());
00390       if (mNeutralEmEnergy) mNeutralEmEnergy->Fill(jet->neutralEmEnergy());
00391       if (mChargedMultiplicity ) mChargedMultiplicity->Fill(jet->chargedMultiplicity());
00392       if (mNeutralMultiplicity ) mNeutralMultiplicity->Fill(jet->neutralMultiplicity());
00393       if (mMuonMultiplicity )mMuonMultiplicity->Fill (jet-> muonMultiplicity());
00394       //_______________________________________________________
00395       if (mNeutralFraction) mNeutralFraction->Fill (jet->neutralMultiplicity()/jet->nConstituents());
00396       
00397       //calculate correctly the dphi
00398       if(numofjets>1) {
00399         diff = fabs(fstPhi - sndPhi);
00400         corr = 2*acos(-1.) - diff;
00401         if(diff < acos(-1.)) { 
00402           dphi = diff; 
00403         } else { 
00404           dphi = corr;
00405         }
00406       } // numofjets>1
00407     } // JetPt>_ptThreshold
00408   } // PF jet loop
00409   if (mNJets)    mNJets->Fill (numofjets);
00410   if (mDPhi)    mDPhi->Fill (dphi);
00411 }