AliRoot Core  edcc906 (edcc906)
TestESD.C
Go to the documentation of this file.
1 
17 #if !defined(__CINT__) || defined(__MAKECINT__)
18 
19 //Root include files
20 #include <Riostream.h>
21 #include <TFile.h>
22 //#include <TSystem.h>
23 #include <TH1F.h>
24 #include <TParticle.h>
25 #include <TRefArray.h>
26 
27 //AliRoot include files
28 #include "AliRunLoader.h"
29 #include "AliStack.h"
30 #include "AliESDEvent.h"
31 #include "AliESDVertex.h"
32 #include "AliESDCaloCluster.h"
33 #include "AliESDCaloCells.h"
34 #include "AliESDCaloTrigger.h"
35 #include "AliPID.h"
36 #include "AliEMCALGeometry.h"
37 
38 #endif
39 
40 // Change the bool depending on what information you want to print
41 // when all FALSE, prints minimum cluster information.
42 Bool_t kPrintKine = kFALSE;
43 Bool_t kPrintCaloCells = kFALSE;
44 Bool_t kPrintCaloTrigger = kFALSE;
45 Bool_t kPrintTrackMatches = kFALSE;
46 Bool_t kPrintCaloCluster = kTRUE ;
47 Bool_t kPrintClusterCells = kFALSE;
48 Bool_t kPrintClusterPID = kFALSE;
49 Bool_t kPrintMisalMatrix = kFALSE;
50 
54 void TestESD()
55 {
56  // Init some example histograms
57  // ESD
58  TH1F * hEta = (TH1F*) new TH1F("hEta","reco #eta",1000, -0.71,0.71);
59  TH1F * hPhi = (TH1F*) new TH1F("hPhi","reco #phi",360, 0,360);
60  TH1F * hE = (TH1F*) new TH1F("hE" ,"reco e",300, 0,30);
61  TH1F * hTime = (TH1F*) new TH1F("hTime" ,"reco time",1000, 0,1000);
62  hEta ->SetXTitle("#eta");
63  hPhi ->SetXTitle("#phi (deg)");
64  hE ->SetXTitle("E (GeV)");
65  hTime->SetXTitle("time (ns)");
66 
67  // Monte Carlo
68  TH1F * hMCEta = (TH1F*) new TH1F("hMCEta","MC #eta",1000, -0.71,0.71);
69  TH1F * hMCPhi = (TH1F*) new TH1F("hMCPhi","MC #phi",360, 0,360);
70  TH1F * hMCE = (TH1F*) new TH1F("hMCE" ,"MC e",300, 0,30);
71  hMCEta->SetXTitle("#eta");
72  hMCPhi->SetXTitle("#phi (deg)");
73  hMCE ->SetXTitle("E (GeV)");
74 
75  // ESD - MonteCarlo
76  TH1F * hMCDEta = (TH1F*) new TH1F("hMCDEta"," #eta cluster - #eta MC",500, -0.05,0.05);
77  TH1F * hMCDPhi = (TH1F*) new TH1F("hMCDPhi"," #phi cluster - #phi MC",500, -0.05,0.05);
78  TH1F * hMCDE = (TH1F*) new TH1F("hMCDE" ,"e cluster - e MC",200, -10,10);
79  hMCDEta->SetXTitle("#eta_{reco}-#eta_{MC}");
80  hMCDPhi->SetXTitle("#phi_{reco}-#phi_{MC} (rad)");
81  hMCDE ->SetXTitle("E_{reco}-E_{MC} (GeV)");
82 
83  // ESD - Track-matching
84  TH1F * hTMDEta = (TH1F*) new TH1F("hTMDEta"," #eta cluster - eta track",500, -0.05,0.05);
85  TH1F * hTMDPhi = (TH1F*) new TH1F("hTMDPhi"," #phi cluster - phi track",500, -0.05,0.05);
86  TH1F * hTMDEtaOut= (TH1F*) new TH1F("hTMDEtaOut"," #eta cluster - eta track-outer",500, -0.05,0.05);
87  TH1F * hTMDPhiOut= (TH1F*) new TH1F("hTMDPhiOut"," #phi cluster - phi track-outer",500, -0.05,0.05);
88  TH1F * hTMResEta = (TH1F*) new TH1F("hTMResEta"," Residual #eta cluster - eta track",500, -0.05,0.05);
89  TH1F * hTMResPhi = (TH1F*) new TH1F("hTMResPhi"," Residual #phi cluster - phi track",500, -0.05,0.05);
90  hTMDEta ->SetXTitle("#eta_{cluster}-#eta_{track}");
91  hTMDPhi ->SetXTitle("#phi_{cluster}-#phi_{track} (rad)");
92  hTMDEtaOut->SetXTitle("#eta_{cluster}-#eta_{track-out}");
93  hTMDPhiOut->SetXTitle("#phi_{cluster}-#phi_{track-out} (rad)");
94 
95  TH1F * hTMEOverP = (TH1F*) new TH1F("hTMEOverP"," E_{cluster}/ p_{Track}",100,0,2);
96  TH1F * hTMEOverPOut = (TH1F*) new TH1F("hTMEOverPOut"," E_{cluster}/ p_{Track-out}",100,0,2);
97  hTMEOverP ->SetXTitle("E_{cluster}/ p_{Track}");
98  hTMEOverPOut ->SetXTitle("E_{cluster}/ p_{Track-out}");
99 
100  // L1 trigger
101  TH2I * hEGAPatch = new TH2I("hEGAPatch","EGA trigger",51,-0.5,50.5,65,-0.5,64.5);
102  TH2I * hEJEPatch = new TH2I("hEJEPatch","EJE trigger",51,-0.5,50.5,65,-0.5,64.5);
103  hEGAPatch->SetXTitle("column (#eta direction)");
104  hEGAPatch->SetYTitle("row (#phi direction)");
105  hEJEPatch->SetXTitle("column (#eta direction)");
106  hEJEPatch->SetYTitle("row (#phi direction)");
107 
108  // Open the ESD file, get the tree with events
109  TFile* f = new TFile("AliESDs.root");
110  TTree* esdTree = (TTree*)f->Get("esdTree");
111  AliESDEvent* esd = new AliESDEvent();
112  esd->ReadFromTree(esdTree);
113 
114  // Init geometry and array that will contain the clusters, Run2
115  AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance("EMCAL_COMPLETE12SMV1_DCAL_8SM") ;
116 
117  // Run1:
118  //AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance("EMCAL_COMPLETE12SMV1") ;
119  //AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance("EMCAL_COMPLETEV1") ;
120  //AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance("EMCAL_FIRSTYEARV1") ;
121 
122  // Pass the geometry transformation matrix from ESDs to geometry
123  for(Int_t mod = 0; mod < (geom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
124  {
125  //printf("matrix %d\n",mod);
126  if(esd->GetEMCALMatrix(mod))
127  {
128  if(kPrintMisalMatrix)
129  {
130  printf("Misalign matrix for: mod %d, matrix %p\n",mod, esd->GetEMCALMatrix(mod));
131 
132  (esd->GetEMCALMatrix(mod))->Print("");
133  }
134 
135  // Fill it in case cell index to global position is requested
136  geom->SetMisalMatrix(esd->GetEMCALMatrix(mod),mod) ;
137  }//matrix
138  }//module
139 
140  // Loop of events
141 
142  Float_t pos[3];
143 
144  TRefArray* caloClusters = new TRefArray();
145 
146  Int_t nEvt = esdTree->GetEntries();
147 
148  for(Int_t iev = 0; iev < nEvt; iev++)
149  {
150  cout << "<<<< Event: " << iev+1 << "/" << nEvt << " >>>>";
151  cout <<" Triggers: "<< esd->GetFiredTriggerClasses();
152  cout<<endl;
153 
154  esdTree->GetEvent(iev);
155 
156  // In case you want to play with MC data, get stack
157  AliStack *stack = 0;
158  if(kPrintKine)
159  {
161  rl->LoadKinematics();
162  rl->GetEvent(iev);
163  stack = rl->Stack();
164  }
165 
166  // Get reconstructed vertex position, only used if
167  // the momentum of the cluster is calculated, see comment below
168  //
169  //Double_t vertex_position[3];
170  //esd->GetVertex()->GetXYZ(vertex_position);
171 
172  //------------------------------------------------------
173  // Get Cells Array, all cells in event, print cells info
174  //------------------------------------------------------
175  AliVCaloCells &cells = *(esd->GetEMCALCells());
176 
177  if(kPrintCaloCells)
178  {
179  Int_t nTotalCells = cells.GetNumberOfCells() ;
180  //Int_t type = cells.GetType();
181  for (Int_t icell= 0; icell < nTotalCells; icell++)
182  {
183  Int_t cellid = cells.GetCellNumber(icell);
184  cout<<"Cell : " <<icell<<"/"<<nTotalCells <<" - ID: "<<cells.GetCellNumber(icell)<<"; High Gain? "<<cells.GetHighGain(icell);
185  cout<<"; Amplitude: "<<cells.GetAmplitude(icell)<<"; Time: "<<cells.GetTime(icell)*1e9;
186  cout<<"; MC label " <<cells.GetMCLabel(icell) <<"; Embeded E fraction "<<cells.GetEFraction(icell);
187  cout<<endl;
188  }// cell loop
189  }
190 
191  //------------------------------------------------------
192  // Calo Trigger, L1
193  //------------------------------------------------------
194 
196  {
197  // Bits to know what L1 trigger,
198  Int_t bitEGA = 6; // 4 for old data
199  Int_t bitEJE = 8; // 5 for old data
200 
201  AliESDCaloTrigger& trg = *(esd->GetCaloTrigger("EMCAL"));
202 
203  trg.Reset();
204  while (trg.Next())
205  {
206  Int_t col, row;
207  trg.GetPosition(col, row);
208 
209  if (col > -1 && row > -1)
210  {
211  Int_t bit = 0;
212  trg.GetTriggerBits(bit);
213 
214  Int_t ts = 0;
215  trg.GetL1TimeSum(ts);
216 
217  // Check if it is EGA or EJE (it does not work??)
218  Bool_t isEGA1 = ((bit>> bitEGA ) & 0x1);
219  Bool_t isEGA2 = ((bit>> (bitEGA+1)) & 0x1);
220  Bool_t isEJE1 = ((bit>> bitEJE ) & 0x1);
221  Bool_t isEJE2 = ((bit>> (bitEJE+1)) & 0x1);
222 
223  if ( ts >= 0 )
224  {
225  printf("Bit %d (G1 %d,G2 %d,J1 %d,J2 %d), cell in patch position (ieta,iphi)=(%d,%d), L1 amplitude %d\n ",
226  bit, isEGA1, isEGA2, isEJE1, isEJE2, col, row, ts);
227  if(bit>10) hEJEPatch->Fill(col,row);
228  else hEGAPatch->Fill(col,row);
229  }
230  }
231  }
232  }
233 
234  //------------------------------------------------------
235  // Calo Clusters
236  //------------------------------------------------------
237 
238  // Get CaloClusters Array
239  esd->GetEMCALClusters(caloClusters);
240 
241  // Loop over clusters
242  Int_t nclus = caloClusters->GetEntries();
243  for (Int_t icl = 0; icl < nclus; icl++)
244  {
245  AliVCluster* clus = (AliVCluster*)caloClusters->At(icl);
246  Float_t energy = clus->E();
247  clus->GetPosition(pos);
248  TVector3 vpos(pos[0],pos[1],pos[2]);
249 
250  // We can get a momentum TLorentzVector per cluster,
251  // corrected by the vertex position, see above
252  //
253  //TLorentzVector p;
254  //clus->GetMomentum(p,vertex_position);
255 
256  Double_t cphi = vpos.Phi();
257  if(cphi < 0) cphi +=TMath::TwoPi();
258  Double_t ceta = vpos.Eta();
259 
260  Int_t nMatched = clus->GetNTracksMatched();
261  Int_t trackIndex = clus->GetTrackMatchedIndex();
262  Int_t nLabels = clus->GetNLabels();
263  Int_t labelIndex = clus->GetLabel();
264  Int_t nCells = clus->GetNCells();
265  Int_t nlm = clus->GetNExMax() ;
266 
267  // Fill some histograms
268  hEta ->Fill(ceta);
269  hPhi ->Fill(cphi*TMath::RadToDeg());
270  hE ->Fill(energy);
271  hTime->Fill(clus->GetTOF()*1e9);
272 
273  // Print basic cluster information
274  if ( kPrintCaloCluster )
275  {
276  cout << "Cluster: " << icl+1 << "/" << nclus << " Energy: " << energy << "; Phi: "
277  << cphi*TMath::RadToDeg() << "; Eta: " << ceta
278  << "; NCells: " << nCells << "; NLM: " << nlm
279  << "; #Labels: " << nLabels << " Index: " << labelIndex
280  << "; Time "<<clus->GetTOF()*1e9<<" ns "<<endl;
281  }
282 
283  if ( nMatched > 0 )
284  {
285  if ( kPrintCaloCluster )
286  printf("\t N matches %d, Residual phi %2.4f, eta %2.4f\n",
287  nMatched,clus->GetTrackDx(),clus->GetTrackDz());
288  hTMResEta->Fill(clus->GetTrackDz());
289  hTMResPhi->Fill(clus->GetTrackDx());
290  }
291 
292  //Print primary info
293  if(stack && kPrintKine)
294  {
295  if(labelIndex >= 0 && labelIndex < stack->GetNtrack())
296  {
297  TParticle * particle = stack->Particle(labelIndex);
298  // Fill histograms with primary info
299  hMCEta ->Fill(particle->Eta());
300  hMCPhi ->Fill(particle->Phi()*TMath::RadToDeg());
301  hMCE ->Fill(particle->Energy());
302  hMCDEta ->Fill(ceta-particle->Eta());
303  hMCDPhi ->Fill(cphi-particle->Phi());
304  hMCDE ->Fill(energy-particle->Energy());
305 
306  // Print primary values
307  cout<<" MC More contributing primary: "<<particle->GetName()<<"; with kinematics: "<<endl;
308  cout<<" \t Energy: "<<particle->Energy()<<"; Phi: "<<particle->Phi()*TMath::RadToDeg()<<"; Eta: "<<particle->Eta()<<endl;
309  cout<<" \t dE: "<<energy-particle->Energy()<<"; dPhi: "<<(cphi-particle->Phi())*TMath::RadToDeg()<<"; dEta: "<<ceta-particle->Eta()<<endl;
310  cout<<" \t deposited energy fraction = "<<clus->GetClusterMCEdepFraction(0)<<endl;
311 
312  //printf("List of e fractions %p; cell edep fraction map %p\n",clus->GetClusterMCEdepFraction(),clus->GetCellsMCEdepFractionMap());
313 
314  for(Int_t i = 1; i < nLabels; i++)
315  {
316  //particle = stack->Particle((((AliESDCaloCluster*)clus)->GetLabelsArray())->At(i));
317  particle = stack->Particle((clus->GetLabels())[i]);
318  //or Int_t *labels = clus->GetLabels();
319  //particle = stack->Particle(labels[i]);
320  cout<<" Other contributing primary: "<<particle->GetName()<< "; Energy "<<particle->Energy()<<
321  "; deposited energy fraction = "<<clus->GetClusterMCEdepFraction(i)<<endl;
322  }
323  }
324  else if( labelIndex >= stack->GetNtrack())
325  cout <<"PROBLEM, label is too large : "<<labelIndex<<" >= particles in stack "<< stack->GetNtrack() <<endl;
326  else
327  cout<<"Negative label!!! : "<<labelIndex<<endl;
328  } // play with stack
329 
330  // Matching results
331  if(kPrintTrackMatches && trackIndex >= 0)
332  {
333  AliESDtrack* track = esd->GetTrack(trackIndex);
334  if(track)
335  {
336  Double_t tphi = track->Phi();
337  Double_t teta = track->Eta();
338  Double_t tmom = track->P();
339 
340  Double_t deta = teta - ceta;
341  Double_t dphi = tphi - cphi;
342 
343  if(dphi > TMath::Pi()) dphi -= 2*TMath::Pi();
344  if(dphi < -TMath::Pi()) dphi += 2*TMath::Pi();
345 
346  Double_t dR = sqrt(dphi*dphi + deta*deta);
347 
348  Double_t pOverE = tmom/energy;
349 
350  cout << "\t Track Momentum : " << tmom << " phi: " << tphi << " eta: " << teta << endl;
351 
352  hTMEOverP->Fill(pOverE);
353  hTMDEta ->Fill(deta);
354  hTMDPhi ->Fill(dphi);
355 
356  if(track->GetOuterParam())
357  {
358  tphi = track->GetOuterParam()->Phi();
359  teta = track->GetOuterParam()->Eta();
360  tmom = track->GetOuterParam()->P();
361 
362  cout << "\t Track Momentum - Outer: " << tmom << " phi: " << tphi << " eta: " << teta << endl;
363 
364  deta = teta - ceta;
365  dphi = tphi - cphi;
366 
367  if(dphi > TMath::Pi()) dphi -= 2*TMath::Pi();
368  if(dphi < -TMath::Pi()) dphi += 2*TMath::Pi();
369 
370  dR = sqrt(dphi*dphi + deta*deta);
371 
372  pOverE = tmom/energy;
373 
374  if(dR < 0.02 && pOverE < 1.8 && nCells > 1)
375  {
376  cout << "\t Excellent MATCH! dR = " << dR << " p/E = " << pOverE << "Residuals phi "<< dphi*TMath::RadToDeg()<< " eta "<< deta << endl;
377  }
378 
379  hTMEOverPOut->Fill(pOverE);
380  hTMDEtaOut->Fill(deta);
381  hTMDPhiOut->Fill(dphi);
382  }
383  else printf("!!! NO OUTER PARAM !!!\n");
384  }
385  else printf("!!! NO TRACK !!!\n");
386  }// matching
387 
388  // Get PID weights and print them
389  if(kPrintClusterPID)
390  {
391  const Double_t *pid = clus->GetPID();
392  printf("PID weights: ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f, hadrons: pion %0.2f, kaon %0.2f, proton %0.2f , neutron %0.2f, kaon %0.2f \n",
397  } // PID
398 
399  // Get CaloCells of cluster and print their info, position.
401  {
402  UShort_t * index = clus->GetCellsAbsId() ;
403  Double_t * fraction = clus->GetCellsAmplitudeFraction() ;
404  Int_t sm = -1;
405 
406  for(Int_t i = 0; i < nCells ; i++)
407  {
408  Int_t absId = index[i]; // or clus->GetCellNumber(i) ;
409  Double_t ampFract = fraction[i];
410  Float_t amp = cells.GetCellAmplitude(absId) ;
411  Double_t time = cells.GetCellTime(absId);
412  Bool_t hg = cells.GetCellHighGain(absId);
413 
414  cout<<"\t Cluster Cell: AbsID : "<< absId << " == "<<clus->GetCellAbsId(i) <<"; Amplitude "<< amp;
415  cout<< "; Fraction "<<ampFract<<"; Time " <<time*1e9<<"; High Gain? "<<hg<<endl;
416 
417  if(kPrintKine)
418  {
419  Float_t eDepFrac[4];
420  clus->GetCellMCEdepFractionArray(i,eDepFrac);
421  Int_t map = 0;
422  if(clus->GetCellsMCEdepFractionMap()) map = clus->GetCellsMCEdepFractionMap()[i];
423  printf("\t \t cell MC label index %d, cell MC deposited energy map %d: p1 %2.2f, p2 %2.2f, p3 %2.2f p4 %2.2f\n",
424  cells.GetCellMCLabel(absId),map,eDepFrac[0],eDepFrac[1],eDepFrac[2],eDepFrac[3]);
425  }
426 
427  // Geometry methods
428  Int_t iSupMod = 0 ;
429  Int_t iTower = 0 ;
430  Int_t iIphi = 0 ;
431  Int_t iIeta = 0 ;
432  Int_t iphi = 0 ;
433  Int_t ieta = 0 ;
434 
435  if(geom)
436  {
437  geom->GetCellIndex(absId,iSupMod,iTower,iIphi,iIeta);
438  //Gives SuperModule and Tower numbers
439  geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,
440  iIphi, iIeta,iphi,ieta);
441  //Gives label of cell in eta-phi position per each supermodule
442  Float_t cellPhi = 0;
443  Float_t cellEta = 0;
444  geom->EtaPhiFromIndex(absId,cellEta,cellPhi);
445  cout<< " SModule "<<iSupMod<<"; Tower "<<iTower <<"; Eta "<<iIeta
446  <<"; Phi "<<iIphi<<"; Index: Cell Eta "<<ieta<<"; Cell Phi "<<iphi
447  <<"; Global: Cell Eta "<<cellEta<<"; Cell Phi "<<cellPhi*TMath::RadToDeg()
448  <<endl;
449 
450  if ( i == 0 ) sm = iSupMod;
451  else if ( sm != iSupMod) printf(" ******CLUSTER SHARED BY 2 SuperModules!!!!\n");
452 
453  }// geometry on
454  }// cluster cell loop
455  }// print cell clusters
456  } //cluster loop
457  } // event loop
458 
459 
460  //Write histograms in a file
461  TFile * fhisto = (TFile*) new TFile("histos.root","recreate");
462  hEta->Write();
463  hPhi->Write();
464  hE ->Write();
465  hTime->Write();
466 
467  if(kPrintKine)
468  {
469  hMCEta->Write();
470  hMCPhi->Write();
471  hMCE ->Write();
472  hMCDEta->Write();
473  hMCDPhi->Write();
474  hMCDE ->Write();
475  }
476 
478  {
479  hTMEOverP->Write();
480  hTMDEta ->Write();
481  hTMDPhi ->Write();
482 
483  hTMEOverPOut->Write();
484  hTMDEtaOut ->Write();
485  hTMDPhiOut ->Write();
486 
487  hTMResEta->Write();
488  hTMResPhi->Write();
489  }
490 
492  {
493  hEJEPatch->Write();
494  hEGAPatch->Write();
495  }
496 
497  fhisto->Close();
498 }
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
static const TString & GetDefaultEventFolderName()
Definition: AliConfig.h:58
Bool_t kPrintMisalMatrix
Print the alignment matrices stored in ESDs.
Definition: TestESD.C:49
virtual void GetCellMCEdepFractionArray(Int_t, Float_t *) const
Definition: AliVCluster.h:134
AliTPCcalibPID * pid
Definition: CalibPID.C:69
virtual UShort_t * GetClusterMCEdepFraction() const
Definition: AliVCluster.h:139
virtual Double_t * GetCellsAmplitudeFraction()
Definition: AliVCluster.h:122
virtual Double_t E() const
Definition: AliVCluster.h:75
virtual Double_t GetCellTime(Short_t cellNumber)=0
virtual Short_t GetCellNumber(Short_t pos) const =0
const AliExternalTrackParam * GetOuterParam() const
Definition: AliESDtrack.h:147
virtual Short_t GetNumberOfCells() const =0
void GetTriggerBits(Int_t &bits) const
Get the trigger bits for a given fastor position.
virtual Int_t GetNTracksMatched() const
Definition: AliVCluster.h:142
Bool_t kPrintCaloTrigger
Print trigger patches information.
Definition: TestESD.C:44
Virtual class for calorimeter cell data handling.
Definition: AliVCaloCells.h:19
void SetMisalMatrix(const TGeoHMatrix *m, Int_t smod) const
void EtaPhiFromIndex(Int_t absId, Double_t &eta, Double_t &phi) const
virtual Bool_t Next()
Forward to next trigger entry (fastor / L0 patch)
virtual Int_t GetCellMCLabel(Short_t cellNumber)=0
virtual Double_t GetTrackDz(void) const
Definition: AliVCluster.h:112
virtual UShort_t * GetCellsAbsId()
Definition: AliVCluster.h:120
AliTPCfastTrack * track
Virtual class for calorimeter cluster data handling.
Definition: AliVCluster.h:20
Bool_t kPrintClusterCells
Print cells in clusters information.
Definition: TestESD.C:47
void GetPosition(Int_t &col, Int_t &row) const
Access to position of the current fastor channel.
static AliRunLoader * Open(const char *filename="galice.root", const char *eventfoldername=AliConfig::GetDefaultEventFolderName(), Option_t *option="READ")
void GetCellPhiEtaIndexInSModule(Int_t nSupMod, Int_t nModule, Int_t nIphi, Int_t nIeta, Int_t &iphi, Int_t &ieta) const
void Print(const char *method, TStopwatch &timer, Int_t n)
virtual Int_t GetMCLabel(Short_t pos) const =0
virtual Int_t GetNtrack() const
Definition: AliStack.h:134
virtual Double_t GetAmplitude(Short_t pos) const =0
virtual Bool_t GetCellHighGain(Short_t cellNumber)=0
AliESDCaloTrigger * GetCaloTrigger(TString calo) const
Definition: AliESDEvent.h:509
Bool_t GetCellIndex(Int_t absId, Int_t &nSupMod, Int_t &nModule, Int_t &nIphi, Int_t &nIeta) const
AliESDCaloCells * GetEMCALCells() const
Definition: AliESDEvent.h:506
Bool_t kPrintCaloCluster
Print cluster parameters.
Definition: TestESD.C:46
Container with calorimeter trigger information in the ESD event.
TString GetFiredTriggerClasses() const
Definition: AliESDEvent.h:210
Int_t GetEMCALClusters(TRefArray *clusters) const
virtual Int_t * GetLabels() const
Definition: AliVCluster.h:128
Bool_t kPrintCaloCells
Print cells parameters.
Definition: TestESD.C:43
virtual Double_t GetCellAmplitude(Short_t cellNumber)=0
void TestESD()
Definition: TestESD.C:54
virtual Double_t GetTime(Short_t pos) const =0
virtual const Double_t * GetPID() const
Definition: AliVCluster.h:86
virtual Double_t GetTrackDx(void) const
Definition: AliVCluster.h:111
TF1 * f
Definition: interpolTest.C:21
virtual UInt_t GetNLabels() const
Definition: AliVCluster.h:129
virtual Int_t GetLabel() const
Definition: AliVCluster.h:126
Int_t GetEvent(Int_t evno)
virtual Double_t GetTOF() const
Definition: AliVCluster.h:106
TParticle * Particle(Int_t id, Bool_t useInEmbedding=kFALSE)
Definition: AliStack.cxx:688
virtual UChar_t GetNExMax() const
Definition: AliVCluster.h:103
void ReadFromTree(TTree *tree, Option_t *opt="")
virtual void Reset()
AliEMCALEMCGeometry * GetEMCGeometry() const
virtual Int_t GetTrackMatchedIndex(Int_t=0) const
Only for ESDs.
Definition: AliVCluster.h:148
virtual Double_t GetEFraction(Short_t pos) const =0
virtual UInt_t * GetCellsMCEdepFractionMap() const
Definition: AliVCluster.h:133
AliESDtrack * GetTrack(Int_t i) const
Definition: AliESDEvent.h:405
const TGeoHMatrix * GetEMCALMatrix(Int_t i) const
Definition: AliESDEvent.h:161
static AliEMCALGeometry * GetInstance()
TEveGeoShape * geom
Definition: tpc_tracks.C:10
virtual Int_t GetNCells() const
Definition: AliVCluster.h:118
virtual Bool_t GetHighGain(Short_t pos) const =0
void GetL1TimeSum(Int_t &timesum) const
Get the L1 time sums (L1 ADC values) for the current fastor.
Bool_t kPrintClusterPID
Print clusters PID (bayesian) weights.
Definition: TestESD.C:48
EMCal geometry, singleton.
Bool_t kPrintTrackMatches
Print cluster-track matching information.
Definition: TestESD.C:45
Bool_t kPrintKine
Print MC related information. Do not use for raw data.
Definition: TestESD.C:42
virtual Int_t GetCellAbsId(Int_t) const
Definition: AliVCluster.h:123
Int_t LoadKinematics(Option_t *option="READ")
AliStack * Stack() const
Definition: AliRunLoader.h:95
virtual void GetPosition(Float_t *) const
Definition: AliVCluster.h:82