CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Validation/RecoJets/plugins/CaloJetTester.cc

Go to the documentation of this file.
00001 // Producer for validation histograms for CaloJet objects
00002 // F. Ratnikov, Sept. 7, 2006
00003 // Modified by J F Novak July 10, 2008
00004 // $Id: CaloJetTester.cc,v 1.24 2010/08/07 14:56:05 wmtan Exp $
00005 
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/ServiceRegistry/interface/Service.h"
00011 
00012 #include "DataFormats/Math/interface/deltaR.h"
00013 
00014 #include "DQMServices/Core/interface/DQMStore.h"
00015 #include "DQMServices/Core/interface/MonitorElement.h"
00016 
00017 #include "DataFormats/JetReco/interface/CaloJet.h"
00018 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00019 #include "DataFormats/JetReco/interface/GenJet.h"
00020 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00021 
00022 //I don't know which of these I actually need yet
00023 #include "DataFormats/METReco/interface/CaloMET.h"
00024 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00025 #include "DataFormats/METReco/interface/GenMET.h"
00026 #include "DataFormats/METReco/interface/GenMETCollection.h"
00027 #include "DataFormats/METReco/interface/MET.h"
00028 #include "DataFormats/METReco/interface/METCollection.h"
00029 
00030 #include "RecoJets/JetProducers/interface/JetMatchingTools.h"
00031 
00032 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00033 
00034 #include "CaloJetTester.h"
00035 
00036 #include <cmath>
00037 
00038 using namespace edm;
00039 using namespace reco;
00040 using namespace std;
00041 
00042 namespace {
00043   bool is_B (const reco::Jet& fJet) {return fabs (fJet.eta()) < 1.3;}
00044   bool is_E (const reco::Jet& fJet) {return fabs (fJet.eta()) >= 1.3 && fabs (fJet.eta()) < 3.;}
00045   bool is_F (const reco::Jet& fJet) {return fabs (fJet.eta()) >= 3.;}
00046 }
00047 
00048 CaloJetTester::CaloJetTester(const edm::ParameterSet& iConfig)
00049   : mInputCollection (iConfig.getParameter<edm::InputTag>( "src" )),
00050     mInputGenCollection (iConfig.getParameter<edm::InputTag>( "srcGen" )),
00051     mOutputFile (iConfig.getUntrackedParameter<std::string>("outputFile", "")),
00052     mMatchGenPtThreshold (iConfig.getParameter<double>("genPtThreshold")),
00053     mGenEnergyFractionThreshold (iConfig.getParameter<double>("genEnergyFractionThreshold")),
00054     mReverseEnergyFractionThreshold (iConfig.getParameter<double>("reverseEnergyFractionThreshold")),
00055     mRThreshold (iConfig.getParameter<double>("RThreshold")),
00056     mTurnOnEverything (iConfig.getUntrackedParameter<std::string>("TurnOnEverything",""))
00057 {
00058     numberofevents
00059     = mEta = mEtaFineBin = mPhi = mPhiFineBin = mE = mE_80 = mE_3000
00060     = mP = mP_80 = mP_3000 = mPt = mPt_80 = mPt_3000
00061     = mMass = mMass_80 = mMass_3000 = mConstituents = mConstituents_80
00062     = mEtaFirst = mPhiFirst = mEFirst = mEFirst_80 = mEFirst_3000 = mPtFirst = mPtFirst_80 = mPtFirst_3000
00063     = mMjj = mMjj_3000 = mDelEta = mDelPhi = mDelPt 
00064     = mMaxEInEmTowers = mMaxEInHadTowers 
00065     = mHadEnergyInHO = mHadEnergyInHB = mHadEnergyInHF = mHadEnergyInHE 
00066     = mHadEnergyInHO_80 = mHadEnergyInHB_80 = mHadEnergyInHE_80 
00067     = mHadEnergyInHO_3000 = mHadEnergyInHB_3000 = mHadEnergyInHE_3000 
00068     = mEmEnergyInEB = mEmEnergyInEE = mEmEnergyInHF 
00069     = mEmEnergyInEB_80 = mEmEnergyInEE_80
00070     = mEmEnergyInEB_3000 = mEmEnergyInEE_3000
00071     = mEnergyFractionHadronic = mEnergyFractionEm 
00072     = mHFLong = mHFTotal = mHFLong_80 = mHFLong_3000 = mHFShort = mHFShort_80 = mHFShort_3000
00073     = mN90
00074     = mCaloMEx = mCaloMEx_3000 = mCaloMEy = mCaloMEy_3000 = mCaloMETSig = mCaloMETSig_3000
00075     = mCaloMET = mCaloMET_3000 =  mCaloMETPhi = mCaloSumET  = mCaloSumET_3000   
00076     = mHadTiming = mEmTiming 
00077     = mNJetsEtaC = mNJetsEtaF = mNJets1 = mNJets2
00078     = mAllGenJetsPt = mMatchedGenJetsPt = mAllGenJetsEta = mMatchedGenJetsEta 
00079     = mGenJetMatchEnergyFraction = mReverseMatchEnergyFraction = mRMatch
00080     = mDeltaEta = mDeltaPhi = mEScale = mlinEScale = mDeltaE
00081     = mHadEnergyProfile = mEmEnergyProfile = mJetEnergyProfile = mHadJetEnergyProfile = mEMJetEnergyProfile
00082     = mEScale_pt10 = mEScaleFineBin
00083     = mpTScaleB_s = mpTScaleE_s = mpTScaleF_s 
00084     = mpTScaleB_d = mpTScaleE_d = mpTScaleF_d
00085     = mpTScale_60_120_s = mpTScale_200_300_s = mpTScale_600_900_s = mpTScale_2700_3500_s
00086     = mpTScale_60_120_d = mpTScale_200_300_d = mpTScale_600_900_d = mpTScale_2700_3500_d
00087     = mpTScale1DB_60_120    = mpTScale1DE_60_120    = mpTScale1DF_60_120 
00088     = mpTScale1DB_200_300   = mpTScale1DE_200_300   = mpTScale1DF_200_300 
00089     = mpTScale1DB_600_900   = mpTScale1DE_600_900   = mpTScale1DF_600_900 
00090     = mpTScale1DB_2700_3500 = mpTScale1DE_2700_3500 = mpTScale1DF_2700_3500
00091     = mpTScale1D_60_120 = mpTScale1D_200_300 = mpTScale1D_600_900 = mpTScale1D_2700_3500
00092     = mHBEne = mHBTime = mHEEne = mHETime = mHFEne = mHFTime = mHOEne = mHOTime
00093     = mEBEne = mEBTime = mEEEne = mEETime 
00094     = mPthat_80 = mPthat_3000
00095     = 0;
00096   
00097   DQMStore* dbe = &*edm::Service<DQMStore>();
00098   if (dbe) {
00099     dbe->setCurrentFolder("JetMET/RecoJetsV/CaloJetTask_" + mInputCollection.label());
00100     //
00101     numberofevents    = dbe->book1D("numberofevents","numberofevents", 3, 0 , 2);
00102     //
00103     mEta              = dbe->book1D("Eta", "Eta", 100, -5, 5); 
00104     mEtaFineBin       = dbe->book1D("EtaFineBin_Pt10", "EtaFineBin_Pt10", 500, -5, 5); 
00105     mEtaFineBin1p     = dbe->book1D("EtaFineBin1p_Pt10", "EtaFineBin1p_Pt10", 100, 0, 1.3); 
00106     mEtaFineBin2p     = dbe->book1D("EtaFineBin2p_Pt10", "EtaFineBin2p_Pt10", 100, 1.3, 3); 
00107     mEtaFineBin3p     = dbe->book1D("EtaFineBin3p_Pt10", "EtaFineBin3p_Pt10", 100, 3, 5); 
00108     mEtaFineBin1m     = dbe->book1D("EtaFineBin1m_Pt10", "EtaFineBin1m_Pt10", 100, -1.3, 0); 
00109     mEtaFineBin2m     = dbe->book1D("EtaFineBin2m_Pt10", "EtaFineBin2m_Pt10", 100, -3, -1.3); 
00110     mEtaFineBin3m     = dbe->book1D("EtaFineBin3m_Pt10", "EtaFineBin3m_Pt10", 100, -5, -3); 
00111     //
00112     mPhi              = dbe->book1D("Phi", "Phi", 70, -3.5, 3.5); 
00113     mPhiFineBin       = dbe->book1D("PhiFineBin_Pt10", "PhiFineBin_Pt10", 350, -3.5, 3.5); 
00114     //
00115     mE                = dbe->book1D("E", "E", 100, 0, 500); 
00116     mE_80             = dbe->book1D("E_80", "E_80", 100, 0, 4500); 
00117     mE_3000           = dbe->book1D("E_3000", "E_3000", 100, 0, 6000); 
00118     //
00119     mP                = dbe->book1D("P", "P", 100, 0, 500); 
00120     mP_80             = dbe->book1D("P_80", "P_80", 100, 0, 4500); 
00121     mP_3000           = dbe->book1D("P_3000", "P_3000", 100, 0, 6000); 
00122     //
00123     mPt               = dbe->book1D("Pt", "Pt", 100, 0, 50); 
00124     mPt_80            = dbe->book1D("Pt_80", "Pt_80", 100, 0, 140); 
00125     mPt_3000          = dbe->book1D("Pt_3000", "Pt_3000", 100, 0, 4000); 
00126     //
00127     mMass             = dbe->book1D("Mass", "Mass", 100, 0, 25); 
00128     mMass_80          = dbe->book1D("Mass_80", "Mass_80", 100, 0, 120); 
00129     mMass_3000        = dbe->book1D("Mass_3000", "Mass_3000", 100, 0, 500); 
00130     //
00131     mConstituents     = dbe->book1D("Constituents", "# of Constituents", 100, 0, 100); 
00132     mConstituents_80  = dbe->book1D("Constituents_80", "# of Constituents_80", 40, 0, 40); 
00133     //
00134     mEtaFirst         = dbe->book1D("EtaFirst", "EtaFirst", 100, -5, 5); 
00135     mPhiFirst         = dbe->book1D("PhiFirst", "PhiFirst", 70, -3.5, 3.5);     
00136     mEFirst           = dbe->book1D("EFirst", "EFirst", 100, 0, 1000); 
00137     mEFirst_80        = dbe->book1D("EFirst_80", "EFirst_80", 100, 0, 180); 
00138     mEFirst_3000      = dbe->book1D("EFirst_3000", "EFirst_3000", 100, 0, 4000); 
00139     mPtFirst          = dbe->book1D("PtFirst", "PtFirst", 100, 0, 50); 
00140     mPtFirst_80       = dbe->book1D("PtFirst_80", "PtFirst_80", 100, 0, 140);
00141     mPtFirst_3000     = dbe->book1D("PtFirst_3000", "PtFirst_3000", 100, 0, 4000);
00142     //
00143     mMjj              = dbe->book1D("Mjj", "Mjj", 100, 0, 2000); 
00144     mMjj_3000         = dbe->book1D("Mjj_3000", "Mjj_3000", 100, 0, 10000); 
00145     mDelEta           = dbe->book1D("DelEta", "DelEta", 100, -.5, .5); 
00146     mDelPhi           = dbe->book1D("DelPhi", "DelPhi", 100, -.5, .5); 
00147     mDelPt            = dbe->book1D("DelPt", "DelPt", 100, -1, 1); 
00148     //
00149     mMaxEInEmTowers   = dbe->book1D("MaxEInEmTowers", "MaxEInEmTowers", 100, 0, 100); 
00150     mMaxEInHadTowers  = dbe->book1D("MaxEInHadTowers", "MaxEInHadTowers", 100, 0, 100); 
00151     mHadEnergyInHO    = dbe->book1D("HadEnergyInHO", "HadEnergyInHO", 100, 0, 10); 
00152     mHadEnergyInHB    = dbe->book1D("HadEnergyInHB", "HadEnergyInHB", 100, 0, 50); 
00153     mHadEnergyInHF    = dbe->book1D("HadEnergyInHF", "HadEnergyInHF", 100, 0, 50); 
00154     mHadEnergyInHE    = dbe->book1D("HadEnergyInHE", "HadEnergyInHE", 100, 0, 100); 
00155     //
00156     mHadEnergyInHO_80    = dbe->book1D("HadEnergyInHO_80", "HadEnergyInHO_80", 100, 0, 50); 
00157     mHadEnergyInHB_80    = dbe->book1D("HadEnergyInHB_80", "HadEnergyInHB_80", 100, 0, 200); 
00158     mHadEnergyInHE_80    = dbe->book1D("HadEnergyInHE_80", "HadEnergyInHE_80", 100, 0, 1000); 
00159     mHadEnergyInHO_3000  = dbe->book1D("HadEnergyInHO_3000", "HadEnergyInHO_3000", 100, 0, 500); 
00160     mHadEnergyInHB_3000  = dbe->book1D("HadEnergyInHB_3000", "HadEnergyInHB_3000", 100, 0, 3000); 
00161     mHadEnergyInHE_3000  = dbe->book1D("HadEnergyInHE_3000", "HadEnergyInHE_3000", 100, 0, 2000); 
00162     //
00163     mEmEnergyInEB     = dbe->book1D("EmEnergyInEB", "EmEnergyInEB", 100, 0, 50); 
00164     mEmEnergyInEE     = dbe->book1D("EmEnergyInEE", "EmEnergyInEE", 100, 0, 50); 
00165     mEmEnergyInHF     = dbe->book1D("EmEnergyInHF", "EmEnergyInHF", 120, -20, 100); 
00166     mEmEnergyInEB_80  = dbe->book1D("EmEnergyInEB_80", "EmEnergyInEB_80", 100, 0, 200); 
00167     mEmEnergyInEE_80  = dbe->book1D("EmEnergyInEE_80", "EmEnergyInEE_80", 100, 0, 1000); 
00168     mEmEnergyInEB_3000= dbe->book1D("EmEnergyInEB_3000", "EmEnergyInEB_3000", 100, 0, 3000); 
00169     mEmEnergyInEE_3000= dbe->book1D("EmEnergyInEE_3000", "EmEnergyInEE_3000", 100, 0, 2000); 
00170     mEnergyFractionHadronic = dbe->book1D("EnergyFractionHadronic", "EnergyFractionHadronic", 120, -0.1, 1.1); 
00171     mEnergyFractionEm = dbe->book1D("EnergyFractionEm", "EnergyFractionEm", 120, -0.1, 1.1); 
00172     //
00173     mHFTotal          = dbe->book1D("HFTotal", "HFTotal", 100, 0, 500);
00174     mHFTotal_80       = dbe->book1D("HFTotal_80", "HFTotal_80", 100, 0, 3000);
00175     mHFTotal_3000     = dbe->book1D("HFTotal_3000", "HFTotal_3000", 100, 0, 6000);
00176     mHFLong           = dbe->book1D("HFLong", "HFLong", 100, 0, 500);
00177     mHFLong_80        = dbe->book1D("HFLong_80", "HFLong_80", 100, 0, 200);
00178     mHFLong_3000      = dbe->book1D("HFLong_3000", "HFLong_3000", 100, 0, 1500);
00179     mHFShort          = dbe->book1D("HFShort", "HFShort", 100, 0, 500);
00180     mHFShort_80       = dbe->book1D("HFShort_80", "HFShort_80", 100, 0, 200);
00181     mHFShort_3000     = dbe->book1D("HFShort_3000", "HFShort_3000", 100, 0, 1500);
00182     //
00183     mN90              = dbe->book1D("N90", "N90", 50, 0, 50); 
00184     //
00185     mGenEta           = dbe->book1D("GenEta", "GenEta", 100, -5, 5);
00186     mGenPhi           = dbe->book1D("GenPhi", "GenPhi", 70, -3.5, 3.5);
00187     mGenPt            = dbe->book1D("GenPt", "GenPt", 100, 0, 50);
00188     mGenPt_80         = dbe->book1D("GenPt_80", "GenPt_80", 100, 0, 140);
00189     mGenPt_3000       = dbe->book1D("GenPt_3000", "GenPt_3000", 100, 0, 1500);
00190     //
00191     mGenEtaFirst      = dbe->book1D("GenEtaFirst", "GenEtaFirst", 100, -5, 5);
00192     mGenPhiFirst      = dbe->book1D("GenPhiFirst", "GenPhiFirst", 70, -3.5, 3.5);
00193     //
00194     mCaloMEx          = dbe->book1D("CaloMEx","CaloMEx",200,-150,150);
00195     mCaloMEx_3000     = dbe->book1D("CaloMEx_3000","CaloMEx_3000",100,-500,500);
00196     mCaloMEy          = dbe->book1D("CaloMEy","CaloMEy",200,-150,150);
00197     mCaloMEy_3000     = dbe->book1D("CaloMEy_3000","CaloMEy_3000",100,-500,500);
00198     mCaloMETSig       = dbe->book1D("CaloMETSig","CaloMETSig",100,0,15);
00199     mCaloMETSig_3000  = dbe->book1D("CaloMETSig_3000","CaloMETSig_3000",100,0,50);
00200     mCaloMET          = dbe->book1D("CaloMET","CaloMET",100,0,150);
00201     mCaloMET_3000     = dbe->book1D("CaloMET_3000","CaloMET_3000",100,0,1000);
00202     mCaloMETPhi       = dbe->book1D("CaloMETPhi","CaloMETPhi",70, -3.5, 3.5);
00203     mCaloSumET        = dbe->book1D("CaloSumET","CaloSumET",100,0,500);
00204     mCaloSumET_3000   = dbe->book1D("CaloSumET_3000","CaloSumET_3000",100,3000,8000);
00205     //
00206     mHadTiming        = dbe->book1D("HadTiming", "HadTiming", 75, -50, 100);
00207     mEmTiming         = dbe->book1D("EMTiming", "EMTiming", 75, -50, 100);
00208     //
00209     mNJetsEtaC        = dbe->book1D("NJetsEtaC_Pt10", "NJetsEtaC_Pt10", 15, 0, 15);
00210     mNJetsEtaF        = dbe->book1D("NJetsEtaF_Pt10", "NJetsEtaF_Pt10", 15, 0, 15);
00211     //
00212     mNJets1           = dbe->bookProfile("NJets1", "NJets1", 100, 0, 200,  100, 0, 50, "s");
00213     mNJets2           = dbe->bookProfile("NJets2", "NJets2", 100, 0, 4000, 100, 0, 50, "s");
00214     //
00215     mHBEne     = dbe->book1D( "HBEne",  "HBEne", 1000, -20, 100 );
00216     mHBTime    = dbe->book1D( "HBTime", "HBTime", 200, -200, 200 );
00217     mHEEne     = dbe->book1D( "HEEne",  "HEEne", 1000, -20, 100 );
00218     mHETime    = dbe->book1D( "HETime", "HETime", 200, -200, 200 );
00219     mHOEne     = dbe->book1D( "HOEne",  "HOEne", 1000, -20, 100 );
00220     mHOTime    = dbe->book1D( "HOTime", "HOTime", 200, -200, 200 );
00221     mHFEne     = dbe->book1D( "HFEne",  "HFEne", 1000, -20, 100 );
00222     mHFTime    = dbe->book1D( "HFTime", "HFTime", 200, -200, 200 );
00223     mEBEne     = dbe->book1D( "EBEne",  "EBEne", 1000, -20, 100 );
00224     mEBTime    = dbe->book1D( "EBTime", "EBTime", 200, -200, 200 );
00225     mEEEne     = dbe->book1D( "EEEne",  "EEEne", 1000, -20, 100 );
00226     mEETime    = dbe->book1D( "EETime", "EETime", 200, -200, 200 );
00227     //
00228     mPthat_80            = dbe->book1D("Pthat_80", "Pthat_80", 100, 40.0, 160.0); 
00229     mPthat_3000          = dbe->book1D("Pthat_3000", "Pthat_3000", 100, 2500.0, 4000.0); 
00230     //
00231     double log10PtMin = 0.5; //=3.1622766
00232     double log10PtMax = 3.75; //=5623.41325
00233     int log10PtBins = 26; 
00234     double etaMin = -5.;
00235     double etaMax = 5.;
00236     int etaBins = 50;
00237 
00238     double linPtMin = 5;
00239     double linPtMax = 155;
00240     int linPtBins = 15;
00241 
00242     int log10PtFineBins = 50;
00243 
00244     mAllGenJetsPt = dbe->book1D("GenJetLOGpT", "GenJet LOG(pT_gen)", 
00245                                 log10PtBins, log10PtMin, log10PtMax);
00246     mMatchedGenJetsPt = dbe->book1D("MatchedGenJetLOGpT", "MatchedGenJet LOG(pT_gen)", 
00247                                     log10PtBins, log10PtMin, log10PtMax);
00248     mAllGenJetsEta = dbe->book2D("GenJetEta", "GenJet Eta vs LOG(pT_gen)", 
00249                                  log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax);
00250     mMatchedGenJetsEta = dbe->book2D("MatchedGenJetEta", "MatchedGenJet Eta vs LOG(pT_gen)", 
00251                                      log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax);
00252     //
00253     if (mTurnOnEverything.compare("yes")==0) {
00254       mHadEnergyProfile = dbe->bookProfile2D("HadEnergyProfile", "HadEnergyProfile", 82, -41, 41, 73, 0, 73, 100, 0, 10000, "s");
00255       mEmEnergyProfile  = dbe->bookProfile2D("EmEnergyProfile", "EmEnergyProfile", 82, -41, 41, 73, 0, 73, 100, 0, 10000, "s");
00256     }
00257     mJetEnergyProfile = dbe->bookProfile2D("JetEnergyProfile", "JetEnergyProfile", 50, -5, 5, 36, -3.1415987, 3.1415987, 100, 0, 10000, "s");
00258     mHadJetEnergyProfile = dbe->bookProfile2D("HadJetEnergyProfile", "HadJetEnergyProfile", 50, -5, 5, 36, -3.1415987, 3.1415987, 100, 0, 10000, "s");
00259     mEMJetEnergyProfile = dbe->bookProfile2D("EMJetEnergyProfile", "EMJetEnergyProfile", 50, -5, 5, 36, -3.1415987, 3.1415987, 100, 0, 10000, "s");
00260     //
00261     if (mTurnOnEverything.compare("yes")==0) {
00262     mGenJetMatchEnergyFraction  = dbe->book3D("GenJetMatchEnergyFraction", "GenJetMatchEnergyFraction vs LOG(pT_gen) vs eta", 
00263                                               log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 101, 0, 1.01);
00264     mReverseMatchEnergyFraction  = dbe->book3D("ReverseMatchEnergyFraction", "ReverseMatchEnergyFraction vs LOG(pT_gen) vs eta", 
00265                                                log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 101, 0, 1.01);
00266     mRMatch  = dbe->book3D("RMatch", "delta(R)(Gen-Calo) vs LOG(pT_gen) vs eta", 
00267                            log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 60, 0, 3);
00268     mDeltaEta = dbe->book3D("DeltaEta", "DeltaEta vs LOG(pT_gen) vs eta", 
00269                               log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 100, -1, 1);
00270     mDeltaPhi = dbe->book3D("DeltaPhi", "DeltaPhi vs LOG(pT_gen) vs eta", 
00271                               log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 100, -1, 1);
00272     mEScale = dbe->book3D("EScale", "EnergyScale vs LOG(pT_gen) vs eta", 
00273                             log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 100, 0, 2);
00274     mlinEScale = dbe->book3D("linEScale", "EnergyScale vs LOG(pT_gen) vs eta", 
00275                             linPtBins, linPtMin, linPtMax, etaBins, etaMin, etaMax, 100, 0, 2);
00276     mDeltaE = dbe->book3D("DeltaE", "DeltaE vs LOG(pT_gen) vs eta", 
00277                             log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 2000, -200, 200);
00278     //
00279     mEScale_pt10 = dbe->book3D("EScale_pt10", "EnergyScale vs LOG(pT_gen) vs eta", 
00280                             log10PtBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 100, 0, 2);
00281     mEScaleFineBin = dbe->book3D("EScaleFineBins", "EnergyScale vs LOG(pT_gen) vs eta", 
00282                             log10PtFineBins, log10PtMin, log10PtMax, etaBins, etaMin, etaMax, 100, 0, 2);
00283     }
00284 
00285     mpTScaleB_s = dbe->bookProfile("pTScaleB_s", "pTScale_s_0<|eta|<1.3",
00286                                     log10PtBins, log10PtMin, log10PtMax, 0, 2, "s");
00287     mpTScaleE_s = dbe->bookProfile("pTScaleE_s", "pTScale_s_1.3<|eta|<3.0",
00288                                     log10PtBins, log10PtMin, log10PtMax, 0, 2, "s");
00289     mpTScaleF_s = dbe->bookProfile("pTScaleF_s", "pTScale_s_3.0<|eta|<5.0",
00290                                     log10PtBins, log10PtMin, log10PtMax, 0, 2, "s");
00291     mpTScaleB_d = dbe->bookProfile("pTScaleB_d", "pTScale_d_0<|eta|<1.3",
00292                                    log10PtBins, log10PtMin, log10PtMax, 0, 2, " ");
00293     mpTScaleE_d = dbe->bookProfile("pTScaleE_d", "pTScale_d_1.3<|eta|<3.0",
00294                                    log10PtBins, log10PtMin, log10PtMax, 0, 2, " ");
00295     mpTScaleF_d = dbe->bookProfile("pTScaleF_d", "pTScale_d_3.0<|eta|<5.0",
00296                                    log10PtBins, log10PtMin, log10PtMax, 0, 2, " ");
00297 
00298     mpTScale_60_120_s    = dbe->bookProfile("pTScale_60_120_s", "pTScale_s_60<pT<120",
00299                                           etaBins, etaMin, etaMax, 0., 2., "s");
00300     mpTScale_200_300_s   = dbe->bookProfile("pTScale_200_300_s", "pTScale_s_200<pT<300",
00301                                           etaBins, etaMin, etaMax, 0., 2., "s");
00302     mpTScale_600_900_s   = dbe->bookProfile("pTScale_600_900_s", "pTScale_s_600<pT<900",
00303                                           etaBins, etaMin, etaMax, 0., 2., "s");
00304     mpTScale_2700_3500_s = dbe->bookProfile("pTScale_2700_3500_s", "pTScale_s_2700<pt<3500",
00305                                           etaBins, etaMin, etaMax, 0., 2., "s");
00306     mpTScale_60_120_d    = dbe->bookProfile("pTScale_60_120_d", "pTScale_d_60<pT<120",
00307                                           etaBins, etaMin, etaMax, 0., 2., " ");
00308     mpTScale_200_300_d   = dbe->bookProfile("pTScale_200_300_d", "pTScale_d_200<pT<300",
00309                                           etaBins, etaMin, etaMax, 0., 2., " ");
00310     mpTScale_600_900_d   = dbe->bookProfile("pTScale_600_900_d", "pTScale_d_600<pT<900",
00311                                           etaBins, etaMin, etaMax, 0., 2., " ");
00312     mpTScale_2700_3500_d = dbe->bookProfile("pTScale_2700_3500_d", "pTScale_d_2700<pt<3500",
00313                                           etaBins, etaMin, etaMax, 0., 2., " ");
00314 
00315     mpTScale1DB_60_120 = dbe->book1D("pTScale1DB_60_120", "pTScale_distribution_for_0<|eta|<1.3_60_120",
00316                                    100, 0, 2);
00317     mpTScale1DE_60_120 = dbe->book1D("pTScale1DE_60_120", "pTScale_distribution_for_1.3<|eta|<3.0_60_120",
00318                                    100, 0, 2);
00319     mpTScale1DF_60_120 = dbe->book1D("pTScale1DF_60_120", "pTScale_distribution_for_3.0<|eta|<5.0_60_120",
00320                                    100, 0, 2);
00321 
00322     mpTScale1DB_200_300 = dbe->book1D("pTScale1DB_200_300", "pTScale_distribution_for_0<|eta|<1.3_200_300",
00323                                    100, 0, 2);
00324     mpTScale1DE_200_300 = dbe->book1D("pTScale1DE_200_300", "pTScale_distribution_for_1.3<|eta|<3.0_200_300",
00325                                    100, 0, 2);
00326     mpTScale1DF_200_300 = dbe->book1D("pTScale1DF_200_300", "pTScale_distribution_for_3.0<|eta|<5.0_200_300",
00327                                    100, 0, 2);
00328 
00329     mpTScale1DB_600_900 = dbe->book1D("pTScale1DB_600_900", "pTScale_distribution_for_0<|eta|<1.3_600_900",
00330                                    100, 0, 2);
00331     mpTScale1DE_600_900 = dbe->book1D("pTScale1DE_600_900", "pTScale_distribution_for_1.3<|eta|<3.0_600_900",
00332                                    100, 0, 2);
00333     mpTScale1DF_600_900 = dbe->book1D("pTScale1DF_600_900", "pTScale_distribution_for_3.0<|eta|<5.0_600_900",
00334                                    100, 0, 2);
00335 
00336     mpTScale1DB_2700_3500 = dbe->book1D("pTScale1DB_2700_3500", "pTScale_distribution_for_0<|eta|<1.3_2700_3500",
00337                                    100, 0, 2);
00338     mpTScale1DE_2700_3500 = dbe->book1D("pTScale1DE_2700_3500", "pTScale_distribution_for_1.3<|eta|<3.0_2700_3500",
00339                                    100, 0, 2);
00340     mpTScale1DF_2700_3500 = dbe->book1D("pTScale1DF_2700_3500", "pTScale_distribution_for_3.0<|eta|<5.0_2700_3500",
00341                                    100, 0, 2);
00342 
00343     mpTScale1D_60_120    = dbe->book1D("pTScale1D_60_120", "pTScale_distribution_for_60<pT<120",
00344                                             100, 0, 2);
00345     mpTScale1D_200_300    = dbe->book1D("pTScale1D_200_300", "pTScale_distribution_for_200<pT<300",
00346                                             100, 0, 2);
00347     mpTScale1D_600_900    = dbe->book1D("pTScale1D_600_900", "pTScale_distribution_for_600<pT<900",
00348                                             100, 0, 2);
00349     mpTScale1D_2700_3500 = dbe->book1D("pTScale1D_2700_3500", "pTScale_distribution_for_2700<pt<3500",
00350                                             100, 0, 2);
00351 
00352   }
00353 
00354   if (mOutputFile.empty ()) {
00355     LogInfo("OutputInfo") << " CaloJet histograms will NOT be saved";
00356   } 
00357   else {
00358     LogInfo("OutputInfo") << " CaloJethistograms will be saved to file:" << mOutputFile;
00359   }
00360 }
00361    
00362 CaloJetTester::~CaloJetTester()
00363 {
00364 }
00365 
00366 void CaloJetTester::beginJob(){
00367 }
00368 
00369 void CaloJetTester::endJob() {
00370  if (!mOutputFile.empty() && &*edm::Service<DQMStore>()) edm::Service<DQMStore>()->save (mOutputFile);
00371 }
00372 
00373 
00374 void CaloJetTester::analyze(const edm::Event& mEvent, const edm::EventSetup& mSetup)
00375 {
00376   double countsfornumberofevents = 1;
00377   numberofevents->Fill(countsfornumberofevents);
00378   // *********************************
00379   // *** Get pThat
00380   // *********************************
00381 
00382   edm::Handle<HepMCProduct> evt;
00383   mEvent.getByLabel("generator", evt);
00384   if (evt.isValid()) {
00385   HepMC::GenEvent * myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
00386   
00387   double pthat = myGenEvent->event_scale();
00388 
00389   mPthat_80->Fill(pthat);
00390   mPthat_3000->Fill(pthat);
00391 
00392   delete myGenEvent; 
00393   }
00394 
00395   // ***********************************
00396   // *** Get CaloMET
00397   // ***********************************
00398 
00399   const CaloMET *calomet;
00400   edm::Handle<CaloMETCollection> calo;
00401   mEvent.getByLabel("met", calo);
00402   if (!calo.isValid()) {
00403     edm::LogInfo("OutputInfo") << " failed to retrieve data required by MET Task";
00404     edm::LogInfo("OutputInfo") << " MET Task cannot continue...!";
00405   } else {
00406     const CaloMETCollection *calometcol = calo.product();
00407     calomet = &(calometcol->front());
00408     
00409     double caloSumET = calomet->sumEt();
00410     double caloMETSig = calomet->mEtSig();
00411     double caloMET = calomet->pt();
00412     double caloMEx = calomet->px();
00413     double caloMEy = calomet->py();
00414     double caloMETPhi = calomet->phi();
00415 
00416     mCaloMEx->Fill(caloMEx);
00417     mCaloMEx_3000->Fill(caloMEx);
00418     mCaloMEy->Fill(caloMEy);
00419     mCaloMEy_3000->Fill(caloMEy);
00420     mCaloMET->Fill(caloMET);
00421     mCaloMET_3000->Fill(caloMET);
00422     mCaloMETPhi->Fill(caloMETPhi);
00423     mCaloSumET->Fill(caloSumET);
00424     mCaloSumET_3000->Fill(caloSumET);
00425     mCaloMETSig->Fill(caloMETSig);
00426     mCaloMETSig_3000->Fill(caloMETSig);
00427   }
00428 
00429   // ***********************************
00430   // *** Get the CaloTower collection
00431   // ***********************************
00432   Handle<CaloTowerCollection> caloTowers;
00433   mEvent.getByLabel( "towerMaker", caloTowers );
00434   if (caloTowers.isValid()) {
00435   for( CaloTowerCollection::const_iterator cal = caloTowers->begin(); cal != caloTowers->end(); ++ cal ){
00436 
00437     //To compensate for the index
00438     if (mTurnOnEverything.compare("yes")==0) {
00439       if (cal->ieta() >> 0 ){mHadEnergyProfile->Fill (cal->ieta()-1, cal->iphi(), cal->hadEnergy());
00440       mEmEnergyProfile->Fill (cal->ieta()-1, cal->iphi(), cal->emEnergy());}
00441       mHadEnergyProfile->Fill (cal->ieta(), cal->iphi(), cal->hadEnergy());
00442       mEmEnergyProfile->Fill (cal->ieta(), cal->iphi(), cal->emEnergy());
00443     }
00444 
00445     mHadTiming->Fill (cal->hcalTime());
00446     mEmTiming->Fill (cal->ecalTime());    
00447   }
00448   }
00449   
00450   // ***********************************
00451   // *** Get the RecHits collection
00452   // ***********************************
00453   try {
00454     std::vector<edm::Handle<HBHERecHitCollection> > colls;
00455     mEvent.getManyByType(colls);
00456     std::vector<edm::Handle<HBHERecHitCollection> >::iterator i;
00457     for (i=colls.begin(); i!=colls.end(); i++) {
00458       for (HBHERecHitCollection::const_iterator j=(*i)->begin(); j!=(*i)->end(); j++) {
00459         //      std::cout << *j << std::endl;
00460         if (j->id().subdet() == HcalBarrel) {
00461           mHBEne->Fill(j->energy()); 
00462           mHBTime->Fill(j->time()); 
00463         }
00464         if (j->id().subdet() == HcalEndcap) {
00465           mHEEne->Fill(j->energy()); 
00466           mHETime->Fill(j->time()); 
00467         }
00468 
00469       }
00470     }
00471   } catch (...) {
00472     edm::LogInfo("OutputInfo") << " No HB/HE RecHits.";
00473   }
00474 
00475   try {
00476     std::vector<edm::Handle<HFRecHitCollection> > colls;
00477     mEvent.getManyByType(colls);
00478     std::vector<edm::Handle<HFRecHitCollection> >::iterator i;
00479     for (i=colls.begin(); i!=colls.end(); i++) {
00480       for (HFRecHitCollection::const_iterator j=(*i)->begin(); j!=(*i)->end(); j++) {
00481         //      std::cout << *j << std::endl;
00482         if (j->id().subdet() == HcalForward) {
00483           mHFEne->Fill(j->energy()); 
00484           mHFTime->Fill(j->time()); 
00485         }
00486       }
00487     }
00488   } catch (...) {
00489     edm::LogInfo("OutputInfo") << " No HF RecHits.";
00490   }
00491 
00492   try {
00493     std::vector<edm::Handle<HORecHitCollection> > colls;
00494     mEvent.getManyByType(colls);
00495     std::vector<edm::Handle<HORecHitCollection> >::iterator i;
00496     for (i=colls.begin(); i!=colls.end(); i++) {
00497       for (HORecHitCollection::const_iterator j=(*i)->begin(); j!=(*i)->end(); j++) {
00498         if (j->id().subdet() == HcalOuter) {
00499           mHOEne->Fill(j->energy()); 
00500           mHOTime->Fill(j->time()); 
00501         }
00502       }
00503     }
00504   } catch (...) {
00505     edm::LogInfo("OutputInfo") << " No HO RecHits.";
00506   }
00507   try {
00508     std::vector<edm::Handle<EBRecHitCollection> > colls;
00509     mEvent.getManyByType(colls);
00510     std::vector<edm::Handle<EBRecHitCollection> >::iterator i;
00511     for (i=colls.begin(); i!=colls.end(); i++) {
00512       for (EBRecHitCollection::const_iterator j=(*i)->begin(); j!=(*i)->end(); j++) {
00513         //      if (j->id() == EcalBarrel) {
00514         mEBEne->Fill(j->energy()); 
00515         mEBTime->Fill(j->time()); 
00516         //    }
00517         //      std::cout << *j << std::endl;
00518         //      std::cout << j->id() << std::endl;
00519       }
00520     }
00521   } catch (...) {
00522     edm::LogInfo("OutputInfo") << " No EB RecHits.";
00523   }
00524 
00525   try {
00526     std::vector<edm::Handle<EERecHitCollection> > colls;
00527     mEvent.getManyByType(colls);
00528     std::vector<edm::Handle<EERecHitCollection> >::iterator i;
00529     for (i=colls.begin(); i!=colls.end(); i++) {
00530       for (EERecHitCollection::const_iterator j=(*i)->begin(); j!=(*i)->end(); j++) {
00531         //      if (j->id().subdet() == EcalEndcap) {
00532         mEEEne->Fill(j->energy()); 
00533         mEETime->Fill(j->time()); 
00534         //    }
00535         //      std::cout << *j << std::endl;
00536       }
00537     }
00538   } catch (...) {
00539     edm::LogInfo("OutputInfo") << " No EE RecHits.";
00540   }
00541 
00542 
00543   //***********************************
00544   //*** Get the Jet collection
00545   //***********************************
00546   math::XYZTLorentzVector p4tmp[2];
00547   Handle<CaloJetCollection> caloJets;
00548   mEvent.getByLabel(mInputCollection, caloJets);
00549   if (!caloJets.isValid()) return;
00550   CaloJetCollection::const_iterator jet = caloJets->begin ();
00551   int jetIndex = 0;
00552   int nJet = 0;
00553   int nJetF = 0;
00554   int nJetC = 0;
00555   for (; jet != caloJets->end (); jet++, jetIndex++) {
00556     if (mEta) mEta->Fill (jet->eta());
00557 
00558     if (jet->pt() > 10.) {
00559       if (fabs(jet->eta()) > 1.3) 
00560         nJetF++;
00561       else 
00562         nJetC++;          
00563     }
00564     if (jet->pt() > 10.) {
00565       if (mEtaFineBin) mEtaFineBin->Fill (jet->eta());
00566       if (mEtaFineBin1p) mEtaFineBin1p->Fill (jet->eta());
00567       if (mEtaFineBin2p) mEtaFineBin2p->Fill (jet->eta());
00568       if (mEtaFineBin3p) mEtaFineBin3p->Fill (jet->eta());
00569       if (mEtaFineBin1m) mEtaFineBin1m->Fill (jet->eta());
00570       if (mEtaFineBin2m) mEtaFineBin2m->Fill (jet->eta());
00571       if (mEtaFineBin3m) mEtaFineBin3m->Fill (jet->eta());
00572       if (mPhiFineBin) mPhiFineBin->Fill (jet->phi());
00573     }
00574     if (mPhi) mPhi->Fill (jet->phi());
00575     if (mE) mE->Fill (jet->energy());
00576     if (mE_80) mE_80->Fill (jet->energy());
00577     if (mE_3000) mE_3000->Fill (jet->energy());
00578     if (mP) mP->Fill (jet->p());
00579     if (mP_80) mP_80->Fill (jet->p());
00580     if (mP_3000) mP_3000->Fill (jet->p());
00581     if (mPt) mPt->Fill (jet->pt());
00582     if (mPt_80) mPt_80->Fill (jet->pt());
00583     if (mPt_3000) mPt_3000->Fill (jet->pt());
00584     if (mMass) mMass->Fill (jet->mass());
00585     if (mMass_80) mMass_80->Fill (jet->mass());
00586     if (mMass_3000) mMass_3000->Fill (jet->mass());
00587     if (mConstituents) mConstituents->Fill (jet->nConstituents());
00588     if (mConstituents_80) mConstituents_80->Fill (jet->nConstituents());
00589     if (jet == caloJets->begin ()) { // first jet
00590       if (mEtaFirst) mEtaFirst->Fill (jet->eta());
00591       if (mPhiFirst) mPhiFirst->Fill (jet->phi());
00592       if (mEFirst) mEFirst->Fill (jet->energy());
00593       if (mEFirst_80) mEFirst_80->Fill (jet->energy());
00594       if (mEFirst_3000) mEFirst_3000->Fill (jet->energy());
00595       if (mPtFirst) mPtFirst->Fill (jet->pt());
00596       if (mPtFirst_80) mPtFirst_80->Fill (jet->pt());
00597       if (mPtFirst_3000) mPtFirst_3000->Fill (jet->pt());
00598     }
00599     if (jetIndex == 0) {
00600       nJet++;
00601       p4tmp[0] = jet->p4();     
00602     }
00603     if (jetIndex == 1) {
00604       nJet++;
00605       p4tmp[1] = jet->p4();     
00606     }
00607 
00608     if (mMaxEInEmTowers) mMaxEInEmTowers->Fill (jet->maxEInEmTowers());
00609     if (mMaxEInHadTowers) mMaxEInHadTowers->Fill (jet->maxEInHadTowers());
00610     if (mHadEnergyInHO) mHadEnergyInHO->Fill (jet->hadEnergyInHO());
00611     if (mHadEnergyInHO_80)   mHadEnergyInHO_80->Fill (jet->hadEnergyInHO());
00612     if (mHadEnergyInHO_3000) mHadEnergyInHO_3000->Fill (jet->hadEnergyInHO());
00613     if (mHadEnergyInHB) mHadEnergyInHB->Fill (jet->hadEnergyInHB());
00614     if (mHadEnergyInHB_80)   mHadEnergyInHB_80->Fill (jet->hadEnergyInHB());
00615     if (mHadEnergyInHB_3000) mHadEnergyInHB_3000->Fill (jet->hadEnergyInHB());
00616     if (mHadEnergyInHF) mHadEnergyInHF->Fill (jet->hadEnergyInHF());
00617     if (mHadEnergyInHE) mHadEnergyInHE->Fill (jet->hadEnergyInHE());
00618     if (mHadEnergyInHE_80)   mHadEnergyInHE_80->Fill (jet->hadEnergyInHE());
00619     if (mHadEnergyInHE_3000) mHadEnergyInHE_3000->Fill (jet->hadEnergyInHE());
00620     if (mEmEnergyInEB) mEmEnergyInEB->Fill (jet->emEnergyInEB());
00621     if (mEmEnergyInEB_80)   mEmEnergyInEB_80->Fill (jet->emEnergyInEB());
00622     if (mEmEnergyInEB_3000) mEmEnergyInEB_3000->Fill (jet->emEnergyInEB());
00623     if (mEmEnergyInEE) mEmEnergyInEE->Fill (jet->emEnergyInEE());
00624     if (mEmEnergyInEE_80)   mEmEnergyInEE_80->Fill (jet->emEnergyInEE());
00625     if (mEmEnergyInEE_3000) mEmEnergyInEE_3000->Fill (jet->emEnergyInEE());
00626     if (mEmEnergyInHF) mEmEnergyInHF->Fill (jet->emEnergyInHF());
00627     if (mEnergyFractionHadronic) mEnergyFractionHadronic->Fill (jet->energyFractionHadronic());
00628     if (mEnergyFractionEm) mEnergyFractionEm->Fill (jet->emEnergyFraction());
00629 
00630     if (mHFTotal)      mHFTotal->Fill (jet->hadEnergyInHF()+jet->emEnergyInHF());
00631     if (mHFTotal_80)   mHFTotal_80->Fill (jet->hadEnergyInHF()+jet->emEnergyInHF());
00632     if (mHFTotal_3000) mHFTotal_3000->Fill (jet->hadEnergyInHF()+jet->emEnergyInHF());
00633     if (mHFLong)       mHFLong->Fill (jet->hadEnergyInHF()*0.5+jet->emEnergyInHF());
00634     if (mHFLong_80)    mHFLong_80->Fill (jet->hadEnergyInHF()*0.5+jet->emEnergyInHF());
00635     if (mHFLong_3000)  mHFLong_3000->Fill (jet->hadEnergyInHF()*0.5+jet->emEnergyInHF());
00636     if (mHFShort)      mHFShort->Fill (jet->hadEnergyInHF()*0.5);
00637     if (mHFShort_80)   mHFShort_80->Fill (jet->hadEnergyInHF()*0.5);
00638     if (mHFShort_3000) mHFShort_3000->Fill (jet->hadEnergyInHF()*0.5);
00639 
00640 
00641 
00642     if (mN90) mN90->Fill (jet->n90());
00643     mJetEnergyProfile->Fill (jet->eta(), jet->phi(), jet->energy());
00644     mHadJetEnergyProfile->Fill (jet->eta(), jet->phi(), jet->hadEnergyInHO()+jet->hadEnergyInHB()+jet->hadEnergyInHF()+jet->hadEnergyInHE());
00645     mEMJetEnergyProfile->Fill (jet->eta(), jet->phi(), jet->emEnergyInEB()+jet->emEnergyInEE()+jet->emEnergyInHF());
00646   }
00647 
00648 
00649 
00650   if (mNJetsEtaC) mNJetsEtaC->Fill( nJetC );
00651   if (mNJetsEtaF) mNJetsEtaF->Fill( nJetF );
00652 
00653   if (nJet == 2) {
00654     if (mMjj) mMjj->Fill( (p4tmp[0]+p4tmp[1]).mass() );
00655     if (mMjj_3000) mMjj_3000->Fill( (p4tmp[0]+p4tmp[1]).mass() );
00656   }
00657 
00658   // Count Jets above Pt cut
00659   for (int istep = 0; istep < 100; ++istep) {
00660     int     njet = 0;
00661     float ptStep = (istep * (200./100.));
00662 
00663     for ( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) {
00664       if ( cal->pt() > ptStep ) njet++;
00665     }
00666     mNJets1->Fill( ptStep, njet );
00667   }
00668 
00669   for (int istep = 0; istep < 100; ++istep) {
00670     int     njet = 0;
00671     float ptStep = (istep * (4000./100.));
00672     for ( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) {
00673       if ( cal->pt() > ptStep ) njet++;
00674     }
00675     mNJets2->Fill( ptStep, njet );
00676   }
00677 
00678   // Gen jet analysis
00679   Handle<GenJetCollection> genJets;
00680   mEvent.getByLabel(mInputGenCollection, genJets);
00681   if (!genJets.isValid()) return;
00682   GenJetCollection::const_iterator gjet = genJets->begin ();
00683   int gjetIndex = 0;
00684   for (; gjet != genJets->end (); gjet++, gjetIndex++) {
00685     if (mGenEta) mGenEta->Fill (gjet->eta());
00686     if (mGenPhi) mGenPhi->Fill (gjet->phi());
00687     if (mGenPt) mGenPt->Fill (gjet->pt());
00688     if (mGenPt_80) mGenPt_80->Fill (gjet->pt());
00689     if (mGenPt_3000) mGenPt_3000->Fill (gjet->pt());
00690     if (gjet == genJets->begin ()) { // first jet
00691       if (mGenEtaFirst) mGenEtaFirst->Fill (gjet->eta());
00692       if (mGenPhiFirst) mGenPhiFirst->Fill (gjet->phi());
00693     }
00694   }
00695 
00696 
00697   // now match CaloJets to GenJets
00698   JetMatchingTools jetMatching (mEvent);
00699   if (!(mInputGenCollection.label().empty())) {
00700     //    Handle<GenJetCollection> genJets;
00701     //    mEvent.getByLabel(mInputGenCollection, genJets);
00702 
00703     std::vector <std::vector <const reco::GenParticle*> > genJetConstituents (genJets->size());
00704     std::vector <std::vector <const reco::GenParticle*> > caloJetConstituents (caloJets->size());
00705     if (mRThreshold > 0) { 
00706     }
00707     else {
00708       for (unsigned iGenJet = 0; iGenJet < genJets->size(); ++iGenJet) {
00709         genJetConstituents [iGenJet] = jetMatching.getGenParticles ((*genJets) [iGenJet]);
00710       }
00711       
00712       for (unsigned iCaloJet = 0; iCaloJet < caloJets->size(); ++iCaloJet) {
00713         caloJetConstituents [iCaloJet] = jetMatching.getGenParticles ((*caloJets) [iCaloJet], false);
00714       }
00715     }
00716 
00717     for (unsigned iGenJet = 0; iGenJet < genJets->size(); ++iGenJet) {               //****************************************************************
00718     //for (unsigned iGenJet = 0; iGenJet < 1; ++iGenJet) {                           // only FIRST Jet !!!!
00719       const GenJet& genJet = (*genJets) [iGenJet];
00720       double genJetPt = genJet.pt();
00721 
00722       //std::cout << iGenJet <<". Genjet: pT = " << genJetPt << "GeV" << std::endl;  //  *****************************************************
00723 
00724       if (fabs(genJet.eta()) > 5.) continue; // out of detector 
00725       if (genJetPt < mMatchGenPtThreshold) continue; // no low momentum 
00726       double logPtGen = log10 (genJetPt);
00727       mAllGenJetsPt->Fill (logPtGen);
00728       mAllGenJetsEta->Fill (logPtGen, genJet.eta());
00729       if (caloJets->size() <= 0) continue; // no CaloJets - nothing to match
00730       if (mRThreshold > 0) {
00731         unsigned iCaloJetBest = 0;
00732         double deltaRBest = 999.;
00733         for (unsigned iCaloJet = 0; iCaloJet < caloJets->size(); ++iCaloJet) {
00734           double dR = deltaR (genJet.eta(), genJet.phi(), (*caloJets) [iCaloJet].eta(), (*caloJets) [iCaloJet].phi());
00735           if (deltaRBest < mRThreshold && dR < mRThreshold && genJet.pt() > 5.) {
00736             /*
00737             std::cout << "Yet another matched jet for GenJet pt=" << genJet.pt()
00738                       << " previous CaloJet pt/dr: " << (*caloJets) [iCaloJetBest].pt() << '/' << deltaRBest
00739                       << " new CaloJet pt/dr: " << (*caloJets) [iCaloJet].pt() << '/' << dR
00740                       << std::endl;
00741             */
00742           }
00743           if (dR < deltaRBest) {
00744             iCaloJetBest = iCaloJet;
00745             deltaRBest = dR;
00746           }
00747         }
00748         if (mTurnOnEverything.compare("yes")==0) {
00749           mRMatch->Fill (logPtGen, genJet.eta(), deltaRBest);
00750         }
00751         if (deltaRBest < mRThreshold) { // Matched
00752           fillMatchHists (genJet, (*caloJets) [iCaloJetBest]);
00753         }
00754       }
00755       else {
00756         unsigned iCaloJetBest = 0;
00757         double energyFractionBest = 0.;
00758         for (unsigned iCaloJet = 0; iCaloJet < caloJets->size(); ++iCaloJet) {
00759           double energyFraction = jetMatching.overlapEnergyFraction (genJetConstituents [iGenJet], 
00760                                                                      caloJetConstituents [iCaloJet]);
00761           if (energyFraction > energyFractionBest) {
00762             iCaloJetBest = iCaloJet;
00763             energyFractionBest = energyFraction;
00764           }
00765         }
00766         if (mTurnOnEverything.compare("yes")==0) {
00767           mGenJetMatchEnergyFraction->Fill (logPtGen, genJet.eta(), energyFractionBest);
00768         }
00769         if (energyFractionBest > mGenEnergyFractionThreshold) { // good enough
00770           double reverseEnergyFraction = jetMatching.overlapEnergyFraction (caloJetConstituents [iCaloJetBest], 
00771                                                                             genJetConstituents [iGenJet]);
00772           if (mTurnOnEverything.compare("yes")==0) {
00773             mReverseMatchEnergyFraction->Fill (logPtGen, genJet.eta(), reverseEnergyFraction);
00774           }
00775           if (reverseEnergyFraction > mReverseEnergyFractionThreshold) { // Matched
00776             fillMatchHists (genJet, (*caloJets) [iCaloJetBest]);
00777           }
00778         }
00779       }
00780     }
00781   }
00782 }
00783 
00784 
00785 void CaloJetTester::fillMatchHists (const reco::GenJet& fGenJet, const reco::CaloJet& fCaloJet) {
00786   double logPtGen = log10 (fGenJet.pt());
00787   double PtGen = fGenJet.pt();
00788   double PtCalo = fCaloJet.pt();
00789   mMatchedGenJetsPt->Fill (logPtGen);
00790   mMatchedGenJetsEta->Fill (logPtGen, fGenJet.eta());
00791 
00792   double PtThreshold = 10.;
00793 
00794   if (mTurnOnEverything.compare("yes")==0) {
00795     mDeltaEta->Fill (logPtGen, fGenJet.eta(), fCaloJet.eta()-fGenJet.eta());
00796     mDeltaPhi->Fill (logPtGen, fGenJet.eta(), fCaloJet.phi()-fGenJet.phi());
00797     mEScale->Fill (logPtGen, fGenJet.eta(), fCaloJet.energy()/fGenJet.energy());
00798     mlinEScale->Fill (fGenJet.pt(), fGenJet.eta(), fCaloJet.energy()/fGenJet.energy());
00799     mDeltaE->Fill (logPtGen, fGenJet.eta(), fCaloJet.energy()-fGenJet.energy());
00800 
00801     mEScaleFineBin->Fill (logPtGen, fGenJet.eta(), fCaloJet.energy()/fGenJet.energy());
00802   
00803     if (fGenJet.pt()>PtThreshold) {
00804       mEScale_pt10->Fill (logPtGen, fGenJet.eta(), fCaloJet.energy()/fGenJet.energy());
00805 
00806     }
00807 
00808   }
00809   if (fCaloJet.pt() > PtThreshold) {
00810     mDelEta->Fill (fGenJet.eta()-fCaloJet.eta());
00811     mDelPhi->Fill (fGenJet.phi()-fCaloJet.phi());
00812     mDelPt->Fill  ((fGenJet.pt()-fCaloJet.pt())/fGenJet.pt());
00813   }
00814 
00815   if (fabs(fGenJet.eta())<1.3) {
00816 
00817     mpTScaleB_s->Fill (log10(PtGen), PtCalo/PtGen);
00818     mpTScaleB_d->Fill (log10(PtGen), PtCalo/PtGen);
00819     if (PtGen>60.0 && PtGen<120.0) {
00820       mpTScale1DB_60_120->Fill (fCaloJet.pt()/fGenJet.pt());
00821     }
00822     if (PtGen>200.0 && PtGen<300.0) {
00823       mpTScale1DB_200_300->Fill (fCaloJet.pt()/fGenJet.pt());
00824     }
00825     if (PtGen>600.0 && PtGen<900.0) {
00826       mpTScale1DB_600_900->Fill (fCaloJet.pt()/fGenJet.pt());
00827     }
00828     if (PtGen>2700.0 && PtGen<3500.0) {
00829       mpTScale1DB_2700_3500->Fill (fCaloJet.pt()/fGenJet.pt());
00830     }
00831   }
00832 
00833   if (fabs(fGenJet.eta())>1.3 && fabs(fGenJet.eta())<3.0) {
00834 
00835     mpTScaleE_s->Fill (log10(PtGen), PtCalo/PtGen);
00836     mpTScaleE_d->Fill (log10(PtGen), PtCalo/PtGen);
00837     if (PtGen>60.0 && PtGen<120.0) {
00838       mpTScale1DE_60_120->Fill (fCaloJet.pt()/fGenJet.pt());
00839     }
00840     if (PtGen>200.0 && PtGen<300.0) {
00841       mpTScale1DE_200_300->Fill (fCaloJet.pt()/fGenJet.pt());
00842     }
00843     if (PtGen>600.0 && PtGen<900.0) {
00844       mpTScale1DE_600_900->Fill (fCaloJet.pt()/fGenJet.pt());
00845     }
00846     if (PtGen>2700.0 && PtGen<3500.0) {
00847       mpTScale1DE_2700_3500->Fill (fCaloJet.pt()/fGenJet.pt());
00848     }
00849   }
00850 
00851   if (fabs(fGenJet.eta())>3.0 && fabs(fGenJet.eta())<5.0) {
00852 
00853     mpTScaleF_s->Fill (log10(PtGen), PtCalo/PtGen);
00854     mpTScaleF_d->Fill (log10(PtGen), PtCalo/PtGen);
00855     if (PtGen>60.0 && PtGen<120.0) {
00856       mpTScale1DF_60_120->Fill (fCaloJet.pt()/fGenJet.pt());
00857     }
00858     if (PtGen>200.0 && PtGen<300.0) {
00859       mpTScale1DF_200_300->Fill (fCaloJet.pt()/fGenJet.pt());
00860     }
00861     if (PtGen>600.0 && PtGen<900.0) {
00862       mpTScale1DF_600_900->Fill (fCaloJet.pt()/fGenJet.pt());
00863     }
00864     if (PtGen>2700.0 && PtGen<3500.0) {
00865       mpTScale1DF_2700_3500->Fill (fCaloJet.pt()/fGenJet.pt());
00866     }
00867   }
00868 
00869   if (fGenJet.pt()>60.0 && fGenJet.pt()<120.0) {
00870     mpTScale_60_120_s->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00871     mpTScale_60_120_d->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00872     mpTScale1D_60_120->Fill (fCaloJet.pt()/fGenJet.pt());
00873   }
00874 
00875   if (fGenJet.pt()>200.0 && fGenJet.pt()<300.0) {
00876     mpTScale_200_300_s->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00877     mpTScale_200_300_d->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00878     mpTScale1D_200_300->Fill (fCaloJet.pt()/fGenJet.pt());
00879   }
00880 
00881   if (fGenJet.pt()>600.0 && fGenJet.pt()<900.0) {
00882     mpTScale_600_900_s->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00883     mpTScale_600_900_d->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00884     mpTScale1D_600_900->Fill (fCaloJet.pt()/fGenJet.pt());
00885   }
00886 
00887   if (fGenJet.pt()>2700.0 && fGenJet.pt()<3500.0) {
00888     mpTScale_2700_3500_s->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00889     mpTScale_2700_3500_d->Fill (fGenJet.eta(),fCaloJet.pt()/fGenJet.pt());
00890     mpTScale1D_2700_3500->Fill (fCaloJet.pt()/fGenJet.pt());
00891   }
00892 
00893 
00894 
00895 }