]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerRec.cxx
Fixes for report #63583: High CPU time spent in TMath::Erf
[u/mrichter/AliRoot.git] / TPC / AliTPCQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 /* $Id: $ */
18
19 /*
20   Based on AliPHOSQADataMaker
21   Produces the data needed to calculate the quality assurance. 
22   All data must be mergeable objects.
23   P. Christiansen, Lund, January 2008
24 */
25
26 /*
27   Implementation:
28
29   We have chosen to have the histograms as non-persistent meber to
30   allow better debugging. In the copy constructor we then have to
31   assign the pointers to the existing histograms in the copied
32   list. This have been implemented but not tested.
33
34   For the QA of the RAW data we use the class, AliTPCdataQA, from the
35   existing TPC Calibration framework (which is more advanced than the
36   standard QA framework) and extract the histograms at the end. This
37   has been tested with zero-suppressed data. The Analyse method of the
38   AliTPCdataQA class is called in the method, EndOfDetectorCycle, and
39   there also: 1d histogram(s) are projected and added to the QA list.
40 */
41
42 /*
43   TODO:
44   Sumw2 for RAW histogram(s)?
45   RecPoints and ESD could have many more histograms
46 */
47
48 #include "AliTPCQADataMakerRec.h"
49
50 // --- ROOT system ---
51 #include <TClonesArray.h>
52 #include <TString.h>
53 #include <TSystem.h>
54 #include <TBox.h>
55 #include <TLine.h>
56
57 // --- Standard library ---
58
59 // --- AliRoot header files ---
60 #include "AliQAChecker.h"
61 #include "AliESDEvent.h"
62 #include "AliESDtrack.h"
63 #include "AliLog.h"
64 #include "AliTPCCalPad.h"
65 #include "AliTPCCalROC.h"
66 #include "AliTPCClustersRow.h"
67 #include "AliTPCclusterMI.h"
68 #include "AliSimDigits.h"
69
70 ClassImp(AliTPCQADataMakerRec)
71
72 //____________________________________________________________________________ 
73 AliTPCQADataMakerRec::AliTPCQADataMakerRec() : 
74 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTPC), 
75                   "TPC Rec Quality Assurance Data Maker"),
76   fTPCdataQA(NULL),
77   fBeautifyOption(1),   // 0:no beautify, !=0:beautify RAW 
78   fOccHighLimit(1e-4),  // high limit for accepting occupancy values
79   fQmaxLowLimit(8),    // low limit for accepting Qmax values
80   fQmaxHighLimit(40)    // high limit for accepting Qmax values
81 {
82   // ctor
83   fTPCdataQA = new AliTPCdataQA*[AliRecoParam::kNSpecies] ;
84   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
85     fTPCdataQA[specie] = NULL ; 
86   
87   for(Int_t i = 0; i < 6; i++)
88     fMapping[i] = 0;
89 }
90
91 //____________________________________________________________________________ 
92 AliTPCQADataMakerRec::AliTPCQADataMakerRec(const AliTPCQADataMakerRec& qadm) :
93   AliQADataMakerRec(),
94   fTPCdataQA(NULL),
95   fBeautifyOption(qadm.GetBeautifyOption()),
96   fOccHighLimit(qadm.GetOccHighLimit()),
97   fQmaxLowLimit(qadm.GetQmaxLowLimit()),
98   fQmaxHighLimit(qadm.GetQmaxHighLimit())
99 {
100   //copy ctor 
101   // Does not copy the calibration object, instead InitRaws have to be
102   // called again
103   SetName((const char*)qadm.GetName()) ; 
104   SetTitle((const char*)qadm.GetTitle()); 
105
106   fTPCdataQA = new AliTPCdataQA*[AliRecoParam::kNSpecies] ;
107   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
108     fTPCdataQA[specie] = NULL ; 
109   
110   for(Int_t i = 0; i < 6; i++)
111     fMapping[i] = 0;
112
113   //
114   // Associate class histogram objects to the copies in the list
115   // Could also be done with the indexes
116   //
117
118 }
119
120 //__________________________________________________________________
121 AliTPCQADataMakerRec& AliTPCQADataMakerRec::operator = (const AliTPCQADataMakerRec& qadm )
122 {
123   // Equal operator.
124   this->~AliTPCQADataMakerRec();
125   new(this) AliTPCQADataMakerRec(qadm);
126   return *this;
127 }
128
129 //__________________________________________________________________
130 AliTPCQADataMakerRec::~AliTPCQADataMakerRec()
131 {
132   // Destructor
133   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
134     if ( fTPCdataQA[specie] != NULL )
135       delete fTPCdataQA[specie] ; 
136   delete[] fTPCdataQA; 
137
138   for(Int_t i = 0; i < 6; i++) 
139     delete fMapping[i];
140 }
141  
142 //____________________________________________________________________________ 
143 void AliTPCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
144 {
145   //Detector specific actions at end of cycle
146   
147   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
148     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
149       continue ; 
150     if(fTPCdataQA[specie] != NULL) { // do the final step of the QA for Raw data
151
152       fTPCdataQA[specie]->Analyse(); // 31/1-08 Analyse is now protected against
153                            //         RAW data files with no TPC data
154
155       SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ; 
156       TH1F * histRawsOccupancy                 = (TH1F*)GetRawsData(kRawsOccupancy) ;
157       TH1F * histRawsOccupancyVsSector         = (TH1F*)GetRawsData(kRawsOccupancyVsSector) ;
158       TH1F * histRawsNClustersPerEventVsSector = (TH1F*)GetRawsData(kRawsNClustersPerEventVsSector) ;
159       TH1F * histRawsQVsSector                 = (TH1F*)GetRawsData(kRawsQVsSector) ;
160       TH1F * histRawsQmaxVsSector              = (TH1F*)GetRawsData(kRawsQmaxVsSector) ;
161       TH1F * histRawsOccupancyVsEvent          = (TH1F*)GetRawsData(kRawsOccupancyVsEvent) ;
162       TH1F * histRawsNclustersVsEvent          = (TH1F*)GetRawsData(kRawsNclustersVsEvent) ;
163       if ( !histRawsOccupancy ||
164            !histRawsOccupancyVsSector ||
165            !histRawsNClustersPerEventVsSector ||
166            !histRawsQVsSector ||
167            !histRawsQmaxVsSector ||
168            !histRawsOccupancyVsEvent ||
169            !histRawsNclustersVsEvent ) {
170         AliError("Something very wrong here, corrupted memory ?????. Please check\n") ; 
171         continue ; 
172       }
173         
174       //Add2RawsList(fTPCdataQA, 0);
175       // get the histograms and add them to the output
176       // 31/8-08 Histogram is only added if the Calibration class 
177       //         receives TPC data 
178       const Int_t eventCounter = fTPCdataQA[specie]->GetEventCounter();
179       if(eventCounter>0) { // some TPC data has been processed
180
181         // Reset histograms and refill them 
182         histRawsOccupancy->Reset();
183         histRawsOccupancyVsSector->Reset();
184         histRawsNClustersPerEventVsSector->Reset();
185         histRawsQVsSector->Reset();
186         histRawsQmaxVsSector->Reset();
187       
188         TH1F* hNormOcc = new TH1F("hNormOcc", 0, 72, 0, 72);
189         hNormOcc->Sumw2();
190         TH1F* hNormNclusters = new TH1F("hNormNclusters", 0, 72, 0, 72);
191         hNormNclusters->Sumw2();
192
193         for (Int_t iSec = 0; iSec < 72; iSec++) {
194         
195           AliTPCCalROC* occupancyROC = 
196           fTPCdataQA[specie]->GetNoThreshold()->GetCalROC(iSec); 
197           AliTPCCalROC* nclusterROC = 
198           fTPCdataQA[specie]->GetNLocalMaxima()->GetCalROC(iSec); 
199           AliTPCCalROC* qROC = 
200           fTPCdataQA[specie]->GetMeanCharge()->GetCalROC(iSec); 
201           AliTPCCalROC* qmaxROC = 
202           fTPCdataQA[specie]->GetMaxCharge()->GetCalROC(iSec); 
203
204           const Int_t nRows = occupancyROC->GetNrows(); 
205           for (Int_t iRow = 0; iRow < nRows; iRow++) {
206
207             const Int_t nPads = occupancyROC->GetNPads(iRow); 
208             for (Int_t iPad = 0; iPad < nPads; iPad++) {
209               
210               histRawsOccupancy->Fill(occupancyROC->GetValue(iRow, iPad));
211               hNormOcc->Fill(iSec);
212               histRawsOccupancyVsSector
213                 ->Fill(iSec, occupancyROC->GetValue(iRow, iPad));
214               
215               const Int_t nClusters = TMath::Nint(nclusterROC->GetValue(iRow, iPad));
216               
217               if(nClusters>0) {
218                 
219                 hNormNclusters->Fill(iSec,nClusters);
220                 histRawsNClustersPerEventVsSector->Fill(iSec, nClusters);
221                 histRawsQVsSector->Fill(iSec, 
222                                         nClusters*qROC->GetValue(iRow, iPad));
223                 histRawsQmaxVsSector->Fill(iSec, 
224                                            nClusters*qmaxROC->GetValue(iRow, iPad));
225               }
226             }
227           }
228         } // end loop over sectors
229       
230         // update event histograms - copy info from TPDdataQA histos
231         TH1F* hQAOccVsEvent = fTPCdataQA[specie]->GetHistOccupancyVsEvent();
232         TH1F* hQANclVsEvent = fTPCdataQA[specie]->GetHistNclustersVsEvent();
233         
234         // the two event histograms should have the same number of bins
235         const Int_t nBins = hQAOccVsEvent->GetXaxis()->GetNbins();
236         for(Int_t bin = 1; bin <= nBins; bin++) {
237           
238           histRawsOccupancyVsEvent->SetBinContent(bin, hQAOccVsEvent->GetBinContent(bin));
239           histRawsNclustersVsEvent->SetBinContent(bin, hQANclVsEvent->GetBinContent(bin));
240         }
241         
242         histRawsOccupancyVsEvent->GetXaxis()->SetRange(hQAOccVsEvent->GetXaxis()->GetFirst(), hQAOccVsEvent->GetXaxis()->GetLast());
243         histRawsNclustersVsEvent->GetXaxis()->SetRange(hQANclVsEvent->GetXaxis()->GetFirst(), hQANclVsEvent->GetXaxis()->GetLast());
244
245         // Normalize histograms
246         histRawsOccupancyVsSector->Divide(hNormOcc);
247         histRawsNClustersPerEventVsSector->Scale(1.0/Float_t(eventCounter));
248         histRawsQVsSector->Divide(hNormNclusters);
249         histRawsQmaxVsSector->Divide(hNormNclusters);
250         delete hNormOcc;
251         delete hNormNclusters;
252
253         if(fBeautifyOption!=0) {
254           // Help make the histogram easier to interpret for the DQM shifter
255           
256           histRawsOccupancyVsSector->ResetBit(AliQAv1::GetQABit());
257           histRawsQmaxVsSector->ResetBit(AliQAv1::GetQABit());
258
259           histRawsOccupancyVsSector->SetMinimum(0.0);
260           if(histRawsOccupancyVsSector->GetMaximum()<1.5*fOccHighLimit)
261             histRawsOccupancyVsSector->SetMaximum(1.5*fOccHighLimit);
262           
263           histRawsQmaxVsSector->SetMinimum(0.0);
264           if(histRawsQmaxVsSector->GetMaximum()<1.5*fQmaxHighLimit)
265             histRawsQmaxVsSector->SetMaximum(1.5*fQmaxHighLimit);
266           
267           Double_t xminOcc = histRawsOccupancyVsSector->GetXaxis()->GetXmin();
268           Double_t xmaxOcc = histRawsOccupancyVsSector->GetXaxis()->GetXmax();
269           Double_t yminOcc = histRawsOccupancyVsSector->GetMinimum();
270           Double_t ymaxOcc = histRawsOccupancyVsSector->GetMaximum();
271           
272           Double_t xminQmax = histRawsQmaxVsSector->GetXaxis()->GetXmin();
273           Double_t xmaxQmax = histRawsQmaxVsSector->GetXaxis()->GetXmax();
274 //        Double_t yminQmax = histRawsQmaxVsSector->GetMinimum();
275 //        Double_t ymaxQmax = histRawsQmaxVsSector->GetMaximum();
276           
277           // For reasons not understood the following stopped working
278           // in the DQM and instead lines were adopted:
279 //        TBox* boxOccOk = new TBox(xminOcc,0,xmaxOcc,fOccHighLimit);
280 //        boxOccOk->SetFillColor(kGreen);
281 //        histRawsOccupancyVsSector->GetListOfFunctions()->Add(boxOccOk);
282
283           TLine* lineOccMin = new TLine(xminOcc,0,xmaxOcc,0);
284           lineOccMin->SetLineColor(kGreen);
285           lineOccMin->SetLineWidth(2);
286           histRawsOccupancyVsSector->GetListOfFunctions()->Add(lineOccMin);
287
288           TLine* lineOccMax = new TLine(xminOcc,fOccHighLimit, xmaxOcc,fOccHighLimit);
289           lineOccMax->SetLineColor(kGreen);
290           lineOccMax->SetLineWidth(2);
291           histRawsOccupancyVsSector->GetListOfFunctions()->Add(lineOccMax);
292
293  
294           // For some reason this dtopped working
295 //        TBox* boxQmaxOk = new TBox(xminQmax,fQmaxLowLimit,xmaxQmax,fQmaxHighLimit);
296 //        boxQmaxOk->SetFillColor(kGreen);
297 //        histRawsQmaxVsSector->GetListOfFunctions()->Add(boxQmaxOk);
298
299           TLine* lineQmaxMin = new TLine(xminQmax,fQmaxLowLimit, xmaxQmax,fQmaxLowLimit);
300           lineQmaxMin->SetLineColor(kGreen);
301           lineQmaxMin->SetLineWidth(2);
302           histRawsQmaxVsSector->GetListOfFunctions()->Add(lineQmaxMin);
303
304           TLine* lineQmaxMax = new TLine(xminQmax,fQmaxHighLimit, xmaxQmax,fQmaxHighLimit);
305           lineQmaxMax->SetLineColor(kGreen);
306           lineQmaxMax->SetLineWidth(2);
307           histRawsQmaxVsSector->GetListOfFunctions()->Add(lineQmaxMax);
308           
309           
310           for(Int_t bin = 1; bin <= 72; bin++) {
311             
312             if(histRawsOccupancyVsSector->GetBinContent(bin)<=0 ||
313                histRawsOccupancyVsSector->GetBinContent(bin)>fOccHighLimit) {
314               
315               histRawsOccupancyVsSector->SetBit(AliQAv1::GetQABit());
316
317               TBox* boxErr = 
318                 new TBox(histRawsOccupancyVsSector->GetXaxis()->GetBinLowEdge(bin), yminOcc,
319                          histRawsOccupancyVsSector->GetXaxis()->GetBinUpEdge(bin), ymaxOcc);
320               boxErr->SetFillColor(kRed);
321               //              histRawsOccupancyVsSector->GetListOfFunctions()->Add(boxErr);
322             }
323             
324             if(histRawsQmaxVsSector->GetBinContent(bin)<fQmaxLowLimit||
325                histRawsQmaxVsSector->GetBinContent(bin)>fQmaxHighLimit) {
326               
327               // Mark that histogram has error
328               histRawsQmaxVsSector->SetBit(AliQAv1::GetQABit());
329
330           // For reasons not understood the following stopped working
331           // in the DQM and instead lines were adopted:
332 //            TBox* boxErr = 
333 //              new TBox(histRawsQmaxVsSector->GetXaxis()->GetBinLowEdge(bin), yminQmax,
334 //                       histRawsQmaxVsSector->GetXaxis()->GetBinUpEdge(bin), ymaxQmax);
335 //            boxErr->SetFillColor(kRed);
336               //              histRawsQmaxVsSector->GetListOfFunctions()->Add(boxErr);
337             }
338           }
339
340           // For reasons not understood the following stopped working
341           // in the DQM and instead lines were adopted:
342           // Now we have to add a copy of the histograms to draw
343           // because the boxes covers the data points
344 //        TH1F* hOccCopy = new TH1F(*histRawsOccupancyVsSector);
345 //        hOccCopy->SetOption("SAME P");
346 //        histRawsOccupancyVsSector->GetListOfFunctions()->Add(hOccCopy);
347
348 //        TH1F* hQmaxCopy = new TH1F(*histRawsQmaxVsSector);
349 //        hQmaxCopy->SetOption("SAME P");
350 //        histRawsQmaxVsSector->GetListOfFunctions()->Add(hQmaxCopy);
351
352         } // end beautify
353       }
354     }
355   }
356   AliQAChecker::Instance()->Run(AliQAv1::kTPC, task, list) ;  
357 }
358
359
360 //____________________________________________________________________________ 
361 void AliTPCQADataMakerRec::InitESDs()
362 {
363   //create ESDs histograms in ESDs subdir  
364   const Bool_t expert   = kTRUE ; 
365   const Bool_t image    = kTRUE ; 
366   
367   TH1F * histESDclusters = 
368     new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
369              160, 0, 160);
370   histESDclusters->Sumw2();
371   Add2ESDsList(histESDclusters, KClusters, !expert, image);
372
373   TH1F * histESDratio = 
374     new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
375              100, 0, 1);
376   histESDratio->Sumw2();
377   Add2ESDsList(histESDratio, kRatio, !expert, image);
378   
379   TH1F * histESDpt = 
380     new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
381              50, 0, 5);
382   histESDpt->Sumw2();
383   Add2ESDsList(histESDpt, kPt, !expert, image);
384
385   // This means we are not running DQM so do not beautify
386   SetBeautifyOption(0);
387 }
388
389 //____________________________________________________________________________ 
390 void AliTPCQADataMakerRec::InitRaws()
391 {
392   //
393   // Adding the raw 
394   //  
395
396   // Modified: 7/7 - 2008
397   // Laurent Aphecetche pointed out that the mapping was read from file
398   // for each event, so now we read in the map here and set if for 
399   // the raw data qa
400   const Bool_t expert   = kTRUE ; 
401   const Bool_t saveCorr = kTRUE ; 
402   const Bool_t image    = kTRUE ; 
403   
404   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
405     
406     // It might happen that we will be in this method a few times because
407     // we create all dataQAs at the first call to this method
408     if(fTPCdataQA[specie]!=0) // data QA already created
409       continue;
410     fTPCdataQA[specie] = new AliTPCdataQA(AliRecoParam::ConvertIndex(specie));
411     LoadMaps(); // Load Altro maps
412     fTPCdataQA[specie]->SetAltroMapping(fMapping); // set Altro mapping
413     fTPCdataQA[specie]->SetRangeTime(100, 920); // set time bin interval 
414 //    Add2RawsList(fTPCdataQA, kTPCdataQ, !expert, image, !saveCorrA); // This is used by the AMORE monitoring <------- THIS WILL FAIL (YS)
415   }
416
417   TH1F * histRawsOccupancy = 
418     new TH1F("hRawsOccupancy", "Occupancy (all pads); Occupancy; Counts",
419              100, 0, 1);
420   histRawsOccupancy->Sumw2();
421   Add2RawsList(histRawsOccupancy, kRawsOccupancy, expert, !image, !saveCorr);
422   
423   TH1F * histRawsOccupancyVsSector = 
424     new TH1F("hRawsOccupancyVsSector", "Occupancy vs sector; Sector; Occupancy",
425              72, 0, 72);
426   histRawsOccupancyVsSector->Sumw2();
427   histRawsOccupancyVsSector->SetMarkerStyle(20);
428   histRawsOccupancyVsSector->SetOption("P");
429   histRawsOccupancyVsSector->SetStats(kFALSE);
430   Add2RawsList(histRawsOccupancyVsSector, kRawsOccupancyVsSector, !expert, image, !saveCorr);
431
432   TH1F * histRawsNClustersPerEventVsSector = 
433     new TH1F("hRawsNClustersPerEventVsSector", "Nclusters per event vs sector; Sector; Nclusters per event",
434              72, 0, 72);
435   histRawsNClustersPerEventVsSector->Sumw2();
436   Add2RawsList(histRawsNClustersPerEventVsSector, kRawsNClustersPerEventVsSector, expert, !image, !saveCorr);
437   
438   TH1F * histRawsQVsSector = 
439     new TH1F("hRawsQVsSector", "<Q> vs sector; Sector; <Q>",
440              72, 0, 72);
441   histRawsQVsSector->Sumw2();
442   Add2RawsList(histRawsQVsSector, kRawsQVsSector, expert, !image, !saveCorr);
443
444   TH1F * histRawsQmaxVsSector = 
445     new TH1F("hRawsQmaxVsSector", "<Qmax> vs sector; Sector; <Qmax>",
446              72, 0, 72);
447   histRawsQmaxVsSector->Sumw2();
448   histRawsQmaxVsSector->SetMarkerStyle(20);
449   histRawsQmaxVsSector->SetOption("P");
450   histRawsQmaxVsSector->SetStats(kFALSE);
451   Add2RawsList(histRawsQmaxVsSector, kRawsQmaxVsSector, !expert, image, !saveCorr);
452
453   // Get histogram information from data QA to build copy
454   TH1F* hOccHelp = fTPCdataQA[0]->GetHistOccupancyVsEvent();
455   TH1F * histRawsOccupancyVsEvent = 
456     new TH1F("hRawsOccupancyVsEvent", hOccHelp->GetTitle(),
457              hOccHelp->GetXaxis()->GetNbins(),
458              hOccHelp->GetXaxis()->GetXmin(), hOccHelp->GetXaxis()->GetXmax());
459   histRawsOccupancyVsEvent->GetXaxis()->SetTitle(hOccHelp->GetXaxis()->GetTitle());
460   histRawsOccupancyVsEvent->GetYaxis()->SetTitle(hOccHelp->GetYaxis()->GetTitle());
461   histRawsOccupancyVsEvent->SetMarkerStyle(20);
462   histRawsOccupancyVsEvent->SetOption("P");
463   histRawsOccupancyVsEvent->SetStats(kFALSE);
464   Add2RawsList(histRawsOccupancyVsEvent, kRawsOccupancyVsEvent, !expert, image, !saveCorr);
465
466   // Get histogram information from data QA to build copy
467   TH1F* hNclHelp = fTPCdataQA[0]->GetHistNclustersVsEvent();
468   TH1F * histRawsNclustersVsEvent = 
469     new TH1F("hRawsNclustersVsEvent", hNclHelp->GetTitle(),
470              hNclHelp->GetXaxis()->GetNbins(),
471              hNclHelp->GetXaxis()->GetXmin(), hNclHelp->GetXaxis()->GetXmax());
472   histRawsNclustersVsEvent->GetXaxis()->SetTitle(hNclHelp->GetXaxis()->GetTitle());
473   histRawsNclustersVsEvent->GetYaxis()->SetTitle(hNclHelp->GetYaxis()->GetTitle());
474   histRawsNclustersVsEvent->SetMarkerStyle(20);
475   histRawsNclustersVsEvent->SetOption("P");
476   histRawsNclustersVsEvent->SetStats(kFALSE);
477   Add2RawsList(histRawsNclustersVsEvent, kRawsNclustersVsEvent, !expert, image, !saveCorr);
478 }
479
480 //____________________________________________________________________________ 
481 void AliTPCQADataMakerRec::InitDigits()
482 {
483   const Bool_t expert   = kTRUE ; 
484   const Bool_t image    = kTRUE ; 
485   TH1F * histDigitsADC = 
486     new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
487              1000, 0, 1000);
488   histDigitsADC->Sumw2();
489   Add2DigitsList(histDigitsADC, kDigitsADC, !expert, image);
490
491   // This means we are not running DQM so do not beautify
492   SetBeautifyOption(0);
493 }
494
495 //____________________________________________________________________________ 
496 void AliTPCQADataMakerRec::InitRecPoints()
497 {
498   const Bool_t expert   = kTRUE ; 
499   const Bool_t image    = kTRUE ; 
500   
501   TH1F * histRecPointsQmaxShort = 
502     new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
503              100, 0, 300);
504   histRecPointsQmaxShort->Sumw2();
505   Add2RecPointsList(histRecPointsQmaxShort, kQmaxShort, !expert, image);
506
507   TH1F * histRecPointsQmaxMedium = 
508     new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
509              100, 0, 300);
510   histRecPointsQmaxMedium->Sumw2();
511   Add2RecPointsList(histRecPointsQmaxMedium, kQmaxMedium, !expert, image);
512
513   TH1F * histRecPointsQmaxLong = 
514     new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
515              100, 0, 300);
516   histRecPointsQmaxLong->Sumw2();
517   Add2RecPointsList(histRecPointsQmaxLong, kQmaxLong, !expert, image);
518
519   TH1F * histRecPointsQShort = 
520     new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
521              100, 0, 2000);
522   histRecPointsQShort->Sumw2();
523   Add2RecPointsList(histRecPointsQShort, kQShort, !expert, image);
524
525   TH1F * histRecPointsQMedium = 
526     new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
527              100, 0, 2000);
528   histRecPointsQMedium->Sumw2();
529   Add2RecPointsList(histRecPointsQMedium, kQMedium, !expert, image);
530
531   TH1F * histRecPointsQLong = 
532     new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
533              100, 0, 2000);
534   histRecPointsQLong->Sumw2();
535   Add2RecPointsList(histRecPointsQLong, kQLong, !expert, image);
536
537   TH1F * histRecPointsRow = 
538     new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
539              159, 0, 159);
540   histRecPointsRow->Sumw2();
541   Add2RecPointsList(histRecPointsRow, kRow, !expert, image);
542
543   // This means we are not running DQM so do not beautify
544   SetBeautifyOption(0);
545 }
546
547 //____________________________________________________________________________
548 void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
549 {
550   // make QA data from ESDs
551  
552   const Int_t nESDTracks = esd->GetNumberOfTracks();
553   Int_t nTPCtracks = 0; 
554   for(Int_t i = 0; i < nESDTracks; i++) {
555     
556     AliESDtrack * track = esd->GetTrack(i);
557     
558     if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
559       continue;
560     
561     nTPCtracks++;
562     
563     Int_t nTPCclusters         = track->GetTPCNcls();
564     Int_t nTPCclustersFindable = track->GetTPCNclsF();
565     if ( nTPCclustersFindable<=0) continue;
566     GetESDsData(KClusters)->Fill(nTPCclusters);
567     GetESDsData(kRatio)->Fill(Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
568     GetESDsData(kPt)->Fill(track->Pt()); 
569   }
570 }
571
572 //____________________________________________________________________________
573 void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
574 {
575   //
576   // To make QA for the RAW data we use the TPC Calibration framework 
577   // to handle the data and then in the end extract the data
578   //
579   
580   GetRawsData(0); // dummy call to init raw data
581   rawReader->Reset() ; 
582   if (! fTPCdataQA[AliRecoParam::AConvert(fEventSpecie)] ) {
583     AliError("Something unexpected here!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!") ; 
584   } else {  
585     fTPCdataQA[AliRecoParam::AConvert(fEventSpecie)]->ProcessEvent(rawReader);  
586   }
587 }
588
589 //____________________________________________________________________________
590 void AliTPCQADataMakerRec::MakeDigits(TTree* digitTree)
591 {
592  
593   TBranch* branch = digitTree->GetBranch("Segment");
594   AliSimDigits* digArray = 0;
595   branch->SetAddress(&digArray);
596   
597   Int_t nEntries = Int_t(digitTree->GetEntries());
598   
599   for (Int_t n = 0; n < nEntries; n++) {
600     
601     digitTree->GetEvent(n);
602     
603     if (digArray->First())
604       do {
605         Float_t dig = digArray->CurrentDigit();
606         
607         GetDigitsData(kDigitsADC)->Fill(dig);
608       } while (digArray->Next());    
609   }
610 }
611
612 //____________________________________________________________________________
613 void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
614 {
615
616   AliTPCClustersRow clrow;
617   clrow.SetClass("AliTPCclusterMI");
618   clrow.SetArray(0);
619   clrow.GetArray()->ExpandCreateFast(10000);
620   AliTPCClustersRow * pclrow = &clrow;
621   TBranch* branch = recTree->GetBranch("Segment");
622   
623   branch->SetAddress(&pclrow);
624
625   const Int_t nEntries = Int_t(recTree->GetEntries());
626   for (Int_t i = 0; i < nEntries; i++) {
627     
628     branch->GetEntry(i);
629     
630     const Int_t nClusters = clrow.GetArray()->GetEntriesFast();
631     for (Int_t icl=0; icl < nClusters; icl++){
632       
633       AliTPCclusterMI* cluster = 
634         (AliTPCclusterMI*)clrow.GetArray()->At(icl);
635       
636       Float_t Qmax = cluster->GetMax();
637       Float_t Q    = cluster->GetQ();
638       Int_t   row  = cluster->GetRow();
639
640       if(cluster->GetDetector()<36) { // IROC (short pads)
641
642         GetRecPointsData(kQmaxShort)->Fill(Qmax);
643         GetRecPointsData(kQShort)->Fill(Q);
644       } else { // OROC (medium and long pads)
645         row += 63;
646         if(cluster->GetRow()<64) { // medium pads
647
648           GetRecPointsData(kQmaxMedium)->Fill(Qmax);
649           GetRecPointsData(kQMedium)->Fill(Q);
650         } else { // long pads
651
652           GetRecPointsData(kQmaxLong)->Fill(Qmax);
653           GetRecPointsData(kQLong)->Fill(Q);
654         }
655       }
656       
657       GetRecPointsData(kRow)->Fill(row);
658     } // end loop over clusters
659   } // end loop over tree
660
661 }
662
663 //____________________________________________________________________________
664 void AliTPCQADataMakerRec::LoadMaps()
665 {
666   TString path = gSystem->Getenv("ALICE_ROOT");
667   path += "/TPC/mapping/Patch";
668
669   for(Int_t i = 0; i < 6; i++) {
670
671     if(fMapping[i]!=0) // mapping already loaded
672       continue;
673     TString path2 = path;
674     path2 += i;
675     path2 += ".data";
676     fMapping[i] = new AliTPCAltroMapping(path2.Data());
677   }
678 }
679
680 //____________________________________________________________________________
681 void AliTPCQADataMakerRec::ResetDetector(AliQAv1::TASKINDEX_t task)
682 {
683   // Overwrites general method for RAW data.
684   // The AliTPCdataQA elements that does the internal processing are
685   // in the case they have processed data deleted and new are created
686
687   if ( task != AliQAv1::kRAWS )
688     AliQADataMakerRec::ResetDetector(task);
689   
690   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
691     
692     if ( fTPCdataQA[specie] != NULL) { // exist
693       
694       if(fTPCdataQA[specie]->GetEventCounter()>0) { // has processed data
695         
696         // old configuration
697         Int_t  firstTime    = fTPCdataQA[specie]->GetFirstTimeBin();
698         Int_t  lastTime     = fTPCdataQA[specie]->GetLastTimeBin();
699         Int_t  minADC       = fTPCdataQA[specie]->GetAdcMin();
700         Int_t  maxADC       = fTPCdataQA[specie]->GetAdcMax();
701         Int_t  maxEvents    = fTPCdataQA[specie]->GetMaxEvents();
702         Int_t  eventsPerBin = fTPCdataQA[specie]->GetEventsPerBin();
703
704         //delete old
705         delete fTPCdataQA[specie]; 
706
707         // create new
708         fTPCdataQA[specie] = new AliTPCdataQA(AliRecoParam::ConvertIndex(specie));
709         // configure new
710         LoadMaps(); // Load Altro maps
711         fTPCdataQA[specie]->SetAltroMapping(fMapping);
712         fTPCdataQA[specie]->SetRangeTime(firstTime, lastTime);
713         fTPCdataQA[specie]->SetRangeAdc(minADC, maxADC);
714         fTPCdataQA[specie]->SetMaxEvents(maxEvents);
715         fTPCdataQA[specie]->SetEventsPerBin(eventsPerBin);
716       }
717     }
718   }
719 }