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