00001
00008
00009
00010
00011
00012
00013
00014 #include "RecoEcal/EgammaClusterProducers/interface/EgammaSimpleAnalyzer.h"
00015
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/Utilities/interface/Exception.h"
00018 #include "DataFormats/Common/interface/Handle.h"
00019
00020 #include "TFile.h"
00021 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00022 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00023 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00024 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00025 #include "DataFormats/EgammaReco/interface/ClusterShape.h"
00026
00027
00028
00029 EgammaSimpleAnalyzer::EgammaSimpleAnalyzer( const edm::ParameterSet& ps )
00030
00031 {
00032
00033 xMinHist_ = ps.getParameter<double>("xMinHist");
00034 xMaxHist_ = ps.getParameter<double>("xMaxHist");
00035 nbinHist_ = ps.getParameter<int>("nbinHist");
00036
00037 islandBarrelBasicClusterCollection_ = ps.getParameter<std::string>("islandBarrelBasicClusterCollection");
00038 islandBarrelBasicClusterProducer_ = ps.getParameter<std::string>("islandBarrelBasicClusterProducer");
00039 islandBarrelBasicClusterShapes_ = ps.getParameter<std::string>("islandBarrelBasicClusterShapes");
00040
00041 islandEndcapBasicClusterCollection_ = ps.getParameter<std::string>("islandEndcapBasicClusterCollection");
00042 islandEndcapBasicClusterProducer_ = ps.getParameter<std::string>("islandEndcapBasicClusterProducer");
00043 islandEndcapBasicClusterShapes_ = ps.getParameter<std::string>("islandEndcapBasicClusterShapes");
00044
00045 islandEndcapSuperClusterCollection_ = ps.getParameter<std::string>("islandEndcapSuperClusterCollection");
00046 islandEndcapSuperClusterProducer_ = ps.getParameter<std::string>("islandEndcapSuperClusterProducer");
00047
00048 correctedIslandEndcapSuperClusterCollection_ = ps.getParameter<std::string>("correctedIslandEndcapSuperClusterCollection");
00049 correctedIslandEndcapSuperClusterProducer_ = ps.getParameter<std::string>("correctedIslandEndcapSuperClusterProducer");
00050
00051 hybridSuperClusterCollection_ = ps.getParameter<std::string>("hybridSuperClusterCollection");
00052 hybridSuperClusterProducer_ = ps.getParameter<std::string>("hybridSuperClusterProducer");
00053
00054 correctedHybridSuperClusterCollection_ = ps.getParameter<std::string>("correctedHybridSuperClusterCollection");
00055 correctedHybridSuperClusterProducer_ = ps.getParameter<std::string>("correctedHybridSuperClusterProducer");
00056
00057 outputFile_ = ps.getParameter<std::string>("outputFile");
00058 rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
00059
00060 }
00061
00062
00063
00064 EgammaSimpleAnalyzer::~EgammaSimpleAnalyzer()
00065
00066 {
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080 delete rootFile_;
00081 }
00082
00083
00084 void
00085 EgammaSimpleAnalyzer::beginJob() {
00086
00087
00088
00089 rootFile_->cd();
00090
00091 h1_nIslandEBBC_ = new TH1F("nIslandEBBC","# basic clusters with island in barrel",11,-0.5,10.5);
00092 h1_nIslandEEBC_ = new TH1F("nIslandEEBC","# basic clusters with island in endcap",11,-0.5,10.5);
00093
00094 h1_nIslandEESC_ = new TH1F("nIslandEESC","# super clusters with island in endcap",11,-0.5,10.5);
00095 h1_nHybridSC_ = new TH1F("nHybridSC","# super clusters with hybrid",11,-0.5,10.5);
00096
00097 h1_islandEBBCEnergy_ = new TH1F("islandEBBCEnergy","Energy of basic clusters with island algo - barrel",nbinHist_,xMinHist_,xMaxHist_);
00098 h1_islandEBBCXtals_ = new TH1F("islandEBBCXtals","#xtals in basic cluster - island barrel",51,-0.5,50.5);
00099
00100 h1_islandEBBCe9over25_= new TH1F("islandEBBCe9over25","e3x3/e5x5 of basic clusters with island algo - barrel",35,0.5,1.2);
00101 h1_islandEBBCe5x5_ = new TH1F("islandEBBCe5x5","e5x5 of basic clusters with island algo - barrel",nbinHist_,xMinHist_,xMaxHist_);
00102 h1_islandEEBCe5x5_ = new TH1F("islandEEBCe5x5","e5x5 of basic clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_);
00103 h1_islandEEBCEnergy_ = new TH1F("islandEEBCEnergy","Energy of basic clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_);
00104 h1_islandEEBCXtals_ = new TH1F("islandEEBCXtals","#xtals in basic cluster - island endcap",51,-0.5,50.5);
00105
00106 h1_islandEESCEnergy_ = new TH1F("islandEESCEnergy","Energy of super clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_);
00107 h1_corrIslandEESCEnergy_ = new TH1F("corrIslandEESCEnergy","Corrected Energy of super clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_);
00108 h1_corrIslandEESCET_ = new TH1F("corrIslandEESCET","Corrected Transverse Energy of super clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_);
00109 h1_islandEESCClusters_ = new TH1F("islandEESCClusters","# basic clusters in super cluster - island endcap",11,-0.5,10.5);
00110
00111 h1_hybridSCEnergy_ = new TH1F("hybridSCEnergy","Energy of super clusters with hybrid algo",nbinHist_,xMinHist_,xMaxHist_);
00112 h1_corrHybridSCEnergy_ = new TH1F("corrHybridSCEnergy","Corrected Energy of super clusters with hybrid algo",nbinHist_,xMinHist_,xMaxHist_);
00113 h1_corrHybridSCET_ = new TH1F("corrHybridSCET","Corrected Transverse Energy of super clusters with hybrid algo",nbinHist_,xMinHist_,xMaxHist_);
00114 h1_corrHybridSCEta_ = new TH1F("corrHybridSCEta","Eta of super clusters with hybrid algo",40,-3.,3.);
00115 h1_corrHybridSCPhi_ = new TH1F("corrHybridSCPhi","Phi of super clusters with hybrid algo",40,0.,6.28);
00116 h1_hybridSCClusters_ = new TH1F("hybridSCClusters","# basic clusters in super cluster - hybrid",11,-0.5,10.5);
00117
00118 }
00119
00120
00121
00122 void
00123 EgammaSimpleAnalyzer::analyze( const edm::Event& evt, const edm::EventSetup& es ) {
00124
00125
00126 using namespace edm;
00127
00128
00129
00130
00131
00132 Handle<reco::BasicClusterCollection> pIslandBarrelBasicClusters;
00133 evt.getByLabel(islandBarrelBasicClusterProducer_, islandBarrelBasicClusterCollection_, pIslandBarrelBasicClusters);
00134 const reco::BasicClusterCollection* islandBarrelBasicClusters = pIslandBarrelBasicClusters.product();
00135 h1_nIslandEBBC_->Fill(islandBarrelBasicClusters->size());
00136
00137
00138 Handle<reco::ClusterShapeCollection> pIslandEBShapes;
00139 evt.getByLabel(islandBarrelBasicClusterProducer_, islandBarrelBasicClusterShapes_, pIslandEBShapes);
00140 const reco::ClusterShapeCollection* islandEBShapes = pIslandEBShapes.product();
00141
00142 std::ostringstream str;
00143 str << "# island basic clusters in barrel: " << islandBarrelBasicClusters->size()
00144 << "\t# associated cluster shapes: " << islandEBShapes->size() << "\n"
00145 << "Loop over island basic clusters in barrel" << "\n";
00146
00147
00148 int iClus=0;
00149 for(reco::BasicClusterCollection::const_iterator aClus = islandBarrelBasicClusters->begin();
00150 aClus != islandBarrelBasicClusters->end(); aClus++) {
00151 h1_islandEBBCEnergy_->Fill( aClus->energy() );
00152 h1_islandEBBCXtals_->Fill( aClus->size() );
00153 str << "energy: " << aClus->energy()
00154 << "\te5x5: " << (*islandEBShapes)[iClus].e5x5()
00155 << "\te2x2: " << (*islandEBShapes)[iClus].e2x2()
00156 << "\n";
00157 h1_islandEBBCe5x5_->Fill( (*islandEBShapes)[iClus].e5x5() );
00158
00159 iClus++;
00160 }
00161 edm::LogInfo("EgammaSimpleAnalyzer") << str.str();
00162
00163
00164
00165
00166
00167
00168 Handle<reco::BasicClusterCollection> pIslandEndcapBasicClusters;
00169 evt.getByLabel(islandEndcapBasicClusterProducer_, islandEndcapBasicClusterCollection_, pIslandEndcapBasicClusters);
00170 const reco::BasicClusterCollection* islandEndcapBasicClusters = pIslandEndcapBasicClusters.product();
00171 h1_nIslandEEBC_->Fill(islandEndcapBasicClusters->size());
00172
00173
00174 Handle<reco::ClusterShapeCollection> pIslandEEShapes;
00175 evt.getByLabel(islandEndcapBasicClusterProducer_, islandEndcapBasicClusterShapes_, pIslandEEShapes);
00176 const reco::ClusterShapeCollection* islandEEShapes = pIslandEEShapes.product();
00177
00178
00179 iClus=0;
00180 for(reco::BasicClusterCollection::const_iterator aClus = islandEndcapBasicClusters->begin();
00181 aClus != islandEndcapBasicClusters->end(); aClus++) {
00182 h1_islandEEBCEnergy_->Fill( aClus->energy() );
00183 h1_islandEEBCXtals_->Fill( aClus->size() );
00184 h1_islandEEBCe5x5_->Fill( (*islandEEShapes)[iClus].e5x5() );
00185 h1_islandEBBCe9over25_->Fill( (*islandEEShapes)[iClus].e3x3()/(*islandEEShapes)[iClus].e5x5() );
00186 iClus++;
00187 }
00188 edm::LogInfo("EgammaSimpleAnalyzer") << str.str();
00189
00190
00191 Handle<reco::SuperClusterCollection> pIslandEndcapSuperClusters;
00192 evt.getByLabel(islandEndcapSuperClusterProducer_, islandEndcapSuperClusterCollection_, pIslandEndcapSuperClusters);
00193 const reco::SuperClusterCollection* islandEndcapSuperClusters = pIslandEndcapSuperClusters.product();
00194
00195
00196 for(reco::SuperClusterCollection::const_iterator aClus = islandEndcapSuperClusters->begin();
00197 aClus != islandEndcapSuperClusters->end(); aClus++) {
00198 h1_islandEESCEnergy_->Fill( aClus->energy() );
00199 }
00200
00201
00202
00203 Handle<reco::SuperClusterCollection> pCorrectedIslandEndcapSuperClusters;
00204 evt.getByLabel(correctedIslandEndcapSuperClusterProducer_, correctedIslandEndcapSuperClusterCollection_, pCorrectedIslandEndcapSuperClusters);
00205 const reco::SuperClusterCollection* correctedIslandEndcapSuperClusters = pCorrectedIslandEndcapSuperClusters.product();
00206 h1_nIslandEESC_->Fill(islandEndcapSuperClusters->size());
00207
00208
00209 for(reco::SuperClusterCollection::const_iterator aClus = correctedIslandEndcapSuperClusters->begin();
00210 aClus != correctedIslandEndcapSuperClusters->end(); aClus++) {
00211 h1_corrIslandEESCEnergy_->Fill( aClus->energy() );
00212 h1_corrIslandEESCET_->Fill( aClus->energy()*sin(aClus->position().theta()) );
00213 h1_islandEESCClusters_->Fill( aClus->clustersSize() );
00214 }
00215
00216
00217
00218
00219
00220
00221
00222 Handle<reco::SuperClusterCollection> pHybridSuperClusters;
00223 evt.getByLabel(hybridSuperClusterProducer_, hybridSuperClusterCollection_, pHybridSuperClusters);
00224 const reco::SuperClusterCollection* hybridSuperClusters = pHybridSuperClusters.product();
00225
00226
00227 for(reco::SuperClusterCollection::const_iterator aClus = hybridSuperClusters->begin();
00228 aClus != hybridSuperClusters->end(); aClus++) {
00229 h1_hybridSCEnergy_->Fill( aClus->energy() );
00230 }
00231
00232
00233
00234 Handle<reco::SuperClusterCollection> pCorrectedHybridSuperClusters;
00235 evt.getByLabel(correctedHybridSuperClusterProducer_, correctedHybridSuperClusterCollection_, pCorrectedHybridSuperClusters);
00236 const reco::SuperClusterCollection* correctedHybridSuperClusters = pCorrectedHybridSuperClusters.product();
00237 h1_nHybridSC_->Fill(correctedHybridSuperClusters->size());
00238
00239
00240
00241 for(reco::SuperClusterCollection::const_iterator aClus = correctedHybridSuperClusters->begin();
00242 aClus != correctedHybridSuperClusters->end(); aClus++) {
00243 h1_hybridSCClusters_->Fill( aClus->clustersSize() );
00244 h1_corrHybridSCEnergy_->Fill( aClus->energy() );
00245 h1_corrHybridSCET_->Fill( aClus->energy()*sin(aClus->position().theta()) );
00246 h1_corrHybridSCEta_->Fill( aClus->position().eta() );
00247 h1_corrHybridSCPhi_->Fill( aClus->position().phi() );
00248 }
00249
00250 }
00251
00252
00253 void
00254 EgammaSimpleAnalyzer::endJob() {
00255
00256
00257
00258 rootFile_->cd();
00259
00260 h1_nIslandEBBC_->Write();
00261 h1_nIslandEEBC_->Write();
00262 h1_nIslandEESC_->Write();
00263 h1_nHybridSC_->Write();
00264
00265 h1_islandEBBCe9over25_->Write();
00266 h1_islandEBBCe5x5_->Write();
00267 h1_islandEBBCEnergy_->Write();
00268 h1_islandEBBCXtals_->Write();
00269
00270 h1_islandEEBCe5x5_->Write();
00271 h1_islandEEBCEnergy_->Write();
00272 h1_islandEEBCXtals_->Write();
00273
00274 h1_islandEESCEnergy_->Write();
00275 h1_corrIslandEESCEnergy_->Write();
00276 h1_corrIslandEESCET_->Write();
00277 h1_islandEESCClusters_->Write();
00278
00279 h1_hybridSCClusters_->Write();
00280 h1_hybridSCEnergy_->Write();
00281 h1_corrHybridSCEnergy_->Write();
00282 h1_corrHybridSCET_->Write();
00283 h1_corrHybridSCEta_->Write();
00284 h1_corrHybridSCPhi_->Write();
00285
00286 rootFile_->Close();
00287 }