![]() |
![]() |
#include <RecoEcal/EgammaClusterProducers/interface/EgammaSimpleAnalyzer.h>
Description: <one line="" class="" summary>="">.
Implementation: \
Implementation: <Notes on="" implementation>="">Definition at line 35 of file EgammaSimpleAnalyzer.h.
EgammaSimpleAnalyzer::EgammaSimpleAnalyzer | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 28 of file EgammaSimpleAnalyzer.cc.
References correctedHybridSuperClusterCollection_, correctedHybridSuperClusterProducer_, correctedIslandEndcapSuperClusterCollection_, correctedIslandEndcapSuperClusterProducer_, edm::ParameterSet::getParameter(), hybridSuperClusterCollection_, hybridSuperClusterProducer_, islandBarrelBasicClusterCollection_, islandBarrelBasicClusterProducer_, islandBarrelBasicClusterShapes_, islandEndcapBasicClusterCollection_, islandEndcapBasicClusterProducer_, islandEndcapBasicClusterShapes_, islandEndcapSuperClusterCollection_, islandEndcapSuperClusterProducer_, nbinHist_, outputFile_, rootFile_, xMaxHist_, and xMinHist_.
00030 { 00031 00032 xMinHist_ = ps.getParameter<double>("xMinHist"); 00033 xMaxHist_ = ps.getParameter<double>("xMaxHist"); 00034 nbinHist_ = ps.getParameter<int>("nbinHist"); 00035 00036 islandBarrelBasicClusterCollection_ = ps.getParameter<std::string>("islandBarrelBasicClusterCollection"); 00037 islandBarrelBasicClusterProducer_ = ps.getParameter<std::string>("islandBarrelBasicClusterProducer"); 00038 islandBarrelBasicClusterShapes_ = ps.getParameter<std::string>("islandBarrelBasicClusterShapes"); 00039 00040 islandEndcapBasicClusterCollection_ = ps.getParameter<std::string>("islandEndcapBasicClusterCollection"); 00041 islandEndcapBasicClusterProducer_ = ps.getParameter<std::string>("islandEndcapBasicClusterProducer"); 00042 islandEndcapBasicClusterShapes_ = ps.getParameter<std::string>("islandEndcapBasicClusterShapes"); 00043 00044 islandEndcapSuperClusterCollection_ = ps.getParameter<std::string>("islandEndcapSuperClusterCollection"); 00045 islandEndcapSuperClusterProducer_ = ps.getParameter<std::string>("islandEndcapSuperClusterProducer"); 00046 00047 correctedIslandEndcapSuperClusterCollection_ = ps.getParameter<std::string>("correctedIslandEndcapSuperClusterCollection"); 00048 correctedIslandEndcapSuperClusterProducer_ = ps.getParameter<std::string>("correctedIslandEndcapSuperClusterProducer"); 00049 00050 hybridSuperClusterCollection_ = ps.getParameter<std::string>("hybridSuperClusterCollection"); 00051 hybridSuperClusterProducer_ = ps.getParameter<std::string>("hybridSuperClusterProducer"); 00052 00053 correctedHybridSuperClusterCollection_ = ps.getParameter<std::string>("correctedHybridSuperClusterCollection"); 00054 correctedHybridSuperClusterProducer_ = ps.getParameter<std::string>("correctedHybridSuperClusterProducer"); 00055 00056 outputFile_ = ps.getParameter<std::string>("outputFile"); 00057 rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE"); // open output file to store histograms 00058 00059 }
EgammaSimpleAnalyzer::~EgammaSimpleAnalyzer | ( | ) |
Definition at line 63 of file EgammaSimpleAnalyzer.cc.
References rootFile_.
00065 { 00066 00067 // apparently ROOT takes ownership of histograms 00068 // created after opening a new TFile... no delete is needed 00069 // ... mysteries of root... 00070 /* 00071 delete h1_islandBCEnergy_; 00072 delete h1_islandSCEnergy_; 00073 delete h1_corrIslandSCEnergy_; 00074 delete h1_hybridSCEnergy_; 00075 delete h1_corrHybridSCEnergy_; 00076 delete h1_corrHybridSCEta_; 00077 delete h1_corrHybridSCPhi_; 00078 */ 00079 delete rootFile_; 00080 }
void EgammaSimpleAnalyzer::analyze | ( | const edm::Event & | evt, | |
const edm::EventSetup & | es | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 122 of file EgammaSimpleAnalyzer.cc.
References correctedHybridSuperClusterCollection_, correctedHybridSuperClusterProducer_, correctedHybridSuperClusters_cfi::correctedHybridSuperClusters, correctedIslandEndcapSuperClusterCollection_, correctedIslandEndcapSuperClusterProducer_, correctedIslandEndcapSuperClusters_cfi::correctedIslandEndcapSuperClusters, edm::Event::getByLabel(), h1_corrHybridSCEnergy_, h1_corrHybridSCET_, h1_corrHybridSCEta_, h1_corrHybridSCPhi_, h1_corrIslandEESCEnergy_, h1_corrIslandEESCET_, h1_hybridSCClusters_, h1_hybridSCEnergy_, h1_islandEBBCe5x5_, h1_islandEBBCe9over25_, h1_islandEBBCEnergy_, h1_islandEBBCXtals_, h1_islandEEBCe5x5_, h1_islandEEBCEnergy_, h1_islandEEBCXtals_, h1_islandEESCClusters_, h1_islandEESCEnergy_, h1_nHybridSC_, h1_nIslandEBBC_, h1_nIslandEEBC_, h1_nIslandEESC_, hybridSuperClusterCollection_, hybridSuperClusterProducer_, hybridSuperClusters_cfi::hybridSuperClusters, islandBarrelBasicClusterCollection_, islandBarrelBasicClusterProducer_, islandBarrelBasicClusterShapes_, islandEndcapBasicClusterCollection_, islandEndcapBasicClusterProducer_, islandEndcapBasicClusterShapes_, islandEndcapSuperClusterCollection_, islandEndcapSuperClusterProducer_, and funct::sin().
00122 { 00123 //======================================================================== 00124 00125 using namespace edm; // needed for all fwk related classes 00126 00127 00128 // ----- barrel with island 00129 00130 // Get island basic clusters 00131 Handle<reco::BasicClusterCollection> pIslandBarrelBasicClusters; 00132 evt.getByLabel(islandBarrelBasicClusterProducer_, islandBarrelBasicClusterCollection_, pIslandBarrelBasicClusters); 00133 const reco::BasicClusterCollection* islandBarrelBasicClusters = pIslandBarrelBasicClusters.product(); 00134 h1_nIslandEBBC_->Fill(islandBarrelBasicClusters->size()); 00135 00136 // fetch cluster shapes of island basic clusters in barrel 00137 Handle<reco::ClusterShapeCollection> pIslandEBShapes; 00138 evt.getByLabel(islandBarrelBasicClusterProducer_, islandBarrelBasicClusterShapes_, pIslandEBShapes); 00139 const reco::ClusterShapeCollection* islandEBShapes = pIslandEBShapes.product(); 00140 00141 std::ostringstream str; 00142 str << "# island basic clusters in barrel: " << islandBarrelBasicClusters->size() 00143 << "\t# associated cluster shapes: " << islandEBShapes->size() << "\n" 00144 << "Loop over island basic clusters in barrel" << "\n"; 00145 00146 // loop over the Basic clusters and fill the histogram 00147 int iClus=0; 00148 for(reco::BasicClusterCollection::const_iterator aClus = islandBarrelBasicClusters->begin(); 00149 aClus != islandBarrelBasicClusters->end(); aClus++) { 00150 h1_islandEBBCEnergy_->Fill( aClus->energy() ); 00151 h1_islandEBBCXtals_->Fill( aClus->getHitsByDetId().size() ); 00152 str << "energy: " << aClus->energy() 00153 << "\te5x5: " << (*islandEBShapes)[iClus].e5x5() 00154 << "\te2x2: " << (*islandEBShapes)[iClus].e2x2() 00155 << "\n"; 00156 h1_islandEBBCe5x5_->Fill( (*islandEBShapes)[iClus].e5x5() ); 00157 00158 iClus++; 00159 } 00160 edm::LogInfo("EgammaSimpleAnalyzer") << str.str(); 00161 00162 // extract energy corrections applied 00163 00164 // ---- island in endcap 00165 00166 // Get island basic clusters 00167 Handle<reco::BasicClusterCollection> pIslandEndcapBasicClusters; 00168 evt.getByLabel(islandEndcapBasicClusterProducer_, islandEndcapBasicClusterCollection_, pIslandEndcapBasicClusters); 00169 const reco::BasicClusterCollection* islandEndcapBasicClusters = pIslandEndcapBasicClusters.product(); 00170 h1_nIslandEEBC_->Fill(islandEndcapBasicClusters->size()); 00171 00172 // fetch cluster shapes of island basic clusters in endcap 00173 Handle<reco::ClusterShapeCollection> pIslandEEShapes; 00174 evt.getByLabel(islandEndcapBasicClusterProducer_, islandEndcapBasicClusterShapes_, pIslandEEShapes); 00175 const reco::ClusterShapeCollection* islandEEShapes = pIslandEEShapes.product(); 00176 00177 // loop over the Basic clusters and fill the histogram 00178 iClus=0; 00179 for(reco::BasicClusterCollection::const_iterator aClus = islandEndcapBasicClusters->begin(); 00180 aClus != islandEndcapBasicClusters->end(); aClus++) { 00181 h1_islandEEBCEnergy_->Fill( aClus->energy() ); 00182 h1_islandEEBCXtals_->Fill( aClus->getHitsByDetId().size() ); 00183 h1_islandEEBCe5x5_->Fill( (*islandEEShapes)[iClus].e5x5() ); 00184 h1_islandEBBCe9over25_->Fill( (*islandEEShapes)[iClus].e3x3()/(*islandEEShapes)[iClus].e5x5() ); 00185 iClus++; 00186 } 00187 edm::LogInfo("EgammaSimpleAnalyzer") << str.str(); 00188 00189 // Get island super clusters 00190 Handle<reco::SuperClusterCollection> pIslandEndcapSuperClusters; 00191 evt.getByLabel(islandEndcapSuperClusterProducer_, islandEndcapSuperClusterCollection_, pIslandEndcapSuperClusters); 00192 const reco::SuperClusterCollection* islandEndcapSuperClusters = pIslandEndcapSuperClusters.product(); 00193 00194 // loop over the super clusters and fill the histogram 00195 for(reco::SuperClusterCollection::const_iterator aClus = islandEndcapSuperClusters->begin(); 00196 aClus != islandEndcapSuperClusters->end(); aClus++) { 00197 h1_islandEESCEnergy_->Fill( aClus->energy() ); 00198 } 00199 00200 00201 // Get island super clusters after energy correction 00202 Handle<reco::SuperClusterCollection> pCorrectedIslandEndcapSuperClusters; 00203 evt.getByLabel(correctedIslandEndcapSuperClusterProducer_, correctedIslandEndcapSuperClusterCollection_, pCorrectedIslandEndcapSuperClusters); 00204 const reco::SuperClusterCollection* correctedIslandEndcapSuperClusters = pCorrectedIslandEndcapSuperClusters.product(); 00205 h1_nIslandEESC_->Fill(islandEndcapSuperClusters->size()); 00206 00207 // loop over the super clusters and fill the histogram 00208 for(reco::SuperClusterCollection::const_iterator aClus = correctedIslandEndcapSuperClusters->begin(); 00209 aClus != correctedIslandEndcapSuperClusters->end(); aClus++) { 00210 h1_corrIslandEESCEnergy_->Fill( aClus->energy() ); 00211 h1_corrIslandEESCET_->Fill( aClus->energy()*sin(aClus->position().theta()) ); 00212 h1_islandEESCClusters_->Fill( aClus->clustersSize() ); 00213 } 00214 00215 // extract energy corrections applied 00216 00217 00218 // ----- hybrid 00219 00220 // Get hybrid super clusters 00221 Handle<reco::SuperClusterCollection> pHybridSuperClusters; 00222 evt.getByLabel(hybridSuperClusterProducer_, hybridSuperClusterCollection_, pHybridSuperClusters); 00223 const reco::SuperClusterCollection* hybridSuperClusters = pHybridSuperClusters.product(); 00224 00225 // loop over the super clusters and fill the histogram 00226 for(reco::SuperClusterCollection::const_iterator aClus = hybridSuperClusters->begin(); 00227 aClus != hybridSuperClusters->end(); aClus++) { 00228 h1_hybridSCEnergy_->Fill( aClus->energy() ); 00229 } 00230 00231 00232 // Get hybrid super clusters after energy correction 00233 Handle<reco::SuperClusterCollection> pCorrectedHybridSuperClusters; 00234 evt.getByLabel(correctedHybridSuperClusterProducer_, correctedHybridSuperClusterCollection_, pCorrectedHybridSuperClusters); 00235 const reco::SuperClusterCollection* correctedHybridSuperClusters = pCorrectedHybridSuperClusters.product(); 00236 h1_nHybridSC_->Fill(correctedHybridSuperClusters->size()); 00237 00238 00239 // loop over the super clusters and fill the histogram 00240 for(reco::SuperClusterCollection::const_iterator aClus = correctedHybridSuperClusters->begin(); 00241 aClus != correctedHybridSuperClusters->end(); aClus++) { 00242 h1_hybridSCClusters_->Fill( aClus->clustersSize() ); 00243 h1_corrHybridSCEnergy_->Fill( aClus->energy() ); 00244 h1_corrHybridSCET_->Fill( aClus->energy()*sin(aClus->position().theta()) ); 00245 h1_corrHybridSCEta_->Fill( aClus->position().eta() ); 00246 h1_corrHybridSCPhi_->Fill( aClus->position().phi() ); 00247 } 00248 00249 }
void EgammaSimpleAnalyzer::beginJob | ( | edm::EventSetup const & | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 84 of file EgammaSimpleAnalyzer.cc.
References h1_corrHybridSCEnergy_, h1_corrHybridSCET_, h1_corrHybridSCEta_, h1_corrHybridSCPhi_, h1_corrIslandEESCEnergy_, h1_corrIslandEESCET_, h1_hybridSCClusters_, h1_hybridSCEnergy_, h1_islandEBBCe5x5_, h1_islandEBBCe9over25_, h1_islandEBBCEnergy_, h1_islandEBBCXtals_, h1_islandEEBCe5x5_, h1_islandEEBCEnergy_, h1_islandEEBCXtals_, h1_islandEESCClusters_, h1_islandEESCEnergy_, h1_nHybridSC_, h1_nIslandEBBC_, h1_nIslandEEBC_, h1_nIslandEESC_, nbinHist_, rootFile_, xMaxHist_, and xMinHist_.
00084 { 00085 //======================================================================== 00086 00087 // go to *OUR* rootfile and book histograms 00088 rootFile_->cd(); 00089 00090 h1_nIslandEBBC_ = new TH1F("nIslandEBBC","# basic clusters with island in barrel",11,-0.5,10.5); 00091 h1_nIslandEEBC_ = new TH1F("nIslandEEBC","# basic clusters with island in endcap",11,-0.5,10.5); 00092 00093 h1_nIslandEESC_ = new TH1F("nIslandEESC","# super clusters with island in endcap",11,-0.5,10.5); 00094 h1_nHybridSC_ = new TH1F("nHybridSC","# super clusters with hybrid",11,-0.5,10.5); 00095 00096 h1_islandEBBCEnergy_ = new TH1F("islandEBBCEnergy","Energy of basic clusters with island algo - barrel",nbinHist_,xMinHist_,xMaxHist_); 00097 h1_islandEBBCXtals_ = new TH1F("islandEBBCXtals","#xtals in basic cluster - island barrel",51,-0.5,50.5); 00098 00099 h1_islandEBBCe9over25_= new TH1F("islandEBBCe9over25","e3x3/e5x5 of basic clusters with island algo - barrel",35,0.5,1.2); 00100 h1_islandEBBCe5x5_ = new TH1F("islandEBBCe5x5","e5x5 of basic clusters with island algo - barrel",nbinHist_,xMinHist_,xMaxHist_); 00101 h1_islandEEBCe5x5_ = new TH1F("islandEEBCe5x5","e5x5 of basic clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_); 00102 h1_islandEEBCEnergy_ = new TH1F("islandEEBCEnergy","Energy of basic clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_); 00103 h1_islandEEBCXtals_ = new TH1F("islandEEBCXtals","#xtals in basic cluster - island endcap",51,-0.5,50.5); 00104 00105 h1_islandEESCEnergy_ = new TH1F("islandEESCEnergy","Energy of super clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_); 00106 h1_corrIslandEESCEnergy_ = new TH1F("corrIslandEESCEnergy","Corrected Energy of super clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_); 00107 h1_corrIslandEESCET_ = new TH1F("corrIslandEESCET","Corrected Transverse Energy of super clusters with island algo - endcap",nbinHist_,xMinHist_,xMaxHist_); 00108 h1_islandEESCClusters_ = new TH1F("islandEESCClusters","# basic clusters in super cluster - island endcap",11,-0.5,10.5); 00109 00110 h1_hybridSCEnergy_ = new TH1F("hybridSCEnergy","Energy of super clusters with hybrid algo",nbinHist_,xMinHist_,xMaxHist_); 00111 h1_corrHybridSCEnergy_ = new TH1F("corrHybridSCEnergy","Corrected Energy of super clusters with hybrid algo",nbinHist_,xMinHist_,xMaxHist_); 00112 h1_corrHybridSCET_ = new TH1F("corrHybridSCET","Corrected Transverse Energy of super clusters with hybrid algo",nbinHist_,xMinHist_,xMaxHist_); 00113 h1_corrHybridSCEta_ = new TH1F("corrHybridSCEta","Eta of super clusters with hybrid algo",40,-3.,3.); 00114 h1_corrHybridSCPhi_ = new TH1F("corrHybridSCPhi","Phi of super clusters with hybrid algo",40,0.,6.28); 00115 h1_hybridSCClusters_ = new TH1F("hybridSCClusters","# basic clusters in super cluster - hybrid",11,-0.5,10.5); 00116 00117 }
Reimplemented from edm::EDAnalyzer.
Definition at line 253 of file EgammaSimpleAnalyzer.cc.
References h1_corrHybridSCEnergy_, h1_corrHybridSCET_, h1_corrHybridSCEta_, h1_corrHybridSCPhi_, h1_corrIslandEESCEnergy_, h1_corrIslandEESCET_, h1_hybridSCClusters_, h1_hybridSCEnergy_, h1_islandEBBCe5x5_, h1_islandEBBCe9over25_, h1_islandEBBCEnergy_, h1_islandEBBCXtals_, h1_islandEEBCe5x5_, h1_islandEEBCEnergy_, h1_islandEEBCXtals_, h1_islandEESCClusters_, h1_islandEESCEnergy_, h1_nHybridSC_, h1_nIslandEBBC_, h1_nIslandEEBC_, h1_nIslandEESC_, and rootFile_.
00253 { 00254 //======================================================================== 00255 00256 // go to *OUR* root file and store histograms 00257 rootFile_->cd(); 00258 00259 h1_nIslandEBBC_->Write(); 00260 h1_nIslandEEBC_->Write(); 00261 h1_nIslandEESC_->Write(); 00262 h1_nHybridSC_->Write(); 00263 00264 h1_islandEBBCe9over25_->Write(); 00265 h1_islandEBBCe5x5_->Write(); 00266 h1_islandEBBCEnergy_->Write(); 00267 h1_islandEBBCXtals_->Write(); 00268 00269 h1_islandEEBCe5x5_->Write(); 00270 h1_islandEEBCEnergy_->Write(); 00271 h1_islandEEBCXtals_->Write(); 00272 00273 h1_islandEESCEnergy_->Write(); 00274 h1_corrIslandEESCEnergy_->Write(); 00275 h1_corrIslandEESCET_->Write(); 00276 h1_islandEESCClusters_->Write(); 00277 00278 h1_hybridSCClusters_->Write(); 00279 h1_hybridSCEnergy_->Write(); 00280 h1_corrHybridSCEnergy_->Write(); 00281 h1_corrHybridSCET_->Write(); 00282 h1_corrHybridSCEta_->Write(); 00283 h1_corrHybridSCPhi_->Write(); 00284 00285 rootFile_->Close(); 00286 }
std::string EgammaSimpleAnalyzer::correctedHybridSuperClusterCollection_ [private] |
Definition at line 65 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::correctedHybridSuperClusterProducer_ [private] |
Definition at line 66 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::correctedIslandEndcapSuperClusterCollection_ [private] |
Definition at line 59 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::correctedIslandEndcapSuperClusterProducer_ [private] |
Definition at line 60 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
TH1F* EgammaSimpleAnalyzer::h1_corrHybridSCEnergy_ [private] |
Definition at line 99 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_corrHybridSCET_ [private] |
Definition at line 100 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_corrHybridSCEta_ [private] |
Definition at line 101 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_corrHybridSCPhi_ [private] |
Definition at line 102 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_corrIslandEESCEnergy_ [private] |
Definition at line 94 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_corrIslandEESCET_ [private] |
Definition at line 95 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_hybridSCClusters_ [private] |
Definition at line 103 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_hybridSCEnergy_ [private] |
Definition at line 98 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEBBCe5x5_ [private] |
Definition at line 85 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEBBCe9over25_ [private] |
Definition at line 84 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEBBCEnergy_ [private] |
Definition at line 86 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEBBCXtals_ [private] |
Definition at line 87 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEEBCe5x5_ [private] |
Definition at line 89 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEEBCEnergy_ [private] |
Definition at line 90 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEEBCXtals_ [private] |
Definition at line 91 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEESCClusters_ [private] |
Definition at line 96 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_islandEESCEnergy_ [private] |
Definition at line 93 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_nHybridSC_ [private] |
Definition at line 82 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_nIslandEBBC_ [private] |
Definition at line 79 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_nIslandEEBC_ [private] |
Definition at line 80 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
TH1F* EgammaSimpleAnalyzer::h1_nIslandEESC_ [private] |
Definition at line 81 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
std::string EgammaSimpleAnalyzer::hybridSuperClusterCollection_ [private] |
Definition at line 62 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::hybridSuperClusterProducer_ [private] |
Definition at line 63 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandBarrelBasicClusterCollection_ [private] |
Definition at line 48 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandBarrelBasicClusterProducer_ [private] |
Definition at line 49 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandBarrelBasicClusterShapes_ [private] |
Definition at line 50 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandEndcapBasicClusterCollection_ [private] |
Definition at line 52 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandEndcapBasicClusterProducer_ [private] |
Definition at line 53 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandEndcapBasicClusterShapes_ [private] |
Definition at line 54 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandEndcapSuperClusterCollection_ [private] |
Definition at line 56 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::islandEndcapSuperClusterProducer_ [private] |
Definition at line 57 of file EgammaSimpleAnalyzer.h.
Referenced by analyze(), and EgammaSimpleAnalyzer().
int EgammaSimpleAnalyzer::nbinHist_ [private] |
Definition at line 74 of file EgammaSimpleAnalyzer.h.
Referenced by beginJob(), and EgammaSimpleAnalyzer().
std::string EgammaSimpleAnalyzer::outputFile_ [private] |
TFile* EgammaSimpleAnalyzer::rootFile_ [private] |
Definition at line 69 of file EgammaSimpleAnalyzer.h.
Referenced by beginJob(), EgammaSimpleAnalyzer(), endJob(), and ~EgammaSimpleAnalyzer().
double EgammaSimpleAnalyzer::xMaxHist_ [private] |
Definition at line 73 of file EgammaSimpleAnalyzer.h.
Referenced by beginJob(), and EgammaSimpleAnalyzer().
double EgammaSimpleAnalyzer::xMinHist_ [private] |
Definition at line 72 of file EgammaSimpleAnalyzer.h.
Referenced by beginJob(), and EgammaSimpleAnalyzer().