]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/TPC/AliTPCClusterHistograms.cxx
Bug in QA class corrected
[u/mrichter/AliRoot.git] / PWG0 / TPC / AliTPCClusterHistograms.cxx
1 /* $Id$ */
2
3 // This class contains a number of histograms for diagnostics of a TPC
4 // read out chamber from the reconstructed clusters.
5 //
6 // TODO:
7 //
8 //
9 //
10
11 #include "AliTPCClusterHistograms.h"
12
13 #include <TStyle.h>
14 #include <TFile.h>
15 #include <TCanvas.h>
16 #include <TH2F.h>
17 #include <TProfile2D.h>
18 #include <TObjArray.h>
19 #include <TLatex.h>
20 #include <TTimeStamp.h>
21 #include <TRandom.h>
22
23 #include <AliTPCclusterMI.h>
24 #include <AliTPCseed.h>
25
26 #include <AliLog.h>
27
28 // uncomment if working with version after release that has new cluster scheme
29 //#define NEWALIROOT
30
31 //____________________________________________________________________
32 ClassImp(AliTPCClusterHistograms)
33
34 //____________________________________________________________________
35 AliTPCClusterHistograms::AliTPCClusterHistograms() 
36   : TNamed(),
37   fTimeStart(0),
38   fTimeStop(0),
39   fhQmaxVsRow(0),
40   fhQtotVsRow(0),          
41   fhQtotProfileVsRow(0),   
42   fhQmaxProfileVsRow(0),
43   fhQtotVsDistanceToEdge(0),
44   fhQtotProfileVsDistanceToEdge(0),
45   fhNClustersYVsRow(0),
46   fhNClustersZVsRow(0),
47   fhSigmaYVsRow(0),        
48   fhSigmaZVsRow(0),                             
49   fhQmaxProfileYVsRow(0), 
50   fhQtotProfileYVsRow(0),
51   fhSigmaYProfileYVsRow(0),
52   fhSigmaZProfileYVsRow(0),
53   fhQmaxProfileZVsRow(0), 
54   fhQtotProfileZVsRow(0),
55   fhSigmaYProfileZVsRow(0),
56   fhSigmaZProfileZVsRow(0),
57   fhMeanQtotVsTime(0),  
58   fhQtotVsTime(0),
59   fhMeanNClustersVsTime(0),    
60   fhNClustersVsTime(0),        
61   fhTrackQtotPerCluster(0),
62   fhTrackQtotPerClusterVsPhi(0),
63   fhTrackQtotPerClusterVsTheta(0),
64   fhTrackMeanQtotPerClusterVsPhi(0),
65   fhTrackMeanQtotPerClusterVsTheta(0),
66   fhMeanNTracksVsTime(),
67   fhNEventsVsTime(),
68   fDetector(0),
69   fIsIROC(kFALSE),
70   fEdgeSuppression(kFALSE),
71   fNClustersInEvent(0),
72   fQtotInEvent(0),
73   fMaxQtotInEvent(0),
74   fKeepEvent(0),
75   fWhyKeepEvent(),
76   fCommentToHistograms()
77 {
78   // default constructor
79 }
80
81 //____________________________________________________________________
82 AliTPCClusterHistograms::AliTPCClusterHistograms(Int_t detector, const Char_t* comment, Int_t timeStart, Int_t timeStop, Bool_t edgeSuppression)
83   : TNamed(),
84   fTimeStart(0),
85   fTimeStop(0),
86   fhQmaxVsRow(0),
87   fhQtotVsRow(0),          
88   fhQtotProfileVsRow(0),   
89   fhQmaxProfileVsRow(0),
90   fhQtotVsDistanceToEdge(0),
91   fhQtotProfileVsDistanceToEdge(0),
92   fhNClustersYVsRow(0),
93   fhNClustersZVsRow(0),
94   fhSigmaYVsRow(0),        
95   fhSigmaZVsRow(0),                             
96   fhQmaxProfileYVsRow(0), 
97   fhQtotProfileYVsRow(0),
98   fhSigmaYProfileYVsRow(0),
99   fhSigmaZProfileYVsRow(0),
100   fhQmaxProfileZVsRow(0), 
101   fhQtotProfileZVsRow(0),
102   fhSigmaYProfileZVsRow(0),
103   fhSigmaZProfileZVsRow(0),
104   fhMeanQtotVsTime(0),  
105   fhQtotVsTime(0),
106   fhMeanNClustersVsTime(0),    
107   fhNClustersVsTime(0),        
108   fhTrackQtotPerCluster(0),
109   fhTrackQtotPerClusterVsPhi(0),
110   fhTrackQtotPerClusterVsTheta(0),
111   fhTrackMeanQtotPerClusterVsPhi(0),
112   fhTrackMeanQtotPerClusterVsTheta(0),
113   fhMeanNTracksVsTime(),
114   fhNEventsVsTime(),
115   fDetector(0),
116   fIsIROC(kFALSE),
117   fEdgeSuppression(edgeSuppression),
118   fNClustersInEvent(0),
119   fQtotInEvent(0),
120   fMaxQtotInEvent(0),
121   fKeepEvent(0),
122   fWhyKeepEvent(),
123   fCommentToHistograms()
124 {
125   // constructor 
126   
127   // make name and title
128   if (detector < 0 || detector >= 72) {
129     AliDebug(AliLog::kError, Form("Detector %d does not exist", detector));
130     return;
131   }
132       
133   TString name(FormDetectorName(detector, edgeSuppression, comment));
134
135   fNClustersInEvent = 0;
136   fQtotInEvent      = 0;
137   fMaxQtotInEvent   = 0;
138
139   fKeepEvent = kFALSE;
140   fWhyKeepEvent = TString("hi");
141
142   fDetector = detector;
143   if (detector < 36)
144     fIsIROC = kTRUE; 
145   
146   SetName(name);
147   SetTitle(Form("%s (detector %d)",name.Data(), detector));
148
149   // rounding down to the closest hour and starting 10 hours before
150   fTimeStart = 3600*UInt_t(timeStart/3600) - 36000;
151   // rounding up to the closest hour
152   fTimeStop  = 3600*UInt_t((3600 + timeStop)/3600);
153   // each time bin covers 5 min
154   Int_t nTimeBins = (fTimeStop-fTimeStart)/300;
155   
156   //  printf(Form(" start time: %d,  stop time: %d \n",fTimeStart, fTimeStop));
157
158   #define BINNING_Z 250, 0, 250
159   
160   Float_t yRange   = 45;
161   Int_t nPadRows   = 96;
162   
163   if (fIsIROC)
164   {
165     yRange   = 25;
166     nPadRows = 63;
167   }
168   
169   // 1 bin for each 0.5 cm
170   Int_t nBinsY = Int_t(4*yRange);
171
172   // do not add this hists to the directory
173   Bool_t oldStatus = TH1::AddDirectoryStatus();
174   TH1::AddDirectory(kFALSE);
175
176   //defining histograms and profile plots
177   fhQmaxVsRow  = new TH2F("QmaxVsPadRow", "Qmax vs. pad row;Pad row;Qmax", nPadRows+2, -1.5, nPadRows+0.5, 500,  0,  500);
178   fhQtotVsRow  = new TH2F("QtotVsPadRow", "Qtot vs. pad row;Pad row;Qtot", nPadRows+2, -1.5, nPadRows+0.5, 400,  0,  4000);
179
180   fhQmaxProfileVsRow = new TProfile("MeanQmaxVsPadRow","Mean Qmax vs. pad row;Pad row;Mean Qmax",nPadRows+2, -1.5, nPadRows+0.5);
181   fhQtotProfileVsRow = new TProfile("MeanQtotVsPadRow","Mean Qtot vs. pad row;Pad row;Mean Qtot",nPadRows+2, -1.5, nPadRows+0.5);
182   
183   fhNClustersYVsRow = new TH2F("NClusters y vs pad row","N clusters y vs pad;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
184   fhNClustersZVsRow = new TH2F("NClusters z vs pad row","N clusters z vs pad;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
185
186   fhSigmaYVsRow = new TH2F("SigmaYVsPadRow", "Sigma Y vs. pad row;Pad row;#sigma_{Y}", nPadRows+2, -1.5, nPadRows+0.5, 100,  0,  0.5);
187   fhSigmaZVsRow = new TH2F("SigmaZVsPadRow", "Sigma Z vs. pad row;Pad row;#sigma_{Z}", nPadRows+2, -1.5, nPadRows+0.5, 100,  0,  0.5);
188   
189   fhQmaxProfileYVsRow = new TProfile2D("MeanQmaxYVsPadRow","Mean Qmax, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
190   fhQtotProfileYVsRow = new TProfile2D("MeanQtotYVsPadRow","Mean Qtot, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
191   fhSigmaYProfileYVsRow = new TProfile2D("MeanSigmaYYVsPadRow","Mean Sigma y, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
192   fhSigmaZProfileYVsRow = new TProfile2D("MeanSigmaZYVsPadRow","Mean Sigma z, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
193
194   fhQmaxProfileZVsRow = new TProfile2D("MeanQmaxZVsPadRow","Mean Qmax, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
195   fhQtotProfileZVsRow = new TProfile2D("MeanQtotZVsPadRow","Mean Qtot, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
196   fhSigmaYProfileZVsRow = new TProfile2D("MeanSigmaYZVsPadRow","Mean Sigma y, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
197   fhSigmaZProfileZVsRow = new TProfile2D("MeanSigmaZZVsPadRow","Mean Sigma z, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
198
199
200   TString start(TTimeStamp(fTimeStart).AsString());
201   //  TString stop(TTimeStamp(fTimeStart).AsString());
202   start.Remove(26);
203   
204   fhMeanQtotVsTime = new TProfile("MeanQtotVsTime",Form("Mean Qtot vs. time (start %s , 1 min bins); time; Qtot",start.Data()),5*nTimeBins, fTimeStart, fTimeStop);
205   fhQtotVsTime     = new TH2F("QtotVsTime",Form("Qtot vs. time (start %s , 1 min bins); time; Qtot",start.Data()),5*nTimeBins, fTimeStart, fTimeStop,400,0,2000);
206
207   fhMeanNClustersVsTime = new TProfile("MeanNClustersVsTime",Form("Mean N Cluster vs. time (start %s , 5 min bins); time; NClusters",start.Data()),nTimeBins, fTimeStart, fTimeStop);
208   fhNClustersVsTime = new TH2F("NClustersVsTime",Form("N Clusters vs. time (start %s , 5 min bins); time; NClusters",start.Data()),nTimeBins, fTimeStart, fTimeStop,400,-0.5,3999.5);
209
210   fhQmaxProfileVsRow->SetLineWidth(2);
211   fhQtotProfileVsRow->SetLineWidth(2);
212
213   fhMeanQtotVsTime->SetLineWidth(2);
214
215   // histograms related to tracks
216
217   fhTrackQtotPerCluster = new TH1F("QtotPerCluster","Qtot per cluster; (Sum Qtot)/clusters",200,0,2000);
218   fhTrackQtotPerCluster->SetMarkerStyle(22);
219   fhTrackQtotPerCluster->SetMarkerSize(1);
220
221   fhTrackQtotPerClusterVsPhi = new TH2F("QtotPerClusterVsPhi","QtotPerCluster vs Phi; Phi; (Sum Qtot)/clusters",40,-2,2,200,0,2000);
222   fhTrackQtotPerClusterVsTheta = new TH2F("QtotPerClusterVsTheta","QtotPerCluster vs Theta; Theta; (Sum Qtot)/clusters",40,-2,2,200,0,2000);
223
224   fhTrackMeanQtotPerClusterVsPhi = new TProfile("MeanQtotPerClusterVsPhi", "QtotPerCluster vs Phi; Phi; Mean (Sum Qtot)/clusters",40,-2,2);
225   fhTrackMeanQtotPerClusterVsTheta = new TProfile("MeanQtotPerClusterVsTheta", "QtotPerCluster vs Theta; Theta; Mean (Sum Qtot)/clusters",40,-2,2);
226
227   fhTrackMeanQtotPerClusterVsPhi->SetLineWidth(2);
228   fhTrackMeanQtotPerClusterVsTheta->SetLineWidth(2);
229
230   fhMeanNTracksVsTime = new TProfile("MeanNTracksVsTime",Form("Mean n tracks vs. time (start %s , 5 min bins); time; N tracks",start.Data()),nTimeBins, fTimeStart, fTimeStop);
231
232   fhNEventsVsTime = new TH1F("NEventsVsTime",Form("N events vs. time (start %s , 5 min bins); time; N events",start.Data()),nTimeBins, fTimeStart, fTimeStop);
233
234   TH1::AddDirectory(oldStatus);
235 }
236
237 //____________________________________________________________________
238 AliTPCClusterHistograms::AliTPCClusterHistograms(const AliTPCClusterHistograms& c) : TNamed(c),
239   fTimeStart(0),
240   fTimeStop(0),
241   fhQmaxVsRow(0),
242   fhQtotVsRow(0),          
243   fhQtotProfileVsRow(0),   
244   fhQmaxProfileVsRow(0),
245   fhQtotVsDistanceToEdge(0),
246   fhQtotProfileVsDistanceToEdge(0),
247   fhNClustersYVsRow(0),
248   fhNClustersZVsRow(0),
249   fhSigmaYVsRow(0),        
250   fhSigmaZVsRow(0),                             
251   fhQmaxProfileYVsRow(0), 
252   fhQtotProfileYVsRow(0),
253   fhSigmaYProfileYVsRow(0),
254   fhSigmaZProfileYVsRow(0),
255   fhQmaxProfileZVsRow(0), 
256   fhQtotProfileZVsRow(0),
257   fhSigmaYProfileZVsRow(0),
258   fhSigmaZProfileZVsRow(0),
259   fhMeanQtotVsTime(0),  
260   fhQtotVsTime(0),
261   fhMeanNClustersVsTime(0),    
262   fhNClustersVsTime(0),        
263   fhTrackQtotPerCluster(0),
264   fhTrackQtotPerClusterVsPhi(0),
265   fhTrackQtotPerClusterVsTheta(0),
266   fhTrackMeanQtotPerClusterVsPhi(0),
267   fhTrackMeanQtotPerClusterVsTheta(0),
268   fhMeanNTracksVsTime(),
269   fhNEventsVsTime(),
270   fDetector(0),
271   fIsIROC(kFALSE),
272   fEdgeSuppression(kFALSE),
273   fNClustersInEvent(0),
274   fQtotInEvent(0),
275   fMaxQtotInEvent(0),
276   fKeepEvent(0),
277   fWhyKeepEvent(),
278   fCommentToHistograms()
279 {
280   // copy constructor
281   ((AliTPCClusterHistograms &)c).Copy(*this);
282 }
283
284 //____________________________________________________________________
285 AliTPCClusterHistograms::~AliTPCClusterHistograms()
286 {
287   //
288   // destructor
289   //
290
291   if (fhQmaxVsRow) {
292     delete fhQmaxVsRow;
293     fhQmaxVsRow = 0;
294   }
295   if (fhQtotVsRow) {
296     delete fhQtotVsRow;
297     fhQtotVsRow = 0; 
298   }
299   if (fhQmaxProfileVsRow) {
300     delete fhQmaxProfileVsRow;
301     fhQmaxProfileVsRow = 0;
302   }
303   if (fhQtotProfileVsRow) {
304     delete fhQtotProfileVsRow;
305     fhQtotProfileVsRow = 0;
306   }
307   if (fhNClustersYVsRow) {
308     delete fhNClustersYVsRow;
309     fhNClustersYVsRow = 0;
310   }
311   if (fhNClustersZVsRow) {
312     delete fhNClustersZVsRow;
313     fhNClustersZVsRow = 0;
314   }
315   if (fhSigmaYVsRow) {
316     delete fhSigmaYVsRow;
317     fhSigmaYVsRow = 0;
318   } 
319   if (fhSigmaZVsRow) {
320     delete fhSigmaZVsRow;
321     fhSigmaZVsRow = 0; 
322   }
323   if (fhQmaxProfileYVsRow) {
324     delete fhQmaxProfileYVsRow;
325     fhQmaxProfileYVsRow = 0;
326   }
327   if (fhQtotProfileYVsRow) {
328     delete fhQtotProfileYVsRow;
329     fhQtotProfileYVsRow = 0;
330   }
331   if (fhSigmaYProfileYVsRow) {
332     delete fhSigmaYProfileYVsRow;
333     fhSigmaYProfileYVsRow = 0;
334   }
335   if (fhSigmaZProfileYVsRow) {
336     delete fhSigmaZProfileYVsRow;
337     fhSigmaZProfileYVsRow = 0;
338   }
339   if (fhQmaxProfileZVsRow) {
340     delete fhQmaxProfileZVsRow;
341     fhQmaxProfileZVsRow = 0;
342   }
343   if (fhQtotProfileZVsRow) {
344     delete fhQtotProfileZVsRow;
345     fhQtotProfileZVsRow = 0;
346   }
347   if (fhSigmaYProfileZVsRow) {
348     delete fhSigmaYProfileZVsRow;
349     fhSigmaYProfileZVsRow = 0;
350   }
351   if (fhSigmaZProfileZVsRow) {
352     delete fhSigmaZProfileZVsRow;
353     fhSigmaZProfileZVsRow = 0;
354   }
355   if (fhMeanQtotVsTime) {
356     delete fhMeanQtotVsTime;
357     fhMeanQtotVsTime = 0;
358   }
359   if (fhQtotVsTime) {
360     delete fhQtotVsTime;
361     fhQtotVsTime = 0;
362   }
363   if (fhMeanNClustersVsTime) {
364     delete fhMeanNClustersVsTime;
365     fhMeanNClustersVsTime = 0;
366   }
367   if (fhNClustersVsTime) {
368     delete fhNClustersVsTime;
369     fhNClustersVsTime = 0;
370   }
371   if (fhTrackQtotPerCluster) {
372     delete fhTrackQtotPerCluster;
373     fhTrackQtotPerCluster = 0;
374   }
375   if (fhTrackQtotPerClusterVsPhi) {
376     delete fhTrackQtotPerClusterVsPhi;
377     fhTrackQtotPerClusterVsPhi = 0;
378   }
379   if (fhTrackQtotPerClusterVsTheta) {
380     delete fhTrackQtotPerClusterVsTheta;
381     fhTrackQtotPerClusterVsTheta = 0;
382   }
383   if (fhTrackMeanQtotPerClusterVsPhi) {
384     delete fhTrackMeanQtotPerClusterVsPhi;
385     fhTrackMeanQtotPerClusterVsPhi = 0;
386   }
387   if (fhTrackMeanQtotPerClusterVsTheta) {
388     delete fhTrackMeanQtotPerClusterVsTheta;
389     fhTrackMeanQtotPerClusterVsTheta = 0;
390   }
391   if (fhMeanNTracksVsTime) {
392     delete fhMeanNTracksVsTime;
393     fhMeanNTracksVsTime = 0;
394   }  
395   if (fhNEventsVsTime) {
396     delete fhNEventsVsTime;
397     fhNEventsVsTime = 0;
398   }  
399 }
400
401 //____________________________________________________________________
402 AliTPCClusterHistograms &AliTPCClusterHistograms::operator=(const AliTPCClusterHistograms &c)
403 {
404   // assigment operator
405
406   if (this != &c)
407     ((AliTPCClusterHistograms &) c).Copy(*this);
408
409   return *this;
410 }
411
412 //____________________________________________________________________
413 TString AliTPCClusterHistograms::FormDetectorName(Int_t detector, Bool_t edgeSuppression, const char* comment)
414 {
415   //
416   // creates a readable name from the detector number
417   //   
418   
419   Int_t sector = detector%18;
420   TString side;
421   TString inout;
422   
423   if (detector<18 || ( detector>=36 && detector<54))
424     side.Form("A");
425   else 
426     side.Form("C");
427   
428   if (detector<36)
429     inout.Form("IROC");
430   else 
431     inout.Form("OROC");
432
433   TString name;
434   name.Form("sector_%s%d_%s", side.Data(), sector, inout.Data());
435
436   if (edgeSuppression)
437     name += "_noedge";
438   
439   if (comment)
440     name += comment;
441
442   return name; 
443 }
444
445 //____________________________________________________________________
446 Long64_t AliTPCClusterHistograms::Merge(TCollection* list)
447 {
448   // Merge a list of AliTPCClusterHistograms objects with this (needed for
449   // PROOF). 
450   // Returns the number of merged objects (including this).
451
452   if (!list)
453     return 0;
454   
455   if (list->IsEmpty())
456     return 1;
457
458   TIterator* iter = list->MakeIterator();
459   TObject* obj;
460
461   // collections of measured and generated histograms
462   TList* collectionQmaxVsRow     = new TList;
463   TList* collectionQtotVsRow     = new TList;
464
465   TList* collectionQmaxProfileVsRow = new TList;
466   TList* collectionQtotProfileVsRow = new TList;
467
468   TList* collectionNClustersYVsRow = new TList;
469   TList* collectionNClustersZVsRow = new TList;
470
471   TList* collectionSigmaYVsRow   = new TList;
472   TList* collectionSigmaZVsRow   = new TList;
473                                         
474   TList* collectionQmaxProfileYVsRow    = new TList;
475   TList* collectionQtotProfileYVsRow    = new TList;
476   TList* collectionSigmaYProfileYVsRow  = new TList;
477   TList* collectionSigmaZProfileYVsRow  = new TList;
478
479   TList* collectionQmaxProfileZVsRow    = new TList;
480   TList* collectionQtotProfileZVsRow    = new TList;
481   TList* collectionSigmaYProfileZVsRow  = new TList;
482   TList* collectionSigmaZProfileZVsRow  = new TList;
483
484   TList* collectionMeanQtotVsTime  = new TList;
485   TList* collectionQtotVsTime      = new TList;
486
487   TList* collectionMeanNClustersVsTime = new TList; 
488   TList* collectionNClustersVsTime     = new TList;
489
490   TList* collectionTrackQtotPerCluster = new TList;
491
492   TList* collectionTrackQtotPerClusterVsPhi = new TList;
493   TList* collectionTrackQtotPerClusterVsTheta = new TList;
494
495   TList* collectionTrackMeanQtotPerClusterVsPhi = new TList;
496   TList* collectionTrackMeanQtotPerClusterVsTheta = new TList;
497
498   TList* collectionMeanNTracksVsTime = new TList;
499   TList* collectionNEventsVsTime = new TList;
500
501   Int_t count = 0;
502   while ((obj = iter->Next())) {
503     
504     AliTPCClusterHistograms* entry = dynamic_cast<AliTPCClusterHistograms*> (obj);
505     if (entry == 0) 
506       continue;
507     
508     collectionQmaxVsRow          ->Add(entry->fhQmaxVsRow          );
509     collectionQtotVsRow   ->Add(entry->fhQtotVsRow         );
510     
511     collectionQmaxProfileVsRow   ->Add(entry->fhQmaxProfileVsRow  );
512     collectionQtotProfileVsRow    ->Add(entry->fhQtotProfileVsRow  );
513     
514     collectionNClustersYVsRow    ->Add(entry->fhNClustersYVsRow);
515     collectionNClustersZVsRow    ->Add(entry->fhNClustersZVsRow);
516     
517     collectionSigmaYVsRow         ->Add(entry->fhSigmaYVsRow       );
518     collectionSigmaZVsRow         ->Add(entry->fhSigmaZVsRow       );
519     
520     collectionQmaxProfileYVsRow  ->Add(entry->fhQmaxProfileYVsRow );
521     collectionQtotProfileYVsRow  ->Add(entry->fhQtotProfileYVsRow );
522     collectionSigmaYProfileYVsRow->Add(entry->fhSigmaYProfileYVsRow);
523     collectionSigmaZProfileYVsRow->Add(entry->fhSigmaZProfileYVsRow);
524     
525     collectionQmaxProfileZVsRow  ->Add(entry->fhQmaxProfileZVsRow );
526     collectionQtotProfileZVsRow  ->Add(entry->fhQtotProfileZVsRow );
527     collectionSigmaYProfileZVsRow->Add(entry->fhSigmaYProfileZVsRow);
528     collectionSigmaZProfileZVsRow->Add(entry->fhSigmaZProfileZVsRow);
529     
530     collectionMeanQtotVsTime     ->Add(entry->fhMeanQtotVsTime);
531     collectionQtotVsTime         ->Add(entry->fhQtotVsTime);
532
533     collectionMeanNClustersVsTime->Add(entry->fhMeanNClustersVsTime);  
534     collectionNClustersVsTime    ->Add(entry->fhNClustersVsTime);
535     
536     collectionTrackQtotPerCluster->Add(entry->fhTrackQtotPerCluster);
537     
538     collectionTrackQtotPerClusterVsPhi->Add(entry->fhTrackQtotPerClusterVsPhi);
539     collectionTrackQtotPerClusterVsTheta->Add(entry->fhTrackQtotPerClusterVsTheta);
540     
541     collectionTrackMeanQtotPerClusterVsPhi->Add(entry->fhTrackMeanQtotPerClusterVsPhi);
542     collectionTrackMeanQtotPerClusterVsTheta->Add(entry->fhTrackMeanQtotPerClusterVsTheta);
543     
544     collectionMeanNTracksVsTime->Add(entry->fhMeanNTracksVsTime);
545     collectionNEventsVsTime->Add(entry->fhNEventsVsTime);
546     
547     count++;
548   }
549   
550   fhQmaxVsRow          ->Merge(collectionQmaxVsRow       );        
551   fhQtotVsRow          ->Merge(collectionQtotVsRow        );       
552   
553   fhQmaxProfileVsRow   ->Merge(collectionQtotProfileVsRow);
554   fhQtotProfileVsRow   ->Merge(collectionQtotProfileVsRow);
555   
556   fhNClustersYVsRow    ->Merge(collectionNClustersYVsRow);
557   fhNClustersZVsRow    ->Merge(collectionNClustersZVsRow);
558   
559   fhSigmaYVsRow        ->Merge(collectionSigmaYVsRow      );       
560   fhSigmaZVsRow        ->Merge(collectionSigmaZVsRow      );       
561   
562   fhQmaxProfileYVsRow  ->Merge(collectionQmaxProfileYVsRow  ); 
563   fhQtotProfileYVsRow  ->Merge(collectionQtotProfileYVsRow  );
564   fhSigmaYProfileYVsRow->Merge(collectionSigmaYProfileYVsRow);
565   fhSigmaZProfileYVsRow->Merge(collectionSigmaZProfileYVsRow);
566   
567   fhQmaxProfileZVsRow  ->Merge(collectionQmaxProfileZVsRow  ); 
568   fhQtotProfileZVsRow  ->Merge(collectionQtotProfileZVsRow  );
569   fhSigmaYProfileZVsRow->Merge(collectionSigmaYProfileZVsRow);
570   fhSigmaZProfileZVsRow->Merge(collectionSigmaZProfileZVsRow);
571   
572   fhMeanQtotVsTime     ->Merge(collectionMeanQtotVsTime);
573   fhQtotVsTime         ->Merge(collectionQtotVsTime);
574
575   fhMeanNClustersVsTime->Merge(collectionMeanNClustersVsTime );
576   fhNClustersVsTime    ->Merge(collectionNClustersVsTime);
577   
578   fhTrackQtotPerCluster->Merge(collectionTrackQtotPerCluster);
579   
580   fhTrackQtotPerClusterVsPhi->Merge(collectionTrackQtotPerClusterVsPhi);
581   fhTrackQtotPerClusterVsTheta->Merge(collectionTrackQtotPerClusterVsTheta);
582   
583   fhTrackMeanQtotPerClusterVsPhi->Merge(collectionTrackMeanQtotPerClusterVsPhi);
584   fhTrackMeanQtotPerClusterVsTheta->Merge(collectionTrackMeanQtotPerClusterVsTheta);
585   
586   fhMeanNTracksVsTime->Merge(collectionMeanNTracksVsTime);
587   fhNEventsVsTime->Merge(collectionNEventsVsTime);
588   
589   delete collectionQmaxVsRow;          
590   delete collectionQtotVsRow;  
591   
592   delete collectionQmaxProfileVsRow;
593   delete collectionQtotProfileVsRow;
594   
595   delete collectionNClustersYVsRow;
596   delete collectionNClustersZVsRow;
597   
598   delete collectionSigmaYVsRow;   
599   delete collectionSigmaZVsRow;   
600   
601   delete collectionQmaxProfileYVsRow;  
602   delete collectionQtotProfileYVsRow;  
603   delete collectionSigmaYProfileYVsRow;
604   delete collectionSigmaZProfileYVsRow;
605   
606   delete collectionQmaxProfileZVsRow;  
607   delete collectionQtotProfileZVsRow;  
608   delete collectionSigmaYProfileZVsRow;
609   delete collectionSigmaZProfileZVsRow;
610   
611   delete collectionMeanQtotVsTime;
612   delete collectionQtotVsTime;
613
614   delete collectionMeanNClustersVsTime; 
615   delete collectionNClustersVsTime;     
616   
617   delete collectionTrackQtotPerCluster;
618   
619   delete collectionTrackQtotPerClusterVsPhi; 
620   delete collectionTrackQtotPerClusterVsTheta;
621   
622   delete collectionTrackMeanQtotPerClusterVsPhi; 
623   delete collectionTrackMeanQtotPerClusterVsTheta;
624   
625   delete collectionMeanNTracksVsTime;
626   delete collectionNEventsVsTime;
627
628   return count+1;
629 }
630
631
632 //____________________________________________________________________
633 void AliTPCClusterHistograms::FillCluster(AliTPCclusterMI* cluster, Int_t time) {
634   //
635   // Fills the different histograms with the information from the cluster.
636   //
637
638   Int_t padRow =   cluster->GetRow(); 
639   Float_t qMax =   cluster->GetMax();
640   Float_t qTot =   cluster->GetQ();
641   Float_t sigmaY = cluster->GetSigmaY2();
642   Float_t sigmaZ = cluster->GetSigmaZ2();
643   Float_t y      = cluster->GetY();
644   Float_t z      = cluster->GetZ();
645
646   // check if this is ok!!!
647   z = TMath::Abs(z);
648
649   if (qMax<=0) {
650     printf(Form("\n WARNING: Hi Marian! How can we have Qmax = %f ??? \n \n", qMax));
651   }
652   if (qTot<=0) {
653     printf(Form("\n WARNING: Hi Marian! How can we have Qtot = %f ??? \n \n ", qTot));
654   }   
655   
656   // check if the cluster is accepted
657   if (fEdgeSuppression)
658     if (IsClusterOnEdge(cluster))
659       return;
660
661   fNClustersInEvent++;
662   fQtotInEvent = fQtotInEvent + qTot;
663   if (qTot > fMaxQtotInEvent)
664     fMaxQtotInEvent = qTot;
665
666   fhQmaxVsRow           ->Fill(padRow, qMax);
667   fhQtotVsRow           ->Fill(padRow, qTot);
668
669   fhQmaxProfileVsRow    ->Fill(padRow, qMax);
670   fhQtotProfileVsRow    ->Fill(padRow, qTot);
671
672   fhNClustersYVsRow     ->Fill(padRow, y, 1);
673   fhNClustersZVsRow     ->Fill(padRow, z, 1);
674                         
675   fhSigmaYVsRow         ->Fill(padRow, sigmaY);
676   fhSigmaZVsRow         ->Fill(padRow, sigmaZ);
677                         
678   fhQmaxProfileYVsRow   ->Fill(padRow, y, qMax);
679   fhQtotProfileYVsRow   ->Fill(padRow, y, qTot); 
680   fhSigmaYProfileYVsRow ->Fill(padRow, y, sigmaY);
681   fhSigmaZProfileYVsRow ->Fill(padRow, y, sigmaZ);
682
683   fhQmaxProfileZVsRow   ->Fill(padRow, z, qMax);
684   fhQtotProfileZVsRow   ->Fill(padRow, z, qTot); 
685   fhSigmaYProfileZVsRow ->Fill(padRow, z, sigmaY);
686   fhSigmaZProfileZVsRow ->Fill(padRow, z, sigmaZ);
687
688   if (time>0 && fTimeStart>0 && fTimeStop>0 && time>fTimeStart) {
689     //Float_t timeFraction = (time - fTimeStart)/(fTimeStop-fTimeStart); 
690
691     fhMeanQtotVsTime->Fill(time,qTot);
692     fhQtotVsTime->Fill(time,qTot);
693   }
694 }
695
696 //____________________________________________________________________
697 void AliTPCClusterHistograms::FillTrack(const AliTPCseed* seed) {
698   //
699   // fill histograms related to tracks
700   //
701
702   Float_t totalQtot = 0;
703   Int_t   nClusters = 0;
704   for (Int_t clusterID = 0; clusterID < 160; clusterID++) {
705     AliTPCclusterMI* cluster = 0;
706 #ifdef NEWALIROOT
707     cluster = seed->GetClusterPointer(clusterID);
708 #endif
709     if (!cluster)
710       continue;
711     
712     // only use clusters within this detector
713     if (cluster->GetDetector()!=fDetector)
714       continue;
715     
716     // check if the cluster is accepted
717     if (fEdgeSuppression)
718       if (IsClusterOnEdge(cluster))
719         return;
720
721     //    Int_t padRow =   cluster->GetRow(); 
722     //Float_t qMax =   cluster->GetMax();
723     Float_t qTot =   cluster->GetQ();    
724
725     nClusters++;
726     totalQtot += qTot;
727     
728   }
729   if (nClusters==0) 
730     return;
731   
732   Float_t meanQtot = totalQtot/nClusters;
733   
734   // azimuthal angle 
735   Float_t phi    =  TMath::ASin(seed->GetSnp() + seed->GetAlpha());
736   // angle with respect to the central membrane
737   Float_t theta  =  TMath::ATan(seed->GetTgl());
738
739   fhTrackQtotPerCluster->Fill(meanQtot);
740
741   fhTrackMeanQtotPerClusterVsPhi->Fill(phi,   meanQtot);
742   fhTrackMeanQtotPerClusterVsTheta->Fill(theta, meanQtot);
743
744   fhTrackQtotPerClusterVsPhi->Fill(phi, meanQtot);
745   fhTrackQtotPerClusterVsTheta->Fill(theta, meanQtot);
746 }
747
748 //____________________________________________________________________
749 void AliTPCClusterHistograms::FillEvent(Int_t time, Int_t nTracks) {
750   //
751   // fill event
752   //
753
754   fhMeanNTracksVsTime->Fill(time, nTracks);
755
756   //  fhNEventsVsTime->Fill(time);
757 }
758
759
760 //____________________________________________________________________
761 Bool_t AliTPCClusterHistograms::IsClusterOnEdge(AliTPCclusterMI* clusterMI) {
762   //
763   // check if the cluster is on the edge
764   //
765
766   Int_t padRow =   clusterMI->GetRow(); 
767   Float_t y      = clusterMI->GetY();
768   
769   Float_t limit = 0;
770   if (fIsIROC)
771     {
772       limit = 12 + padRow * (20.0 - 12.0) / 63; 
773     }
774   else
775     limit = 16 + padRow * (36.0 - 16.0) / 96;
776   
777   if (TMath::Abs(y) > limit)
778     return kTRUE;
779   
780   return kFALSE;
781 }
782
783 //____________________________________________________________________
784 Float_t AliTPCClusterHistograms::DistanceToEdge(AliTPCclusterMI* clusterMI) {
785   //
786   // get the y-distance to closest edge 
787   //
788   
789   Int_t detector = clusterMI->GetDetector();
790   Int_t padRow   = clusterMI->GetRow(); 
791   Float_t y      = clusterMI->GetY();
792   
793   Float_t yEdge = -9999;
794   Float_t d     = 0;
795   
796   // IROC
797   if (detector < 36) {    
798     yEdge = 14 + padRow * 0.1333;
799     
800   }    
801   else { // OROC
802     if (padRow<64) // small pads
803       yEdge = 22.5 + padRow * 0.1746;
804     else           // large pads
805       yEdge = 34.0 + (padRow-64) * 0.2581;      
806   }
807   if (y<=0) yEdge = -yEdge;
808   
809   d = yEdge - y;
810   
811   return d;
812 }
813
814
815 //____________________________________________________________________
816 Bool_t AliTPCClusterHistograms::KeepThisEvent(TString& why) {
817   //
818   // is this event interesting? 
819   //
820   // the criteria are ...
821   // 
822   
823   if (fKeepEvent) {
824     why = TString(fWhyKeepEvent);
825     return kTRUE;
826   }
827
828   if (fNClustersInEvent>20000) {
829     why.Append("_moreThan20000clusters");
830     fWhyKeepEvent = TString(why);
831     fKeepEvent = kTRUE;
832     return kTRUE;
833   }
834   
835   if (fMaxQtotInEvent>10000) {
836     why.Append("_clusterWithQtot20000plus");
837     fWhyKeepEvent = TString(why);
838     fKeepEvent = kTRUE;
839     return kTRUE;
840   }
841
842   if (gRandom->Uniform()<0.001) {
843     why.Append("_random");
844     fWhyKeepEvent = TString(why);
845     fKeepEvent = kTRUE;
846     return kTRUE;
847   }
848
849   return kFALSE;
850 }
851
852 //____________________________________________________________________
853 void AliTPCClusterHistograms::StartEvent() {
854   //
855   // reset counters 
856   // 
857  
858   fNClustersInEvent  = 0; 
859   fQtotInEvent       = 0; 
860   fMaxQtotInEvent    = 0; 
861   fKeepEvent         = kFALSE; 
862   fWhyKeepEvent      = TString("");
863  
864 }
865
866
867 //____________________________________________________________________
868 void AliTPCClusterHistograms::FinishEvent(Int_t timeStamp) {
869   //
870   // fill histograms related to the event
871   //
872  
873   fhMeanNClustersVsTime->Fill(timeStamp, fNClustersInEvent); 
874   fhNClustersVsTime    ->Fill(timeStamp, fNClustersInEvent); 
875
876   fhNEventsVsTime->Fill(timeStamp);
877  
878 }
879
880
881 //____________________________________________________________________
882 void AliTPCClusterHistograms::SaveHistograms()
883 {
884   //
885   // saves the histograms
886   //
887
888   gDirectory->mkdir(fName.Data());
889   gDirectory->cd(fName.Data());
890
891   fhQmaxVsRow           ->Write();
892   fhQtotVsRow           ->Write();
893
894   fhQmaxProfileVsRow    ->Write();
895   fhQtotProfileVsRow    ->Write();
896
897   fhNClustersYVsRow     ->Write();
898   fhNClustersZVsRow     ->Write();
899                         
900   fhSigmaYVsRow         ->Write();
901   fhSigmaZVsRow         ->Write();
902                         
903   fhQmaxProfileYVsRow   ->Write();
904   fhQtotProfileYVsRow   ->Write();
905   fhSigmaYProfileYVsRow ->Write();
906   fhSigmaZProfileYVsRow ->Write();
907
908   fhQmaxProfileZVsRow   ->Write();
909   fhQtotProfileZVsRow   ->Write();
910   fhSigmaYProfileZVsRow ->Write();
911   fhSigmaZProfileZVsRow ->Write();
912
913   TNamed* comment = new TNamed("comment", fCommentToHistograms.Data());
914   comment->Write();
915
916   if (fhMeanQtotVsTime->GetEntries()>0)
917     fhMeanQtotVsTime->Write();
918
919   if (fhQtotVsTime->GetEntries()>0)
920     fhQtotVsTime->Write();
921
922   if (fhMeanNClustersVsTime->GetEntries()>0)
923     fhMeanNClustersVsTime->Write();
924
925   if (fhNClustersVsTime->GetEntries()>0)
926     fhNClustersVsTime->Write();
927
928   if (fhNEventsVsTime->GetEntries()>0)
929     fhNEventsVsTime->Write();
930
931   gDirectory->mkdir("track_hists");
932   gDirectory->cd("track_hists");
933
934   fhTrackQtotPerCluster->Write();
935
936   fhTrackQtotPerClusterVsPhi->Write();
937   fhTrackQtotPerClusterVsTheta->Write();
938
939   fhTrackMeanQtotPerClusterVsPhi->Write();
940   fhTrackMeanQtotPerClusterVsTheta->Write();
941
942   fhMeanNTracksVsTime->Write();
943
944   gDirectory->cd("../");
945
946   gDirectory->cd("../");
947
948 }
949
950
951 //____________________________________________________________________
952 TCanvas* AliTPCClusterHistograms::DrawHistograms(const Char_t* opt) {
953   //
954   // Draws some histograms and save the canvas as eps and gif file.
955   //  
956
957   TCanvas* c = new TCanvas(fName.Data(), fName.Data(), 1200, 800);
958
959   gStyle->SetOptStat(0);
960   gStyle->SetOptFit(0);
961
962   gStyle->SetPadLeftMargin(0.1);
963
964   c->Divide(3, 2);
965
966   c->Draw();
967
968   c->cd(1);
969
970   // this is not really a nice way to do it...
971   c->GetPad(1)->Delete();
972
973   TLatex* tName = new TLatex(0.05,0.9,fName.Data());
974   tName->SetTextSize(0.02);
975   tName->DrawClone();
976
977   TLatex* tEdge;
978   if (fEdgeSuppression)
979     tEdge = new TLatex(0.05,0.85,"(edges cut)");
980   else
981     tEdge = new TLatex(0.05,0.85,"(no edge cut)");
982
983   tEdge->SetTextSize(0.015);
984   tEdge->DrawClone();
985
986   tName = new TLatex(0.05,0.7, Form("Run: %s", opt));
987   tName->SetTextSize(0.02);
988   tName->DrawClone();
989
990   c->cd(2);
991   fhNClustersYVsRow->Draw("colz");
992
993   c->cd(3);
994   fhQtotVsRow->Draw("colz");
995   fhQtotProfileVsRow->SetMarkerStyle(3);
996   fhQtotProfileVsRow->Draw("same");
997
998   c->cd(4);
999   fhQtotProfileYVsRow   ->Draw("colz");
1000
1001   c->cd(5);
1002   fhQtotProfileZVsRow   ->Draw("colz");
1003
1004   c->cd(6);
1005   fhQtotVsTime->Draw("COLZ");
1006   fhMeanQtotVsTime->SetMarkerStyle(3);
1007   fhMeanQtotVsTime->Draw("SAME");
1008
1009   return c;
1010 }