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