AliRoot Core  edcc906 (edcc906)
TestEMCALData.C
Go to the documentation of this file.
1 
13 #if !defined(__CINT__) || defined(__MAKECINT__)
14 
15 //Root include files
16 #include <TH1I.h>
17 #include <TH1F.h>
18 #include <TCanvas.h>
19 #include <TClonesArray.h>
20 #include <Riostream.h>
21 
22 //AliRoot include files
23 #include "AliStack.h"
24 #include "AliRunLoader.h"
25 #include "AliEMCALLoader.h"
26 #include "AliEMCAL.h"
27 #include "AliEMCALHit.h"
28 #include "AliEMCALDigit.h"
29 #include "AliEMCALRecPoint.h"
30 
31 #endif
32 
36 void TestEMCALData() {
37 
38  TH1F* hE = new TH1F("hEmcalHits", "Hits energy distribution in EMCAL", 100, 0., 10.) ;
39  hE->Sumw2() ;
40  TH1I* hM = new TH1I("hEmcalHitsMul", "Hits multiplicity distribution in EMCAL", 500, 0., 10000) ;
41  hM->Sumw2() ;
42 
43  TH1I * dA = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL", 500, 0, 5000) ;
44  dA->Sumw2() ;
45  TH1I * dM = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL", 500, 0, 1000) ;
46  dM->Sumw2() ;
47 
48  TH1F * sE = new TH1F("hEmcalSDigits", "SDigits energy distribution in EMCAL", 200, 0, 1000) ;
49  sE->Sumw2() ;
50  TH1I * sM = new TH1I("hEmcalSDigitsMul", "SDigits multiplicity distribution in EMCAL", 500, 0, 1000) ;
51  sM->Sumw2() ;
52 
53  TH1F * cE = new TH1F("hEmcalRecPoints", "RecPoints energy distribution in EMCAL", 100, 0, 20) ;
54  cE->Sumw2() ;
55  TH1I * cM = new TH1I("hEmcalRecPointsMul", "RecPoints multiplicity distribution in EMCAL", 500, 0, 1000) ;
56  cM->Sumw2() ;
57 
59  if (rl == 0x0)
60  cout<<"Can not instatiate the Run Loader"<<endl;
61 
62  AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>
63  (rl->GetDetectorLoader("EMCAL"));
64 
65  // Load Hits
66  rl->LoadHits("EMCAL");
67  //Load Digits
68  rl->LoadDigits("EMCAL");
69  //Load SDigits
70  rl->LoadSDigits("EMCAL");
71  //Load RecPoints
72  rl->LoadRecPoints("EMCAL");
73 
74  // One way to do it that works for all three
75  AliEMCALHit* hit;
76  AliEMCALDigit* dig;
77  AliEMCALRecPoint* rp;
78 
79  rl->GetEvent(0);
80 
81  TClonesArray *hits = 0;
82  // Fill array of digits
83  TClonesArray *sdigits = emcalLoader->SDigits();
84  // Fill array of digits
85  TClonesArray *digits = emcalLoader->Digits();
86  // Fill array of clusters/rec-points
87  // TObjArray *clusters = emcalLoader->RecPoints(); //It should work, need to FIX
88  TObjArray *clusters = 0;
89 
90  //
91  // Get hits from the list
92  // Hits are stored in different branches in the hits Tree,
93  // first get the branch and then access the hits in the branch
94  //
95  Int_t nHit = 0;
96  TTree *treeH = emcalLoader->TreeH();
97  if (treeH)
98  {
99  // TreeH exists, get the branch
100  Int_t nTrack = treeH->GetEntries(); // TreeH has array of hits for every primary
101  TBranch * branchH = treeH->GetBranch("EMCAL");
102  branchH->SetAddress(&hits);
103 
104  // Now get the hits in this branch
105  for (Int_t iTrack = 0; iTrack < nTrack; iTrack++)
106  {
107  branchH->GetEntry(iTrack);
108  Int_t nHit = hits->GetEntriesFast();
109  for(Int_t ihit = 0; ihit< nHit;ihit++)
110  {
111  hit = static_cast<AliEMCALHit *> (hits->At(ihit));//hits->At(ihit);
112  if(!hit) continue;
113 
114  nHit++;
115  hE->Fill(hit->GetEnergy());
116  cout<<"Hit Info "<<hit->GetId()<<" ELoss "<<hit->GetEnergy()<<endl;
117  }//hit loop
118  }// track loop
119  }//treeH?
120 
121  hM->Fill(nHit);
122 
123  //
124  // Get sdigits from the list
125  //
126  if(sdigits)
127  {
128  sM->Fill(sdigits->GetEntries());
129  for(Int_t isdig = 0; isdig< sdigits->GetEntries();isdig++)
130  {
131  //cout<<">> idig "<<idig<<endl;
132  dig = static_cast<AliEMCALDigit *>(sdigits->At(isdig)) ;
133 
134  if(!dig) continue;
135 
136  sE->Fill(dig->GetAmplitude());
137  cout << "SDigit info " << dig->GetId() << " " << dig->GetAmplitude() << endl;
138  }
139  }
140  else printf("No Sdigits available\n");
141 
142  //
143  // Get digits from the list
144  //
145  if(digits)
146  {
147  dM->Fill(digits->GetEntries());
148  for(Int_t idig = 0; idig< digits->GetEntries();idig++)
149  {
150  //cout<<">> idig "<<idig<<endl;
151  dig = static_cast<AliEMCALDigit *>(digits->At(idig)) ;
152 
153  if(!dig) continue;
154  dA->Fill(dig->GetAmplitude());
155  cout << "Digit info " << dig->GetId() << " " << dig->GetAmplitude() << endl;
156  }
157  }
158  else printf("No digits available\n");
159 
160  // Get clusters from the list
161  TTree *treeR = emcalLoader->TreeR();
162  TBranch * branchR = treeR->GetBranch("EMCALECARP");
163  branchR->SetAddress(&clusters);
164  branchR->GetEntry(0);
165 
166  if(clusters)
167  {
168  cM->Fill(clusters->GetEntries());
169  for(Int_t iclu = 0; iclu< clusters->GetEntries();iclu++)
170  {
171  //cout<<">> idig "<<idig<<endl;
172  rp = static_cast<AliEMCALRecPoint *>(clusters->At(iclu)) ;
173 
174  if(!rp) continue;
175  cE->Fill(rp->GetEnergy());
176  cout << "RecPoint info " << rp->GetAbsId(0) << " " << rp->GetEnergy() << endl;
177  }
178  }
179  else printf("No recpoints available\n");
180 
181  /*
182  //another way to do it
183  TTree *hTree = rl->GetTreeH("EMCAL",false);
184  TTree *dTree = rl->GetTreeD("EMCAL",false);
185  TTree *sTree = rl->GetTreeS("EMCAL",false);
186  TTree *cTree = rl->GetTreeR("EMCAL",false);
187 
188  TObjArray *harr=NULL;
189  TBranch *hbranch=hTree->GetBranch("EMCAL");
190  hbranch->SetAddress(&harr);
191 
192  TObjArray *darr=NULL;
193  TBranch *dbranch=dTree->GetBranch("EMCAL");
194  dbranch->SetAddress(&darr);
195 
196  TObjArray *sarr=NULL;
197  TBranch *sbranch=sTree->GetBranch("EMCAL");
198  sbranch->SetAddress(&sarr);
199 
200  TObjArray *carr=NULL;
201  TBranch *cbranch=cTree->GetBranch("EMCALECARP");
202  cbranch->SetAddress(&carr);
203 
204 
205  if(hbranch->GetEvent(0)) {
206  for(Int_t ih = 0; ih < harr->GetEntriesFast(); ih++) {
207  hM->Fill(harr->GetEntriesFast());
208  AliEMCALHit* hit =(AliEMCALHit*)harr->UncheckedAt(ih);
209  if(hit != 0){
210  hE->Fill(hit->GetEnergy());
211  cout << "Hit info " << hit->GetId() << " " << hit->GetEnergy()*10.5 << endl;
212  }
213  }
214  }
215 
216  if(dbranch->GetEvent(0)) {
217  for(Int_t id = 0; id < darr->GetEntriesFast(); id++) {
218  dM->Fill(darr->GetEntriesFast());
219  AliEMCALDigit* dig =(AliEMCALDigit*)darr->UncheckedAt(id);
220  if(dig != 0){
221  dA->Fill(dig->GetAmp());
222  cout << "Digit info " << dig->GetId() << " " << dig->GetAmp() << endl;
223  }
224  }
225  }
226 
227  if(sbranch->GetEvent(0)) {
228  for(Int_t id = 0; id < sarr->GetEntriesFast(); id++) {
229  sM->Fill(sarr->GetEntriesFast());
230  AliEMCALDigit* sdig =(AliEMCALDigit*)sarr->UncheckedAt(id);
231  if(sdig != 0){
232  sE->Fill(sdig->GetAmp()/1.e+6);
233  cout << "SDigit info " << sdig->GetId() << " " << sdig->GetAmp()/1.e+6 << endl;
234  }
235  }
236  }
237 
238  if(cbranch->GetEvent(0)) {
239  for(Int_t ic = 0; ic < carr->GetEntriesFast(); ic++) {
240  cE->Fill(carr->GetEntriesFast());
241  AliEMCALRecPoint* rp =(AliEMCALRecPoint*)carr->UncheckedAt(ic);
242  if(rp != 0){
243  cE->Fill(rp->GetEnergy());
244  cout << "RecPoint info " << rp->GetAbsId() << " " << rp->GetEnergy() << endl;
245  }
246  }
247  }
248 
249  */
250 
251  TCanvas *chits = new TCanvas("chits","Hits",20,20,800,400);
252  chits->Divide(2,1);
253  chits->cd(1);
254  hE->Draw();
255  chits->cd(2);
256  hM->Draw();
257 
258  TCanvas *cdig = new TCanvas("cdig","Digits",20,40,800,400);
259  cdig->Divide(2,1);
260  cdig->cd(1);
261  dA->Draw();
262  cdig->cd(2);
263  dM->Draw();
264 
265 
266  TCanvas *csdig = new TCanvas("csdig","SDigits",20,60,800,400);
267  csdig->Divide(2,1);
268  csdig->cd(1);
269  sE->Draw();
270  csdig->cd(2);
271  sM->Draw();
272 
273  TCanvas *cclu = new TCanvas("cclu","Clusters",20,60,800,400);
274  cclu->Divide(2,1);
275  cclu->cd(1);
276  cE->Draw();
277  cclu->cd(2);
278  cM->Draw();
279 }
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
static const TString & GetDefaultEventFolderName()
Definition: AliConfig.h:58
#define TObjArray
Int_t LoadSDigits(Option_t *detectors="all", Option_t *opt="READ")
TClonesArray * SDigits()
Int_t LoadDigits(Option_t *detectors="all", Option_t *opt="READ")
AliLoader * GetDetectorLoader(const char *detname)
virtual Float_t GetEnergy() const
void TestEMCALData()
Definition: TestEMCALData.C:36
EMCal digits object.
Definition: AliEMCALDigit.h:30
Int_t GetId(void) const
Definition: AliEMCALHit.h:44
TClonesArray * Digits()
Give access to hits, digits, recpoints arrays and OCDB.
static AliRunLoader * Open(const char *filename="galice.root", const char *eventfoldername=AliConfig::GetDefaultEventFolderName(), Option_t *option="READ")
TTree * TreeR() const
Definition: AliLoader.h:89
EMCal hits object.
Definition: AliEMCALHit.h:24
Int_t * GetAbsId() const
Int_t LoadRecPoints(Option_t *detectors="all", Option_t *opt="READ")
EMCal rec_points object.
Int_t GetEvent(Int_t evno)
TTree * TreeH() const
Definition: AliLoader.h:83
Int_t LoadHits(Option_t *detectors="all", Option_t *opt="READ")
Int_t GetId() const
Definition: AliDigitNew.h:23
Float_t GetEnergy(void) const
Definition: AliEMCALHit.h:41
Float_t GetAmplitude() const
Definition: AliEMCALDigit.h:55