]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCdataQA.cxx
Changing the class version after modifications of the class layout in previous revisions.
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCdataQA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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   July 2011:
21
22   Changes to accomodate updates of general DQM/QA changes to have per trigger
23   histograms (for a given event specie).
24
25   AliTPCdataQA has a new flag for only keeping DQM info event by
26   event!
27   The expert/DA functionality has been kept exactly the same!
28   
29
30   June 2010
31
32   This update should solve two problems mainly:
33   * The vs event histograms have been limited to a fixed size for the
34   DQM. The 500k seemed to be a big size but is no longer so, so we
35   need to dynamically expand the range. The non-trivial point is that
36   we also have to do it for the copy owned by AliTPCQADataMakerRec.
37   * The amoreGui now remembers the x-range of the first visualization so
38   the trick of setting the relevant event range as the histogram is
39   filled no longer works.
40
41   The fix is a bit crude but avoids creating a new histogram. Instead
42   the range is expanded (max events and events per bin is doubled) but
43   the number of bins is kept constant! In this way we can change just
44   double the max of the X-axis of the hist and rebin the data. The
45   same can easily be done for the copy owned by AliTPCQADataMakerRec.
46
47   CAUTION:
48   If we change the number of bins we could crash the whole system
49   because ROOT does not create space for extra bins! (but we do not do
50   this). In that way it is a crude solution.
51   The rebinning in the code only works for an even number of bins.
52
53   In addition to the above a bug in the reading of the config file was
54   also found and corrected. fAdcMax was set instead of fEventsPerBin.
55
56   Finally cout was changes to AliInfo.
57
58   February 2008
59
60   The code has been heavily modified so that now the RAW data is
61   "expanded" for each sector and stored in a big signal array. Then a
62   simple version of the code in AliTPCclusterer is used to identify
63   the local maxima and these are then used for QA. This gives a better
64   estimate of the charge (both max and total) and also limits the
65   effect of noise.
66
67   Implementation:
68
69   In Update the RAW signals >= 3 ADC channels are stored in the arrays.
70   
71   There are 3 arrays:
72   Float_t** fAllBins       2d array [row][bin(pad, time)] ADC signal
73   Int_t**   fAllSigBins    2d array [row][signal#] bin(with signal)
74   Int_t*    fAllNSigBins;  1d array [row] Nsignals
75
76   This is done sector by sector.
77
78   When data from a new sector is encountered, the method
79   FindLocalMaxima is called on the data from the previous sector, and
80   the calibration/data objects are updated with the "cluster"
81   info. Finally the arrays are cleared.
82
83   The requirements for a local maxima is:
84   Charge in bin is >= 5 ADC channels.
85   Charge in bin is larger than all the 8 neighboring bins.
86   At least one of the two pad neighbors has a signal.
87   At least one of the two time neighbors has a signal.
88
89   Before accessing the data it is expected that the Analyse method is
90   called. This normalizes some of the data objects to per event or per
91   cluster. 
92   If more data is passed to the class after Analyse has been called
93   the normalization is reversed and Analyse has to be called again.
94 */
95
96
97 //Root includes
98 #include <TH1F.h>
99 #include <TString.h>
100 #include <TMath.h>
101 #include <TDirectory.h>
102 #include <TFile.h>
103 #include <TError.h>
104 #include <TMap.h>
105 #include <TProfile.h>
106 //AliRoot includes
107 #include "AliRawReader.h"
108 #include "AliRawReaderRoot.h"
109 #include "AliRawReaderDate.h"
110 #include "AliTPCRawStreamV3.h"
111 #include "AliTPCCalROC.h"
112 #include "AliTPCROC.h"
113 #include "AliMathBase.h"
114 #include "TTreeStream.h"
115
116 //date
117 #include "event.h"
118 #include "AliTPCCalPad.h"
119 #include "AliTPCPreprocessorOnline.h"
120
121 //header file
122 #include "AliTPCdataQA.h"
123 #include "AliLog.h"
124
125
126 ClassImp(AliTPCdataQA)
127
128 AliTPCdataQA::AliTPCdataQA() : /*FOLD00*/  
129   fFirstTimeBin(60),
130   fLastTimeBin(1000),
131   fAdcMin(1),
132   fAdcMax(100),
133   fMapping(NULL),
134   fPedestal(0),
135   fNoise(0),
136   fNLocalMaxima(0),
137   fMaxCharge(0),
138   fMeanCharge(0),
139   fNoThreshold(0),
140   fNTimeBins(0),
141   fNPads(0),
142   fTimePosition(0),
143   fOverThreshold10(0),
144   fOverThreshold20(0),
145   fOverThreshold30(0),
146   fHistQVsTimeSideA(0),
147   fHistQVsTimeSideC(0),
148   fHistQMaxVsTimeSideA(0),
149   fHistQMaxVsTimeSideC(0),
150   fHistOccupancyVsEvent(0),
151   fHistNclustersVsEvent(0),
152   fEventCounter(0),
153   fIsAnalysed(kFALSE),
154   fMaxEvents(500000),           // Max events for event histograms
155   fEventsPerBin(1000),          // Events per bin for event histograms
156   fSignalCounter(0),            // Signal counter
157   fClusterCounter(0),           // Cluster counter
158   fAllBins(0),
159   fAllSigBins(0),
160   fAllNSigBins(0),
161   fRowsMax(0),
162   fPadsMax(0),
163   fTimeBinsMax(0),
164   fIsDQM(kFALSE),
165   fHistOccVsSector(0x0),
166   fHistOcc2dVsSector(0x0),
167   fHistQVsSector(0x0),
168   fHistQmaxVsSector(0x0),
169   fOccVec(0x0),
170   fOccMaxVec(0x0),
171   fOccVecFine(0x0),
172   fOccMaxVecFine(0x0)
173 {
174   //
175   // default constructor
176   //
177 }
178
179 //_____________________________________________________________________
180 AliTPCdataQA::AliTPCdataQA(const AliTPCdataQA &ped) : /*FOLD00*/
181   TH1F(ped),
182   fFirstTimeBin(ped.GetFirstTimeBin()),
183   fLastTimeBin(ped.GetLastTimeBin()),
184   fAdcMin(ped.GetAdcMin()),
185   fAdcMax(ped.GetAdcMax()),
186   fMapping(NULL),
187   fPedestal(0),
188   fNoise(0),
189   fNLocalMaxima(0),
190   fMaxCharge(0),
191   fMeanCharge(0),
192   fNoThreshold(0),
193   fNTimeBins(0),
194   fNPads(0),
195   fTimePosition(0),
196   fOverThreshold10(0),
197   fOverThreshold20(0),
198   fOverThreshold30(0),
199   fHistQVsTimeSideA(0),
200   fHistQVsTimeSideC(0),
201   fHistQMaxVsTimeSideA(0),
202   fHistQMaxVsTimeSideC(0),
203   fHistOccupancyVsEvent(0),
204   fHistNclustersVsEvent(0),
205   fEventCounter(ped.GetEventCounter()),
206   fIsAnalysed(ped.GetIsAnalysed()),
207   fMaxEvents(ped.GetMaxEvents()),
208   fEventsPerBin(ped.GetEventsPerBin()),
209   fSignalCounter(ped.GetSignalCounter()),
210   fClusterCounter(ped.GetClusterCounter()),
211   fAllBins(0),
212   fAllSigBins(0),
213   fAllNSigBins(0),
214   fRowsMax(0),
215   fPadsMax(0),
216   fTimeBinsMax(0),
217   fIsDQM(ped.GetIsDQM()),
218   fHistOccVsSector(0x0),
219   fHistOcc2dVsSector(0x0),
220   fHistQVsSector(0x0),
221   fHistQmaxVsSector(0x0),
222   fOccVec(0x0),
223   fOccMaxVec(0x0),
224   fOccVecFine(0x0),
225   fOccMaxVecFine(0x0)
226 {
227   //
228   // copy constructor
229   //
230   if(ped.GetNLocalMaxima())
231     fNLocalMaxima  = new AliTPCCalPad(*ped.GetNLocalMaxima());
232   if(ped.GetMaxCharge())
233     fMaxCharge      = new AliTPCCalPad(*ped.GetMaxCharge());
234   if(ped.GetMeanCharge())
235     fMeanCharge     = new AliTPCCalPad(*ped.GetMeanCharge());
236   if(ped.GetNoThreshold())
237     fNoThreshold  = new AliTPCCalPad(*ped.GetNoThreshold());
238   if(ped.GetNTimeBins())
239     fNTimeBins  = new AliTPCCalPad(*ped.GetNTimeBins());
240   if(ped.GetNPads())
241     fNPads  = new AliTPCCalPad(*ped.GetNPads());
242   if(ped.GetTimePosition())
243     fTimePosition  = new AliTPCCalPad(*ped.GetTimePosition());
244   if(ped.GetOverThreshold10())
245     fOverThreshold10  = new AliTPCCalPad(*ped.GetOverThreshold10());
246   if(ped.GetOverThreshold20())
247     fOverThreshold20  = new AliTPCCalPad(*ped.GetOverThreshold20());
248   if(ped.GetOverThreshold30())
249     fOverThreshold30  = new AliTPCCalPad(*ped.GetOverThreshold30());
250   if(ped.GetHistQVsTimeSideA()) {
251     fHistQVsTimeSideA = new TProfile(*ped.GetHistQVsTimeSideA());
252     fHistQVsTimeSideA->SetDirectory(0);
253   }
254   if(ped.GetHistQVsTimeSideC()) {
255     fHistQVsTimeSideC = new TProfile(*ped.GetHistQVsTimeSideC());
256     fHistQVsTimeSideC->SetDirectory(0);
257   }
258   if(ped.GetHistQMaxVsTimeSideA()) {
259     fHistQMaxVsTimeSideA = new TProfile(*ped.GetHistQMaxVsTimeSideA());
260     fHistQMaxVsTimeSideA->SetDirectory(0);
261   }
262   if(ped.GetHistQMaxVsTimeSideC()) {
263     fHistQMaxVsTimeSideC = new TProfile(*ped.GetHistQMaxVsTimeSideC());
264     fHistQMaxVsTimeSideC->SetDirectory(0);
265   }
266   if(ped.GetHistOccupancyVsEventConst()) {
267     fHistOccupancyVsEvent  = new TH1F(*ped.GetHistOccupancyVsEventConst());
268     fHistOccupancyVsEvent->SetDirectory(0);
269   }
270   if(ped.GetHistNclustersVsEventConst()) {
271     fHistNclustersVsEvent  = new TH1F(*ped.GetHistNclustersVsEventConst());
272     fHistNclustersVsEvent->SetDirectory(0);
273   }
274 }
275
276 //_____________________________________________________________________
277 AliTPCdataQA::AliTPCdataQA(const TMap *config) : /*FOLD00*/  
278   TH1F("TPCRAW","TPCRAW",100,0,100),
279   fFirstTimeBin(60),
280   fLastTimeBin(1000),
281   fAdcMin(1),
282   fAdcMax(100),
283   fMapping(NULL),
284   fPedestal(0),
285   fNoise(0),
286   fNLocalMaxima(0),
287   fMaxCharge(0),
288   fMeanCharge(0),
289   fNoThreshold(0),
290   fNTimeBins(0),
291   fNPads(0),
292   fTimePosition(0),
293   fOverThreshold10(0),
294   fOverThreshold20(0),
295   fOverThreshold30(0),
296   fHistQVsTimeSideA(0),
297   fHistQVsTimeSideC(0),
298   fHistQMaxVsTimeSideA(0),
299   fHistQMaxVsTimeSideC(0),
300   fHistOccupancyVsEvent(0),
301   fHistNclustersVsEvent(0),
302   fEventCounter(0),
303   fIsAnalysed(kFALSE),
304   fMaxEvents(500000),
305   fEventsPerBin(1000),
306   fSignalCounter(0),
307   fClusterCounter(0),
308   fAllBins(0),
309   fAllSigBins(0),
310   fAllNSigBins(0),
311   fRowsMax(0),
312   fPadsMax(0),
313   fTimeBinsMax(0),
314   fIsDQM(kFALSE),
315   fHistOccVsSector(0x0),
316   fHistOcc2dVsSector(0x0),
317   fHistQVsSector(0x0),
318   fHistQmaxVsSector(0x0),
319   fOccVec(0x0),
320   fOccMaxVec(0x0),
321   fOccVecFine(0x0),
322   fOccMaxVecFine(0x0)
323 {
324   //
325   // default constructor
326   //
327   if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
328   if (config->GetValue("LastTimeBin"))  fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
329   if (config->GetValue("AdcMin"))       fAdcMin = ((TObjString*)config->GetValue("AdcMin"))->GetString().Atoi();
330   if (config->GetValue("AdcMax"))       fAdcMax = ((TObjString*)config->GetValue("AdcMax"))->GetString().Atoi();
331   if (config->GetValue("MaxEvents"))    fMaxEvents = ((TObjString*)config->GetValue("MaxEvents"))->GetString().Atoi();
332   if (config->GetValue("EventsPerBin")) fEventsPerBin = ((TObjString*)config->GetValue("EventsPerBin"))->GetString().Atoi();
333 }
334
335 //_____________________________________________________________________
336 AliTPCdataQA& AliTPCdataQA::operator = (const  AliTPCdataQA &source)
337 {
338   //
339   // assignment operator
340   //
341   if (&source == this) return *this;
342   new (this) AliTPCdataQA(source);
343
344   return *this;
345 }
346
347
348 //_____________________________________________________________________
349 AliTPCdataQA::~AliTPCdataQA() /*FOLD00*/
350 {
351   //
352   // destructor
353   //
354
355   // do not delete fMapping, because we do not own it.
356   // do not delete fMapping, because we do not own it.
357   // do not delete fNoise and fPedestal, because we do not own them.
358
359   delete fNLocalMaxima;
360   delete fMaxCharge;
361   delete fMeanCharge;
362   delete fNoThreshold;
363   delete fNTimeBins;
364   delete fNPads;
365   delete fTimePosition;
366   delete fOverThreshold10;
367   delete fOverThreshold20;
368   delete fOverThreshold30;
369   delete fHistQVsTimeSideA;
370   delete fHistQVsTimeSideC;
371   delete fHistQMaxVsTimeSideA;
372   delete fHistQMaxVsTimeSideC;
373   delete fHistOccupancyVsEvent;
374   delete fHistNclustersVsEvent;
375
376   // DQM
377   delete fHistOccVsSector;
378   delete fHistOcc2dVsSector;
379   delete fHistQVsSector;
380   delete fHistQmaxVsSector;
381   delete fOccVec;
382   delete fOccMaxVec;
383   delete fOccVecFine;
384   delete fOccMaxVecFine;
385   
386   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
387     delete [] fAllBins[iRow];
388     delete [] fAllSigBins[iRow];
389   }  
390   delete [] fAllBins;
391   delete [] fAllSigBins;
392   delete [] fAllNSigBins;
393 }
394
395 //_____________________________________________________________________
396 TH1F* AliTPCdataQA::GetHistOccupancyVsEvent()
397 {
398   //
399   // Create Occupancy vs event histogram
400   // (we create this histogram differently then the other histograms
401   //  because this we want to be able to access and copy
402   //  from AliTPCQAMakerRec before it normally would be created)
403   //
404   if(!fHistOccupancyVsEvent) {
405
406     Int_t nBins = fMaxEvents/fEventsPerBin;
407     fHistOccupancyVsEvent = new TH1F("hOccupancyVsEvent", "Occupancy vs event number (~time); Event number; Occupancy", nBins, 0, nBins*fEventsPerBin);
408     fHistOccupancyVsEvent->SetDirectory(0);
409   }
410   
411   return fHistOccupancyVsEvent;
412 }
413
414 //_____________________________________________________________________
415 TH1F* AliTPCdataQA::GetHistNclustersVsEvent()
416 {
417   //
418   // Create Nclusters vs event histogram
419   // (we create this histogram differently then the other histograms
420   //  because this we want to be able to access and copy
421   //  from AliTPCQAMakerRec before it normally would be created)
422   //
423   if(!fHistNclustersVsEvent) {
424
425     Int_t nBins = fMaxEvents/fEventsPerBin;
426     fHistNclustersVsEvent = new TH1F("hNclustersVsEvent", "Nclusters vs event number (~time); Event number; Nclusters per event", nBins, 0, nBins*fEventsPerBin);
427     fHistNclustersVsEvent->SetDirectory(0);
428   }
429   
430   return fHistNclustersVsEvent;
431 }
432
433 //_____________________________________________________________________
434 void AliTPCdataQA::UpdateEventHistograms()
435 {
436   // Update histograms that display occupancy and 
437   // number of clusters as a function of number of 
438   // events
439   if (!fHistOccupancyVsEvent)
440     GetHistOccupancyVsEvent();
441   if (!fHistNclustersVsEvent)
442     GetHistNclustersVsEvent();
443   
444   if(fEventCounter > fMaxEvents) {
445     
446     // we have to expand the histogram to handle the larger number of
447     // events. The way it is done now is to double the range and the
448     // number of events per bin (so the number of histogram bins stays
449     // constant)
450     fEventsPerBin *= 2;
451     fMaxEvents *= 2;
452
453     // Change histogram limits
454     const Int_t nBins = fHistOccupancyVsEvent->GetXaxis()->GetNbins();
455     fHistOccupancyVsEvent->GetXaxis()->Set(nBins, fHistOccupancyVsEvent->GetXaxis()->GetNbins(), fMaxEvents);
456     fHistNclustersVsEvent->GetXaxis()->Set(nBins, fHistNclustersVsEvent->GetXaxis()->GetNbins(), fMaxEvents);
457
458     // Rebin the histogram
459     for(Int_t bin = 1; bin <= nBins; bin+=2) {
460
461       Int_t newBin = TMath::Nint(Float_t(bin+1)/2.0);
462       Float_t newContent = (fHistOccupancyVsEvent->GetBinContent(bin)+
463                             fHistOccupancyVsEvent->GetBinContent(bin+1))/2.0;
464       fHistOccupancyVsEvent->SetBinContent(newBin, newContent); 
465
466       newContent = (fHistNclustersVsEvent->GetBinContent(bin)+
467                     fHistNclustersVsEvent->GetBinContent(bin+1))/2.0;
468       fHistNclustersVsEvent->SetBinContent(newBin, newContent); 
469     }
470
471     // Set the remaining bins to 0
472     Int_t lastHalf = nBins/2 +1;
473     for(Int_t bin = lastHalf; bin <= nBins; bin++) {
474
475       fHistOccupancyVsEvent->SetBinContent(bin, 0); 
476       fHistNclustersVsEvent->SetBinContent(bin, 0); 
477     }
478
479     // In this case we should nut update but wait untill the new
480     // number of events per bin is reached!
481     return;
482   }
483
484   const Int_t bin = TMath::Nint(Float_t(fEventCounter)/fEventsPerBin);
485
486   Float_t averageOccupancy =
487     Float_t(fSignalCounter)/fEventsPerBin/(fLastTimeBin - fFirstTimeBin +1.0)
488     / 570132.0; // 570,132 is number of pads
489   fHistOccupancyVsEvent->SetBinContent(bin, averageOccupancy);
490   fSignalCounter = 0;
491   
492   Float_t averageNclusters =
493     Float_t(fClusterCounter)/fEventsPerBin;
494   fHistNclustersVsEvent->SetBinContent(bin, averageNclusters);
495   fClusterCounter = 0;
496 }
497
498 //_____________________________________________________________________
499 Bool_t AliTPCdataQA::ProcessEvent(AliTPCRawStreamV3 *const rawStreamV3)
500 {
501   //
502   // Event Processing loop - AliTPCRawStreamV3
503   //
504   Bool_t withInput = kFALSE;
505   Int_t nSignals = 0;
506   Int_t lastSector = -1;
507   
508   while ( rawStreamV3->NextDDL() ){
509
510     while ( rawStreamV3->NextChannel() ){
511     
512       Int_t iSector = rawStreamV3->GetSector(); //  current sector
513       Int_t iRow    = rawStreamV3->GetRow();    //  current row
514       Int_t iPad    = rawStreamV3->GetPad();    //  current pad
515       Int_t iPatch  = rawStreamV3->GetPatchIndex(); //  current patch
516       Int_t iBranch = rawStreamV3->GetBranch();    //  current branch
517       if (iRow<0 || iPad<0) continue;
518       // Call local maxima finder if the data is in a new sector
519       if(iSector != lastSector) {
520         
521         if(nSignals>0)
522           FindLocalMaxima(lastSector);
523         
524         CleanArrays();
525         lastSector = iSector;
526         nSignals = 0;
527       }
528       
529       while ( rawStreamV3->NextBunch() ){
530
531         Int_t  startTbin    = (Int_t)rawStreamV3->GetStartTimeBin();
532         Int_t  bunchlength  = (Int_t)rawStreamV3->GetBunchLength();
533         const UShort_t *sig = rawStreamV3->GetSignals();
534         
535         for (Int_t iTimeBin = 0; iTimeBin<bunchlength; iTimeBin++){
536           Float_t signal=(Float_t)sig[iTimeBin];
537           nSignals += Update(iSector,iRow,iPad,startTbin--,signal, iPatch, iBranch);
538           withInput = kTRUE;
539         }
540       }
541     }
542   }
543     
544   if (lastSector>=0&&nSignals>0)  
545     FindLocalMaxima(lastSector);
546   
547   return withInput;
548 }
549
550 //_____________________________________________________________________
551 Bool_t AliTPCdataQA::ProcessEvent(AliRawReader *const rawReader)
552 {
553   //
554   //  Event processing loop - AliRawReader
555   //
556   AliTPCRawStreamV3 rawStreamV3(rawReader,(AliAltroMapping**)fMapping);
557   Bool_t res=ProcessEvent(&rawStreamV3);
558   if(res) {
559     fEventCounter++; // only increment event counter if there is TPC data
560
561     if(fEventCounter%fEventsPerBin==0) 
562       UpdateEventHistograms();
563   }
564   return res;
565 }
566
567 //_____________________________________________________________________
568 Bool_t AliTPCdataQA::ProcessEvent(eventHeaderStruct *const event)
569 {
570   //
571   //  process date event
572   //
573
574   AliRawReaderDate rawReader((void*)event);
575   Bool_t result=ProcessEvent(&rawReader);
576   return result;
577 }
578
579
580
581 //_____________________________________________________________________
582 void AliTPCdataQA::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
583 {
584   //
585   //  Write class to file
586   //
587
588   TString sDir(dir);
589   TString option;
590
591   if ( append )
592     option = "update";
593   else
594     option = "recreate";
595
596   TDirectory *backup = gDirectory;
597   TFile f(filename,option.Data());
598   f.cd();
599   if ( !sDir.IsNull() ){
600     f.mkdir(sDir.Data());
601     f.cd(sDir);
602   }
603   this->Write();
604   f.Close();
605
606   if ( backup ) backup->cd();
607 }
608
609
610 //_____________________________________________________________________
611 Int_t AliTPCdataQA::Update(const Int_t iSector, /*FOLD00*/
612                            const Int_t iRow,
613                            const Int_t iPad,
614                            const Int_t iTimeBin,
615                            Float_t signal,
616                            const Int_t iPatch,
617                            const Int_t iBranch)
618 {
619   //
620   // Signal filling method
621   //
622   
623   //
624   // Define the calibration objects the first time Update is called
625   // NB! This has to be done first even if the data is rejected by the time 
626   // cut to make sure that the objects are available in Analyse
627   //
628   if(!fIsDQM) {
629
630     if (!fNLocalMaxima) fNLocalMaxima = new AliTPCCalPad("NLocalMaxima","NLocalMaxima");
631     if (!fMaxCharge) fMaxCharge = new AliTPCCalPad("MaxCharge","MaxCharge");
632     if (!fMeanCharge) fMeanCharge = new AliTPCCalPad("MeanCharge","MeanCharge");
633     if (!fNoThreshold) fNoThreshold = new AliTPCCalPad("NoThreshold","NoThreshold");
634     if (!fNTimeBins) fNTimeBins = new AliTPCCalPad("NTimeBins","NTimeBins");
635     if (!fNPads) fNPads = new AliTPCCalPad("NPads","NPads");
636     if (!fTimePosition) fTimePosition = new AliTPCCalPad("TimePosition","TimePosition");
637     if (!fOverThreshold10) fOverThreshold10 = new AliTPCCalPad("OverThreshold10","OverThreshold10");
638     if (!fOverThreshold20) fOverThreshold20 = new AliTPCCalPad("OverThreshold20","OverThreshold20");
639     if (!fOverThreshold30) fOverThreshold30 = new AliTPCCalPad("OverThreshold30","OverThreshold30");
640     if (!fHistQVsTimeSideA) {
641       fHistQVsTimeSideA  = new TProfile("hQVsTimeSideA", "Q vs time (side A); Time [Timebin]; Q [ADC ch]", 100, 0, 1000);
642       fHistQVsTimeSideA->SetDirectory(0);
643     }
644     if (!fHistQVsTimeSideC) {
645       fHistQVsTimeSideC  = new TProfile("hQVsTimeSideC", "Q vs time (side C); Time [Timebin]; Q [ADC ch]", 100, 0, 1000);
646       fHistQVsTimeSideC->SetDirectory(0);
647   }
648     if (!fHistQMaxVsTimeSideA) {
649       fHistQMaxVsTimeSideA  = new TProfile("hQMaxVsTimeSideA", "Q_{MAX} vs time (side A); Time [Timebin]; Q_{MAX} [ADC ch]", 100, 0, 1000);
650       fHistQMaxVsTimeSideA->SetDirectory(0);
651     }
652     if (!fHistQMaxVsTimeSideC) {
653       fHistQMaxVsTimeSideC  = new TProfile("hQMaxVsTimeSideC", "Q_{MAX} vs time (side C); Time [Timebin]; Q_{MAX} [ADC ch]", 100, 0, 1000);
654       fHistQMaxVsTimeSideC->SetDirectory(0);
655     }
656   } else { // DQM histograms and array
657     
658     if (!fHistOccVsSector) {
659       fHistOccVsSector  = new TProfile("hOccVsSector", "Occupancy vs sector; Sector; Occupancy", 72, 0, 72);
660       fHistOccVsSector->SetDirectory(0);
661
662       fHistOcc2dVsSector  = new TProfile2D("hOcc2dVsSector", "Occupancy vs sector and patch; Sector; Patch", 72, 0, 36, 6, 0, 6);
663       fHistOcc2dVsSector->SetDirectory(0);
664
665       fHistQVsSector  = new TProfile("hQVsSector", "Q vs sector; Sector; Q [ADC ch]", 72, 0, 72);
666       fHistQVsSector->SetDirectory(0);
667
668       fHistQmaxVsSector  = new TProfile("hQmaxVsSector", "Qmax vs sector; Sector; Qmax [ADC ch]", 72, 0, 72);
669       fHistQmaxVsSector->SetDirectory(0);
670
671       fOccVec = new TArrayD(72);
672       for(Int_t i = 0; i < 72; i++)
673         fOccVec->GetArray()[i] = 0;
674
675       fOccMaxVec = new TArrayD(72);
676       const Double_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
677       for(Int_t i = 0; i < 72; i++)
678         
679         if(i<36) // IROCs (5504 pads)
680           fOccMaxVec->GetArray()[i] = nTimeBins*5504;
681         else     // OROCs (9984 pads)
682           fOccMaxVec->GetArray()[i] = nTimeBins*9984;
683
684       // 12 branches for each full sector
685       const Int_t nBranches = 36*12;
686       fOccVecFine = new TArrayD(nBranches);
687       for(Int_t i = 0; i < nBranches; i++)
688         fOccVecFine->GetArray()[i] = 0;
689
690       // Pads per patch same for all sectors
691       Int_t nPads0[6] = {1152, 1536, 1152, 1280, 1280, 1280};
692       Int_t nPads1[6] = {1152, 1664, 1152, 1280, 1280, 1280};
693
694       fOccMaxVecFine = new TArrayD(nBranches);
695       for(Int_t i = 0; i < nBranches; i++) {
696         
697         const Int_t fullSector = Int_t(i/12);
698         Int_t branch = i - fullSector*12;
699         R__ASSERT(branch>=0 && branch<12);
700         
701         const Int_t patch = Int_t(branch/2);
702         branch -= patch*2;
703         
704         R__ASSERT(branch>=0 && branch<2);
705         if(branch == 0)
706           fOccMaxVecFine->GetArray()[i] = nTimeBins*nPads0[patch];
707         else     // OROCs (9984 pads)
708           fOccMaxVecFine->GetArray()[i] = nTimeBins*nPads1[patch];
709       }
710     }
711   }
712   // Make the arrays for expanding the data
713
714   if (!fAllBins)
715     MakeArrays();
716
717   //
718   // If Analyse has been previously called we need now to denormalize the data
719   // as more data is coming
720   //
721   if(fIsAnalysed == kTRUE && !fIsDQM) {
722     
723     const Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
724     const Float_t denormalization = Float_t(fEventCounter * nTimeBins);
725     fNoThreshold->Multiply(denormalization);  
726     
727     fMeanCharge->Multiply(fNLocalMaxima);
728     fMaxCharge->Multiply(fNLocalMaxima);
729     fNTimeBins->Multiply(fNLocalMaxima);
730     fNPads->Multiply(fNLocalMaxima);
731     fTimePosition->Multiply(fNLocalMaxima);
732     fIsAnalysed = kFALSE;
733   }
734
735   //
736   // TimeBin cut
737   //
738   if (iTimeBin<fFirstTimeBin) return 0;
739   if (iTimeBin>fLastTimeBin) return 0;
740   
741   // if pedestal calibrations are loaded subtract pedestals
742   if(fPedestal) {
743
744     Float_t ped = fPedestal->GetCalROC(iSector)->GetValue(iRow, iPad);
745     // Don't use data from pads where pedestals are abnormally small or large
746     if(ped<10 || ped>90)
747       return 0;
748     signal -= ped;
749   }
750   
751   if(fIsDQM) {
752
753     fOccVec->GetArray()[iSector] += 1.0;
754     // To change before committing
755     if(iPatch>=0 && iBranch>=0 && iPatch<=5 && iBranch <= 1)
756       fOccVecFine->GetArray()[(iSector%36)*12+iPatch*2+iBranch] += 1.0;
757   } else {
758     // In fNoThreshold we fill all data to estimate the ZS volume
759     Float_t count = fNoThreshold->GetCalROC(iSector)->GetValue(iRow, iPad);
760     fNoThreshold->GetCalROC(iSector)->SetValue(iRow, iPad,count+1);
761   }
762
763   // Require at least 3 ADC channels
764   if (signal < 3.0)
765     return 0;
766
767   // if noise calibrations are loaded require at least 3*sigmaNoise
768   if(fNoise) {
769     
770     Float_t noise = fNoise->GetCalROC(iSector)->GetValue(iRow, iPad);
771     
772     if(signal < noise*3.0)
773       return 0;
774   }
775
776   //
777   // This signal is ok and we store it in the cluster map
778   //
779
780   SetExpandDigit(iRow, iPad, iTimeBin, signal);
781
782   fSignalCounter++;
783   
784   return 1; // signal was accepted
785 }
786
787 //_____________________________________________________________________
788 void AliTPCdataQA::FindLocalMaxima(const Int_t iSector)
789 {
790   //
791   // This method is called after the data from each sector has been
792   // exapanded into an array
793   // Loop over the signals and identify local maxima and fill the
794   // calibration objects with the information
795   //
796
797   Int_t nLocalMaxima = 0;
798   const Int_t maxTimeBin = fTimeBinsMax+4; // Used to step between neighboring pads 
799                                            // Because we have tha pad-time data in a 
800                                            // 1d array
801
802   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
803
804     Float_t* allBins = fAllBins[iRow];
805     Int_t* sigBins   = fAllSigBins[iRow];
806     const Int_t nSigBins   = fAllNSigBins[iRow];
807     
808     for (Int_t iSig = 0; iSig < nSigBins; iSig++) {
809
810       Int_t bin  = sigBins[iSig];
811       Float_t *b = &allBins[bin];
812
813       //
814       // Now we check if this is a local maximum
815       //
816
817       Float_t qMax = b[0];
818
819       // First check that the charge is bigger than the threshold
820       if (qMax<5) 
821         continue;
822       
823       // Require at least one neighboring pad with signal
824       if (b[-maxTimeBin]+b[maxTimeBin]<=0) continue;
825
826       // Require at least one neighboring time bin with signal
827       if (b[-1]+b[1]<=0) continue;
828       
829       //
830       // Check that this is a local maximum
831       // Note that the checking is done so that if 2 charges has the same
832       // qMax then only 1 cluster is generated 
833       // (that is why there is BOTH > and >=)
834       //
835       if (b[-maxTimeBin]   >= qMax) continue;
836       if (b[-1  ]          >= qMax) continue; 
837       if (b[+maxTimeBin]   > qMax)  continue; 
838       if (b[+1  ]          > qMax)  continue; 
839       if (b[-maxTimeBin-1] >= qMax) continue;
840       if (b[+maxTimeBin-1] >= qMax) continue; 
841       if (b[+maxTimeBin+1] > qMax)  continue; 
842       if (b[-maxTimeBin+1] >= qMax) continue;
843       
844       //
845       // Now we accept the local maximum and fill the calibration/data objects
846       //
847       nLocalMaxima++;
848
849       Int_t iPad, iTimeBin;
850       GetPadAndTimeBin(bin, iPad, iTimeBin);
851       
852       if(!fIsDQM) {
853         Float_t count = fNLocalMaxima->GetCalROC(iSector)->GetValue(iRow, iPad);
854         fNLocalMaxima->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
855         
856         count = fTimePosition->GetCalROC(iSector)->GetValue(iRow, iPad);
857         fTimePosition->GetCalROC(iSector)->SetValue(iRow, iPad, count+iTimeBin);
858       
859         Float_t charge = fMaxCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
860         fMaxCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qMax);
861         
862         if(qMax>=10) {
863           count = fOverThreshold10->GetCalROC(iSector)->GetValue(iRow, iPad);
864           fOverThreshold10->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
865         }
866         if(qMax>=20) {
867           count = fOverThreshold20->GetCalROC(iSector)->GetValue(iRow, iPad);
868           fOverThreshold20->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
869         }
870         if(qMax>=30) {
871           count = fOverThreshold30->GetCalROC(iSector)->GetValue(iRow, iPad);
872           fOverThreshold30->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
873         }
874       }
875
876       //
877       // Calculate the total charge as the sum over the region:
878       //
879       //    o o o o o
880       //    o i i i o
881       //    o i C i o
882       //    o i i i o
883       //    o o o o o
884       //
885       // with qmax at the center C.
886       //
887       // The inner charge (i) we always add, but we only add the outer
888       // charge (o) if the neighboring inner bin (i) has a signal.
889       //
890       Int_t minP = 0, maxP = 0, minT = 0, maxT = 0;
891       Float_t qTot = qMax;
892       for(Int_t i = -1; i<=1; i++) {
893         for(Int_t j = -1; j<=1; j++) {
894           
895           if(i==0 && j==0)
896             continue;
897           
898           Float_t charge1 = GetQ(b, i, j, maxTimeBin, minT, maxT, minP, maxP);
899           qTot += charge1;
900           if(charge1>0) {
901             // see if the next neighbor is also above threshold
902             if(i*j==0) {
903               qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP); 
904             } else {
905               // we are in a diagonal corner
906               qTot += GetQ(b,   i, 2*j, maxTimeBin, minT, maxT, minP, maxP); 
907               qTot += GetQ(b, 2*i,   j, maxTimeBin, minT, maxT, minP, maxP); 
908               qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP); 
909             }
910           }
911         }
912       }
913       
914       if(fIsDQM) {
915         fHistQVsSector->Fill(iSector, qTot);
916         fHistQmaxVsSector->Fill(iSector, qMax);
917       } else {
918         Float_t charge = fMeanCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
919         fMeanCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qTot);
920         
921         Float_t count = fNTimeBins->GetCalROC(iSector)->GetValue(iRow, iPad);
922         fNTimeBins->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxT-minT+1);
923         
924         count = fNPads->GetCalROC(iSector)->GetValue(iRow, iPad);
925         fNPads->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxP-minP+1);
926       
927         if((iSector%36)<18) { // A side
928           fHistQVsTimeSideA->Fill(iTimeBin, qTot);
929           fHistQMaxVsTimeSideA->Fill(iTimeBin, qMax);
930         } else {
931           fHistQVsTimeSideC->Fill(iTimeBin, qTot);
932           fHistQMaxVsTimeSideC->Fill(iTimeBin, qMax);      
933         }
934       }
935     } // end loop over signals
936   } // end loop over rows
937   
938   fClusterCounter += nLocalMaxima;
939 }
940
941 //_____________________________________________________________________
942 void AliTPCdataQA::Analyse()
943 {
944   //
945   //  Calculate calibration constants
946   //
947   
948   AliInfo("Analyse called");
949
950   if(fIsDQM == kTRUE) {
951     
952     AliInfo("DQM flas is set -> No 2d information to analyze");
953     return;
954   }
955
956   if(fIsAnalysed == kTRUE) {
957     
958     AliInfo("No new data since Analyse was called last time");
959     return;
960   }
961
962   if(fEventCounter==0) {
963     
964       AliInfo("EventCounter == 0, Cannot analyse");
965     return;
966   }
967   
968   Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
969   AliInfo(Form("EventCounter: %d , TimeBins: %d", fEventCounter, nTimeBins));
970
971   Float_t normalization = 1.0 / Float_t(fEventCounter * nTimeBins);
972   fNoThreshold->Multiply(normalization);  
973   
974   fMeanCharge->Divide(fNLocalMaxima);
975   fMaxCharge->Divide(fNLocalMaxima);
976   fNTimeBins->Divide(fNLocalMaxima);
977   fNPads->Divide(fNLocalMaxima);
978   fTimePosition->Divide(fNLocalMaxima);
979
980   fIsAnalysed = kTRUE;
981 }
982
983
984 //_____________________________________________________________________
985 void AliTPCdataQA::MakeTree(const char *fname) const {
986   //
987   // Export result to the tree -located in the file
988   // This file can be analyzed using AliTPCCalibViewer
989   // 
990   AliTPCPreprocessorOnline preprocesor;
991
992   if (fNLocalMaxima) preprocesor.AddComponent(fNLocalMaxima);
993   if (fMaxCharge) preprocesor.AddComponent(fMaxCharge);  
994   if (fMeanCharge) preprocesor.AddComponent(fMeanCharge);  
995   if (fNoThreshold) preprocesor.AddComponent(fNoThreshold);
996   if (fNTimeBins) preprocesor.AddComponent(fNTimeBins);
997   if (fNPads) preprocesor.AddComponent(fNPads);
998   if (fTimePosition) preprocesor.AddComponent(fTimePosition);
999   if (fOverThreshold10) preprocesor.AddComponent(fOverThreshold10);
1000   if (fOverThreshold20) preprocesor.AddComponent(fOverThreshold20);
1001   if (fOverThreshold30) preprocesor.AddComponent(fOverThreshold30);
1002
1003   preprocesor.DumpToFile(fname);  
1004 }
1005
1006
1007 //_____________________________________________________________________
1008 void AliTPCdataQA::MakeArrays(){
1009   //
1010   // The arrays for expanding the raw data are defined and 
1011   // som parameters are intialised
1012   //
1013   AliTPCROC * roc = AliTPCROC::Instance();
1014   //
1015   // To make the array big enough for all sectors we take 
1016   // the dimensions from the outer row of an OROC (the last sector)
1017   //
1018   fRowsMax     = roc->GetNRows(roc->GetNSector()-1);
1019   fPadsMax     = roc->GetNPads(roc->GetNSector()-1,fRowsMax-1);
1020   fTimeBinsMax = fLastTimeBin - fFirstTimeBin +1; 
1021
1022   //
1023   // Since we have added 2 pads (TimeBins) before and after the real pads (TimeBins) 
1024   // to make sure that we can always query the exanded table even when the 
1025   // max is on the edge
1026   //
1027
1028  
1029   fAllBins = new Float_t*[fRowsMax];
1030   fAllSigBins = new Int_t*[fRowsMax];
1031   fAllNSigBins = new Int_t[fRowsMax];
1032
1033   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
1034     //
1035     Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);  
1036     fAllBins[iRow] = new Float_t[maxBin];
1037     memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin); // set all values to 0
1038     fAllSigBins[iRow] = new Int_t[maxBin];
1039     fAllNSigBins[iRow] = 0;
1040   }
1041 }
1042
1043
1044 //_____________________________________________________________________
1045 void AliTPCdataQA::CleanArrays(){
1046   //
1047   //
1048   //
1049
1050   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
1051
1052     // To speed up the performance by a factor 2 on cosmic data (and
1053     // presumably pp data as well) where the ocupancy is low, the
1054     // memset is only called if there is more than 1000 signals for a
1055     // row (of the order 1% occupancy)
1056     if(fAllNSigBins[iRow]<1000) {
1057       
1058       Float_t* allBins = fAllBins[iRow];
1059       Int_t* sigBins   = fAllSigBins[iRow];
1060       const Int_t nSignals = fAllNSigBins[iRow];
1061       for(Int_t i = 0; i < nSignals; i++)
1062         allBins[sigBins[i]]=0;      
1063     } else {
1064
1065       Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4); 
1066       memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1067     }
1068
1069     fAllNSigBins[iRow]=0;
1070   }
1071 }
1072
1073 //_____________________________________________________________________
1074 void AliTPCdataQA::GetPadAndTimeBin(Int_t bin, Int_t& iPad, Int_t& iTimeBin){
1075   //
1076   // Return pad and timebin for a given bin
1077   //
1078   
1079   //  Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
1080   iTimeBin  = bin%(fTimeBinsMax+4);
1081   iPad      = (bin-iTimeBin)/(fTimeBinsMax+4);
1082
1083   iPad     -= 2;
1084   iTimeBin -= 2;
1085   iTimeBin += fFirstTimeBin;
1086
1087   R__ASSERT(iPad>=0 && iPad<=fPadsMax);
1088   R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
1089 }
1090
1091 //_____________________________________________________________________
1092 void AliTPCdataQA::SetExpandDigit(const Int_t iRow, Int_t iPad, 
1093                                   Int_t iTimeBin, const Float_t signal) 
1094 {
1095   //
1096   // 
1097   //
1098   R__ASSERT(iRow>=0 && iRow<fRowsMax);
1099   R__ASSERT(iPad>=0 && iPad<=fPadsMax);
1100   R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
1101   
1102   iTimeBin -= fFirstTimeBin;
1103   iPad     += 2;
1104   iTimeBin += 2;
1105   
1106   Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
1107   fAllBins[iRow][bin] = signal;
1108   fAllSigBins[iRow][fAllNSigBins[iRow]] = bin;
1109   fAllNSigBins[iRow]++;
1110 }
1111
1112 //______________________________________________________________________________
1113 Float_t AliTPCdataQA::GetQ(const Float_t* adcArray, const Int_t time, 
1114                            const Int_t pad, const Int_t maxTimeBins, 
1115                            Int_t& timeMin, Int_t& timeMax, 
1116                            Int_t& padMin,  Int_t& padMax) const
1117 {
1118   //
1119   // This methods return the charge in the bin time+pad*maxTimeBins
1120   // If the charge is above 0 it also updates the padMin, padMax, timeMin
1121   // and timeMax if necessary
1122   //
1123   Float_t charge = adcArray[time + pad*maxTimeBins];
1124   if(charge > 0) {
1125     timeMin = TMath::Min(time, timeMin); timeMax = TMath::Max(time, timeMax);
1126     padMin = TMath::Min(pad, padMin); padMax = TMath::Max(pad, padMax);
1127   }
1128   return charge; 
1129 }
1130
1131 //______________________________________________________________________________
1132 void AliTPCdataQA::Streamer(TBuffer &xRuub)
1133 {
1134   // Automatic schema evolution was first used from revision 4
1135   // Code based on:
1136   // http://root.cern.ch/root/roottalk/roottalk02/3207.html
1137
1138    UInt_t xRuus, xRuuc;
1139    if (xRuub.IsReading()) {
1140       Version_t xRuuv = xRuub.ReadVersion(&xRuus, &xRuuc);
1141       //we use the automatic algorithm for class version > 3
1142       if (xRuuv > 3) {
1143         AliTPCdataQA::Class()->ReadBuffer(xRuub, this, xRuuv, xRuus,
1144                                           xRuuc);
1145         return;
1146       }
1147       TH1F::Streamer(xRuub);
1148       xRuub >> fFirstTimeBin;
1149       xRuub >> fLastTimeBin;
1150       xRuub >> fAdcMin;
1151       xRuub >> fAdcMax;
1152       xRuub >> fNLocalMaxima;
1153       xRuub >> fMaxCharge;
1154       xRuub >> fMeanCharge;
1155       xRuub >> fNoThreshold;
1156       xRuub >> fNTimeBins;
1157       xRuub >> fNPads;
1158       xRuub >> fTimePosition;
1159       xRuub >> fEventCounter;
1160       xRuub >> fIsAnalysed;
1161       xRuub.CheckByteCount(xRuus, xRuuc, AliTPCdataQA::IsA());
1162    } else {
1163      AliTPCdataQA::Class()->WriteBuffer(xRuub,this);
1164    }
1165 }
1166
1167 //____________________________________________________________________________________________
1168 void AliTPCdataQA::FillOccupancyProfile()
1169 {
1170   // This has to be filled at the end of the loop over data
1171   if(!fIsDQM) 
1172     AliInfo("Method only meaningful for DQM");
1173   
1174   for(Int_t i = 0; i < 72; i++) {
1175
1176     fOccVec->GetArray()[i] /= fOccMaxVec->GetArray()[i];
1177     fHistOccVsSector->Fill(i, fOccVec->GetArray()[i]);
1178   }
1179
1180   const Int_t nBranches = 36*12;
1181   for(Int_t i = 0; i < nBranches; i++) {
1182
1183     fOccVecFine->GetArray()[i] /= fOccMaxVecFine->GetArray()[i];
1184
1185     const Int_t fullSector = Int_t(i/12);
1186
1187     Int_t branch = i - fullSector*12;
1188     const Int_t patch = Int_t(branch/2);
1189
1190     branch -= patch*2;
1191
1192     fHistOcc2dVsSector->Fill(fullSector+0.5*branch+0.1, patch+0.5, fOccVecFine->GetArray()[i]);
1193   }
1194 }
1195
1196 //____________________________________________________________________________________________
1197 void AliTPCdataQA::ResetProfiles()
1198 {
1199   if(!fIsDQM) 
1200     AliInfo("Method only meaningful for DQM");
1201   
1202   if(fHistQVsSector)
1203     fHistQVsSector->Reset();
1204   if(fHistQmaxVsSector)
1205     fHistQmaxVsSector->Reset();
1206   if(fHistOccVsSector)
1207     fHistOccVsSector->Reset();
1208   if(fHistOcc2dVsSector)
1209     fHistOcc2dVsSector->Reset();
1210
1211   if(fOccVec)
1212     for(Int_t i = 0; i < 72; i++)
1213       fOccVec->GetArray()[i] = 0.0;
1214   if(fOccVecFine)
1215     for(Int_t i = 0; i < 36*12; i++)
1216       fOccVecFine->GetArray()[i] = 0.0;
1217 }