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");
66 #include <TObjString.h>
69 #include "AliTPCCalROC.h"
70 #include "AliAltroRawStream.h"
73 #include "AliTPCCalibRaw.h"
75 ClassImp(AliTPCCalibRaw)
77 AliTPCCalibRaw::AliTPCCalibRaw() :
83 fNSecTime(600), //default 10 minutes
84 fNBinsTime(60), //default 60*10 minutes = 10 hours
85 fPadProcessed(kFALSE),
96 fArrCurrentPhaseDist(4),
97 fArrALTROL1Phase(1000),
98 fArrALTROL1PhaseEvent(216),
99 fArrALTROL1PhaseFailEvent(216),
105 SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
110 //_____________________________________________________________________
111 AliTPCCalibRaw::AliTPCCalibRaw(const TMap *config) :
112 AliTPCCalibRawBase(),
117 fNSecTime(600), //default 10 minutes
118 fNBinsTime(60), //default 60*10 minutes = 10 hours
119 fPadProcessed(kFALSE),
130 fArrCurrentPhaseDist(4),
131 fArrALTROL1Phase(1000),
132 fArrALTROL1PhaseEvent(216),
133 fArrALTROL1PhaseFailEvent(216),
139 SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
143 if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
144 if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
147 //_____________________________________________________________________
148 AliTPCCalibRaw::~AliTPCCalibRaw()
155 //_____________________________________________________________________
156 // AliTPCCalibRaw& AliTPCCalibRaw::operator = (const AliTPCCalibRaw &source)
159 // // assignment operator
161 // if (&source == this) return *this;
162 // new (this) AliTPCCalibRaw(source);
167 //_____________________________________________________________________
168 Int_t AliTPCCalibRaw::Update(const Int_t isector, const Int_t iRow, const Int_t iPad,
169 const Int_t iTimeBin, const Float_t signal)
172 // Data filling method
174 if (iRow<0) return 0;
175 if (iPad<0) return 0;
176 if (iTimeBin<0) return 0;
177 if (!fFirstTimeStamp) fFirstTimeStamp=GetTimeStamp();
178 if ( (iTimeBin>fLastTimeBin) || (iTimeBin<fFirstTimeBin) ) return 0;
179 //don't process edge pads
180 if (IsEdgePad(isector,iRow,iPad)) return 0;
181 // Double_t x[kHnBinsDV]={1,isector,0};
182 // fHnDrift->Fill(x);
183 Int_t iChannel = fROC->GetRowIndexes(isector)[iRow]+iPad; // global pad position in sector
184 if (fCurrentChannel==iChannel){
185 if (fPadProcessed) return 0;
187 fPadProcessed=kFALSE;
193 // Double_t x2[kHnBinsDV]={2,isector,0};
194 // fHnDrift->Fill(x2);
197 if (signal>fLastSignal) ++fNOkPlus;
198 else if(signal<fLastSignal && fNOkPlus>=fPeakDetPlus){
200 if (!fPeakTimeBin) fPeakTimeBin=fLastTimeBinProc;
201 if ( fNOkMinus>=fPeakDetMinus ) {
202 Double_t x[kHnBinsDV]={fPeakTimeBin,isector,(fTimeStamp-fFirstTimeStamp)/fNSecTime};
212 fLastTimeBinProc=iTimeBin;
213 fLastSignal=TMath::Nint(signal);
214 fCurrentChannel = iChannel;
217 //_____________________________________________________________________
218 void AliTPCCalibRaw::UpdateDDL(){
220 // fill ALTRO L1 information
222 // if (fCurrDDLNum!=fPrevDDLNum){
223 TVectorF *arr=MakeArrL1PhaseRCU(fCurrDDLNum,kTRUE);
224 if (arr->GetNrows()<=fNevents) arr->ResizeTo(arr->GetNrows()+1000);
225 // phase as a position of a quarter time bin
226 Int_t phase=(Int_t)(GetL1PhaseTB()*4.);
227 // printf("DDL: %03d, phase: %d (%f))\n",fCurrDDLNum,phase,GetL1PhaseTB());
228 //Fill pahse information of current rcu and event
229 (arr->GetMatrixArray())[fNevents]=phase;
230 //increase phase counter
231 ++((fArrCurrentPhaseDist.GetMatrixArray())[phase]);
232 // printf("RCUId: %03d (%03d), DDL: %03d, sector: %02d\n",fCurrRCUId, fPrevRCUId, fCurrDDLNum, isector);
236 //_____________________________________________________________________
237 void AliTPCCalibRaw::ResetEvent()
240 // Reset event counters
244 fArrCurrentPhaseDist.Zero();
246 //_____________________________________________________________________
247 void AliTPCCalibRaw::EndEvent()
250 // End event analysis
254 //find phase of the current event
255 Int_t phaseMaxEntries=-1;
257 for (Int_t i=0;i<fArrCurrentPhaseDist.GetNrows();++i){
258 Int_t entries=(Int_t)fArrCurrentPhaseDist[i];
259 if (maxEntries<entries) {
264 // store phase of current event
265 if (fArrALTROL1Phase.GetNrows()<=GetNevents())
266 fArrALTROL1Phase.ResizeTo(GetNevents()+1000);
267 (fArrALTROL1Phase.GetMatrixArray())[GetNevents()]=phaseMaxEntries;
269 //loop over RCUs and test failures
270 for (Int_t ircu=0;ircu<216;++ircu){
271 const TVectorF *arr=GetALTROL1PhaseEventsRCU(ircu);//MakeArrL1PhaseRCU(ircu);
273 TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu,kTRUE);
274 if (arrF->GetNrows()<=fNevents) arrF->ResizeTo(arrF->GetNrows()+1000);
275 if ((arr->GetMatrixArray())[fNevents]!=phaseMaxEntries){
276 (arrF->GetMatrixArray())[fNevents]=1;
282 //_____________________________________________________________________
283 TH2C *AliTPCCalibRaw::MakeHistL1RCUEvents(Int_t type)
285 // Create a 2D histo RCU:Events indicating the there was a deviation
286 // from the mean L1 phase of the event
288 //type: 0=Failures, 1=Phases
289 TH2C *h2 = new TH2C("hL1FailRCUEvents","L1 Failures;RCU;Event",216,0,216,GetNevents(),0,GetNevents());
291 for (Int_t ircu=0;ircu<216;++ircu) {
294 v=GetALTROL1PhaseFailEventsRCU(ircu);
298 } else if (type==1) {
299 v=GetALTROL1PhaseEventsRCU(ircu);
305 for (Int_t iev=0;iev<GetNevents();++iev) {
306 h2->SetBinContent(ircu+1,iev+1,(*v)(iev)+add);
311 //_____________________________________________________________________
312 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsIROC(Int_t type)
315 // Create a 2D histo RCU:Events indicating the there was a deviation
316 // from the mean L1 phase of the event
318 TH2C *h2 = new TH2C("hL1FailRCUEventsIROC","L1 Failures IROCs;RCU;Event",72,0,36,GetNevents(),0,GetNevents());
319 for (Int_t ircu=0;ircu<72;++ircu) {
321 if (type==0) v=GetALTROL1PhaseFailEventsRCU(ircu);
322 else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
324 for (Int_t iev=0;iev<GetNevents();++iev) {
325 h2->SetBinContent(ircu+1,iev+1,(*v)(iev));
330 //_____________________________________________________________________
331 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsOROC(Int_t type)
334 // Create a 2D histo RCU:Events indicating the there was a deviation
335 // from the mean L1 phase of the event
337 TH2C *h2 = new TH2C("hL1FailRCUEventsOROC","L1 Failures OROCs;RCU;Event",144,0,36,GetNevents(),0,GetNevents());
338 for (Int_t ircu=72;ircu<216;++ircu) {
340 if (type==0) v=GetALTROL1PhaseFailEventsRCU(ircu);
341 else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
343 for (Int_t iev=0;iev<GetNevents();++iev) {
344 h2->SetBinContent(ircu-72+1,iev+1,(*v)(iev));
349 //_____________________________________________________________________
350 void AliTPCCalibRaw::CreateDVhist()
353 // Setup the HnSparse for the drift velocity determination
355 if (fHnDrift) return;
357 //time bin, roc, time
358 Int_t bins[kHnBinsDV] = {fLastTimeBin-fFirstTimeBin, 72, fNBinsTime};
359 Double_t xmin[kHnBinsDV] = {fFirstTimeBin,0,0};
360 Double_t xmax[kHnBinsDV] = {fLastTimeBin,72,fNBinsTime};
361 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);
364 //_____________________________________________________________________
365 void AliTPCCalibRaw::Analyse()
372 fArrALTROL1Phase.ResizeTo(GetNevents());
373 for (Int_t ircu=0;ircu<216;++ircu){
374 TVectorF *arr=MakeArrL1PhaseRCU(ircu);//MakeArrL1PhaseRCU(ircu);
376 TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu);
377 arr->ResizeTo(GetNevents());
378 arrF->ResizeTo(GetNevents());
381 //Analyse drift velocity