]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibRaw.cxx
QA for the on-the-fly V0 (ana marin)
[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 <TH2C.h>
65
66 //AliRoot includes
67 #include "AliTPCCalROC.h"
68 #include "AliAltroRawStream.h"
69 #include "AliLog.h"
70
71 //class header
72 #include "AliTPCCalibRaw.h"
73
74 ClassImp(AliTPCCalibRaw)
75
76 AliTPCCalibRaw::AliTPCCalibRaw() :
77   AliTPCCalibRawBase(),
78   fPeakDetMinus(1),
79   fPeakDetPlus(2),
80   fNFailL1Phase(0),
81   fFirstTimeStamp(0),
82   fNSecTime(600), //default 10 minutes
83   fNBinsTime(60), //default 60*10 minutes = 10 hours
84   fPadProcessed(kFALSE),
85   fCurrentChannel(-1),
86   fCurrentSector(-1),
87   fLastSector(-2),
88   fCurrentRow(-1),
89   fCurrentPad(-1),
90   fLastTimeBinProc(0),
91   fPeakTimeBin(0),
92   fLastSignal(0),
93   fNOkPlus(0),
94   fNOkMinus(0),
95   fArrCurrentPhaseDist(4),
96   fArrALTROL1Phase(1000),
97   fArrALTROL1PhaseEvent(216),
98   fArrALTROL1PhaseFailEvent(216),
99   fHnDrift(0x0)
100 {
101   //
102   // Default ctor
103   //
104   SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
105   CreateDVhist();
106   fFirstTimeBin=850;
107   fLastTimeBin=1000;
108 }
109 //_____________________________________________________________________
110 AliTPCCalibRaw::~AliTPCCalibRaw()
111 {
112   //
113   // dtor
114   //
115   delete fHnDrift;  
116 }
117 //_____________________________________________________________________
118 // AliTPCCalibRaw& AliTPCCalibRaw::operator = (const  AliTPCCalibRaw &source)
119 // {
120 //   //
121 //   // assignment operator
122 //   //
123 //   if (&source == this) return *this;
124 //   new (this) AliTPCCalibRaw(source);
125 //   
126 //   return *this;
127 // }
128
129 //_____________________________________________________________________
130 Int_t AliTPCCalibRaw::Update(const Int_t isector, const Int_t iRow, const Int_t iPad,
131              const Int_t iTimeBin, const Float_t signal)
132 {
133   //
134   // Data filling method
135   //
136   if (iRow<0) return 0;
137   if (iPad<0) return 0;
138   if (iTimeBin<0) return 0;
139   if (!fFirstTimeStamp) fFirstTimeStamp=GetTimeStamp();
140   if (fCurrDDLNum!=fPrevDDLNum){
141     TVectorF *arr=MakeArrL1PhaseRCU(fCurrDDLNum,kTRUE);
142     if (arr->GetNrows()<=fNevents) arr->ResizeTo(arr->GetNrows()+1000);
143     // phase as a position of a quarter time bin
144     Int_t phase=(Int_t)(GetL1PhaseTB()*4.);
145 //     printf("DDL: %03d, phase: %d (%f))\n",fCurrDDLNum,phase,GetL1PhaseTB());
146     //Fill pahse information of current rcu and event
147     (arr->GetMatrixArray())[fNevents]=phase;
148     //increase phase counter 
149     ++((fArrCurrentPhaseDist.GetMatrixArray())[phase]);
150 //     printf("RCUId: %03d (%03d), DDL: %03d, sector: %02d\n",fCurrRCUId, fPrevRCUId, fCurrDDLNum, isector);
151   }
152
153   if ( (iTimeBin>fLastTimeBin) || (iTimeBin<fFirstTimeBin)   ) return 0;
154   //don't process edge pads
155   if (IsEdgePad(isector,iRow,iPad)) return 0;
156 //   Double_t x[kHnBinsDV]={1,isector,0};
157 //   fHnDrift->Fill(x);
158   Int_t iChannel  = fROC->GetRowIndexes(isector)[iRow]+iPad; //  global pad position in sector
159   if (fCurrentChannel==iChannel){
160     if (fPadProcessed) return 0;
161   } else {
162     fPadProcessed=kFALSE;
163     fNOkPlus=0;
164     fNOkMinus=0;
165     fPeakTimeBin=0;
166     fLastSignal=0;
167   }
168 //   Double_t x2[kHnBinsDV]={2,isector,0};
169 //   fHnDrift->Fill(x2);
170   
171
172   if (signal>fLastSignal) ++fNOkPlus;
173   else if(signal<fLastSignal && fNOkPlus>=fPeakDetPlus){
174     ++fNOkMinus;
175     if (!fPeakTimeBin) fPeakTimeBin=fLastTimeBinProc;
176     if ( fNOkMinus>=fPeakDetMinus ) {
177       Double_t x[kHnBinsDV]={fPeakTimeBin,isector,(fTimeStamp-fFirstTimeStamp)/fNSecTime};
178       fHnDrift->Fill(x);
179     } 
180   } else {
181     fNOkPlus=0;
182     fNOkMinus=0;
183     fPeakTimeBin=0;
184     fLastSignal=0;
185   }
186
187   fLastTimeBinProc=iTimeBin;
188   fLastSignal=TMath::Nint(signal);
189   fCurrentChannel = iChannel;
190   return 0;
191 }
192 //_____________________________________________________________________
193 void AliTPCCalibRaw::ResetEvent()
194 {
195   //
196   // Reset event counters
197   //
198
199   fCurrentChannel=-1;
200   fArrCurrentPhaseDist.Zero();
201 }
202 //_____________________________________________________________________
203 void AliTPCCalibRaw::EndEvent()
204 {
205   //
206   // End event analysis
207   //
208
209   
210   //find phase of the current event
211   Int_t phaseMaxEntries=-1;
212   Int_t maxEntries=0;
213   for (Int_t i=0;i<fArrCurrentPhaseDist.GetNrows();++i){
214     Int_t entries=(Int_t)fArrCurrentPhaseDist[i];
215     if (maxEntries<entries) {
216       maxEntries=entries;
217       phaseMaxEntries=i;
218     }
219   }
220   // store phase of current event
221   if (fArrALTROL1Phase.GetNrows()<=GetNevents())
222     fArrALTROL1Phase.ResizeTo(GetNevents()+1000);
223   (fArrALTROL1Phase.GetMatrixArray())[GetNevents()]=phaseMaxEntries;
224   
225   //loop over RCUs and test failures
226   for (Int_t ircu=0;ircu<216;++ircu){
227     const TVectorF *arr=GetALTROL1PhaseEventsRCU(ircu);//MakeArrL1PhaseRCU(ircu);
228     if (!arr) continue;
229     TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu,kTRUE);
230     if (arrF->GetNrows()<=fNevents) arrF->ResizeTo(arrF->GetNrows()+1000);
231     if ((arr->GetMatrixArray())[fNevents]!=phaseMaxEntries){
232       (arrF->GetMatrixArray())[fNevents]=1;
233       ++fNFailL1Phase;
234       }
235   }
236   IncrementNevents();
237 }
238 //_____________________________________________________________________
239 TH2C *AliTPCCalibRaw::MakeHistL1RCUEvents(Int_t type)
240 {
241   // Create a 2D histo RCU:Events indicating the there was a deviation
242   // from the mean L1 phase of the event
243   //
244   //type: 0=Failures, 1=Phases
245   TH2C *h2 = new TH2C("hL1FailRCUEvents","L1 Failures;RCU;Event",216,0,216,GetNevents(),0,GetNevents());
246   for (Int_t ircu=0;ircu<216;++ircu) {
247     const TVectorF *v=0;
248     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
249     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
250     if (!v) continue;
251     for (Int_t iev=0;iev<GetNevents();++iev) {
252       h2->SetBinContent(ircu+1,iev+1,(*v)(iev)+1);
253     }
254   }
255   return h2;
256 }
257 //_____________________________________________________________________
258 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsIROC(Int_t type)
259 {
260   //
261   // Create a 2D histo RCU:Events indicating the there was a deviation
262   // from the mean L1 phase of the event
263   //
264   TH2C *h2 = new TH2C("hL1FailRCUEventsIROC","L1 Failures IROCs;RCU;Event",72,0,36,GetNevents(),0,GetNevents());
265   for (Int_t ircu=0;ircu<72;++ircu) {
266     const TVectorF *v=0;
267     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
268     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
269     if (!v) continue;
270     for (Int_t iev=0;iev<GetNevents();++iev) {
271       h2->SetBinContent(ircu+1,iev+1,(*v)(iev)+1);
272     }
273   }
274   return h2;
275 }
276 //_____________________________________________________________________
277 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsOROC(Int_t type)
278 {
279   //
280   // Create a 2D histo RCU:Events indicating the there was a deviation
281   // from the mean L1 phase of the event
282   //
283   TH2C *h2 = new TH2C("hL1FailRCUEventsOROC","L1 Failures OROCs;RCU;Event",144,0,36,GetNevents(),0,GetNevents());
284   for (Int_t ircu=72;ircu<216;++ircu) {
285     const TVectorF *v=0;
286     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
287     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
288     if (!v) continue;
289     for (Int_t iev=0;iev<GetNevents();++iev) {
290       h2->SetBinContent(ircu-72+1,iev+1,(*v)(iev)+1);
291     }
292   }
293   return h2;
294 }
295 //_____________________________________________________________________
296 void AliTPCCalibRaw::CreateDVhist()
297 {
298   //
299   // Setup the HnSparse for the drift velocity determination
300   //
301   if (fHnDrift) return;
302   //HnSparse bins
303   //time bin, roc, time
304   Int_t    bins[kHnBinsDV] = {1000, 72, fNBinsTime};
305   Double_t xmin[kHnBinsDV] = {0,0,0};
306   Double_t xmax[kHnBinsDV] = {1000,72,fNBinsTime};
307   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);
308     
309 }
310 //_____________________________________________________________________
311 void AliTPCCalibRaw::Analyse()
312 {
313   //
314   // Analyse Data
315   //
316
317   //resize arrays
318   fArrALTROL1Phase.ResizeTo(GetNevents());
319   for (Int_t ircu=0;ircu<216;++ircu){
320     TVectorF *arr=MakeArrL1PhaseRCU(ircu);//MakeArrL1PhaseRCU(ircu);
321     if (!arr) continue;
322     TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu);
323     arr->ResizeTo(GetNevents());
324     arrF->ResizeTo(GetNevents());
325   }
326
327   //Analyse drift velocity
328   
329 }
330