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