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