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