1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////////////////
20 // Implementation of the TPC Raw drift velocity and Altro L1 Phase calibration //
22 // Origin: Jens Wiechula, J.Wiechula@gsi.de //
24 ////////////////////////////////////////////////////////////////////////////////////////
27 // *************************************************************************************
28 // * Class Description *
29 // *************************************************************************************
33 TFile f("CalibAltro.root");
34 AliTPCCalibRaw *al=(AliTPCCalibRaw*)f.Get(f.GetListOfKeys()->At(0)->GetName())
36 TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
37 if (!c1) c1=new TCanvas("c1","c1");
40 TH2F h2f("h2","h2;RCU;fail",216,0,216,al->GetNevents(),0,al->GetNevents());
44 TVectorF *v=al->GetALTROL1PhaseFailEventsRCU(i);
46 for (iev=0;iev<al->GetNevents();++iev) {
47 h2f->SetBinContent(i+1,iev+1,(*v)(iev));
50 // h.SetLineColor(i/216.*50+50);
51 // ((TH1F*)h.Clone(Form("h%d",i)))->Draw(first?"":"same");
67 #include <TObjString.h>
70 #include "AliTPCCalROC.h"
71 #include "AliAltroRawStream.h"
74 #include "AliTPCCalibRaw.h"
76 ClassImp(AliTPCCalibRaw)
78 AliTPCCalibRaw::AliTPCCalibRaw() :
83 fNFailL1PhaseEvent(0),
85 fNSecTime(600), //default 10 minutes
86 fNBinsTime(60), //default 60*10 minutes = 10 hours
87 fPadProcessed(kFALSE),
98 fArrCurrentPhaseDist(4),
99 fArrCurrentPhase(kNRCU),
100 fArrFailEventNumber(100),
101 fArrALTROL1Phase(1000),
102 fArrALTROL1PhaseEvent(kNRCU),
103 fArrALTROL1PhaseFailEvent(kNRCU),
109 SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
111 for (Int_t ircu=0;ircu<kNRCU;++ircu) fArrCurrentPhase.GetMatrixArray()[ircu]=-1;
115 //_____________________________________________________________________
116 AliTPCCalibRaw::AliTPCCalibRaw(const TMap *config) :
117 AliTPCCalibRawBase(),
121 fNFailL1PhaseEvent(0),
123 fNSecTime(600), //default 10 minutes
124 fNBinsTime(60), //default 60*10 minutes = 10 hours
125 fPadProcessed(kFALSE),
136 fArrCurrentPhaseDist(4),
137 fArrCurrentPhase(kNRCU),
138 fArrFailEventNumber(100),
139 fArrALTROL1Phase(1000),
140 fArrALTROL1PhaseEvent(kNRCU),
141 fArrALTROL1PhaseFailEvent(kNRCU),
147 SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
149 for (Int_t ircu=0;ircu<kNRCU;++ircu) fArrCurrentPhase.GetMatrixArray()[ircu]=-1;
152 if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
153 if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
156 //_____________________________________________________________________
157 AliTPCCalibRaw::~AliTPCCalibRaw()
164 //_____________________________________________________________________
165 // AliTPCCalibRaw& AliTPCCalibRaw::operator = (const AliTPCCalibRaw &source)
168 // // assignment operator
170 // if (&source == this) return *this;
171 // new (this) AliTPCCalibRaw(source);
176 //_____________________________________________________________________
177 Int_t AliTPCCalibRaw::Update(const Int_t isector, const Int_t iRow, const Int_t iPad,
178 const Int_t iTimeBin, const Float_t signal)
181 // Data filling method
183 if (iRow<0) return 0;
184 if (iPad<0) return 0;
185 if (iTimeBin<0) return 0;
186 if (!fFirstTimeStamp) fFirstTimeStamp=GetTimeStamp();
187 if ( (iTimeBin>fLastTimeBin) || (iTimeBin<fFirstTimeBin) ) return 0;
188 //don't process edge pads
189 if (IsEdgePad(isector,iRow,iPad)) return 0;
190 // Double_t x[kHnBinsDV]={1,isector,0};
191 // fHnDrift->Fill(x);
192 Int_t iChannel = fROC->GetRowIndexes(isector)[iRow]+iPad; // global pad position in sector
193 if (fCurrentChannel==iChannel){
194 if (fPadProcessed) return 0;
196 fPadProcessed=kFALSE;
202 // Double_t x2[kHnBinsDV]={2,isector,0};
203 // fHnDrift->Fill(x2);
206 if (signal>fLastSignal) ++fNOkPlus;
207 else if(signal<fLastSignal && fNOkPlus>=fPeakDetPlus){
209 if (!fPeakTimeBin) fPeakTimeBin=fLastTimeBinProc;
210 if ( fNOkMinus>=fPeakDetMinus ) {
211 Double_t x[kHnBinsDV]={fPeakTimeBin,isector,(fTimeStamp-fFirstTimeStamp)/fNSecTime};
221 fLastTimeBinProc=iTimeBin;
222 fLastSignal=TMath::Nint(signal);
223 fCurrentChannel = iChannel;
226 //_____________________________________________________________________
227 void AliTPCCalibRaw::UpdateDDL(){
229 // fill ALTRO L1 information
233 Int_t phase=(Int_t)(GetL1PhaseTB()*4.);
234 //Fill pahse information of current rcu and event
235 fArrCurrentPhase.GetMatrixArray()[fCurrDDLNum]=phase;
236 //increase phase counter
237 ++((fArrCurrentPhaseDist.GetMatrixArray())[phase]);
240 //_____________________________________________________________________
241 void AliTPCCalibRaw::ResetEvent()
244 // Reset event counters
248 fArrCurrentPhaseDist.Zero();
250 //_____________________________________________________________________
251 void AliTPCCalibRaw::EndEvent()
254 // End event analysis
258 //find phase of the current event
259 Int_t phaseMaxEntries=-1;
261 for (Int_t i=0;i<fArrCurrentPhaseDist.GetNrows();++i){
262 Int_t entries=(Int_t)fArrCurrentPhaseDist[i];
263 if (maxEntries<entries) {
268 // store phase of current event
269 if (fArrALTROL1Phase.GetNrows()<=GetNevents())
270 fArrALTROL1Phase.ResizeTo(GetNevents()+1000);
271 (fArrALTROL1Phase.GetMatrixArray())[GetNevents()]=phaseMaxEntries;
273 //loop over RCUs and test failures
275 for (Int_t ircu=0;ircu<kNRCU;++ircu){
276 Int_t phase=(Int_t)fArrCurrentPhase[ircu];
277 if (phase<0) continue;
278 if (phase!=phaseMaxEntries){
279 TVectorF *arr=MakeArrL1PhaseRCU(fCurrDDLNum,kTRUE);
280 if (arr->GetNrows()<=(Int_t)fNFailL1PhaseEvent) arr->ResizeTo(arr->GetNrows()+100);
281 (arr->GetMatrixArray())[fNFailL1PhaseEvent]=phase;
285 //reset current phase information
286 fArrCurrentPhase[ircu]=-1;
289 if (fArrFailEventNumber.GetNrows()<=(Int_t)fNFailL1PhaseEvent) fArrFailEventNumber.ResizeTo(fArrFailEventNumber.GetNrows()+100);
290 fArrFailEventNumber.GetMatrixArray()[fNFailL1PhaseEvent]=GetNevents();
292 fNFailL1PhaseEvent+=fail;
295 //_____________________________________________________________________
296 TH2C *AliTPCCalibRaw::MakeHistL1RCUEvents(Int_t type)
298 // Create a 2D histo RCU:Events indicating the there was a deviation
299 // from the mean L1 phase of the event
301 //type: 0=Failures, 1=Phases
303 //number of relavant events, depending on version
304 Int_t nevents=GetNevents();
306 Bool_t newVersion=kFALSE;
307 for (Int_t ircu=0; ircu<kNRCU; ++ircu){
308 const TVectorF *v=GetALTROL1PhaseEventsRCU(ircu);
310 if ((UInt_t)(v->GetNrows())==fNFailL1PhaseEvent){
312 nevents=fNFailL1PhaseEvent;
316 TH2C *h2 = new TH2C("hL1FailRCUEvents","L1 Failures;RCU;Event",kNRCU,0,kNRCU,nevents,0,nevents);
318 for (Int_t ircu=0;ircu<kNRCU;++ircu) {
319 const TVectorF *v=GetALTROL1PhaseEventsRCU(ircu);
324 } else if (type==1) {
330 for (Int_t iev=0;iev<nevents;++iev) {
331 Float_t val=(*v)(iev);
332 Float_t phase=fArrALTROL1Phase.GetMatrixArray()[iev];
334 Int_t event=(Int_t)fArrFailEventNumber.GetMatrixArray()[iev];
335 phase=fArrALTROL1Phase.GetMatrixArray()[event];
337 if (type==0) val=(val!=phase);
338 h2->SetBinContent(ircu+1,iev+1,val+add);
343 //_____________________________________________________________________
344 TH1F *AliTPCCalibRaw::MakeHistL1PhaseDist()
347 // L1 phase distribution. Should be flat in ideal case
349 TH1F *h=new TH1F("L1phaseDist","Normalized L1 phase distribution;phase;fraction of events",4,0,4);
351 for (Int_t iev=0;iev<GetNevents();++iev) h->Fill(fArrALTROL1Phase.GetMatrixArray()[iev]);
352 if (GetNevents()>0) h->Scale(1./GetNevents());
357 //_____________________________________________________________________
358 TVectorF *AliTPCCalibRaw::MakeVectL1PhaseDist()
361 // L1 phase distribution. Should be flat in ideal case
363 TVectorF *v=new TVectorF(4);
364 for (Int_t iev=0;iev<GetNevents();++iev) {
365 Int_t phase=(Int_t)fArrALTROL1Phase.GetMatrixArray()[iev];
366 ((v->GetMatrixArray())[phase])+=1./GetNevents();
370 //_____________________________________________________________________
371 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsIROC(Int_t type)
374 // Create a 2D histo RCU:Events indicating the there was a deviation
375 // from the mean L1 phase of the event
377 TH2C *h2 = new TH2C("hL1FailRCUEventsIROC","L1 Failures IROCs;RCU;Event",72,0,36,GetNevents(),0,GetNevents());
378 for (Int_t ircu=0;ircu<72;++ircu) {
380 if (type==0) v=GetALTROL1PhaseFailEventsRCU(ircu);
381 else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
383 for (Int_t iev=0;iev<GetNevents();++iev) {
384 h2->SetBinContent(ircu+1,iev+1,(*v)(iev));
389 //_____________________________________________________________________
390 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsOROC(Int_t type)
393 // Create a 2D histo RCU:Events indicating the there was a deviation
394 // from the mean L1 phase of the event
396 TH2C *h2 = new TH2C("hL1FailRCUEventsOROC","L1 Failures OROCs;RCU;Event",144,0,36,GetNevents(),0,GetNevents());
397 for (Int_t ircu=72;ircu<kNRCU;++ircu) {
399 if (type==0) v=GetALTROL1PhaseFailEventsRCU(ircu);
400 else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
402 for (Int_t iev=0;iev<GetNevents();++iev) {
403 h2->SetBinContent(ircu-72+1,iev+1,(*v)(iev));
408 //_____________________________________________________________________
409 void AliTPCCalibRaw::CreateDVhist()
412 // Setup the HnSparse for the drift velocity determination
414 if (fHnDrift) return;
416 //time bin, roc, time
417 Int_t bins[kHnBinsDV] = {fLastTimeBin-fFirstTimeBin, 72, fNBinsTime};
418 Double_t xmin[kHnBinsDV] = {fFirstTimeBin,0,0};
419 Double_t xmax[kHnBinsDV] = {fLastTimeBin,72,fNBinsTime};
420 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);
423 //_____________________________________________________________________
424 void AliTPCCalibRaw::Analyse()
431 fArrALTROL1Phase.ResizeTo(GetNevents());
432 for (Int_t ircu=0;ircu<kNRCU;++ircu){
433 TVectorF *arr=MakeArrL1PhaseRCU(ircu);//MakeArrL1PhaseRCU(ircu);
435 arr->ResizeTo(fNFailL1PhaseEvent);
436 fArrFailEventNumber.ResizeTo(fNFailL1PhaseEvent);
437 // TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu);
438 // arrF->ResizeTo(1);
441 //Analyse drift velocity