]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCCalibRaw.cxx
Extracting Branch and Revision from Git.
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCCalibRaw.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 /* $Id:  $ */
17
18 ////////////////////////////////////////////////////////////////////////////////////////
19 //                                                                                    //
20 //     Implementation of the TPC Raw drift velocity and Altro L1 Phase  calibration   //
21 //                                                                                    //
22 //               Origin: Jens Wiechula, J.Wiechula@gsi.de                             //
23 //                                                                                    //
24 ////////////////////////////////////////////////////////////////////////////////////////
25 //
26 //
27 // *************************************************************************************
28 // *                                Class Description                                  *
29 // *************************************************************************************
30 /*
31
32 ----example---
33 TFile f("CalibAltro.root");
34 AliTPCCalibRaw *al=(AliTPCCalibRaw*)f.Get(f.GetListOfKeys()->At(0)->GetName())
35 {
36 TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
37 if (!c1) c1=new TCanvas("c1","c1");
38 c1->Clear();
39
40 TH2F h2f("h2","h2;RCU;fail",216,0,216,al->GetNevents(),0,al->GetNevents());
41 Bool_t first=kTRUE;
42 Int_t i,iev;
43 for (i=0;i<216;++i) {
44   TVectorF *v=al->GetALTROL1PhaseFailEventsRCU(i);
45   if (!v) continue;
46   for (iev=0;iev<al->GetNevents();++iev) {
47     h2f->SetBinContent(i+1,iev+1,(*v)(iev));
48   }
49 //   TH1F h(*v);
50 //   h.SetLineColor(i/216.*50+50);
51 //   ((TH1F*)h.Clone(Form("h%d",i)))->Draw(first?"":"same");
52 //   c1->Modified();
53 //   c1->Update();
54   first=kFALSE;
55 }
56 h2f->Draw("col");
57 }
58
59 */
60
61
62
63 //Root includes
64 #include <TROOT.h>
65 #include <TH2C.h>
66 #include <TH1F.h>
67 #include <TMap.h>
68 #include <TGraph.h>
69 #include <TObjString.h>
70 #include <TTimeStamp.h>
71 #include <TCanvas.h>
72 //AliRoot includes
73 #include "AliTPCCalROC.h"
74 #include "AliAltroRawStream.h"
75 #include "AliLog.h"
76 //class header
77 #include "AliTPCCalibRaw.h"
78
79 ClassImp(AliTPCCalibRaw)
80
81 AliTPCCalibRaw::AliTPCCalibRaw() :
82   AliTPCCalibRawBase(),
83   fPeakDetMinus(1),
84   fPeakDetPlus(2),
85   fNFailL1Phase(0),
86   fNFailL1PhaseEvent(0),
87   fNSecTime(600), //default 10 minutes
88   fNBinsTime(60), //default 60*10 minutes = 10 hours
89   fPadProcessed(kFALSE),
90   fCurrentChannel(-1),
91   fCurrentSector(-1),
92   fLastSector(-2),
93   fCurrentRow(-1),
94   fCurrentPad(-1),
95   fLastTimeBinProc(0),
96   fPeakTimeBin(0),
97   fLastSignal(0),
98   fNOkPlus(0),
99   fNOkMinus(0),
100   fNanoSec(0),
101   fArrCurrentPhaseDist(4),
102   fArrCurrentPhase(kNRCU),
103   fArrFailEventNumber(100),
104   fArrALTROL1Phase(100000),
105   fArrALTROL1PhaseEvent(kNRCU),
106   fArrALTROL1PhaseFailEvent(kNRCU),
107   fHnDrift(0x0),
108   fVOccupancyEvent(100000),
109   fVSignalSumEvent(100000),
110   fVOccupancySenEvent(100000),
111   fVSignalSumSenEvent(100000),
112   fVNfiredPadsSenEvent(100000),
113   fVTimeStampEvent(100000)
114 {
115   //
116   // Default ctor
117   //
118   SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
119   CreateDVhist();
120   for (Int_t ircu=0;ircu<kNRCU;++ircu) fArrCurrentPhase.GetMatrixArray()[ircu]=-1;
121   fFirstTimeBin=850;
122   fLastTimeBin=1020;
123 }
124 //_____________________________________________________________________
125 AliTPCCalibRaw::AliTPCCalibRaw(const TMap *config) :
126 AliTPCCalibRawBase(),
127 fPeakDetMinus(1),
128 fPeakDetPlus(2),
129 fNFailL1Phase(0),
130 fNFailL1PhaseEvent(0),
131 fNSecTime(600), //default 10 minutes
132 fNBinsTime(60), //default 60*10 minutes = 10 hours
133 fPadProcessed(kFALSE),
134 fCurrentChannel(-1),
135 fCurrentSector(-1),
136 fLastSector(-2),
137 fCurrentRow(-1),
138 fCurrentPad(-1),
139 fLastTimeBinProc(0),
140 fPeakTimeBin(0),
141 fLastSignal(0),
142 fNOkPlus(0),
143 fNOkMinus(0),
144 fNanoSec(0),
145 fArrCurrentPhaseDist(4),
146 fArrCurrentPhase(kNRCU),
147 fArrFailEventNumber(100),
148 fArrALTROL1Phase(100000),
149 fArrALTROL1PhaseEvent(kNRCU),
150 fArrALTROL1PhaseFailEvent(kNRCU),
151 fHnDrift(0x0),
152 fVOccupancyEvent(100000),
153 fVSignalSumEvent(100000),
154 fVOccupancySenEvent(100000),
155 fVSignalSumSenEvent(100000),
156 fVNfiredPadsSenEvent(100000),
157 fVTimeStampEvent(100000)
158 {
159   //
160   // Default ctor
161   //
162   SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
163   CreateDVhist();
164   for (Int_t ircu=0;ircu<kNRCU;++ircu) fArrCurrentPhase.GetMatrixArray()[ircu]=-1;
165   fFirstTimeBin=850;
166   fLastTimeBin=1020;
167   if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
168   if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
169   if (config->GetValue("DebugLevel")) fDebugLevel = ((TObjString*)config->GetValue("DebugLevel"))->GetString().Atoi();
170 }
171
172 //_____________________________________________________________________
173 AliTPCCalibRaw::~AliTPCCalibRaw()
174 {
175   //
176   // dtor
177   //
178   delete fHnDrift;  
179 }
180 //_____________________________________________________________________
181 // AliTPCCalibRaw& AliTPCCalibRaw::operator = (const  AliTPCCalibRaw &source)
182 // {
183 //   //
184 //   // assignment operator
185 //   //
186 //   if (&source == this) return *this;
187 //   new (this) AliTPCCalibRaw(source);
188 //   
189 //   return *this;
190 // }
191
192 //_____________________________________________________________________
193 Int_t AliTPCCalibRaw::Update(const Int_t isector, const Int_t iRow, const Int_t iPad,
194              const Int_t iTimeBin, const Float_t signal)
195 {
196   //
197   // Data filling method
198   //
199   if (iRow<0) return 0;
200   if (iPad<0) return 0;
201   if (iTimeBin<0) return 0;
202   //
203   Int_t iChannel  = fROC->GetRowIndexes(isector)[iRow]+iPad; //  global pad position in sector
204   //occupancy
205   fVOccupancyEvent.GetMatrixArray()[GetNevents()]++;
206   fVSignalSumEvent.GetMatrixArray()[GetNevents()]+=signal;
207   //occupancy in sensitive regions
208   Int_t npads=(Int_t)fROC->GetNPads(isector,iRow);
209   Int_t cpad=iPad-npads/2;
210   if (isector<(Int_t)fROC->GetNInnerSector()){
211     //IROC case (spot)
212     if ( iRow>19 && iRow<46 ){
213       if ( TMath::Abs(cpad)<7 ){
214         fVOccupancySenEvent.GetMatrixArray()[GetNevents()]++;
215         fVSignalSumSenEvent.GetMatrixArray()[GetNevents()]+=signal;
216         if (iChannel!=fCurrentChannel) fVNfiredPadsSenEvent.GetMatrixArray()[GetNevents()]++;
217       }
218     }
219   } else if ( iRow>75 ){
220     //OROC case (outer corners and last three rows are sensitive)
221     Int_t padEdge=(Int_t)TMath::Min(iPad,npads-iPad);
222     Int_t nrows=(Int_t)fROC->GetNRows(isector);
223     if ((nrows-iRow-1)<3 || padEdge<((((Int_t)iRow-76)/4+1))*2){
224       fVOccupancySenEvent.GetMatrixArray()[GetNevents()]++;
225       fVSignalSumSenEvent.GetMatrixArray()[GetNevents()]+=signal;
226       if (iChannel!=fCurrentChannel) fVNfiredPadsSenEvent.GetMatrixArray()[GetNevents()]++;
227     }
228   }
229   //
230   if ( (iTimeBin>fLastTimeBin) || (iTimeBin<fFirstTimeBin)   ) return 0;
231   //don't process edge pads
232   if (IsEdgePad(isector,iRow,iPad)) return 0;
233 //   Double_t x[kHnBinsDV]={1,isector,0};
234 //   fHnDrift->Fill(x);
235   if (fCurrentChannel==iChannel){
236     if (fPadProcessed) return 0;
237   } else {
238     fPadProcessed=kFALSE;
239     fNOkPlus=0;
240     fNOkMinus=0;
241     fPeakTimeBin=0;
242     fLastSignal=0;
243   }
244 //   Double_t x2[kHnBinsDV]={2,isector,0};
245 //   fHnDrift->Fill(x2);
246   
247
248   if (signal>fLastSignal) ++fNOkPlus;
249   else if(signal<fLastSignal && fNOkPlus>=fPeakDetPlus){
250     ++fNOkMinus;
251     if (!fPeakTimeBin) fPeakTimeBin=fLastTimeBinProc;
252     if ( fNOkMinus>=fPeakDetMinus ) {
253       Double_t x[kHnBinsDV]={fPeakTimeBin,isector,(fTimeStamp-fFirstTimeStamp)/fNSecTime};
254       fHnDrift->Fill(x);
255     } 
256   } else {
257     fNOkPlus=0;
258     fNOkMinus=0;
259     fPeakTimeBin=0;
260     fLastSignal=0;
261   }
262
263   fLastTimeBinProc=iTimeBin;
264   fLastSignal=TMath::Nint(signal);
265   fCurrentChannel = iChannel;
266   return 0;
267 }
268 //_____________________________________________________________________
269 void AliTPCCalibRaw::UpdateDDL(){
270   //
271   // fill ALTRO L1 information
272   //
273   
274   //set nanoseconds
275   if (!fNanoSec) {
276     TTimeStamp s;
277     fNanoSec=s.GetNanoSec();
278   }
279   // current phase
280   Int_t phase=(Int_t)(GetL1PhaseTB()*4.);
281   //Fill pahse information of current rcu and event
282   fArrCurrentPhase.GetMatrixArray()[fCurrDDLNum]=phase;
283   //increase phase counter
284   ++((fArrCurrentPhaseDist.GetMatrixArray())[phase]);
285   
286 }
287 //_____________________________________________________________________
288 void AliTPCCalibRaw::ResetEvent()
289 {
290   //
291   // Reset event counters
292   //
293
294   fCurrentChannel=-1;
295   fCurrentRow=-1;
296   fCurrentPad=-1;
297   fArrCurrentPhaseDist.Zero();
298 }
299 //_____________________________________________________________________
300 void AliTPCCalibRaw::EndEvent()
301 {
302   //
303   // End event analysis
304   //
305
306   
307   //find phase of the current event
308   Int_t phaseMaxEntries=-1;
309   Int_t maxEntries=0;
310   for (Int_t i=0;i<fArrCurrentPhaseDist.GetNrows();++i){
311     Int_t entries=(Int_t)fArrCurrentPhaseDist[i];
312     if (maxEntries<entries) {
313       maxEntries=entries;
314       phaseMaxEntries=i;
315     }
316   }
317   // store phase of current event
318   if (fArrALTROL1Phase.GetNrows()-1<=GetNevents())
319     fArrALTROL1Phase.ResizeTo(GetNevents()+10000);
320   (fArrALTROL1Phase.GetMatrixArray())[GetNevents()]=phaseMaxEntries;
321   
322   //loop over RCUs and test failures
323   UInt_t fail=0;
324   for (Int_t ircu=0;ircu<kNRCU;++ircu){
325     Int_t phase=(Int_t)fArrCurrentPhase[ircu];
326     if (phase<0) continue;
327     if (phase!=phaseMaxEntries){
328       TVectorF *arr=MakeArrL1PhaseRCU(fCurrDDLNum,kTRUE);
329       if (arr->GetNrows()-1<=(Int_t)fNFailL1PhaseEvent) arr->ResizeTo(arr->GetNrows()+100);
330       (arr->GetMatrixArray())[fNFailL1PhaseEvent]=phase;
331       ++fNFailL1Phase;
332       fail=1;
333       }
334     //reset current phase information
335     fArrCurrentPhase[ircu]=-1;
336   }
337   if (fail){
338     if (fArrFailEventNumber.GetNrows()-1<=(Int_t)fNFailL1PhaseEvent) fArrFailEventNumber.ResizeTo(fArrFailEventNumber.GetNrows()+100);
339     fArrFailEventNumber.GetMatrixArray()[fNFailL1PhaseEvent]=GetNevents();
340   }
341   fNFailL1PhaseEvent+=fail;
342   //time stamps
343   fVTimeStampEvent.GetMatrixArray()[GetNevents()]=GetTimeStamp()-fFirstTimeStamp+fNanoSec*1e-9;
344   fNanoSec=0;
345   //occupance related
346   if (fVOccupancyEvent.GetNrows()-1<=GetNevents()){
347     fVOccupancyEvent.ResizeTo(GetNevents()+10000);
348     fVSignalSumEvent.ResizeTo(GetNevents()+10000);
349     fVOccupancySenEvent.ResizeTo(GetNevents()+10000);
350     fVSignalSumSenEvent.ResizeTo(GetNevents()+10000);
351     fVTimeStampEvent.ResizeTo(GetNevents()+10000);
352     fVNfiredPadsSenEvent.ResizeTo(GetNevents()+10000);
353   }
354   IncrementNevents();
355 }
356 //_____________________________________________________________________
357 TH2C *AliTPCCalibRaw::MakeHistL1RCUEvents(Int_t type)
358 {
359   // Create a 2D histo RCU:Events indicating the there was a deviation
360   // from the mean L1 phase of the event
361   //
362   //type: 0=Failures, 1=Phases
363
364   //number of relavant events, depending on version
365   Int_t nevents=GetNevents();
366   //check version
367   Bool_t newVersion=kFALSE;
368   for (Int_t ircu=0; ircu<kNRCU; ++ircu){
369     const TVectorF *v=GetALTROL1PhaseEventsRCU(ircu);
370     if (!v) continue;
371     if ((UInt_t)(v->GetNrows())==fNFailL1PhaseEvent){
372       newVersion=kTRUE;
373       nevents=fNFailL1PhaseEvent;
374     }
375     break;
376   }
377   TH2C *h2 = new TH2C("hL1FailRCUEvents","L1 Failures;RCU;Event",kNRCU,0,kNRCU,nevents,0,nevents);
378   Int_t add=0;
379   for (Int_t ircu=0;ircu<kNRCU;++ircu) {
380     const TVectorF *v=GetALTROL1PhaseEventsRCU(ircu);
381     if (type==0){
382       add=1;
383       h2->SetMinimum(0);
384       h2->SetMaximum(2);
385     } else if (type==1) {
386       add=0;
387       h2->SetMinimum(0);
388       h2->SetMaximum(4);
389     }
390     if (!v) continue;
391     for (Int_t iev=0;iev<nevents;++iev) {
392       Float_t val=(*v)(iev);
393       Float_t phase=fArrALTROL1Phase.GetMatrixArray()[iev];
394       if (newVersion) {
395         Int_t event=(Int_t)fArrFailEventNumber.GetMatrixArray()[iev];
396         phase=fArrALTROL1Phase.GetMatrixArray()[event];
397       }
398       if (type==0) val=(val!=phase);
399       h2->SetBinContent(ircu+1,iev+1,val+add);
400     }
401   }
402   return h2;
403 }
404 //_____________________________________________________________________
405 TH1F *AliTPCCalibRaw::MakeHistL1PhaseDist()
406 {
407   //
408   // L1 phase distribution. Should be flat in ideal case
409   //
410   TH1F *h=new TH1F("L1phaseDist","Normalized L1 phase distribution;phase;fraction of events",4,0,4);
411   h->Sumw2();
412   for (Int_t iev=0;iev<GetNevents();++iev) h->Fill(fArrALTROL1Phase.GetMatrixArray()[iev]);
413   if (GetNevents()>0) h->Scale(1./GetNevents());
414   h->SetMinimum(0);
415   h->SetMaximum(1);
416   return h;
417 }
418 //_____________________________________________________________________
419 TVectorF *AliTPCCalibRaw::MakeVectL1PhaseDist()
420 {
421   //
422   // L1 phase distribution. Should be flat in ideal case
423   //
424   TVectorF *v=new TVectorF(4);
425   for (Int_t iev=0;iev<GetNevents();++iev) {
426     Int_t phase=(Int_t)fArrALTROL1Phase.GetMatrixArray()[iev];
427     ((v->GetMatrixArray())[phase])+=1./GetNevents();
428   }
429   return v;
430 }
431 //_____________________________________________________________________
432 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsIROC(Int_t type)
433 {
434   //
435   // Create a 2D histo RCU:Events indicating the there was a deviation
436   // from the mean L1 phase of the event
437   //
438   TH2C *h2 = new TH2C("hL1FailRCUEventsIROC","L1 Failures IROCs;RCU;Event",72,0,36,GetNevents(),0,GetNevents());
439   for (Int_t ircu=0;ircu<72;++ircu) {
440     const TVectorF *v=0;
441     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
442     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
443     if (!v) continue;
444     for (Int_t iev=0;iev<GetNevents();++iev) {
445       h2->SetBinContent(ircu+1,iev+1,(*v)(iev));
446     }
447   }
448   return h2;
449 }
450 //_____________________________________________________________________
451 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsOROC(Int_t type)
452 {
453   //
454   // Create a 2D histo RCU:Events indicating the there was a deviation
455   // from the mean L1 phase of the event
456   //
457   TH2C *h2 = new TH2C("hL1FailRCUEventsOROC","L1 Failures OROCs;RCU;Event",144,0,36,GetNevents(),0,GetNevents());
458   for (Int_t ircu=72;ircu<kNRCU;++ircu) {
459     const TVectorF *v=0;
460     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
461     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
462     if (!v) continue;
463     for (Int_t iev=0;iev<GetNevents();++iev) {
464       h2->SetBinContent(ircu-72+1,iev+1,(*v)(iev));
465     }
466   }
467   return h2;
468 }
469 //_____________________________________________________________________
470 void AliTPCCalibRaw::CreateDVhist()
471 {
472   //
473   // Setup the HnSparse for the drift velocity determination
474   //
475   if (fHnDrift) return;
476   //HnSparse bins
477   //time bin, roc, time
478   Int_t    bins[kHnBinsDV] = {fLastTimeBin-fFirstTimeBin, 72, fNBinsTime};
479   Double_t xmin[kHnBinsDV] = {fFirstTimeBin,0,0};
480   Double_t xmax[kHnBinsDV] = {fLastTimeBin,72,fNBinsTime};
481   fHnDrift=new THnSparseI("fHnDrift",Form("Drift velocity using last time bin;time bin[#times 100ns];ROC;Time bin [#times %us]",fNSecTime),kHnBinsDV, bins, xmin, xmax);
482     
483 }
484 //_____________________________________________________________________
485 void AliTPCCalibRaw::Analyse()
486 {
487   //
488   // Analyse Data
489   //
490
491   //resize arrays
492   fArrALTROL1Phase.ResizeTo(GetNevents());
493   for (Int_t ircu=0;ircu<kNRCU;++ircu){
494     TVectorF *arr=MakeArrL1PhaseRCU(ircu);//MakeArrL1PhaseRCU(ircu);
495     if (!arr) continue;
496     arr->ResizeTo(fNFailL1PhaseEvent);
497     fArrFailEventNumber.ResizeTo(fNFailL1PhaseEvent);
498 //    TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu);
499 //     arrF->ResizeTo(1);
500   }
501   //resize occupancy arrays only save event occupancy in sensitive regions by default
502   //save the rest in debub mode
503   fVOccupancySenEvent.ResizeTo(GetNevents());
504   if (fDebugLevel>0){
505     fVOccupancyEvent.ResizeTo(GetNevents());
506     fVSignalSumEvent.ResizeTo(GetNevents());
507     fVSignalSumSenEvent.ResizeTo(GetNevents());
508     fVNfiredPadsSenEvent.ResizeTo(GetNevents());
509     fVTimeStampEvent.ResizeTo(GetNevents());
510   } else {
511     fVOccupancyEvent.ResizeTo(0);
512     fVSignalSumEvent.ResizeTo(0);
513     fVSignalSumSenEvent.ResizeTo(0);
514     fVNfiredPadsSenEvent.ResizeTo(0);
515     fVTimeStampEvent.ResizeTo(0);
516   }
517   //Analyse drift velocity TODO
518   
519 }
520 //_____________________________________________________________________
521 TGraph* AliTPCCalibRaw::MakeGraphOccupancy(const Int_t type, const Int_t xType)
522 {
523   //
524   // create occupancy graph (of samples abouve threshold)
525   // type=0:   number of samples
526   // type=1:   mean data volume (ADC counts/sample)
527   // type=2:   data volume (ADC counts)
528   // type=3:   samples per ADC count
529   // type=4:   sample occupancy
530   //
531   // type=5: number of sample sensitive / number of samples
532   //
533   // same in sensitive regions:
534   // type=10:  number of samples
535   // type=11:  mean data volume (ADC counts/sample)
536   // type=12:  data volume (ADC counts)
537   // type=13:  samples per ADC count
538   // type=14:  sample occupancy
539   //
540   // type=16: number of samples sensitive / number of pads sensitive
541   // type=17: pad occupancy in sensitive regions
542   // xType=0:  vs. time stamp
543   // xType=1:  vs. event counter
544   //
545
546   TString title("Event occupancy");
547   TString xTitle("Time");
548   TString yTitle("number of samples");
549   TGraph *gr=new TGraph(GetNevents());
550   if (fVSignalSumEvent.GetNrows()==0&&!(type==10||type==14)) return 0;
551   TVectorF *vOcc=&fVOccupancyEvent;
552   TVectorF *vSum=&fVSignalSumEvent;
553   TVectorF *vPads=&fVNfiredPadsSenEvent;
554   Double_t norm=557568.;
555   if (type!=14&&fVOccupancyEvent.GetNrows()==0){
556     AliWarning("In non debug mode only occupancy in sensitive regions vs. event awailable!!!");
557     return 0;
558   }
559   if (type>=10){
560     vOcc=&fVOccupancySenEvent;
561     vSum=&fVSignalSumSenEvent;
562     vPads=&fVNfiredPadsSenEvent;
563     norm=33012.;
564   }
565   for (Int_t i=0;i<GetNevents(); ++i){
566     Double_t nAboveThreshold=vOcc->GetMatrixArray()[i];
567     
568     Double_t nSumADC        =1;
569     Double_t timestamp      =1;
570     Double_t nPads          =1;
571
572     if (fVOccupancyEvent.GetNrows()>0){
573       nSumADC        =vSum->GetMatrixArray()[i];
574       timestamp      =fVTimeStampEvent.GetMatrixArray()[i]+fFirstTimeStamp;
575       nPads          =vPads->GetMatrixArray()[i];
576     }
577     Double_t x=timestamp;
578     Double_t y=0;
579     //
580     if (xType==1)     x=i;
581     //
582     if (type%10==0)            y=nAboveThreshold;
583     if (type%10==1&&nAboveThreshold>0)      y=nSumADC/nAboveThreshold;
584     if (type%10==2)            y=nSumADC;
585     if (type%10==3&&nSumADC>0) y=nAboveThreshold/nSumADC;
586     if (type%10==4)            y=nAboveThreshold/(norm*(fLastTimeBin-fFirstTimeBin));
587     if (type==5)               y=fVOccupancySenEvent.GetMatrixArray()[i]/fVOccupancyEvent.GetMatrixArray()[i];
588     if (type==16&&nPads>0)     y=nAboveThreshold/nPads;
589     if (type==17)              y=nPads/norm;
590     //
591     gr->SetPoint(i,x,y);
592   }
593   if (xType==1) xTitle="Event";
594   if (type%10==1)      yTitle="Mean ADC counts/sample";
595   else if (type%10==2) yTitle="Data volume [ADC counts]";
596   else if (type%10==3) yTitle="samples per ADC count";
597   else if (type%10==4) yTitle="sample occupancy";
598   if (type==5)         yTitle="N samples (sensitive) / N samples";
599   if (type%10==6)      yTitle="N samples / N pads";
600   if (type==17)        yTitle="Pad Occupancy";
601   if (type>=10)        yTitle+=" (sensitive)";
602   title=yTitle+":"+xTitle;
603   title+=";"+xTitle+";"+yTitle;
604   gr->SetTitle(title.Data());
605   gr->SetEditable(kFALSE);
606   return gr;
607 }
608 //_____________________________________________________________________
609 // TGraph* AliTPCCalibRaw::MakeGraphNoiseEvents()
610 // {
611   //
612   // Not implemented for the moment
613   //
614 //   return 0;  
615 // }
616 //_____________________________________________________________________
617 TCanvas* AliTPCCalibRaw::MakeCanvasOccupancy(const Int_t xType, Bool_t sen)
618 {
619   //
620   // Create a canvas with occupancy information of all 'type's (see MakeGraphOccupancy)
621   // xType=0: vs. timestamp
622   // xType=1: vs. event number
623   //
624   // sen=kTRUE: for sensitive regions
625   //
626
627   TString name("RawOccupancy_");
628   TString title("Raw Occupancy vs. ");
629   if (xType==0){
630     name+="Time";
631     title+="time";
632   } else if (xType==1){
633     name+="Event";
634     title+="event";
635   }
636   if (sen){
637     name+="Sen";
638     title+=" (sensitive)";
639   }
640   TCanvas *c=(TCanvas*)gROOT->GetListOfCanvases()->FindObject(name.Data());
641   if (!c) c=new TCanvas(name.Data(),title.Data());
642   c->Clear();
643   c->Divide(2,2);
644   for (Int_t i=0;i<4;++i){
645     c->cd(i+1);
646     TGraph *gr=MakeGraphOccupancy(i+10*(Int_t)sen,xType);
647     gr->Draw("alp");
648   }
649   return c;
650 }
651
652 //_____________________________________________________________________
653 void AliTPCCalibRaw::Merge(AliTPCCalibRaw * const sig)
654 {
655   //
656   // Merge sig with this instance
657   //
658
659   if (!sig) return;
660   MergeBase(sig);
661   //Add last time bin distribution histogram
662   fHnDrift->Add(sig->fHnDrift);
663
664   //Add occupancy data
665   
666 }
667
668 //_____________________________________________________________________
669 Long64_t AliTPCCalibRaw::Merge(TCollection * const list)
670 {
671   //
672   // Merge all objects of this type in list
673   //
674   
675   Long64_t nmerged=1;
676   
677   TIter next(list);
678   AliTPCCalibRaw *ce=0;
679   TObject *o=0;
680   
681   while ( (o=next()) ){
682     ce=dynamic_cast<AliTPCCalibRaw*>(o);
683     if (ce){
684       Merge(ce);
685       ++nmerged;
686     }
687   }
688   
689   return nmerged;
690 }
691