CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/Validation/RecoEgamma/plugins/EgammaObjects.cc

Go to the documentation of this file.
00001 #include "Validation/RecoEgamma/plugins/EgammaObjects.h"
00002 
00003 #include "FWCore/Utilities/interface/Exception.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 #include "DataFormats/Common/interface/Handle.h"
00007 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00008 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00009 #include "DataFormats/EgammaReco/interface/ClusterShape.h"
00010 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00011 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00012 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00013 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00014 
00015 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00016 
00017 EgammaObjects::EgammaObjects( const edm::ParameterSet& ps )
00018 {
00019   particleID = ps.getParameter<int>("particleID");
00020   EtCut = ps.getParameter<int>("EtCut");
00021 
00022   if( particleID == 22 ) particleString = "Photon";
00023   else if( particleID == 11 ) particleString = "Electron";
00024   else
00025     throw(std::runtime_error("\n\nEgammaObjects: Only 11 = Photon and 22 = Electron are acceptable parictleIDs! Exiting...\n\n"));
00026 
00027   loadCMSSWObjects(ps);
00028   loadHistoParameters(ps);
00029 
00030   rootFile_ = TFile::Open(ps.getParameter<std::string>("outputFile").c_str(),"RECREATE");
00031 
00032   hist_EtaEfficiency_ = 0 ;
00033   hist_EtaNumRecoOverNumTrue_ = 0 ;
00034   hist_deltaEtaVsEt_ = 0 ;
00035   hist_deltaEtaVsE_ = 0 ;
00036   hist_deltaEtaVsEta_ = 0 ;
00037   hist_deltaEtaVsPhi_ = 0 ;
00038   hist_resolutionEtaVsEt_ = 0 ;
00039   hist_resolutionEtaVsE_ = 0 ;
00040   hist_resolutionEtaVsEta_ = 0 ;
00041   hist_resolutionEtaVsPhi_ = 0 ;
00042 
00043   hist_Phi_ = 0 ;
00044   hist_PhiOverTruth_ = 0 ;
00045   hist_PhiEfficiency_ = 0 ;
00046   hist_PhiNumRecoOverNumTrue_ = 0 ;
00047   hist_deltaPhiVsEt_ = 0 ;
00048   hist_deltaPhiVsE_ = 0 ;
00049   hist_deltaPhiVsEta_ = 0 ;
00050   hist_deltaPhiVsPhi_ = 0 ;
00051   hist_resolutionPhiVsEt_ = 0 ;
00052   hist_resolutionPhiVsE_ = 0 ;
00053   hist_resolutionPhiVsEta_ = 0 ;
00054   hist_resolutionPhiVsPhi_ = 0 ;
00055 
00056   hist_All_recoMass_ = 0 ;
00057   hist_BarrelOnly_recoMass_ = 0 ;
00058   hist_EndcapOnly_recoMass_ = 0 ;
00059   hist_Mixed_recoMass_ = 0 ;
00060 
00061   hist_recoMass_withBackgroud_NoEtCut_ = 0 ;
00062   hist_recoMass_withBackgroud_5EtCut_ = 0 ;
00063   hist_recoMass_withBackgroud_10EtCut_ = 0 ;
00064   hist_recoMass_withBackgroud_20EtCut_ = 0 ;
00065 
00066   _TEMP_scatterPlot_EtOverTruthVsEt_ = 0 ;
00067   _TEMP_scatterPlot_EtOverTruthVsE_ = 0 ;
00068   _TEMP_scatterPlot_EtOverTruthVsEta_ = 0 ;
00069   _TEMP_scatterPlot_EtOverTruthVsPhi_ = 0 ;
00070 
00071   _TEMP_scatterPlot_EOverTruthVsEt_ = 0 ;
00072   _TEMP_scatterPlot_EOverTruthVsE_ = 0 ;
00073   _TEMP_scatterPlot_EOverTruthVsEta_ = 0 ;
00074   _TEMP_scatterPlot_EOverTruthVsPhi_ = 0 ;
00075 
00076   _TEMP_scatterPlot_deltaEtaVsEt_ = 0 ;
00077   _TEMP_scatterPlot_deltaEtaVsE_ = 0 ;
00078   _TEMP_scatterPlot_deltaEtaVsEta_ = 0 ;
00079   _TEMP_scatterPlot_deltaEtaVsPhi_ = 0 ;
00080 
00081   _TEMP_scatterPlot_deltaPhiVsEt_ = 0 ;
00082   _TEMP_scatterPlot_deltaPhiVsE_ = 0 ;
00083   _TEMP_scatterPlot_deltaPhiVsEta_ = 0 ;
00084   _TEMP_scatterPlot_deltaPhiVsPhi_ = 0 ;
00085 
00086 }
00087 
00088 void EgammaObjects::loadCMSSWObjects(const edm::ParameterSet& ps)
00089 {
00090   MCTruthCollection_ = ps.getParameter<edm::InputTag>("MCTruthCollection");
00091   RecoCollection_  = ps.getParameter<edm::InputTag>("RecoCollection");
00092 }
00093 
00094 void EgammaObjects::loadHistoParameters(const edm::ParameterSet& ps)
00095 {
00096   hist_min_Et_  = ps.getParameter<double>("hist_min_Et");
00097   hist_max_Et_  = ps.getParameter<double>("hist_max_Et");
00098   hist_bins_Et_ = ps.getParameter<int>   ("hist_bins_Et");
00099 
00100   hist_min_E_  = ps.getParameter<double>("hist_min_E");
00101   hist_max_E_  = ps.getParameter<double>("hist_max_E");
00102   hist_bins_E_ = ps.getParameter<int>   ("hist_bins_E");
00103 
00104   hist_min_Eta_  = ps.getParameter<double>("hist_min_Eta");
00105   hist_max_Eta_  = ps.getParameter<double>("hist_max_Eta");
00106   hist_bins_Eta_ = ps.getParameter<int>   ("hist_bins_Eta");
00107 
00108   hist_min_Phi_  = ps.getParameter<double>("hist_min_Phi");
00109   hist_max_Phi_  = ps.getParameter<double>("hist_max_Phi");
00110   hist_bins_Phi_ = ps.getParameter<int>   ("hist_bins_Phi");
00111 
00112   hist_min_EtOverTruth_  = ps.getParameter<double>("hist_min_EtOverTruth");
00113   hist_max_EtOverTruth_  = ps.getParameter<double>("hist_max_EtOverTruth");
00114   hist_bins_EtOverTruth_ = ps.getParameter<int>   ("hist_bins_EtOverTruth");
00115 
00116   hist_min_EOverTruth_  = ps.getParameter<double>("hist_min_EOverTruth");
00117   hist_max_EOverTruth_  = ps.getParameter<double>("hist_max_EOverTruth");
00118   hist_bins_EOverTruth_ = ps.getParameter<int>   ("hist_bins_EOverTruth");
00119 
00120   hist_min_EtaOverTruth_  = ps.getParameter<double>("hist_min_EtaOverTruth");
00121   hist_max_EtaOverTruth_  = ps.getParameter<double>("hist_max_EtaOverTruth");
00122   hist_bins_EtaOverTruth_ = ps.getParameter<int>   ("hist_bins_EtaOverTruth");
00123 
00124   hist_min_PhiOverTruth_  = ps.getParameter<double>("hist_min_PhiOverTruth");
00125   hist_max_PhiOverTruth_  = ps.getParameter<double>("hist_max_PhiOverTruth");
00126   hist_bins_PhiOverTruth_ = ps.getParameter<int>   ("hist_bins_PhiOverTruth");
00127 
00128   hist_min_deltaEta_  = ps.getParameter<double>("hist_min_deltaEta");
00129   hist_max_deltaEta_  = ps.getParameter<double>("hist_max_deltaEta");
00130   hist_bins_deltaEta_ = ps.getParameter<int>   ("hist_bins_deltaEta");
00131 
00132   hist_min_deltaPhi_  = ps.getParameter<double>("hist_min_deltaPhi");
00133   hist_max_deltaPhi_  = ps.getParameter<double>("hist_max_deltaPhi");
00134   hist_bins_deltaPhi_ = ps.getParameter<int>   ("hist_bins_deltaPhi");
00135 
00136   hist_min_recoMass_ = ps.getParameter<double>("hist_min_recoMass");
00137   hist_max_recoMass_ = ps.getParameter<double>("hist_max_recoMass");
00138   hist_bins_recoMass_ = ps.getParameter<int>   ("hist_bins_recoMass");
00139 }
00140 
00141 EgammaObjects::~EgammaObjects()
00142 {
00143   delete rootFile_;
00144 }
00145 
00146 void EgammaObjects::beginJob()
00147 {
00148   TH1::SetDefaultSumw2(true);
00149 
00150   createBookedHistoObjects();
00151   createTempHistoObjects();
00152 }
00153 
00154 void EgammaObjects::createBookedHistoObjects()
00155 {
00156   hist_Et_
00157     = new TH1D("hist_Et_",("Et Distribution of "+particleString).c_str(),
00158       hist_bins_Et_,hist_min_Et_,hist_max_Et_);
00159   hist_EtOverTruth_
00160     = new TH1D("hist_EtOverTruth_",("Reco Et over True Et of "+particleString).c_str(),
00161       hist_bins_EtOverTruth_,hist_min_EtOverTruth_,hist_max_EtOverTruth_);
00162   hist_EtEfficiency_
00163     = new TH1D("hist_EtEfficiency_",("# of True "+particleString+" Reconstructed over # of True "+particleString+" VS Et of "+particleString).c_str(),
00164       hist_bins_Et_,hist_min_Et_,hist_max_Et_);
00165   hist_EtNumRecoOverNumTrue_
00166     = new TH1D("hist_EtNumRecoOverNumTrue_",("# of Reco "+particleString+" over # of True "+particleString+" VS Et of "+particleString).c_str(),
00167       hist_bins_Et_,hist_min_Et_,hist_max_Et_);
00168 
00169   hist_E_
00170     = new TH1D("hist_E_",("E Distribution of "+particleString).c_str(),
00171       hist_bins_E_,hist_min_E_,hist_max_E_);
00172   hist_EOverTruth_
00173     = new TH1D("hist_EOverTruth_",("Reco E over True E of "+particleString).c_str(),
00174       hist_bins_EOverTruth_,hist_min_EOverTruth_,hist_max_EOverTruth_);
00175   hist_EEfficiency_
00176     = new TH1D("hist_EEfficiency_",("# of True "+particleString+" Reconstructed over # of True "+particleString+" VS E of "+particleString).c_str(),
00177       hist_bins_E_,hist_min_E_,hist_max_E_);
00178   hist_ENumRecoOverNumTrue_
00179     = new TH1D("hist_ENumRecoOverNumTrue_",("# of Reco "+particleString+" over # of True "+particleString+" VS E of "+particleString).c_str(),
00180       hist_bins_E_,hist_min_E_,hist_max_E_);
00181 
00182   hist_Eta_
00183     = new TH1D("hist_Eta_",("Eta Distribution of "+particleString).c_str(),
00184       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_);
00185   hist_EtaOverTruth_
00186     = new TH1D("hist_EtaOverTruth_",("Reco Eta over True Eta of "+particleString).c_str(),
00187       hist_bins_EtaOverTruth_,hist_min_EtaOverTruth_,hist_max_EtaOverTruth_);
00188   hist_EtaEfficiency_
00189     = new TH1D("hist_EtaEfficiency_",("# of True "+particleString+" Reconstructed over # of True "+particleString+" VS Eta of "+particleString).c_str(),
00190       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_);
00191   hist_EtaNumRecoOverNumTrue_
00192     = new TH1D("hist_EtaNumRecoOverNumTrue_",("# of Reco "+particleString+" over # of True "+particleString+" VS Eta of "+particleString).c_str(),
00193       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_);
00194 
00195   hist_Phi_
00196     = new TH1D("hist_Phi_",("Phi Distribution of "+particleString).c_str(),
00197       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_);
00198   hist_PhiOverTruth_
00199     = new TH1D("hist_PhiOverTruth_",("Reco Phi over True Phi of "+particleString).c_str(),
00200       hist_bins_PhiOverTruth_,hist_min_PhiOverTruth_,hist_max_PhiOverTruth_);
00201   hist_PhiEfficiency_
00202     = new TH1D("hist_PhiEfficiency_",("# of True "+particleString+" Reconstructed over # of True "+particleString+" VS Phi of "+particleString).c_str(),
00203       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_);
00204    hist_PhiNumRecoOverNumTrue_
00205     = new TH1D("hist_PhiNumRecoOverNumTrue_",("# of Reco "+particleString+" over # of True "+particleString+" VS Phi of "+particleString).c_str(),
00206       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_);
00207 
00208   std::string recoParticleName;
00209 
00210   if( particleID == 22 ) recoParticleName = "Higgs";
00211   else if( particleID == 11 ) recoParticleName = "Z";
00212 
00213   hist_All_recoMass_
00214     = new TH1D("hist_All_recoMass_",(recoParticleName+" Mass from "+particleString+" in All Regions").c_str(),
00215       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00216   hist_BarrelOnly_recoMass_
00217     = new TH1D("hist_BarrelOnly_recoMass_",(recoParticleName+" Mass from "+particleString+" in Barrel").c_str(),
00218       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00219   hist_EndcapOnly_recoMass_
00220     = new TH1D("hist_EndcapOnly_recoMass_",(recoParticleName+" Mass from "+particleString+" in EndCap").c_str(),
00221       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00222   hist_Mixed_recoMass_
00223     = new TH1D("hist_Mixed_recoMass_",(recoParticleName+" Mass from "+particleString+" in Split Detectors").c_str(),
00224       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00225 
00226   hist_recoMass_withBackgroud_NoEtCut_
00227     = new TH1D("hist_recoMass_withBackgroud_NoEtCut_",(recoParticleName+" Mass from "+particleString+" with Background, No Et Cut").c_str(),
00228       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00229   hist_recoMass_withBackgroud_5EtCut_
00230     = new TH1D("hist_recoMass_withBackgroud_5EtCut_",(recoParticleName+" Mass from "+particleString+" with Background, 5 Et Cut").c_str(),
00231       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00232   hist_recoMass_withBackgroud_10EtCut_
00233     = new TH1D("hist_recoMass_withBackgroud_10EtCut_",(recoParticleName+" Mass from "+particleString+" with Background, 10 Et Cut").c_str(),
00234       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00235   hist_recoMass_withBackgroud_20EtCut_
00236     = new TH1D("hist_recoMass_withBackgroud_20EtCut_",(recoParticleName+" Mass from "+particleString+" with Background, 20 Et Cut").c_str(),
00237       hist_bins_recoMass_,hist_min_recoMass_,hist_max_recoMass_);
00238 }
00239 
00240 void EgammaObjects::createTempHistoObjects()
00241 {
00242   _TEMP_scatterPlot_EtOverTruthVsEt_
00243     = new TH2D("_TEMP_scatterPlot_EtOverTruthVsEt_","_TEMP_scatterPlot_EtOverTruthVsEt_",
00244       hist_bins_Et_,hist_min_Et_,hist_max_Et_,hist_bins_EtOverTruth_,hist_min_EtOverTruth_,hist_max_EtOverTruth_);
00245   _TEMP_scatterPlot_EtOverTruthVsE_
00246     = new TH2D("_TEMP_scatterPlot_EtOverTruthVsE_","_TEMP_scatterPlot_EtOverTruthVsE_",
00247       hist_bins_E_,hist_min_E_,hist_max_E_,hist_bins_EtOverTruth_,hist_min_EtOverTruth_,hist_max_EtOverTruth_);
00248   _TEMP_scatterPlot_EtOverTruthVsEta_
00249     = new TH2D("_TEMP_scatterPlot_EtOverTruthVsEta_","_TEMP_scatterPlot_EtOverTruthVsEta_",
00250       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_,hist_bins_EtOverTruth_,hist_min_EtOverTruth_,hist_max_EtOverTruth_);
00251   _TEMP_scatterPlot_EtOverTruthVsPhi_
00252     = new TH2D("_TEMP_scatterPlot_EtOverTruthVsPhi_","_TEMP_scatterPlot_EtOverTruthVsPhi_",
00253       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_,hist_bins_EtOverTruth_,hist_min_EtOverTruth_,hist_max_EtOverTruth_);
00254 
00255   _TEMP_scatterPlot_EOverTruthVsEt_
00256     = new TH2D("_TEMP_scatterPlot_EOverTruthVsEt_","_TEMP_scatterPlot_EOverTruthVsEt_",
00257       hist_bins_Et_,hist_min_Et_,hist_max_Et_,hist_bins_EOverTruth_,hist_min_EOverTruth_,hist_max_EOverTruth_);
00258   _TEMP_scatterPlot_EOverTruthVsE_
00259     = new TH2D("_TEMP_scatterPlot_EOverTruthVsE_","_TEMP_scatterPlot_EOverTruthVsE_",
00260       hist_bins_E_,hist_min_E_,hist_max_E_,hist_bins_EOverTruth_,hist_min_EOverTruth_,hist_max_EOverTruth_);
00261   _TEMP_scatterPlot_EOverTruthVsEta_
00262     = new TH2D("_TEMP_scatterPlot_EOverTruthVsEta_","_TEMP_scatterPlot_EOverTruthVsEta_",
00263       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_,hist_bins_EOverTruth_,hist_min_EOverTruth_,hist_max_EOverTruth_);
00264   _TEMP_scatterPlot_EOverTruthVsPhi_
00265     = new TH2D("_TEMP_scatterPlot_EOverTruthVsPhi_","_TEMP_scatterPlot_EOverTruthVsPhi_",
00266       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_,hist_bins_EOverTruth_,hist_min_EOverTruth_,hist_max_EOverTruth_);
00267 
00268   _TEMP_scatterPlot_deltaEtaVsEt_
00269     = new TH2D("_TEMP_scatterPlot_deltaEtaVsEt_","_TEMP_scatterPlot_deltaEtaVsEt_",
00270       hist_bins_Et_,hist_min_Et_,hist_max_Et_,hist_bins_deltaEta_,hist_min_deltaEta_,hist_max_deltaEta_);
00271   _TEMP_scatterPlot_deltaEtaVsE_
00272     = new TH2D("_TEMP_scatterPlot_deltaEtaVsE_","_TEMP_scatterPlot_deltaEtaVsE_",
00273       hist_bins_E_,hist_min_E_,hist_max_E_,hist_bins_deltaEta_,hist_min_deltaEta_,hist_max_deltaEta_);
00274   _TEMP_scatterPlot_deltaEtaVsEta_
00275     = new TH2D("_TEMP_scatterPlot_deltaEtaVsEta_","_TEMP_scatterPlot_deltaEtaVsEta_",
00276       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_,hist_bins_deltaEta_,hist_min_deltaEta_,hist_max_deltaEta_);
00277   _TEMP_scatterPlot_deltaEtaVsPhi_
00278     = new TH2D("_TEMP_scatterPlot_deltaEtaVsPhi_","_TEMP_scatterPlot_deltaEtaVsPhi_",
00279       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_,hist_bins_deltaEta_,hist_min_deltaEta_,hist_max_deltaEta_);
00280 
00281   _TEMP_scatterPlot_deltaPhiVsEt_
00282     = new TH2D("_TEMP_scatterPlot_deltaPhiVsEt_","_TEMP_scatterPlot_deltaPhiVsEt_",
00283       hist_bins_Et_,hist_min_Et_,hist_max_Et_,hist_bins_deltaPhi_,hist_min_deltaPhi_,hist_max_deltaPhi_);
00284   _TEMP_scatterPlot_deltaPhiVsE_
00285     = new TH2D("_TEMP_scatterPlot_deltaPhiVsE_","_TEMP_scatterPlot_deltaPhiVsE_",
00286       hist_bins_E_,hist_min_E_,hist_max_E_,hist_bins_deltaPhi_,hist_min_deltaPhi_,hist_max_deltaPhi_);
00287   _TEMP_scatterPlot_deltaPhiVsEta_
00288     = new TH2D("_TEMP_scatterPlot_deltaPhiVsEta_","_TEMP_scatterPlot_deltaPhiVsEta_",
00289       hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_,hist_bins_deltaPhi_,hist_min_deltaPhi_,hist_max_deltaPhi_);
00290   _TEMP_scatterPlot_deltaPhiVsPhi_
00291     = new TH2D("_TEMP_scatterPlot_deltaPhiVsPhi_","_TEMP_scatterPlot_deltaPhiVsPhi_",
00292       hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_,hist_bins_deltaPhi_,hist_min_deltaPhi_,hist_max_deltaPhi_);
00293 }
00294 
00295 void EgammaObjects::analyze( const edm::Event& evt, const edm::EventSetup& es )
00296 {
00297   if( particleID == 22 ) analyzePhotons(evt, es);
00298   else if( particleID == 11 ) analyzeElectrons(evt, es);
00299 }
00300 
00301 void EgammaObjects::analyzePhotons( const edm::Event& evt, const edm::EventSetup& es )
00302 {
00303   edm::Handle<reco::PhotonCollection> pPhotons;
00304   evt.getByLabel(RecoCollection_, pPhotons);
00305   if (!pPhotons.isValid()) {
00306     edm::LogError("EgammaObjects") << "Error! can't get collection with label " << RecoCollection_.label();
00307   }
00308 
00309   const reco::PhotonCollection* photons = pPhotons.product();
00310   std::vector<reco::Photon> photonsMCMatched;
00311 
00312   for(reco::PhotonCollection::const_iterator aClus = photons->begin(); aClus != photons->end(); aClus++)
00313   {
00314     if(aClus->et() >= EtCut)
00315     {
00316       hist_Et_->Fill(aClus->et());
00317       hist_E_->Fill(aClus->energy());
00318       hist_Eta_->Fill(aClus->eta());
00319       hist_Phi_->Fill(aClus->phi());
00320      }
00321   }
00322 
00323   for(int firstPhoton = 0, numPhotons = photons->size(); firstPhoton < numPhotons - 1; firstPhoton++)
00324     for(int secondPhoton = firstPhoton + 1; secondPhoton < numPhotons; secondPhoton++)
00325     {
00326       reco::Photon pOne = (*photons)[firstPhoton];
00327       reco::Photon pTwo = (*photons)[secondPhoton];
00328 
00329       double recoMass = findRecoMass(pOne, pTwo);
00330 
00331       hist_recoMass_withBackgroud_NoEtCut_->Fill(recoMass);
00332 
00333       if(pOne.et() >= 5 && pTwo.et() >= 5)
00334         hist_recoMass_withBackgroud_5EtCut_->Fill(recoMass);
00335 
00336       if(pOne.et() >= 10 && pTwo.et() >= 10)
00337         hist_recoMass_withBackgroud_10EtCut_->Fill(recoMass);
00338 
00339       if(pOne.et() >= 20 && pTwo.et() >= 20)
00340         hist_recoMass_withBackgroud_20EtCut_->Fill(recoMass);
00341     }
00342 
00343   edm::Handle<edm::HepMCProduct> pMCTruth ;
00344   evt.getByLabel(MCTruthCollection_, pMCTruth);
00345   if (!pMCTruth.isValid()) {
00346     edm::LogError("EgammaObjects") << "Error! can't get collection with label " << MCTruthCollection_.label();
00347   }
00348 
00349   const HepMC::GenEvent* genEvent = pMCTruth->GetEvent();
00350 
00351   for(HepMC::GenEvent::particle_const_iterator currentParticle = genEvent->particles_begin();
00352       currentParticle != genEvent->particles_end(); currentParticle++ )
00353   {
00354     if(abs((*currentParticle)->pdg_id())==22 && (*currentParticle)->status()==1
00355       && (*currentParticle)->momentum().e()/cosh(ecalEta((*currentParticle)->momentum().eta(), (*currentParticle)->production_vertex()->position().z()/10.,
00356       (*currentParticle)->production_vertex()->position().perp()/10.)) >= EtCut)
00357     {
00358                         HepMC::FourVector vtx = (*currentParticle)->production_vertex()->position();
00359                         double phiTrue = (*currentParticle)->momentum().phi();
00360                         double etaTrue = ecalEta((*currentParticle)->momentum().eta(), vtx.z()/10., vtx.perp()/10.);
00361       double eTrue  = (*currentParticle)->momentum().e();
00362                         double etTrue  = (*currentParticle)->momentum().e()/cosh(etaTrue);
00363 
00364       double etaCurrent, etaFound = -999;
00365       double phiCurrent, phiFound = -999;
00366       double etCurrent,  etFound  = -999;
00367       double eCurrent,   eFound  = -999;
00368 
00369       reco::Photon bestMatchPhoton;
00370 
00371       double closestParticleDistance = 999;
00372 
00373       for(reco::PhotonCollection::const_iterator aClus = photons->begin(); aClus != photons->end(); aClus++)
00374       {
00375         if(aClus->et() > EtCut)
00376         {
00377           etaCurrent = aClus->eta();
00378           phiCurrent = aClus->phi();
00379           etCurrent  = aClus->et();
00380           eCurrent  = aClus->energy();
00381 
00382           double deltaPhi = phiCurrent-phiTrue;
00383           if(deltaPhi > Geom::pi()) deltaPhi -= 2.*Geom::pi();
00384           if(deltaPhi < -Geom::pi()) deltaPhi += 2.*Geom::pi();
00385           double deltaR = std::sqrt(std::pow(etaCurrent-etaTrue,2)+std::pow(deltaPhi,2));
00386 
00387               if(deltaR < closestParticleDistance)
00388           {
00389             etFound  = etCurrent;
00390             eFound   = eCurrent;
00391             etaFound = etaCurrent;
00392             phiFound = phiCurrent;
00393             closestParticleDistance = deltaR;
00394             bestMatchPhoton = *aClus;
00395           }
00396         }
00397       }
00398 
00399       if(closestParticleDistance < 0.05 && etFound/etTrue > .5 && etFound/etTrue < 1.5)
00400       {
00401         hist_EtOverTruth_->Fill(etFound/etTrue);
00402         hist_EOverTruth_->Fill(eFound/eTrue);
00403         hist_EtaOverTruth_->Fill(etaFound/etaTrue);
00404         hist_PhiOverTruth_->Fill(phiFound/phiTrue);
00405 
00406         hist_EtEfficiency_->Fill(etTrue);
00407         hist_EEfficiency_->Fill(eTrue);
00408         hist_EtaEfficiency_->Fill(etaTrue);
00409         hist_PhiEfficiency_->Fill(phiTrue);
00410 
00411         double deltaPhi = phiFound-phiTrue;
00412         if(deltaPhi > Geom::pi()) deltaPhi -= 2.*Geom::pi();
00413         if(deltaPhi < -Geom::pi()) deltaPhi += 2.*Geom::pi();
00414 
00415         _TEMP_scatterPlot_EtOverTruthVsEt_->Fill(etFound,etFound/etTrue);
00416         _TEMP_scatterPlot_EtOverTruthVsE_->Fill(eFound,etFound/etTrue);
00417         _TEMP_scatterPlot_EtOverTruthVsEta_->Fill(etaFound,etFound/etTrue);
00418         _TEMP_scatterPlot_EtOverTruthVsPhi_->Fill(phiFound,etFound/etTrue);
00419 
00420         _TEMP_scatterPlot_EOverTruthVsEt_->Fill(etFound,eFound/eTrue);
00421         _TEMP_scatterPlot_EOverTruthVsE_->Fill(eFound,eFound/eTrue);
00422         _TEMP_scatterPlot_EOverTruthVsEta_->Fill(etaFound,eFound/eTrue);
00423         _TEMP_scatterPlot_EOverTruthVsPhi_->Fill(phiFound,eFound/eTrue);
00424 
00425         _TEMP_scatterPlot_deltaEtaVsEt_->Fill(etFound,etaFound-etaTrue);
00426         _TEMP_scatterPlot_deltaEtaVsE_->Fill(eFound,etaFound-etaTrue);
00427         _TEMP_scatterPlot_deltaEtaVsEta_->Fill(etaFound,etaFound-etaTrue);
00428         _TEMP_scatterPlot_deltaEtaVsPhi_->Fill(phiFound,etaFound-etaTrue);
00429 
00430         _TEMP_scatterPlot_deltaPhiVsEt_->Fill(etFound,deltaPhi);
00431         _TEMP_scatterPlot_deltaPhiVsE_->Fill(eFound,deltaPhi);
00432         _TEMP_scatterPlot_deltaPhiVsEta_->Fill(etaFound,deltaPhi);
00433         _TEMP_scatterPlot_deltaPhiVsPhi_->Fill(phiFound,deltaPhi);
00434 
00435         photonsMCMatched.push_back(bestMatchPhoton);
00436       }
00437 
00438       hist_EtNumRecoOverNumTrue_->Fill(etTrue);
00439       hist_ENumRecoOverNumTrue_->Fill(eTrue);
00440       hist_EtaNumRecoOverNumTrue_->Fill(etaTrue);
00441       hist_PhiNumRecoOverNumTrue_->Fill(phiTrue);
00442     }
00443   }
00444 
00445   if(photonsMCMatched.size() == 2)
00446   {
00447     reco::Photon pOne = photonsMCMatched[0];
00448     reco::Photon pTwo = photonsMCMatched[1];
00449 
00450     double recoMass = findRecoMass(pOne, pTwo);
00451 
00452     hist_All_recoMass_->Fill(recoMass);
00453 
00454     if(pOne.superCluster()->seed()->algo() == 1 && pTwo.superCluster()->seed()->algo() == 1)
00455       hist_BarrelOnly_recoMass_->Fill(recoMass);
00456     else if(pOne.superCluster()->seed()->algo() == 0 && pTwo.superCluster()->seed()->algo() == 0)
00457       hist_EndcapOnly_recoMass_->Fill(recoMass);
00458     else
00459       hist_Mixed_recoMass_->Fill(recoMass);
00460   }
00461 }
00462 
00463 void EgammaObjects::analyzeElectrons( const edm::Event& evt, const edm::EventSetup& es )
00464 {
00465   edm::Handle<reco::GsfElectronCollection> pElectrons;
00466   evt.getByLabel(RecoCollection_, pElectrons);
00467   if (!pElectrons.isValid()) {
00468     edm::LogError("DOEPlotsProducerElectrons") << "Error! can't get collection with label " << RecoCollection_.label();
00469   }
00470 
00471   const reco::GsfElectronCollection* electrons = pElectrons.product();
00472   std::vector<reco::GsfElectron> electronsMCMatched;
00473 
00474   for(reco::GsfElectronCollection::const_iterator aClus = electrons->begin(); aClus != electrons->end(); aClus++)
00475   {
00476     if(aClus->et() >= EtCut)
00477     {
00478       hist_Et_->Fill(aClus->et());
00479       hist_E_->Fill(aClus->energy());
00480       hist_Eta_->Fill(aClus->eta());
00481       hist_Phi_->Fill(aClus->phi());
00482      }
00483   }
00484 
00485   for(int firstElectron = 0, numElectrons = electrons->size(); firstElectron < numElectrons - 1; firstElectron++)
00486     for(int secondElectron = firstElectron + 1; secondElectron < numElectrons; secondElectron++)
00487     {
00488       reco::GsfElectron eOne = (*electrons)[firstElectron];
00489       reco::GsfElectron eTwo = (*electrons)[secondElectron];
00490 
00491       double recoMass = findRecoMass(eOne, eTwo);
00492 
00493       hist_recoMass_withBackgroud_NoEtCut_->Fill(recoMass);
00494 
00495       if(eOne.et() >= 5 && eTwo.et() >= 5)
00496         hist_recoMass_withBackgroud_5EtCut_->Fill(recoMass);
00497 
00498       if(eOne.et() >= 10 && eTwo.et() >= 10)
00499         hist_recoMass_withBackgroud_10EtCut_->Fill(recoMass);
00500 
00501       if(eOne.et() >= 20 && eTwo.et() >= 20)
00502         hist_recoMass_withBackgroud_20EtCut_->Fill(recoMass);
00503     }
00504 
00505   edm::Handle<edm::HepMCProduct> pMCTruth ;
00506   evt.getByLabel(MCTruthCollection_, pMCTruth);
00507   if (!pMCTruth.isValid()) {
00508     edm::LogError("DOEPlotsProducerElectrons") << "Error! can't get collection with label " << MCTruthCollection_.label();
00509   }
00510 
00511   const HepMC::GenEvent* genEvent = pMCTruth->GetEvent();
00512   for(HepMC::GenEvent::particle_const_iterator currentParticle = genEvent->particles_begin();
00513       currentParticle != genEvent->particles_end(); currentParticle++ )
00514   {
00515     if(abs((*currentParticle)->pdg_id())==11 && (*currentParticle)->status()==1
00516       && (*currentParticle)->momentum().e()/cosh((*currentParticle)->momentum().eta()) >= EtCut)
00517     {
00518                         double phiTrue = (*currentParticle)->momentum().phi();
00519                         double etaTrue = (*currentParticle)->momentum().eta();
00520       double eTrue  = (*currentParticle)->momentum().e();
00521                         double etTrue  = (*currentParticle)->momentum().e()/cosh(etaTrue);
00522 
00523       double etaCurrent, etaFound = -999;
00524       double phiCurrent, phiFound = -999;
00525       double etCurrent,  etFound  = -999;
00526       double eCurrent,   eFound  = -999;
00527 
00528       reco::GsfElectron bestMatchElectron;
00529 
00530       double closestParticleDistance = 999;
00531 
00532       for(reco::GsfElectronCollection::const_iterator aClus = electrons->begin(); aClus != electrons->end(); aClus++)
00533       {
00534         if(aClus->et() > EtCut)
00535         {
00536           etaCurrent = aClus->eta();
00537           phiCurrent = aClus->phi();
00538           etCurrent  = aClus->et();
00539           eCurrent  = aClus->energy();
00540 
00541           double deltaPhi = phiCurrent-phiTrue;
00542           if(deltaPhi > Geom::pi()) deltaPhi -= 2.*Geom::pi();
00543           if(deltaPhi < -Geom::pi()) deltaPhi += 2.*Geom::pi();
00544           double deltaR = std::sqrt(std::pow(etaCurrent-etaTrue,2)+std::pow(deltaPhi,2));
00545 
00546               if(deltaR < closestParticleDistance)
00547           {
00548             etFound  = etCurrent;
00549             eFound   = eCurrent;
00550             etaFound = etaCurrent;
00551             phiFound = phiCurrent;
00552             closestParticleDistance = deltaR;
00553             bestMatchElectron = *aClus;
00554           }
00555         }
00556       }
00557 
00558       if(closestParticleDistance < 0.05 && etFound/etTrue > .5 && etFound/etTrue < 1.5)
00559       {
00560         hist_EtOverTruth_->Fill(etFound/etTrue);
00561         hist_EOverTruth_->Fill(eFound/eTrue);
00562         hist_EtaOverTruth_->Fill(etaFound/etaTrue);
00563         hist_PhiOverTruth_->Fill(phiFound/phiTrue);
00564 
00565         hist_EtEfficiency_->Fill(etTrue);
00566         hist_EEfficiency_->Fill(eTrue);
00567         hist_EtaEfficiency_->Fill(etaTrue);
00568         hist_PhiEfficiency_->Fill(phiTrue);
00569 
00570         double deltaPhi = phiFound-phiTrue;
00571         if(deltaPhi > Geom::pi()) deltaPhi -= 2.*Geom::pi();
00572         if(deltaPhi < -Geom::pi()) deltaPhi += 2.*Geom::pi();
00573 
00574         _TEMP_scatterPlot_EtOverTruthVsEt_->Fill(etFound,etFound/etTrue);
00575         _TEMP_scatterPlot_EtOverTruthVsE_->Fill(eFound,etFound/etTrue);
00576         _TEMP_scatterPlot_EtOverTruthVsEta_->Fill(etaFound,etFound/etTrue);
00577         _TEMP_scatterPlot_EtOverTruthVsPhi_->Fill(phiFound,etFound/etTrue);
00578 
00579         _TEMP_scatterPlot_EOverTruthVsEt_->Fill(etFound,eFound/eTrue);
00580         _TEMP_scatterPlot_EOverTruthVsE_->Fill(eFound,eFound/eTrue);
00581         _TEMP_scatterPlot_EOverTruthVsEta_->Fill(etaFound,eFound/eTrue);
00582         _TEMP_scatterPlot_EOverTruthVsPhi_->Fill(phiFound,eFound/eTrue);
00583 
00584         _TEMP_scatterPlot_deltaEtaVsEt_->Fill(etFound,etaFound-etaTrue);
00585         _TEMP_scatterPlot_deltaEtaVsE_->Fill(eFound,etaFound-etaTrue);
00586         _TEMP_scatterPlot_deltaEtaVsEta_->Fill(etaFound,etaFound-etaTrue);
00587         _TEMP_scatterPlot_deltaEtaVsPhi_->Fill(phiFound,etaFound-etaTrue);
00588 
00589         _TEMP_scatterPlot_deltaPhiVsEt_->Fill(etFound,deltaPhi);
00590         _TEMP_scatterPlot_deltaPhiVsE_->Fill(eFound,deltaPhi);
00591         _TEMP_scatterPlot_deltaPhiVsEta_->Fill(etaFound,deltaPhi);
00592         _TEMP_scatterPlot_deltaPhiVsPhi_->Fill(phiFound,deltaPhi);
00593 
00594         electronsMCMatched.push_back(bestMatchElectron);
00595       }
00596 
00597       hist_EtNumRecoOverNumTrue_->Fill(etTrue);
00598       hist_ENumRecoOverNumTrue_->Fill(eTrue);
00599       hist_EtaNumRecoOverNumTrue_->Fill(etaTrue);
00600       hist_PhiNumRecoOverNumTrue_->Fill(phiTrue);
00601     }
00602   }
00603 
00604   if(electronsMCMatched.size() == 2)
00605   {
00606     reco::GsfElectron eOne = electronsMCMatched[0];
00607     reco::GsfElectron eTwo = electronsMCMatched[1];
00608 
00609     double recoMass = findRecoMass(eOne, eTwo);
00610 
00611     hist_All_recoMass_->Fill(recoMass);
00612 
00613     if(eOne.superCluster()->seed()->algo() == 1 && eTwo.superCluster()->seed()->algo() == 1)
00614       hist_BarrelOnly_recoMass_->Fill(recoMass);
00615     else if(eOne.superCluster()->seed()->algo() == 0 && eTwo.superCluster()->seed()->algo() == 0)
00616       hist_EndcapOnly_recoMass_->Fill(recoMass);
00617     else
00618       hist_Mixed_recoMass_->Fill(recoMass);
00619   }
00620 }
00621 
00622 double EgammaObjects::findRecoMass(reco::Photon pOne, reco::Photon pTwo)
00623 {
00624   double cosTheta
00625     = (cos(pOne.superCluster()->phi() - pTwo.superCluster()->phi()) + sinh(pOne.superCluster()->eta()) * sinh(pTwo.superCluster()->eta())) /
00626       (cosh(pOne.superCluster()->eta()) * cosh(pTwo.superCluster()->eta()));
00627 
00628   double recoMass = sqrt(2 * (pOne.superCluster())->energy() * (pTwo.superCluster())->energy() * (1 - cosTheta));
00629 
00630   return recoMass;
00631 }
00632 
00633 double EgammaObjects::findRecoMass(reco::GsfElectron eOne, reco::GsfElectron eTwo)
00634 {
00635   double cosTheta
00636     = (cos(eOne.caloPosition().phi() - eTwo.caloPosition().phi()) + sinh(eOne.caloPosition().eta()) * sinh(eTwo.caloPosition().eta())) /
00637       (cosh(eOne.caloPosition().eta()) * cosh(eTwo.caloPosition().eta()));
00638 
00639   double recoMass = sqrt(2 * eOne.caloEnergy() * eTwo.caloEnergy() * (1 - cosTheta));
00640 
00641   return recoMass;
00642 }
00643 
00644 float EgammaObjects::ecalEta(float EtaParticle , float Zvertex, float plane_Radius)
00645 {
00646         const float R_ECAL           = 136.5;
00647         const float Z_Endcap         = 328.0;
00648         const float etaBarrelEndcap  = 1.479;
00649 
00650         if(EtaParticle != 0.)
00651         {
00652                 float Theta = 0.0  ;
00653                 float ZEcal = (R_ECAL-plane_Radius)*sinh(EtaParticle)+Zvertex;
00654 
00655                 if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
00656                 if(Theta<0.0) Theta = Theta+Geom::pi() ;
00657 
00658                 float ETA = - log(tan(0.5*Theta));
00659 
00660                 if( std::abs(ETA) > etaBarrelEndcap )
00661                 {
00662                         float Zend = Z_Endcap ;
00663                         if(EtaParticle<0.0 )  Zend = -Zend ;
00664                         float Zlen = Zend - Zvertex ;
00665                         float RR = Zlen/sinh(EtaParticle);
00666                         Theta = atan((RR+plane_Radius)/Zend);
00667                         if(Theta<0.0) Theta = Theta+Geom::pi() ;
00668                         ETA = - log(tan(0.5*Theta));
00669                 }
00670 
00671                 return ETA;
00672         }
00673         else
00674         {
00675                 edm::LogWarning("")  << "[EgammaObjects::ecalEta] Warning: Eta equals to zero, not correcting" ;
00676                 return EtaParticle;
00677         }
00678 }
00679 
00680 
00681 void EgammaObjects::endJob()
00682 {
00683   rootFile_->cd();
00684   rootFile_->mkdir(particleString.c_str());
00685 
00686   getDeltaResHistosViaSlicing();
00687   getEfficiencyHistosViaDividing();
00688   fitHistos();
00689 
00690   applyLabels();
00691   setDrawOptions();
00692   saveHistos();
00693   rootFile_->Close();
00694 }
00695 
00696 void EgammaObjects::getDeltaResHistosViaSlicing()
00697 {
00698   _TEMP_scatterPlot_EtOverTruthVsEt_->FitSlicesY(0,1,hist_bins_Et_,10,"QRG3");
00699   _TEMP_scatterPlot_EtOverTruthVsE_->FitSlicesY(0,1,hist_bins_E_,10,"QRG3");
00700   _TEMP_scatterPlot_EtOverTruthVsEta_->FitSlicesY(0,1,hist_bins_Eta_,10,"QRG2");
00701   _TEMP_scatterPlot_EtOverTruthVsPhi_->FitSlicesY(0,1,hist_bins_Phi_,10,"QRG2");
00702 
00703   _TEMP_scatterPlot_EOverTruthVsEt_->FitSlicesY(0,1,hist_bins_Et_,10,"QRG3");
00704   _TEMP_scatterPlot_EOverTruthVsE_->FitSlicesY(0,1,hist_bins_E_,10,"QRG3");
00705   _TEMP_scatterPlot_EOverTruthVsEta_->FitSlicesY(0,1,hist_bins_Eta_,10,"QRG2");
00706   _TEMP_scatterPlot_EOverTruthVsPhi_->FitSlicesY(0,1,hist_bins_Phi_,10,"QRG2");
00707 
00708   _TEMP_scatterPlot_deltaEtaVsEt_->FitSlicesY(0,1,hist_bins_Et_,10,"QRG3");
00709   _TEMP_scatterPlot_deltaEtaVsE_->FitSlicesY(0,1,hist_bins_E_,10,"QRG3");
00710   _TEMP_scatterPlot_deltaEtaVsEta_->FitSlicesY(0,1,hist_bins_Eta_,10,"QRG2");
00711   _TEMP_scatterPlot_deltaEtaVsPhi_->FitSlicesY(0,1,hist_bins_Phi_,10,"QRG2");
00712 
00713   _TEMP_scatterPlot_deltaPhiVsEt_->FitSlicesY(0,1,hist_bins_Et_,10,"QRG3");
00714   _TEMP_scatterPlot_deltaPhiVsE_->FitSlicesY(0,1,hist_bins_E_,10,"QRG3");
00715   _TEMP_scatterPlot_deltaPhiVsEta_->FitSlicesY(0,1,hist_bins_Eta_,10,"QRG2");
00716   _TEMP_scatterPlot_deltaPhiVsPhi_->FitSlicesY(0,1,hist_bins_Phi_,10,"QRG2");
00717 
00718   hist_EtOverTruthVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsEt__1");
00719   hist_EtOverTruthVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsE__1");
00720   hist_EtOverTruthVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsEta__1");
00721   hist_EtOverTruthVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsPhi__1");
00722 
00723   hist_EOverTruthVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsEt__1");
00724   hist_EOverTruthVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsE__1");
00725   hist_EOverTruthVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsEta__1");
00726   hist_EOverTruthVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsPhi__1");
00727 
00728   hist_deltaEtaVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsEt__1");
00729   hist_deltaEtaVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsE__1");
00730   hist_deltaEtaVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsEta__1");
00731   hist_deltaEtaVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsPhi__1");
00732 
00733   hist_deltaPhiVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsEt__1");
00734   hist_deltaPhiVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsE__1");
00735   hist_deltaPhiVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsEta__1");
00736   hist_deltaPhiVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsPhi__1");
00737 
00738   hist_EtOverTruthVsEt_->SetNameTitle("hist_EtOverTruthVsEt_",("Reco Et over True Et VS Et of "+particleString).c_str());
00739   hist_EtOverTruthVsE_->SetNameTitle("hist_EtOverTruthVsE_",("Reco Et over True Et VS E of "+particleString).c_str());
00740   hist_EtOverTruthVsEta_->SetNameTitle("hist_EtOverTruthVsEta_",("Reco Et over True Et VS Eta of "+particleString).c_str());
00741   hist_EtOverTruthVsPhi_->SetNameTitle("hist_EtOverTruthVsPhi_",("Reco Et over True Et VS Phi of "+particleString).c_str());
00742 
00743   hist_EOverTruthVsEt_->SetNameTitle("hist_EOverTruthVsEt_",("Reco E over True E VS Et of "+particleString).c_str());
00744   hist_EOverTruthVsE_->SetNameTitle("hist_EOverTruthVsE_",("Reco E over True E VS E of "+particleString).c_str());
00745   hist_EOverTruthVsEta_->SetNameTitle("hist_EOverTruthVsEta_",("Reco E over True E VS Eta of "+particleString).c_str());
00746   hist_EOverTruthVsPhi_->SetNameTitle("hist_EOverTruthVsPhi_",("Reco E over True E VS Phi of "+particleString).c_str());
00747 
00748   hist_deltaEtaVsEt_->SetNameTitle("hist_deltaEtaVsEt_",("delta Eta VS Et of "+particleString).c_str());
00749   hist_deltaEtaVsE_->SetNameTitle("hist_deltaEtaVsE_",("delta Eta VS E of "+particleString).c_str());
00750   hist_deltaEtaVsEta_->SetNameTitle("hist_deltaEtaVsEta_",("delta Eta VS Eta of "+particleString).c_str());
00751   hist_deltaEtaVsPhi_->SetNameTitle("hist_deltaEtaVsPhi_",("delta Eta VS Phi of "+particleString).c_str());
00752 
00753   hist_deltaPhiVsEt_->SetNameTitle("hist_deltaPhiVsEt_",("delta Phi VS Et of "+particleString).c_str());
00754   hist_deltaPhiVsE_->SetNameTitle("hist_deltaPhiVsE_",("delta Phi VS E of "+particleString).c_str());
00755   hist_deltaPhiVsEta_->SetNameTitle("hist_deltaPhiVsEta_",("delta Phi VS Eta of "+particleString).c_str());
00756   hist_deltaPhiVsPhi_->SetNameTitle("hist_deltaPhiVsPhi_",("delta Phi VS Phi of "+particleString).c_str());
00757 
00758   hist_resolutionEtVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsEt__2");
00759   hist_resolutionEtVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsE__2");
00760   hist_resolutionEtVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsEta__2");
00761   hist_resolutionEtVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EtOverTruthVsPhi__2");
00762 
00763   hist_resolutionEVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsEt__2");
00764   hist_resolutionEVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsE__2");
00765   hist_resolutionEVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsEta__2");
00766   hist_resolutionEVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_EOverTruthVsPhi__2");
00767 
00768   hist_resolutionEtaVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsEt__2");
00769   hist_resolutionEtaVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsE__2");
00770   hist_resolutionEtaVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsEta__2");
00771   hist_resolutionEtaVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaEtaVsPhi__2");
00772 
00773   hist_resolutionPhiVsEt_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsEt__2");
00774   hist_resolutionPhiVsE_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsE__2");
00775   hist_resolutionPhiVsEta_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsEta__2");
00776   hist_resolutionPhiVsPhi_ = (TH1D*)gDirectory->Get("_TEMP_scatterPlot_deltaPhiVsPhi__2");
00777 
00778   hist_resolutionEtVsEt_->SetNameTitle("hist_resolutionEtVsEt_",("#sigma of Reco Et over True Et VS Et of "+particleString).c_str());
00779   hist_resolutionEtVsE_->SetNameTitle("hist_resolutionEtVsE_",("#sigma of Reco Et over True Et VS E of "+particleString).c_str());
00780   hist_resolutionEtVsEta_->SetNameTitle("hist_resolutionEtVsEta_",("#sigma of Reco Et over True Et VS Eta of "+particleString).c_str());
00781   hist_resolutionEtVsPhi_->SetNameTitle("hist_resolutionEtVsPhi_",("#sigma of Reco Et over True Et VS Phi of "+particleString).c_str());
00782 
00783   hist_resolutionEVsEt_->SetNameTitle("hist_resolutionEVsEt_",("#sigma of Reco E over True E VS Et of "+particleString).c_str());
00784   hist_resolutionEVsE_->SetNameTitle("hist_resolutionEVsE_",("#sigma of Reco E over True E VS E of "+particleString).c_str());
00785   hist_resolutionEVsEta_->SetNameTitle("hist_resolutionEVsEta_",("#sigma of Reco E over True E VS Eta of "+particleString).c_str());
00786   hist_resolutionEVsPhi_->SetNameTitle("hist_resolutionEVsPhi_",("#sigma of Reco E over True E VS Phi of "+particleString).c_str());
00787 
00788   hist_resolutionEtaVsEt_->SetNameTitle("hist_resolutionEtaVsEt_",("#sigma of delta Eta VS Et of "+particleString).c_str());
00789   hist_resolutionEtaVsE_->SetNameTitle("hist_resolutionEtaVsE_",("#sigma of delta Eta VS E of "+particleString).c_str());
00790   hist_resolutionEtaVsEta_->SetNameTitle("hist_resolutionEtaVsEta_",("#sigma of delta Eta VS Eta of "+particleString).c_str());
00791   hist_resolutionEtaVsPhi_->SetNameTitle("hist_resolutionEtaVsPhi_",("#sigma of delta Eta VS Phi of "+particleString).c_str());
00792 
00793   hist_resolutionPhiVsEt_->SetNameTitle("hist_resolutionPhiVsEt_",("#sigma of delta Phi VS Et of "+particleString).c_str());
00794   hist_resolutionPhiVsE_->SetNameTitle("hist_resolutionPhiVsE_",("#sigma of delta Phi VS E of "+particleString).c_str());
00795   hist_resolutionPhiVsEta_->SetNameTitle("hist_resolutionPhiVsEta_",("#sigma of delta Phi VS Eta of "+particleString).c_str());
00796   hist_resolutionPhiVsPhi_->SetNameTitle("hist_resolutionPhiVsPhi_",("#sigma of delta Phi VS Phi of "+particleString).c_str());
00797 }
00798 
00799 void EgammaObjects::getEfficiencyHistosViaDividing()
00800 {
00801   hist_EtEfficiency_->Divide(hist_EtEfficiency_,hist_EtNumRecoOverNumTrue_,1,1);
00802   hist_EEfficiency_->Divide(hist_EEfficiency_,hist_ENumRecoOverNumTrue_,1,1);
00803   hist_EtaEfficiency_->Divide(hist_EtaEfficiency_,hist_EtaNumRecoOverNumTrue_,1,1);
00804   hist_PhiEfficiency_->Divide(hist_PhiEfficiency_,hist_PhiNumRecoOverNumTrue_,1,1);
00805 
00806   hist_EtNumRecoOverNumTrue_->Divide(hist_Et_,hist_EtNumRecoOverNumTrue_,1,1);
00807   hist_ENumRecoOverNumTrue_->Divide(hist_E_,hist_ENumRecoOverNumTrue_,1,1);
00808   hist_EtaNumRecoOverNumTrue_->Divide(hist_Eta_,hist_EtaNumRecoOverNumTrue_,1,1);
00809   hist_PhiNumRecoOverNumTrue_->Divide(hist_Phi_,hist_PhiNumRecoOverNumTrue_,1,1);
00810 }
00811 
00812 void EgammaObjects::fitHistos()
00813 {
00814   hist_EtOverTruth_->Fit("gaus","QEM");
00815 //  hist_EtNumRecoOverNumTrue_->Fit("pol1","QEM");
00816 
00817   hist_EOverTruth_->Fit("gaus","QEM");
00818 //  hist_ENumRecoOverNumTrue_->Fit("pol1","QEM");
00819 
00820   hist_EtaOverTruth_->Fit("gaus","QEM");
00821 //  hist_EtaNumRecoOverNumTrue_->Fit("pol1","QEM");
00822 
00823   hist_PhiOverTruth_->Fit("gaus","QEM");
00824 //  hist_PhiNumRecoOverNumTrue_->Fit("pol1","QEM");
00825 
00826   /*
00827   hist_EtOverTruthVsEt_->Fit("pol1","QEM");
00828   hist_EtOverTruthVsEta_->Fit("pol1","QEM");
00829   hist_EtOverTruthVsPhi_->Fit("pol1","QEM");
00830   hist_resolutionEtVsEt_->Fit("pol1","QEM");
00831   hist_resolutionEtVsEta_->Fit("pol1","QEM");
00832   hist_resolutionEtVsPhi_->Fit("pol1","QEM");
00833 
00834   hist_EOverTruthVsEt_->Fit("pol1","QEM");
00835   hist_EOverTruthVsEta_->Fit("pol1","QEM");
00836   hist_EOverTruthVsPhi_->Fit("pol1","QEM");
00837   hist_resolutionEVsEt_->Fit("pol1","QEM");
00838   hist_resolutionEVsEta_->Fit("pol1","QEM");
00839   hist_resolutionEVsPhi_->Fit("pol1","QEM");
00840 
00841   hist_deltaEtaVsEt_->Fit("pol1","QEM");
00842   hist_deltaEtaVsEta_->Fit("pol1","QEM");
00843   hist_deltaEtaVsPhi_->Fit("pol1","QEM");
00844   hist_resolutionEtaVsEt_->Fit("pol1","QEM");
00845   hist_resolutionEtaVsEta_->Fit("pol1","QEM");
00846   hist_resolutionEtaVsPhi_->Fit("pol1","QEM");
00847 
00848   hist_deltaPhiVsEt_->Fit("pol1","QEM");
00849   hist_deltaPhiVsEta_->Fit("pol1","QEM");
00850   hist_deltaPhiVsPhi_->Fit("pol1","QEM");
00851   hist_resolutionPhiVsEt_->Fit("pol1","QEM");
00852   hist_resolutionPhiVsEta_->Fit("pol1","QEM");
00853   hist_resolutionPhiVsPhi_->Fit("pol1","QEM");
00854   */
00855 }
00856 
00857 void EgammaObjects::applyLabels()
00858 {
00859   hist_Et_->GetXaxis()->SetTitle("Et (GeV)");
00860   hist_Et_->GetYaxis()->SetTitle("# per Et Bin");
00861   hist_EtOverTruth_->GetXaxis()->SetTitle("Reco Et/True Et");
00862   hist_EtOverTruth_->GetYaxis()->SetTitle("# per Ratio Bin");
00863   hist_EtEfficiency_->GetXaxis()->SetTitle("Et (GeV)");
00864   hist_EtEfficiency_->GetYaxis()->SetTitle("# True Reconstructed/# True per Et Bin");
00865   hist_EtNumRecoOverNumTrue_->GetXaxis()->SetTitle("Et (GeV)");
00866   hist_EtNumRecoOverNumTrue_->GetYaxis()->SetTitle("# Reco/# True per Et Bin");
00867   hist_EtOverTruthVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00868   hist_EtOverTruthVsEt_->GetYaxis()->SetTitle("Reco Et/True Et per Et Bin");
00869   hist_EtOverTruthVsE_->GetXaxis()->SetTitle("E (GeV)");
00870   hist_EtOverTruthVsE_->GetYaxis()->SetTitle("Reco Et/True Et per E Bin");
00871   hist_EtOverTruthVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00872   hist_EtOverTruthVsEta_->GetYaxis()->SetTitle("Reco Et/True Et per Eta Bin");
00873   hist_EtOverTruthVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00874   hist_EtOverTruthVsPhi_->GetYaxis()->SetTitle("Reco Et/True Et per Phi Bin");
00875   hist_resolutionEtVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00876   hist_resolutionEtVsEt_->GetYaxis()->SetTitle("#sigma of Reco Et/True Et per Et Bin");
00877   hist_resolutionEtVsE_->GetXaxis()->SetTitle("E (GeV)");
00878   hist_resolutionEtVsE_->GetYaxis()->SetTitle("#sigma of Reco Et/True Et per E Bin");
00879   hist_resolutionEtVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00880   hist_resolutionEtVsEta_->GetYaxis()->SetTitle("#sigma of Reco Et/True Et per Eta Bin");
00881   hist_resolutionEtVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00882   hist_resolutionEtVsPhi_->GetYaxis()->SetTitle("#sigma of Reco Et/True Et per Phi Bin");
00883 
00884   hist_E_->GetXaxis()->SetTitle("E (GeV)");
00885   hist_E_->GetYaxis()->SetTitle("# per E Bin");
00886   hist_EOverTruth_->GetXaxis()->SetTitle("Reco E/True E");
00887   hist_EOverTruth_->GetYaxis()->SetTitle("# per Ratio Bin");
00888   hist_EEfficiency_->GetXaxis()->SetTitle("E (GeV)");
00889   hist_EEfficiency_->GetYaxis()->SetTitle("# True Reconstructed/# True per E Bin");
00890   hist_ENumRecoOverNumTrue_->GetXaxis()->SetTitle("E (GeV)");
00891   hist_ENumRecoOverNumTrue_->GetYaxis()->SetTitle("# Reco/# True per E Bin");
00892   hist_EOverTruthVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00893   hist_EOverTruthVsEt_->GetYaxis()->SetTitle("Reco E/True E per Et Bin");
00894   hist_EOverTruthVsE_->GetXaxis()->SetTitle("E (GeV)");
00895   hist_EOverTruthVsE_->GetYaxis()->SetTitle("Reco E/True E per E Bin");
00896   hist_EOverTruthVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00897   hist_EOverTruthVsEta_->GetYaxis()->SetTitle("Reco E/True E per Eta Bin");
00898   hist_EOverTruthVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00899   hist_EOverTruthVsPhi_->GetYaxis()->SetTitle("Reco E/True E per Phi Bin");
00900   hist_resolutionEVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00901   hist_resolutionEVsEt_->GetYaxis()->SetTitle("#sigma of Reco E/True E per Et Bin");
00902   hist_resolutionEVsE_->GetXaxis()->SetTitle("E (GeV)");
00903   hist_resolutionEVsE_->GetYaxis()->SetTitle("#sigma of Reco E/True E per E Bin");
00904   hist_resolutionEVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00905   hist_resolutionEVsEta_->GetYaxis()->SetTitle("#sigma of Reco E/True E per Eta Bin");
00906   hist_resolutionEVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00907   hist_resolutionEVsPhi_->GetYaxis()->SetTitle("#sigma of Reco E/True E per Phi Bin");
00908 
00909   hist_Eta_->GetXaxis()->SetTitle("#eta (Radians)");
00910   hist_Eta_->GetYaxis()->SetTitle("# per Eta Bin");
00911   hist_EtaOverTruth_->GetXaxis()->SetTitle("Reco Eta/True Eta");
00912   hist_EtaOverTruth_->GetYaxis()->SetTitle("# per Ratio Bin");
00913   hist_EtaEfficiency_->GetXaxis()->SetTitle("#eta (Radians)");
00914   hist_EtaEfficiency_->GetYaxis()->SetTitle("# True Reconstructed/# True per Eta Bin");
00915   hist_EtaNumRecoOverNumTrue_->GetXaxis()->SetTitle("#eta (Radians)");
00916   hist_EtaNumRecoOverNumTrue_->GetYaxis()->SetTitle("# Reco/# True per Eta Bin");
00917   hist_deltaEtaVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00918   hist_deltaEtaVsEt_->GetYaxis()->SetTitle("Reco Eta - True Eta per Et Bin");
00919   hist_deltaEtaVsE_->GetXaxis()->SetTitle("E (GeV)");
00920   hist_deltaEtaVsE_->GetYaxis()->SetTitle("Reco Eta - True Eta per E Bin");
00921   hist_deltaEtaVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00922   hist_deltaEtaVsEta_->GetYaxis()->SetTitle("Reco Eta - True Eta per Eta Bin");
00923   hist_deltaEtaVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00924   hist_deltaEtaVsPhi_->GetYaxis()->SetTitle("Reco Eta - True Eta per Phi Bin");
00925   hist_resolutionEtaVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00926   hist_resolutionEtaVsEt_->GetYaxis()->SetTitle("#sigma of Reco Eta - True Eta per Et Bin");
00927   hist_resolutionEtaVsE_->GetXaxis()->SetTitle("E (GeV)");
00928   hist_resolutionEtaVsE_->GetYaxis()->SetTitle("#sigma of Reco Eta - True Eta per E Bin");
00929   hist_resolutionEtaVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00930   hist_resolutionEtaVsEta_->GetYaxis()->SetTitle("#sigma of Reco Eta - True Eta per Eta Bin");
00931   hist_resolutionEtaVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00932   hist_resolutionEtaVsPhi_->GetYaxis()->SetTitle("#sigma of Reco Eta - True Eta per Phi Bin");
00933 
00934   hist_Phi_->GetXaxis()->SetTitle("#phi (Radians)");
00935   hist_Phi_->GetYaxis()->SetTitle("# per Phi Bin");
00936   hist_PhiOverTruth_->GetXaxis()->SetTitle("Reco Phi/True Phi");
00937   hist_PhiOverTruth_->GetYaxis()->SetTitle("# per Ratio Bin");
00938   hist_PhiEfficiency_->GetXaxis()->SetTitle("#phi (Radians)");
00939   hist_PhiEfficiency_->GetYaxis()->SetTitle("# True Reconstructed/# True per Phi Bin");
00940   hist_PhiNumRecoOverNumTrue_->GetXaxis()->SetTitle("#Phi (Radians)");
00941   hist_PhiNumRecoOverNumTrue_->GetYaxis()->SetTitle("# Reco/# True per Phi Bin");
00942   hist_deltaPhiVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00943   hist_deltaPhiVsEt_->GetYaxis()->SetTitle("Reco Phi - True Phi per Et Bin");
00944   hist_deltaPhiVsE_->GetXaxis()->SetTitle("E (GeV)");
00945   hist_deltaPhiVsE_->GetYaxis()->SetTitle("Reco Phi - True Phi per E Bin");
00946   hist_deltaPhiVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00947   hist_deltaPhiVsEta_->GetYaxis()->SetTitle("Reco Phi - True Phi per Eta Bin");
00948   hist_deltaPhiVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00949   hist_deltaPhiVsPhi_->GetYaxis()->SetTitle("Reco Phi - True Phi per Phi Bin");
00950   hist_resolutionPhiVsEt_->GetXaxis()->SetTitle("Et (GeV)");
00951   hist_resolutionPhiVsEt_->GetYaxis()->SetTitle("#sigma of Reco Phi - True Phi per Et Bin");
00952   hist_resolutionPhiVsE_->GetXaxis()->SetTitle("E (GeV)");
00953   hist_resolutionPhiVsE_->GetYaxis()->SetTitle("#sigma of Reco Phi - True Phi per E Bin");
00954   hist_resolutionPhiVsEta_->GetXaxis()->SetTitle("#eta (Radians)");
00955   hist_resolutionPhiVsEta_->GetYaxis()->SetTitle("#sigma of Reco Phi - True Phi per Eta Bin");
00956   hist_resolutionPhiVsPhi_->GetXaxis()->SetTitle("#phi (Radians)");
00957   hist_resolutionPhiVsPhi_->GetYaxis()->SetTitle("#sigma of Reco Phi - True Phi per Phi Bin");
00958 
00959   std::string recoParticleName;
00960 
00961   if( particleID == 22 ) recoParticleName = "Reco Higgs";
00962   else if( particleID == 11 ) recoParticleName = "Reco Z";
00963 
00964   hist_All_recoMass_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00965   hist_All_recoMass_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00966   hist_BarrelOnly_recoMass_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00967   hist_BarrelOnly_recoMass_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00968   hist_EndcapOnly_recoMass_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00969   hist_EndcapOnly_recoMass_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00970   hist_Mixed_recoMass_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00971   hist_Mixed_recoMass_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00972   hist_recoMass_withBackgroud_NoEtCut_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00973   hist_recoMass_withBackgroud_NoEtCut_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00974   hist_recoMass_withBackgroud_5EtCut_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00975   hist_recoMass_withBackgroud_5EtCut_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00976   hist_recoMass_withBackgroud_10EtCut_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00977   hist_recoMass_withBackgroud_10EtCut_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00978   hist_recoMass_withBackgroud_20EtCut_->GetXaxis()->SetTitle((recoParticleName+" Mass (GeV)").c_str());
00979   hist_recoMass_withBackgroud_20EtCut_->GetYaxis()->SetTitle("# of Reco Masses per Mass Bin");
00980 }
00981 
00982 void EgammaObjects::setDrawOptions()
00983 {
00984   hist_Et_->SetOption("e");
00985   hist_EtOverTruth_->SetOption("e");
00986   hist_EtEfficiency_->SetOption("e");
00987   hist_EtNumRecoOverNumTrue_->SetOption("e");
00988   hist_EtOverTruthVsEt_->SetOption("e");
00989   hist_EtOverTruthVsE_->SetOption("e");
00990   hist_EtOverTruthVsEta_->SetOption("e");
00991   hist_EtOverTruthVsPhi_->SetOption("e");
00992   hist_resolutionEtVsEt_->SetOption("e");
00993   hist_resolutionEtVsE_->SetOption("e");
00994   hist_resolutionEtVsEta_->SetOption("e");
00995   hist_resolutionEtVsPhi_->SetOption("e");
00996 
00997   hist_E_->SetOption("e");
00998   hist_EOverTruth_->SetOption("e");
00999   hist_EEfficiency_->SetOption("e");
01000   hist_ENumRecoOverNumTrue_->SetOption("e");
01001   hist_EOverTruthVsEt_->SetOption("e");
01002   hist_EOverTruthVsE_->SetOption("e");
01003   hist_EOverTruthVsEta_->SetOption("e");
01004   hist_EOverTruthVsPhi_->SetOption("e");
01005   hist_resolutionEVsEt_->SetOption("e");
01006   hist_resolutionEVsE_->SetOption("e");
01007   hist_resolutionEVsEta_->SetOption("e");
01008   hist_resolutionEVsPhi_->SetOption("e");
01009 
01010   hist_Eta_->SetOption("e");
01011   hist_EtaOverTruth_->SetOption("e");
01012   hist_EtaEfficiency_->SetOption("e");
01013   hist_EtaNumRecoOverNumTrue_->SetOption("e");
01014   hist_deltaEtaVsEt_->SetOption("e");
01015   hist_deltaEtaVsE_->SetOption("e");
01016   hist_deltaEtaVsEta_->SetOption("e");
01017   hist_deltaEtaVsPhi_->SetOption("e");
01018   hist_resolutionEtaVsEt_->SetOption("e");
01019   hist_resolutionEtaVsE_->SetOption("e");
01020   hist_resolutionEtaVsEta_->SetOption("e");
01021   hist_resolutionEtaVsPhi_->SetOption("e");
01022 
01023   hist_Phi_->SetOption("e");
01024   hist_PhiOverTruth_->SetOption("e");
01025   hist_PhiEfficiency_->SetOption("e");
01026   hist_PhiNumRecoOverNumTrue_->SetOption("e");
01027   hist_deltaPhiVsEt_->SetOption("e");
01028   hist_deltaPhiVsE_->SetOption("e");
01029   hist_deltaPhiVsEta_->SetOption("e");
01030   hist_deltaPhiVsPhi_->SetOption("e");
01031   hist_resolutionPhiVsEt_->SetOption("e");
01032   hist_resolutionPhiVsE_->SetOption("e");
01033   hist_resolutionPhiVsEta_->SetOption("e");
01034   hist_resolutionPhiVsPhi_->SetOption("e");
01035 
01036   hist_All_recoMass_->SetOption("e");
01037   hist_BarrelOnly_recoMass_->SetOption("e");
01038   hist_EndcapOnly_recoMass_->SetOption("e");
01039   hist_Mixed_recoMass_->SetOption("e");
01040   hist_recoMass_withBackgroud_NoEtCut_->SetOption("e");
01041   hist_recoMass_withBackgroud_5EtCut_->SetOption("e");
01042   hist_recoMass_withBackgroud_10EtCut_->SetOption("e");
01043   hist_recoMass_withBackgroud_20EtCut_->SetOption("e");
01044 }
01045 
01046 void EgammaObjects::saveHistos()
01047 {
01048   rootFile_->cd();
01049   rootFile_->GetDirectory(particleString.c_str())->mkdir("ET");
01050   rootFile_->cd(("/"+particleString+"/ET").c_str());
01051 
01052   hist_Et_->Write();
01053   hist_EtOverTruth_->Write();
01054   hist_EtEfficiency_->Write();
01055   hist_EtNumRecoOverNumTrue_->Write();
01056   hist_EtOverTruthVsEt_->Write();
01057   hist_EtOverTruthVsE_->Write();
01058   hist_EtOverTruthVsEta_->Write();
01059   hist_EtOverTruthVsPhi_->Write();
01060   hist_resolutionEtVsEt_->Write();
01061   hist_resolutionEtVsE_->Write();
01062   hist_resolutionEtVsEta_->Write();
01063   hist_resolutionEtVsPhi_->Write();
01064 
01065   rootFile_->cd();
01066   rootFile_->GetDirectory(particleString.c_str())->mkdir("E");
01067   rootFile_->cd(("/"+particleString+"/E").c_str());
01068 
01069   hist_E_->Write();
01070   hist_EOverTruth_->Write();
01071   hist_EEfficiency_->Write();
01072   hist_ENumRecoOverNumTrue_->Write();
01073   hist_EOverTruthVsEt_->Write();
01074   hist_EOverTruthVsE_->Write();
01075   hist_EOverTruthVsEta_->Write();
01076   hist_EOverTruthVsPhi_->Write();
01077   hist_resolutionEVsEt_->Write();
01078   hist_resolutionEVsE_->Write();
01079   hist_resolutionEVsEta_->Write();
01080   hist_resolutionEVsPhi_->Write();
01081 
01082   rootFile_->cd();
01083   rootFile_->GetDirectory(particleString.c_str())->mkdir("Eta");
01084   rootFile_->cd(("/"+particleString+"/Eta").c_str());
01085 
01086   hist_Eta_->Write();
01087   hist_EtaOverTruth_->Write();
01088   hist_EtaEfficiency_->Write();
01089   hist_EtaNumRecoOverNumTrue_->Write();
01090   hist_deltaEtaVsEt_->Write();
01091   hist_deltaEtaVsE_->Write();
01092   hist_deltaEtaVsEta_->Write();
01093   hist_deltaEtaVsPhi_->Write();
01094   hist_resolutionEtaVsEt_->Write();
01095   hist_resolutionEtaVsE_->Write();
01096   hist_resolutionEtaVsEta_->Write();
01097   hist_resolutionEtaVsPhi_->Write();
01098 
01099   rootFile_->cd();
01100   rootFile_->GetDirectory(particleString.c_str())->mkdir("Phi");
01101   rootFile_->cd(("/"+particleString+"/Phi").c_str());
01102 
01103   hist_Phi_->Write();
01104   hist_PhiOverTruth_->Write();
01105   hist_PhiEfficiency_->Write();
01106   hist_PhiNumRecoOverNumTrue_->Write();
01107   hist_deltaPhiVsEt_->Write();
01108   hist_deltaPhiVsE_->Write();
01109   hist_deltaPhiVsEta_->Write();
01110   hist_deltaPhiVsPhi_->Write();
01111   hist_resolutionPhiVsEt_->Write();
01112   hist_resolutionPhiVsE_->Write();
01113   hist_resolutionPhiVsEta_->Write();
01114   hist_resolutionPhiVsPhi_->Write();
01115 
01116   std::string recoParticleName;
01117 
01118   if( particleID == 22 ) recoParticleName = "HiggsRecoMass";
01119   else if( particleID == 11 ) recoParticleName = "ZRecoMass";
01120 
01121   rootFile_->cd();
01122   rootFile_->GetDirectory(particleString.c_str())->mkdir(recoParticleName.c_str());
01123   rootFile_->cd(("/"+particleString+"/"+recoParticleName).c_str());
01124 
01125   hist_All_recoMass_->Write();
01126   hist_BarrelOnly_recoMass_->Write();
01127   hist_EndcapOnly_recoMass_->Write();
01128   hist_Mixed_recoMass_->Write();
01129   hist_recoMass_withBackgroud_NoEtCut_->Write();
01130   hist_recoMass_withBackgroud_5EtCut_->Write();
01131   hist_recoMass_withBackgroud_10EtCut_->Write();
01132   hist_recoMass_withBackgroud_20EtCut_->Write();
01133 
01134   rootFile_->cd();
01135   rootFile_->GetDirectory(particleString.c_str())->mkdir("_TempScatterPlots");
01136   rootFile_->cd(("/"+particleString+"/_TempScatterPlots").c_str());
01137 
01138   _TEMP_scatterPlot_EtOverTruthVsEt_->Write();
01139   _TEMP_scatterPlot_EtOverTruthVsE_->Write();
01140   _TEMP_scatterPlot_EtOverTruthVsEta_->Write();
01141   _TEMP_scatterPlot_EtOverTruthVsPhi_->Write();
01142 
01143   _TEMP_scatterPlot_EOverTruthVsEt_->Write();
01144   _TEMP_scatterPlot_EOverTruthVsE_->Write();
01145   _TEMP_scatterPlot_EOverTruthVsEta_->Write();
01146   _TEMP_scatterPlot_EOverTruthVsPhi_->Write();
01147 
01148   _TEMP_scatterPlot_deltaEtaVsEt_->Write();
01149   _TEMP_scatterPlot_deltaEtaVsE_->Write();
01150   _TEMP_scatterPlot_deltaEtaVsEta_->Write();
01151   _TEMP_scatterPlot_deltaEtaVsPhi_->Write();
01152 
01153   _TEMP_scatterPlot_deltaPhiVsEt_->Write();
01154   _TEMP_scatterPlot_deltaPhiVsE_->Write();
01155   _TEMP_scatterPlot_deltaPhiVsEta_->Write();
01156   _TEMP_scatterPlot_deltaPhiVsPhi_->Write();
01157 
01158   rootFile_->cd();
01159 }
01160