]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/TPC/AliTPCClusterHistograms.cxx
Added track cut method, warning for q<0, profiles for qmax vs pad row and qtot vs...
[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 <TLatex.h>
19
20 #include <../TPC/AliTPCclusterMI.h>
21 #include <AliLog.h>
22
23
24 //____________________________________________________________________
25 ClassImp(AliTPCClusterHistograms)
26
27 //____________________________________________________________________
28 AliTPCClusterHistograms::AliTPCClusterHistograms() 
29   : TNamed(),
30   fhQmaxVsRow(0),          
31   fhQtotVsRow(0),          
32   fhSigmaYVsRow(0),        
33   fhSigmaZVsRow(0),                             
34   fhQmaxProfileYVsRow(0), 
35   fhQtotProfileYVsRow(0),
36   fhSigmaYProfileYVsRow(0),
37   fhSigmaZProfileYVsRow(0),
38   fhQmaxProfileZVsRow(0), 
39   fhQtotProfileZVsRow(0),
40   fhSigmaYProfileZVsRow(0),
41   fhSigmaZProfileZVsRow(0),
42   fhQtotVsTime(0),  
43   fhQmaxVsTime(0),
44   fIsIROC(kFALSE),
45   fEdgeSuppression(kFALSE)
46 {
47   // default constructor
48 }
49
50 //____________________________________________________________________
51 AliTPCClusterHistograms::AliTPCClusterHistograms(Int_t detector, const Char_t* comment, Int_t timeStart, Int_t timeStop, Bool_t edgeSuppression)
52   : TNamed(),
53   fhQmaxVsRow(0),          
54   fhQtotVsRow(0),          
55   fhSigmaYVsRow(0),        
56   fhSigmaZVsRow(0),                             
57   fhQmaxProfileYVsRow(0), 
58   fhQtotProfileYVsRow(0),
59   fhSigmaYProfileYVsRow(0),
60   fhSigmaZProfileYVsRow(0),
61   fhQmaxProfileZVsRow(0), 
62   fhQtotProfileZVsRow(0),
63   fhSigmaYProfileZVsRow(0),
64   fhSigmaZProfileZVsRow(0),
65   fhQtotVsTime(0),  
66   fhQmaxVsTime(0),
67   fIsIROC(kFALSE),
68   fEdgeSuppression(edgeSuppression)
69 {
70   // constructor 
71   
72   // make name and title
73   if (detector < 0 || detector >= 72) {
74     AliDebug(AliLog::kError, Form("Detector %d does not exist", detector));
75     return;
76   }
77       
78   TString name(FormDetectorName(detector, edgeSuppression, comment));
79
80   if (detector < 36)
81     fIsIROC = kTRUE; 
82   
83   SetName(name);
84   SetTitle(Form("%s (detector %d)",name.Data(), detector));
85
86   fTimeStart = timeStart;
87   fTimeStop  = timeStop;
88   
89   #define BINNING_Z 250, 0, 250
90   
91   Float_t yRange   = 45;
92   Int_t nPadRows   = 96;
93   
94   if (fIsIROC)
95   {
96     yRange = 25;
97     nPadRows = 63;
98   }
99   
100   // 1 bin for each 0.5 cm
101   Int_t nBinsY = Int_t(4*yRange);
102
103   // do not add this hists to the directory
104   Bool_t oldStatus = TH1::AddDirectoryStatus();
105   TH1::AddDirectory(kFALSE);
106
107   //defining histograms and profile plots
108   fhQmaxVsRow  = new TH2F("QmaxVsPadRow", "Qmax vs. pad row;Pad row;Qmax", nPadRows+2, -1.5, nPadRows+0.5, 500,  0,  500);
109   fhQtotVsRow  = new TH2F("QtotVsPadRow", "Qtot vs. pad row;Pad row;Qtot", nPadRows+2, -1.5, nPadRows+0.5, 400,  0,  4000);
110   
111   fhSigmaYVsRow = new TH2F("SigmaYVsPadRow", "Sigma Y vs. pad row;Pad row;#sigma_{Y}", nPadRows+2, -1.5, nPadRows+0.5, 100,  0,  0.5);
112   fhSigmaZVsRow = new TH2F("SigmaZVsPadRow", "Sigma Z vs. pad row;Pad row;#sigma_{Z}", nPadRows+2, -1.5, nPadRows+0.5, 100,  0,  0.5);
113   
114   fhQmaxProfileYVsRow = new TProfile2D("MeanQmaxYVsPadRow","Mean Qmax, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
115   fhQtotProfileYVsRow = new TProfile2D("MeanQtotYVsPadRow","Mean Qtot, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
116   fhSigmaYProfileYVsRow = new TProfile2D("MeanSigmaYYVsPadRow","Mean Sigma y, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
117   fhSigmaZProfileYVsRow = new TProfile2D("MeanSigmaZYVsPadRow","Mean Sigma z, y vs pad row;Pad row;y",nPadRows+2, -1.5, nPadRows+0.5, nBinsY, -yRange, yRange);
118
119   fhQmaxProfileZVsRow = new TProfile2D("MeanQmaxZVsPadRow","Mean Qmax, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
120   fhQtotProfileZVsRow = new TProfile2D("MeanQtotZVsPadRow","Mean Qtot, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
121   fhSigmaYProfileZVsRow = new TProfile2D("MeanSigmaYZVsPadRow","Mean Sigma y, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
122   fhSigmaZProfileZVsRow = new TProfile2D("MeanSigmaZZVsPadRow","Mean Sigma z, z vs pad row;Pad row;z",nPadRows+2, -1.5, nPadRows+0.5, BINNING_Z);
123
124   Int_t nTimeBins  = 100;
125   
126   fhQtotVsTime = new TProfile("MeanQtotVsTime", "Mean Qtot vs. event time stamp; time; Qtot",nTimeBins, fTimeStart, fTimeStop);
127   fhQmaxVsTime = new TProfile("MeanQmaxVsTime", "Mean Qmax vs. event time stamp; time; Qmax",nTimeBins, fTimeStart, fTimeStop);
128   
129   TH1::AddDirectory(oldStatus);
130 }
131
132 //____________________________________________________________________
133 AliTPCClusterHistograms::AliTPCClusterHistograms(const AliTPCClusterHistograms& c) : TNamed(c)
134 {
135   // copy constructor
136   ((AliTPCClusterHistograms &)c).Copy(*this);
137 }
138
139 //____________________________________________________________________
140 AliTPCClusterHistograms::~AliTPCClusterHistograms()
141 {
142   //
143   // destructor
144   //
145
146   if (fhQmaxVsRow) {
147     delete fhQmaxVsRow;
148     fhQmaxVsRow = 0;
149   }
150   if (fhQtotVsRow) {
151     delete fhQtotVsRow;
152     fhQtotVsRow = 0; 
153   }
154   if (fhSigmaYVsRow) {
155     delete fhSigmaYVsRow;
156     fhSigmaYVsRow = 0;
157   } 
158   if (fhSigmaZVsRow) {
159     delete fhSigmaZVsRow;
160     fhSigmaZVsRow = 0; 
161   }
162   if (fhQmaxProfileYVsRow) {
163     delete fhQmaxProfileYVsRow;
164     fhQmaxProfileYVsRow = 0;
165   }
166   if (fhQtotProfileYVsRow) {
167     delete fhQtotProfileYVsRow;
168     fhQtotProfileYVsRow = 0;
169   }
170   if (fhSigmaYProfileYVsRow) {
171     delete fhSigmaYProfileYVsRow;
172     fhSigmaYProfileYVsRow = 0;
173   }
174   if (fhSigmaZProfileYVsRow) {
175     delete fhSigmaZProfileYVsRow;
176     fhSigmaZProfileYVsRow = 0;
177   }
178   if (fhQmaxProfileZVsRow) {
179     delete fhQmaxProfileZVsRow;
180     fhQmaxProfileZVsRow = 0;
181   }
182   if (fhQtotProfileZVsRow) {
183     delete fhQtotProfileZVsRow;
184     fhQtotProfileZVsRow = 0;
185   }
186   if (fhSigmaYProfileZVsRow) {
187     delete fhSigmaYProfileZVsRow;
188     fhSigmaYProfileZVsRow = 0;
189   }
190   if (fhSigmaZProfileZVsRow) {
191     delete fhSigmaZProfileZVsRow;
192     fhSigmaZProfileZVsRow = 0;
193   }
194
195   if (fhQtotVsTime) {
196     delete fhQtotVsTime;
197     fhQtotVsTime = 0;
198   }
199   if (fhQmaxVsTime) {
200     delete fhQmaxVsTime;
201     fhQmaxVsTime = 0;
202   }
203 }
204
205 //____________________________________________________________________
206 AliTPCClusterHistograms &AliTPCClusterHistograms::operator=(const AliTPCClusterHistograms &c)
207 {
208   // assigment operator
209
210   if (this != &c)
211     ((AliTPCClusterHistograms &) c).Copy(*this);
212
213   return *this;
214 }
215
216 //____________________________________________________________________
217 const char* AliTPCClusterHistograms::FormDetectorName(Int_t detector, Bool_t edgeSuppression, const char* comment)
218 {
219   //
220   // creates a readable name from the detector number
221   //   
222   
223   Int_t sector = detector%18;
224   TString side;
225   TString inout;
226   
227   if (detector<18 || ( detector>=36 && detector<54))
228     side.Form("A");
229   else 
230     side.Form("B");
231   
232   if (detector<36)
233     inout.Form("IROC");
234   else 
235     inout.Form("OROC");
236
237   TString name;
238   name.Form("sector_%s%d_%s", side.Data(), sector, inout.Data());
239
240   if (edgeSuppression)
241     name += "_noedge";
242   
243   if (comment)
244     name += comment;
245
246   return name; 
247 }
248
249 //____________________________________________________________________
250 Long64_t AliTPCClusterHistograms::Merge(TCollection* list)
251 {
252   // Merge a list of AliTPCClusterHistograms objects with this (needed for
253   // PROOF). 
254   // Returns the number of merged objects (including this).
255
256   if (!list)
257     return 0;
258   
259   if (list->IsEmpty())
260     return 1;
261
262   TIterator* iter = list->MakeIterator();
263   TObject* obj;
264
265   // collections of measured and generated histograms
266   TList* collectionQmaxVsRow     = new TList;
267   TList* collectionQtotVsRow     = new TList;
268   TList* collectionSigmaYVsRow   = new TList;
269   TList* collectionSigmaZVsRow   = new TList;
270                                         
271   TList* collectionQmaxProfileYVsRow    = new TList;
272   TList* collectionQtotProfileYVsRow    = new TList;
273   TList* collectionSigmaYProfileYVsRow  = new TList;
274   TList* collectionSigmaZProfileYVsRow  = new TList;
275
276   TList* collectionQmaxProfileZVsRow    = new TList;
277   TList* collectionQtotProfileZVsRow    = new TList;
278   TList* collectionSigmaYProfileZVsRow  = new TList;
279   TList* collectionSigmaZProfileZVsRow  = new TList;
280
281   TList* collectionQtotVsTime  = new TList;
282   TList* collectionQmaxVsTime  = new TList;
283
284    Int_t count = 0;
285    while ((obj = iter->Next())) {
286     
287      AliTPCClusterHistograms* entry = dynamic_cast<AliTPCClusterHistograms*> (obj);
288      if (entry == 0) 
289        continue;
290
291      collectionQmaxVsRow          ->Add(entry->fhQmaxVsRow         );
292      collectionQtotVsRow          ->Add(entry->fhQtotVsRow         );
293      collectionSigmaYVsRow        ->Add(entry->fhSigmaYVsRow       );
294      collectionSigmaZVsRow        ->Add(entry->fhSigmaZVsRow       );
295                                                                            
296      collectionQmaxProfileYVsRow  ->Add(entry->fhQmaxProfileYVsRow );
297      collectionQtotProfileYVsRow  ->Add(entry->fhQtotProfileYVsRow );
298      collectionSigmaYProfileYVsRow->Add(entry->fhSigmaYProfileYVsRow);
299      collectionSigmaZProfileYVsRow->Add(entry->fhSigmaZProfileYVsRow);
300
301      collectionQmaxProfileZVsRow  ->Add(entry->fhQmaxProfileZVsRow );
302      collectionQtotProfileZVsRow  ->Add(entry->fhQtotProfileZVsRow );
303      collectionSigmaYProfileZVsRow->Add(entry->fhSigmaYProfileZVsRow);
304      collectionSigmaZProfileZVsRow->Add(entry->fhSigmaZProfileZVsRow);
305
306      collectionQtotVsTime->Add(entry->fhQtotVsTime);
307      collectionQmaxVsTime->Add(entry->fhQmaxVsTime);
308
309
310      count++;
311    }
312
313    fhQmaxVsRow          ->Merge(collectionQmaxVsRow       );       
314    fhQtotVsRow          ->Merge(collectionQtotVsRow       );       
315    fhSigmaYVsRow        ->Merge(collectionSigmaYVsRow     );       
316    fhSigmaZVsRow        ->Merge(collectionSigmaZVsRow     );       
317                                                                      
318    fhQmaxProfileYVsRow  ->Merge(collectionQmaxProfileYVsRow  ); 
319    fhQtotProfileYVsRow  ->Merge(collectionQtotProfileYVsRow  );
320    fhSigmaYProfileYVsRow->Merge(collectionSigmaYProfileYVsRow);
321    fhSigmaZProfileYVsRow->Merge(collectionSigmaZProfileYVsRow);
322
323    fhQmaxProfileZVsRow  ->Merge(collectionQmaxProfileZVsRow  ); 
324    fhQtotProfileZVsRow  ->Merge(collectionQtotProfileZVsRow  );
325    fhSigmaYProfileZVsRow->Merge(collectionSigmaYProfileZVsRow);
326    fhSigmaZProfileZVsRow->Merge(collectionSigmaZProfileZVsRow);
327
328    fhQtotVsTime->Merge(collectionQtotVsTime);
329    fhQmaxVsTime->Merge(collectionQmaxVsTime);
330
331    delete collectionQmaxVsRow;          
332    delete collectionQtotVsRow;  
333    delete collectionSigmaYVsRow;          
334    delete collectionSigmaZVsRow;          
335                                   
336    delete collectionQmaxProfileYVsRow;  
337    delete collectionQtotProfileYVsRow;  
338    delete collectionSigmaYProfileYVsRow;
339    delete collectionSigmaZProfileYVsRow;
340
341    delete collectionQmaxProfileZVsRow;  
342    delete collectionQtotProfileZVsRow;  
343    delete collectionSigmaYProfileZVsRow;
344    delete collectionSigmaZProfileZVsRow;
345
346    delete collectionQtotVsTime;
347    delete collectionQmaxVsTime;
348
349   return count+1;
350 }
351
352 //____________________________________________________________________
353 void AliTPCClusterHistograms::FillCluster(AliTPCclusterMI* cluster, Int_t time) {
354   //
355   // Fills the different histograms with the information from the cluster.
356   //
357
358   Int_t padRow =   cluster->GetRow(); 
359   Float_t qMax =   cluster->GetMax();
360   Float_t qTot =   cluster->GetQ();
361   Float_t sigmaY = cluster->GetSigmaY2();
362   Float_t sigmaZ = cluster->GetSigmaZ2();
363   Float_t y      = cluster->GetY();
364   Float_t z      = cluster->GetZ();
365
366   if (qMax<=0) {
367     printf(Form("\n WARNING: Hi Marian! How can we have Qmax = %f ??? \n \n", qMax));
368     return;
369   }
370   if (qTot<=0) {
371     printf(Form("\n WARNING: Hi Marian! How can we have Qtot = %f ??? \n \n ", qTot));
372     return;
373   } 
374   
375   if (fEdgeSuppression)
376   {
377     Float_t limit = 0;
378     if (fIsIROC)
379     {
380       limit = 12 + padRow * (20.0 - 12.0) / 63; 
381     }
382     else
383       limit = 16 + padRow * (36.0 - 16.0) / 96;
384       
385     if (TMath::Abs(y) > limit)
386       return;  
387   }
388
389   fhQmaxVsRow           ->Fill(padRow, qMax);
390   fhQtotVsRow           ->Fill(padRow, qTot);
391                         
392   fhSigmaYVsRow         ->Fill(padRow, sigmaY);
393   fhSigmaZVsRow         ->Fill(padRow, sigmaZ);
394                         
395   fhQmaxProfileYVsRow   ->Fill(padRow, y, qMax);
396   fhQtotProfileYVsRow   ->Fill(padRow, y, qTot); 
397   fhSigmaYProfileYVsRow ->Fill(padRow, y, sigmaY);
398   fhSigmaZProfileYVsRow ->Fill(padRow, y, sigmaZ);
399
400   fhQmaxProfileZVsRow   ->Fill(padRow, z, qMax);
401   fhQtotProfileZVsRow   ->Fill(padRow, z, qTot); 
402   fhSigmaYProfileZVsRow ->Fill(padRow, z, sigmaY);
403   fhSigmaZProfileZVsRow ->Fill(padRow, z, sigmaZ);
404
405   if (time>0 & fTimeStart>0 & fTimeStop>0 & time>fTimeStart) {
406     //Float_t timeFraction = (time - fTimeStart)/(fTimeStop-fTimeStart); 
407
408     fhQtotVsTime->Fill(time,qTot);
409     fhQmaxVsTime->Fill(time,qMax);
410   }
411
412 }
413
414
415 //____________________________________________________________________
416 void AliTPCClusterHistograms::SaveHistograms()
417 {
418   //
419   // saves the histograms
420   //
421
422   gDirectory->mkdir(fName.Data());
423   gDirectory->cd(fName.Data());
424
425   fhQmaxVsRow           ->Write();
426   fhQtotVsRow           ->Write();
427                         
428   fhSigmaYVsRow         ->Write();
429   fhSigmaZVsRow         ->Write();
430                         
431   fhQmaxProfileYVsRow   ->Write();
432   fhQtotProfileYVsRow   ->Write();
433   fhSigmaYProfileYVsRow ->Write();
434   fhSigmaZProfileYVsRow ->Write();
435
436   fhQmaxProfileZVsRow   ->Write();
437   fhQtotProfileZVsRow   ->Write();
438   fhSigmaYProfileZVsRow ->Write();
439   fhSigmaZProfileZVsRow ->Write();
440
441   TProfile* profileQmaxVsRow = fhQmaxVsRow->ProfileX("MeanQmaxVsRow");
442   TProfile* profileQtotVsRow = fhQtotVsRow->ProfileX("MeanQtotVsRow");
443
444   profileQmaxVsRow->SetTitle("Mean Qmax vs. pad row; Pad row; Mean Qmax");
445   profileQtotVsRow->SetTitle("Mean Qtot vs. pad row; Pad row; Mean Qmax");
446
447   profileQmaxVsRow->Write();
448   profileQtotVsRow->Write();
449
450   if (fhQtotVsTime->GetEntries()>0)
451     fhQtotVsTime->Write();
452
453   if (fhQmaxVsTime->GetEntries()>0)
454     fhQmaxVsTime->Write();
455
456   gDirectory->cd("../");
457
458 }
459
460 //____________________________________________________________________
461 TCanvas* AliTPCClusterHistograms::DrawHistograms(const Char_t* /*opt*/) {
462   //
463   // Draws some histograms and save the canvas as eps and gif file.
464   //  
465
466   TCanvas* c = new TCanvas(fName.Data(), fName.Data(), 1200, 1000);
467
468   gStyle->SetOptStat(0);
469   gStyle->SetOptFit(0);
470
471   gStyle->SetPadLeftMargin(0.1);
472
473   c->Divide(3,3);
474
475   c->Draw();  
476
477   c->cd(1);
478   
479   // this is not really a nice way to do it...
480   c->GetPad(1)->Delete();
481   
482   TLatex* name = new TLatex(0.1,0.8,fName.Data());
483   name->SetTextSize(0.02);
484   name->DrawClone();
485
486   c->cd(2);
487   fhQmaxVsRow->Draw("colz");
488
489   c->cd(3);
490   fhQtotVsRow->Draw("colz");         
491                         
492   c->cd(4);
493   fhQmaxProfileYVsRow   ->Draw("colz");
494   c->cd(5);
495   fhQtotProfileYVsRow   ->Draw("colz");
496   c->cd(6);
497   fhQmaxProfileZVsRow   ->Draw("colz");
498   c->cd(7);
499   fhQtotProfileZVsRow   ->Draw("colz");
500   c->cd(8);
501   fhSigmaYVsRow         ->Draw("colz");
502
503   c->cd(9);
504   fhSigmaZVsRow         ->Draw("colz");
505                         
506   //fhSigmaYProfileYVsRow ->Draw("colz");
507   //fhSigmaZProfileYVsRow ->Draw("colz");
508
509   //fhSigmaYProfileZVsRow ->Draw("colz");
510   //fhSigmaZProfileZVsRow ->Draw("colz");
511   return c;
512 }