#include <XANADOO/domain/interface/XANAEsdEvent.h>
#include <XANADOO/XANAClusters/interface/XANACluster.h>
#include <XANADOO/XANAClusters/interface/XANACaloRecHit.h>
#include <XANADOO/XANAClusters/interface/XANAEmCluster.h>
#include <XANADOO/XANAClusters/interface/XANASuperCluster.h>
#include <XANADOO/XANAClusters/interface/XANAHadCluster.h>
#include <XANADOO/XANATracks/interface/XANATrack.h>
#include <XANADOO/XANAMcInfo/interface/XANAGeneratorEvent.h>
#include <XANADOO/XANAMcInfo/interface/XANAGeantEvent.h>
#include <XANADOO/XANAMcInfo/interface/XANAGeneratorParticle.h>
#include <XANADOO/XANAElectronCandidate/interface/XANAElectronCandidate.h>
#include "TTree.h"
#include "TFile.h"
#include "TBranch.h"
#include "TClonesArray.h"
#include "TObjArray.h"
#include "TH1.h"
#include "TH2.h"
#include "TSystem.h"
#include <string>
Include dependency graph for readESDHist.cpp:
Go to the source code of this file.
Defines | |
#define | SINGLE_ELECTRON |
Functions | |
int | main (int argc, char **argv) |
|
Definition at line 26 of file readESDHist.cpp. |
|
Definition at line 31 of file readESDHist.cpp. References XANAGeantEvent::clear, XANAGeneratorEvent::clear, XANAEsdEvent::clear, XANASuperCluster::getBremClusters, XANAElectronCandidate::getDeltaEtaSuperClusterAtVtx, XANAEsdEvent::getElectronCandidates, XANAEsdEvent::getEmClusterHits, XANAEsdEvent::getEmClusters, XANACaloRecHit::getEnergy, XANASuperCluster::getEnergy, XANACluster::getEnergy, XANAElectronCandidate::getESuperClusterOverP, XANAEsdEvent::getHadClusters, XANAElectronCandidate::getHadronicOverEm, XANAGeneratorParticle::getMomentum, XANATrack::getMomentumAtFirstPoint, XANATrack::getMomentumAtLastPoint, XANATrack::getMomentumAtVertex, XANASuperCluster::getNumberOfBrems, XANAEsdEvent::getNumberOfElectronCandidates, XANAEsdEvent::getNumberOfEmClusterHits, XANAEsdEvent::getNumberOfEmClusters, XANAGeantEvent::getNumberOfGeantTracks, XANAGeantEvent::getNumberOfGeantVertices, XANAEsdEvent::getNumberOfHadClusters, XANAGeneratorEvent::getNumberOfParticles, XANAGeneratorEvent::getNumberOfPileupParticles, XANAEsdEvent::getNumberOfPreshClusters, XANAGeneratorEvent::getNumberOfSignalParticles, XANAGeneratorEvent::getNumberOfStablePileupParticles, XANAGeneratorEvent::getNumberOfStableSignalParticles, XANAEsdEvent::getNumberOfSuperClusters, XANAEsdEvent::getNumberOfTracks, XANAEsdEvent::getNumberOfVertexTracks, XANAEsdEvent::getNumberOfVertices, XANACluster::getPosition, XANAEsdEvent::getPreshClusters, XANASuperCluster::getSeedCluster, XANAGeneratorEvent::getSignalParticles, XANAEsdEvent::getSuperClusters, XANAEsdEvent::getTracks, and XANAEsdEvent::getVertices.
00032 { 00033 string inputFileName = "h300eemm.root"; 00034 string histoFileName = "hist_"+inputFileName; 00035 00036 // Connect to file generated 00037 TFile f(inputFileName.c_str()); 00038 00039 // Booking the histograms 00040 00041 // Clusters 00042 TH1F hClustersNumber("hClustersNumber", "Number of clusters", 20, 0., 20.); 00043 TH1F hClustersEnergy("hClustersEnergy", "Energy of cluster", 100, 0., 500.); 00044 TH1F hSuperClustersNumber("hSuperClustersNumber", "Number of super clusters", 50, 0., 50.); 00045 TH1F hSuperClustersEnergy("hSuperClustersEnergy", "Energy of super cluster", 100, 0., 500.); 00046 TH1F hBremsRecovered("hBremsRecovered","Energy of recovered bremsstrahlung photons",100,0.,200.); 00047 TH2F hBremsDetaDphi("hBremsDetaDphi","Delta phi vs Delta eta for brems clusters",100,-0.3,0.3,100,-0.3,0.3); 00048 00049 TH1F hHadClustersNumber("hHadClustersNumber", "Number of hadron clusters", 20, 0., 20.); 00050 TH1F hHadClustersEnergy("hHadClustersEnergy", "Energy of hadron cluster", 100, 0., 500.); 00051 00052 TH1F hPreshClustersNumber("hPreshClustersNumber", "Number of preshower clusters", 20, 0., 20.); 00053 00054 // Hits 00055 TH1F hHitsNumber("hHitsNumber", "Number of hits", 100, 0., 100.); 00056 TH1F hHitsEnergy("hHitsEnergy", "Energy of hit", 200, 0., 200.); 00057 00058 // Tracks 00059 TH1F hTracksNumber("hTracksNumber", "Number of tracks", 10, 0., 10.); 00060 TH1F hTracksMomentumAtVertex("hTracksMomentumAtVertex", "Momentum of track at vertex", 100, 0., 500.); 00061 TH1F hTracksMomentumAtFirstPoint("hTracksMomentumAtFirstPoint", "Momentum of track at first point", 100, 0., 500.); 00062 TH1F hTracksMomentumAtLastPoint("hTracksMomentumAtLastPoint", "Momentum of track at last point", 100, 0., 500.); 00063 00064 //Vertex 00065 TH1F hVertexNumber("hVertexNumber", "Number of vertices", 10, 0., 10.); 00066 TH1F hVertexTracksNumber("hVertexTracksNumber", "Number of vertex tracks", 10, 0., 10.); 00067 00068 // Monte Carlo information 00069 00070 // generator event 00071 TH1F hGenNumberOfParticles("hGenNumberOfParticles","Generator info: number of particles", 100, 0., 100.); 00072 00073 #ifdef SINGLE_ELECTRON 00074 TH1F hGenElectronMomentum("hGenElectronMomentum","Generator info: Electron momentum", 100, 0., 500.); 00075 TH1F hGenElectronPt("hGenElectronPt","MGeneratorC info: Electron pT", 200, 0., 200.); 00076 TH1F hGenElectronEta("hGenElectronEta","Generator info: Electron eta", 100, -3.,3.); 00077 TH1F hGenElectronPhi("hGenElectronPhi","Generator info: Electron phi", 100, -3.2,3.2); 00078 XANAGeneratorParticle *tmpGenElectron = 0; 00079 #endif 00080 00081 // GEANT event 00082 TH1F hGeantNumberOfTracks("hGeantNumberOfTracks","GEANT info: Number of tracks", 50, 0., 50.); 00083 TH1F hGeantNumberOfVertices("hGeantNumberOfVertices","GEANT info: Number of vertices", 50, 0., 50.); 00084 00085 // electron candidates 00086 TH1F hElectronNumber("hElectronNumber", "Number of electrons", 10, 0., 10.); 00087 TH1F hElectronEscOverP("hElectronEscOverP", "Electron: E supercluster over P", 100, 0.5, 1.5); 00088 TH1F hElectronDetaSc("hElectronDetaSc", "Electron: Delta eta super cluster", 100, -0.2, 0.2); 00089 TH1F hElectronHoE("hElectronHoE", "Electron: Hadronic over Em", 100, 0.5, 1.5); 00090 00091 // Read Tree named "T" in memory. Tree pointer is assigned the same name 00092 TTree *T = (TTree*)f.Get("ESD"); 00093 00094 // Reconstructed event 00095 XANAEsdEvent *event = 0; 00096 T->SetBranchAddress("Event",&event); 00097 00098 Short_t nc = 0, nh = 0, nhc = 0, ns = 0, npc = 0, nt = 0, nv = 0, nvt = 0; 00099 Short_t nbrems = 0, ne = 0; 00100 00101 TObjArray *clusters = 0; 00102 TObjArray *hadClusters = 0; 00103 TObjArray *superClusters = 0; 00104 TObjArray *preshClusters = 0; 00105 TRefArray *hits = 0; 00106 TObjArray *tracks = 0; 00107 TObjArray *vertices = 0; 00108 TObjArray *electrons = 0; 00109 00110 XANAEmCluster *tmpEmCluster = 0; 00111 XANASuperCluster *tmpSuperCluster = 0; 00112 XANAEmCluster *tmpSeedCluster = 0; 00113 XANAEmCluster *tmpBremCluster = 0; 00114 XANAHadCluster *tmpHadCluster = 0; 00115 XANAElectronCandidate *tmpElectronCandidate = 0; 00116 00117 // MC information (generator and GEANT event) 00118 XANAGeneratorEvent *generatorEvent = new XANAGeneratorEvent(); 00119 TBranch *branchGeneratorEvent = T->GetBranch("GeneratorEvent"); 00120 branchGeneratorEvent->SetAddress(&generatorEvent); 00121 00122 XANAGeantEvent *geantEvent = new XANAGeantEvent(); 00123 TBranch *branchGeantEvent = T->GetBranch("GeantEvent"); 00124 branchGeantEvent->SetAddress(&geantEvent); 00125 00126 Short_t mc_np = 0, mc_nsp = 0, mc_nssp = 0, mc_npp = 0, mc_nspp = 0, mc_ngt = 0, mc_ngv = 0; 00127 00128 // Start main loop on all events 00129 Int_t nb = 0; 00130 Int_t nevent = (Int_t)T->GetEntries(); 00131 for (Int_t i=0; i<nevent; i++) { 00132 nb += T->GetEntry(i); 00133 00134 cout << "================ Event: " << i << " ================" << endl; 00135 00136 // Reconstructed event 00137 00138 cout << "Reconstructed:" ; 00139 00140 nc = event->getNumberOfEmClusters(); clusters = event->getEmClusters(); 00141 nh = event->getNumberOfEmClusterHits(); hits = event->getEmClusterHits(); 00142 nhc = event->getNumberOfHadClusters(); hadClusters = event->getHadClusters(); 00143 ns = event->getNumberOfSuperClusters(); superClusters = event->getSuperClusters(); 00144 npc = event->getNumberOfPreshClusters(); preshClusters = event->getPreshClusters(); 00145 nt = event->getNumberOfTracks(); tracks = event->getTracks(); 00146 nv = event->getNumberOfVertices(); vertices = event->getVertices(); 00147 nvt = event->getNumberOfVertexTracks(); 00148 ne = event->getNumberOfElectronCandidates(); electrons = event->getElectronCandidates(); 00149 00150 cout << " clust, hits, hadClu, supers, presho, tracks, vertex, vTrack, nElec" << endl; 00151 printf(" \t \t %d\t%d\t%d\t%d\t%d\t%d\t%d\t%d\t%d\n", nc, nh, nhc, ns, npc, nt, nv, nvt, ne); 00152 00153 // Monte Carlo information 00154 00155 // generator event 00156 mc_np = generatorEvent->getNumberOfParticles(); 00157 mc_nsp = generatorEvent->getNumberOfSignalParticles(); 00158 mc_nssp = generatorEvent->getNumberOfStableSignalParticles(); 00159 mc_npp = generatorEvent->getNumberOfPileupParticles(); 00160 mc_nspp = generatorEvent->getNumberOfStablePileupParticles(); 00161 00162 // GEANT event 00163 mc_ngt = geantEvent->getNumberOfGeantTracks(); 00164 mc_ngv = geantEvent->getNumberOfGeantVertices(); 00165 00166 cout << "MC Info: "; 00167 cout << "parti, signl, stSig, pileup, stPil, geantTr, geantVe" << endl; 00168 printf("\t %d\t %d\t%d\t%d\t%d\t%d\t%d\n", mc_np, mc_nsp, mc_nssp, mc_npp, mc_nspp, mc_ngt, mc_ngv); 00169 00170 00171 // Filling the histograms 00172 00173 // Clusters 00174 00175 hClustersNumber.Fill((Axis_t) nc); 00176 for (int ic = 0; ic < nc; ic++) 00177 { 00178 tmpEmCluster = (XANAEmCluster *) clusters->At(ic); 00179 hClustersEnergy.Fill((Axis_t) tmpEmCluster->getEnergy()); 00180 } 00181 00182 hSuperClustersNumber.Fill((Axis_t) ns); 00183 for (int is = 0; is < ns; is++) 00184 { 00185 tmpSuperCluster = (XANASuperCluster *) superClusters->At(is); 00186 hSuperClustersEnergy.Fill((Axis_t) tmpSuperCluster->getEnergy()); 00187 tmpSeedCluster = (XANAEmCluster *) tmpSuperCluster->getSeedCluster(); 00188 hBremsRecovered.Fill((Axis_t) tmpSuperCluster->getEnergy() - tmpSeedCluster->getEnergy()); 00189 nbrems = tmpSuperCluster->getNumberOfBrems(); 00190 for (int ib = 0; ib < nbrems; ib++) 00191 { 00192 tmpBremCluster = (XANAEmCluster *) tmpSuperCluster->getBremClusters()->At(ib); 00193 hBremsDetaDphi.Fill(tmpSeedCluster->getPosition().getEta()-tmpBremCluster->getPosition().getEta(), 00194 tmpSeedCluster->getPosition().getPhi()-tmpBremCluster->getPosition().getPhi()); 00195 } 00196 } 00197 00198 hHadClustersNumber.Fill((Axis_t) nh); 00199 for (int ihc = 0; ihc < nhc; ihc++) 00200 { 00201 tmpHadCluster = (XANAHadCluster *) hadClusters->At(ihc); 00202 hHadClustersEnergy.Fill((Axis_t) tmpHadCluster->getEnergy()); 00203 } 00204 00205 hPreshClustersNumber.Fill((Axis_t) npc); 00206 00207 // Hits 00208 XANACaloRecHit *tmpCaloRecHit; 00209 00210 hHitsNumber.Fill((Axis_t) nh); 00211 for (int ih = 0; ih < nh; ih++) 00212 { 00213 tmpCaloRecHit = (XANACaloRecHit *) hits->At(ih); 00214 hHitsEnergy.Fill((Axis_t) tmpCaloRecHit->getEnergy()); 00215 } 00216 00217 // Tracks 00218 XANATrack *tmpTrack; 00219 00220 hTracksNumber.Fill((Axis_t) nt); 00221 for (int it = 0; it < nt; it++) 00222 { 00223 tmpTrack = (XANATrack *) tracks->At(it); 00224 hTracksMomentumAtVertex.Fill((Axis_t) tmpTrack->getMomentumAtVertex().mag()); 00225 hTracksMomentumAtFirstPoint.Fill((Axis_t) tmpTrack->getMomentumAtFirstPoint().mag()); 00226 hTracksMomentumAtLastPoint.Fill((Axis_t) tmpTrack->getMomentumAtLastPoint().mag()); 00227 } 00228 00229 // Vertex 00230 hVertexNumber.Fill((Axis_t) nv); 00231 hVertexTracksNumber.Fill((Axis_t) nvt); 00232 00233 // Monte Carlo information 00234 00235 // generator event 00236 hGenNumberOfParticles.Fill((Axis_t) mc_np); 00237 00238 #ifdef SINGLE_ELECTRON 00239 tmpGenElectron = (XANAGeneratorParticle *) generatorEvent->getSignalParticles()->At(0); 00240 hGenElectronMomentum.Fill((Axis_t) tmpGenElectron->getMomentum().vect().mag()); 00241 hGenElectronPt.Fill((Axis_t) tmpGenElectron->getMomentum().perp()); 00242 hGenElectronEta.Fill((Axis_t) tmpGenElectron->getMomentum().pseudoRapidity()); 00243 hGenElectronPhi.Fill((Axis_t) tmpGenElectron->getMomentum().phi()); 00244 #endif 00245 00246 // GEANT event 00247 hGeantNumberOfTracks.Fill((Axis_t) mc_ngt); 00248 hGeantNumberOfVertices.Fill((Axis_t) mc_ngv); 00249 00250 00251 // electron candidates 00252 hElectronNumber.Fill((Axis_t) ne); 00253 if (ne != 0) { 00254 cout << "Electrons: Esc/p, DetaSC, HoE" << endl; 00255 for (int ie = 0; ie < ne ; ie++) { 00256 tmpElectronCandidate = (XANAElectronCandidate *) electrons->At(ie); 00257 00258 hElectronEscOverP.Fill((Axis_t) tmpElectronCandidate->getESuperClusterOverP()); 00259 hElectronDetaSc.Fill((Axis_t) tmpElectronCandidate->getDeltaEtaSuperClusterAtVtx()); 00260 hElectronHoE.Fill((Axis_t) tmpElectronCandidate->getHadronicOverEm()); 00261 00262 printf(" %d, \t%f, %f, %f \n", ie+1, 00263 tmpElectronCandidate->getESuperClusterOverP(), 00264 tmpElectronCandidate->getDeltaEtaSuperClusterAtVtx(), 00265 tmpElectronCandidate->getHadronicOverEm()); 00266 } 00267 } 00268 00269 00270 event->clear(); 00271 generatorEvent->clear(); 00272 geantEvent->clear(); 00273 } 00274 00275 cout << endl; 00276 cout << nevent << " events and " << nb << " bytes read " << endl; 00277 cout << endl; 00278 00279 f.Close(); 00280 00281 00282 // Exporting histograms to file 00283 TFile hFile(histoFileName.c_str(),"RECREATE"); 00284 00285 hClustersNumber.Write(); 00286 hClustersEnergy.Write(); 00287 hSuperClustersNumber.Write(); 00288 hSuperClustersEnergy.Write(); 00289 hHadClustersNumber.Write(); 00290 hHadClustersEnergy.Write(); 00291 hPreshClustersNumber.Write(); 00292 hHitsNumber.Write(); 00293 hHitsEnergy.Write(); 00294 hTracksNumber.Write(); 00295 hTracksMomentumAtVertex.Write(); 00296 hTracksMomentumAtFirstPoint.Write(); 00297 hTracksMomentumAtLastPoint.Write(); 00298 hBremsRecovered.Write(); 00299 hBremsDetaDphi.Write(); 00300 hVertexNumber.Write(); 00301 hVertexTracksNumber.Write(); 00302 hGenNumberOfParticles.Write(); 00303 #ifdef SINGLE_ELECTRON 00304 hGenElectronMomentum.Write(); 00305 hGenElectronPt.Write(); 00306 hGenElectronEta.Write(); 00307 hGenElectronPhi.Write(); 00308 #endif 00309 hGeantNumberOfTracks.Write(); 00310 hGeantNumberOfVertices.Write(); 00311 hElectronNumber.Write(); 00312 hElectronEscOverP.Write(); 00313 hElectronDetaSc.Write(); 00314 hElectronHoE.Write(); 00315 00316 hFile.Close(); 00317 00318 00319 } |