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