]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFClusterFinder.cxx
Removing semaphore .done files.
[u/mrichter/AliRoot.git] / TOF / AliTOFClusterFinder.cxx
CommitLineData
d08a92dd 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/*
17
18Revision 0.03 2005/07/28 A. De Caro
19 Implement public method
20 Raw2Digits(Int_t, AliRawReader *)
21 to convert digits from raw data in MC digits
22 (temporary solution)
23
24Revision 0.02 2005/07/27 A. De Caro
25 Implement public method
26 Digits2RecPoint(Int_t)
27 to convert digits in clusters
28
29Revision 0.02 2005/07/26 A. De Caro
30 Implement private methods
31 InsertCluster(AliTOFcluster *)
32 FindClusterIndex(Double_t)
33 originally implemented in AliTOFtracker
34 by S. Arcelli and C. Zampolli
35
36Revision 0.01 2005/07/25 A. De Caro
37 Implement public methods
38 Digits2RecPoint(AliRawReader *, TTree *)
39 Digits2RecPoint(Int_t, AliRawReader *)
40 to convert raw data in clusters
41 */
42
43////////////////////////////////////////////////////////////////
44// //
45// Class for TOF cluster finder //
46// //
47// Starting from Raw Data, create rec points, //
48// fill TreeR for TOF, //
49// write TOF.RecPoints.root file //
50// //
51////////////////////////////////////////////////////////////////
52
0e46b9ae 53#include "TClonesArray.h"
54#include "TFile.h"
55#include "TTree.h"
d08a92dd 56
d0eb8f39 57#include "AliDAQ.h"
0e46b9ae 58#include "AliLoader.h"
d08a92dd 59#include "AliLog.h"
0e46b9ae 60#include "AliRawReader.h"
d08a92dd 61#include "AliRunLoader.h"
d08a92dd 62
37879eed 63#include "AliTOFcalib.h"
64#include "AliTOFCal.h"
65#include "AliTOFChannel.h"
d08a92dd 66#include "AliTOFClusterFinder.h"
0e46b9ae 67#include "AliTOFcluster.h"
68#include "AliTOFdigit.h"
69#include "AliTOFGeometryV5.h"
70#include "AliTOFGeometry.h"
71#include "AliTOFRawStream.h"
72
73extern TFile *gFile;
d08a92dd 74
75ClassImp(AliTOFClusterFinder)
76
77AliTOFClusterFinder::AliTOFClusterFinder():
78 fRunLoader(0),
79 fTOFLoader(0),
80 fTreeD(0),
81 fTreeR(0),
58d8d9a3 82 fTOFGeometry(new AliTOFGeometryV5()),
d08a92dd 83 fDigits(new TClonesArray("AliTOFdigit", 4000)),
84 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
85 fNumberOfTofClusters(0)
86{
87//
88// Constructor
89//
90
3bbdf45d 91 AliInfo("V5 TOF Geometry is taken as the default");
37879eed 92
d08a92dd 93}
94//______________________________________________________________________________
95
96AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
97 fRunLoader(runLoader),
98 fTOFLoader(runLoader->GetLoader("TOFLoader")),
99 fTreeD(0),
100 fTreeR(0),
58d8d9a3 101 fTOFGeometry(new AliTOFGeometryV5()),
d08a92dd 102 fDigits(new TClonesArray("AliTOFdigit", 4000)),
103 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
104 fNumberOfTofClusters(0)
105{
106//
107// Constructor
108//
109
d3c7bfac 110 runLoader->CdGAFile();
111 TFile *in=(TFile*)gFile;
112 in->cd();
113 fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
d08a92dd 114
7aeeaf38 115}
116
117//------------------------------------------------------------------------
118AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
58d8d9a3 119 :TObject(),
120 fRunLoader(0),
121 fTOFLoader(0),
122 fTreeD(0),
123 fTreeR(0),
124 fTOFGeometry(new AliTOFGeometryV5()),
125 fDigits(new TClonesArray("AliTOFdigit", 4000)),
126 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
127 fNumberOfTofClusters(0)
7aeeaf38 128{
129 // copy constructor
130 this->fDigits=source.fDigits;
131 this->fRecPoints=source.fRecPoints;
132 this->fTOFGeometry=source.fTOFGeometry;
133
134}
135
136//------------------------------------------------------------------------
d0eb8f39 137AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
7aeeaf38 138{
139 // ass. op.
140 this->fDigits=source.fDigits;
141 this->fRecPoints=source.fRecPoints;
142 this->fTOFGeometry=source.fTOFGeometry;
143 return *this;
144
d08a92dd 145}
146//______________________________________________________________________________
147
148AliTOFClusterFinder::~AliTOFClusterFinder()
149{
150
151 //
152 // Destructor
153 //
154
155 if (fDigits)
156 {
157 fDigits->Delete();
158 delete fDigits;
159 fDigits=0;
160 }
161 if (fRecPoints)
162 {
163 fRecPoints->Delete();
164 delete fRecPoints;
165 fRecPoints=0;
166 }
167
168 delete fTOFGeometry;
169
170}
171//______________________________________________________________________________
172
173void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
174{
175 //
176 // Converts digits to recpoints for TOF
177 //
178
179 fRunLoader->GetEvent(iEvent);
180
181 fTreeD = fTOFLoader->TreeD();
182 if (fTreeD == 0x0)
183 {
184 AliFatal("AliTOFClusterFinder: Can not get TreeD");
185 }
186
187 TBranch *branch = fTreeD->GetBranch("TOF");
188 if (!branch) {
189 AliError("can't get the branch with the TOF digits !");
190 return;
191 }
192
193 TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
194 branch->SetAddress(&digits);
195
196 ResetRecpoint();
197
198 fTreeR = fTOFLoader->TreeR();
199 if (fTreeR == 0x0)
200 {
201 fTOFLoader->MakeTree("R");
202 fTreeR = fTOFLoader->TreeR();
203 }
204
205 Int_t bufsize = 32000;
206 fTreeR->Branch("TOF", &fRecPoints, bufsize);
207
208 fTreeD->GetEvent(0);
209 Int_t nDigits = digits->GetEntriesFast();
210 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
211
212 Int_t ii, jj;
213 Int_t dig[5];
214 Float_t g[3];
215 Double_t h[5];
340693af 216 Float_t tToT;
217 Double_t tTdcND;
d08a92dd 218 for (ii=0; ii<nDigits; ii++) {
219 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
220 dig[0]=d->GetSector();
221 dig[1]=d->GetPlate();
222 dig[2]=d->GetStrip();
223 dig[3]=d->GetPadz();
224 dig[4]=d->GetPadx();
225
d0eb8f39 226 //AliInfo(Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
227
d08a92dd 228 for (jj=0; jj<3; jj++) g[jj] = 0.;
229 fTOFGeometry->GetPos(dig,g);
230
231 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
232 h[1] = TMath::ATan2(g[1],g[0]);
233 h[2] = g[2];
234 h[3] = d->GetTdc();
235 h[4] = d->GetAdc();
340693af 236 tToT = d->GetToT();
237 tTdcND = d->GetTdcND();
d08a92dd 238
340693af 239 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
d08a92dd 240 InsertCluster(tofCluster);
241
242 }
243
244 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
245
37879eed 246 CalibrateRecPoint();
d08a92dd 247 FillRecPoint();
248
249 fTreeR->Fill();
250 ResetRecpoint();
251
252 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
253 fTOFLoader->WriteRecPoints("OVERWRITE");
254
255}
256//______________________________________________________________________________
257
258void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
259 TTree *clustersTree)
260{
261 //
262 // Converts RAW data to recpoints for TOF
263 //
264
d0eb8f39 265 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
266 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
d08a92dd 267
268 ResetRecpoint();
269
270 Int_t bufsize = 32000;
271 clustersTree->Branch("TOF", &fRecPoints, bufsize);
272
273 Int_t ii = 0;
274 Int_t indexDDL = 0;
275
276 Int_t detectorIndex[5];
277 Float_t position[3];
278 Double_t cylindricalPosition[5];
340693af 279 Float_t tToT;
280 Double_t tTdcND;
d08a92dd 281
282 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
283
284 rawReader->Reset();
285 AliTOFRawStream tofInput(rawReader);
362c9d61 286 rawReader->Select("TOF", indexDDL, indexDDL);
d08a92dd 287
288 while(tofInput.Next()) {
289
d0eb8f39 290 for (ii=0; ii<5; ii++) detectorIndex[ii] = -1;
291
d08a92dd 292 detectorIndex[0] = tofInput.GetSector();
293 detectorIndex[1] = tofInput.GetPlate();
294 detectorIndex[2] = tofInput.GetStrip();
295 detectorIndex[3] = tofInput.GetPadZ();
296 detectorIndex[4] = tofInput.GetPadX();
297
d0eb8f39 298 //AliInfo(Form(" %2i %1i %2i %1i %2i ",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
299
300 if (detectorIndex[0]==-1 ||
301 detectorIndex[1]==-1 ||
302 detectorIndex[2]==-1 ||
303 detectorIndex[3]==-1 ||
304 detectorIndex[4]==-1) continue;
305
d08a92dd 306 for (ii=0; ii<3; ii++) position[ii] = 0.;
307
308 fTOFGeometry->GetPos(detectorIndex, position);
309
310 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
311 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
312 cylindricalPosition[2] = position[2];
313 cylindricalPosition[3] = tofInput.GetTofBin();
d0eb8f39 314 cylindricalPosition[4] = tofInput.GetToTbin();
315 tToT = tofInput.GetToTbin();
340693af 316 tTdcND = -1.;
d08a92dd 317 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
340693af 318 tofCluster->SetToT(tToT);
319 tofCluster->SetTDCND(tTdcND);
d08a92dd 320 InsertCluster(tofCluster);
321
322 } // while loop
323
324 } // loop on DDL files
325
326 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
327
37879eed 328 CalibrateRecPoint();
d08a92dd 329 FillRecPoint();
330
331 clustersTree->Fill();
37879eed 332
d08a92dd 333 ResetRecpoint();
334
335}
336//______________________________________________________________________________
337
338void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
339{
340 //
341 // Converts RAW data to recpoints for TOF
342 //
343
d0eb8f39 344 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
345 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
d08a92dd 346
347 fRunLoader->GetEvent(iEvent);
348
349 AliDebug(2,Form(" Event number %2i ", iEvent));
350
351 fTreeR = fTOFLoader->TreeR();
352
353 if (fTreeR == 0x0){
354 fTOFLoader->MakeTree("R");
355 fTreeR = fTOFLoader->TreeR();
356 }
357
358 Int_t bufsize = 32000;
359 fTreeR->Branch("TOF", &fRecPoints, bufsize);
360
361 Int_t ii = 0;
362 Int_t indexDDL = 0;
363
d0eb8f39 364 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
d08a92dd 365 Float_t position[3];
366 Double_t cylindricalPosition[5];
340693af 367 Float_t tToT;
368 Double_t tTdcND;
d08a92dd 369
370 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
371
372 rawReader->Reset();
373 AliTOFRawStream tofInput(rawReader);
362c9d61 374 rawReader->Select("TOF", indexDDL, indexDDL);
d08a92dd 375
376 while(tofInput.Next()) {
377
d0eb8f39 378 for (ii=0; ii<5; ii++) detectorIndex[ii] = -1;
379
d08a92dd 380 detectorIndex[0] = (Int_t)tofInput.GetSector();
381 detectorIndex[1] = (Int_t)tofInput.GetPlate();
382 detectorIndex[2] = (Int_t)tofInput.GetStrip();
383 detectorIndex[3] = (Int_t)tofInput.GetPadZ();
384 detectorIndex[4] = (Int_t)tofInput.GetPadX();
385
d0eb8f39 386 if (detectorIndex[0]==-1 ||
387 detectorIndex[1]==-1 ||
388 detectorIndex[2]==-1 ||
389 detectorIndex[3]==-1 ||
390 detectorIndex[4]==-1) continue;
391
392 //AliInfo(Form(" %2i %1i %2i %1i %2i ",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
393
d08a92dd 394 for (ii=0; ii<3; ii++) position[ii] = 0.;
395
396 fTOFGeometry->GetPos(detectorIndex, position);
397
398 cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
399 cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
400 cylindricalPosition[2] = (Double_t)position[2];
401 cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
d0eb8f39 402 cylindricalPosition[4] = (Double_t)tofInput.GetToTbin();
403 tToT = tofInput.GetToTbin();
340693af 404 tTdcND = -1.;
37879eed 405
d08a92dd 406 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
340693af 407 tofCluster->SetToT(tToT);
408 tofCluster->SetTDCND(tTdcND);
d08a92dd 409 InsertCluster(tofCluster);
410
411 } // while loop
412
413 } // DDL Loop
414
415 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
416
37879eed 417 CalibrateRecPoint();
d08a92dd 418 FillRecPoint();
419
420 fTreeR->Fill();
421 ResetRecpoint();
422
423 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
424 fTOFLoader->WriteRecPoints("OVERWRITE");
425
426}
427//______________________________________________________________________________
428
429void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
430{
431 //
432 // Converts RAW data to MC digits for TOF
433 //
434 // (temporary solution)
435 //
436
437 fRunLoader->GetEvent(iEvent);
438
439 fTreeD = fTOFLoader->TreeD();
440 if (fTreeD)
441 {
d0eb8f39 442 AliInfo("TreeD re-creation");
d08a92dd 443 fTreeD = 0x0;
444 fTOFLoader->MakeTree("D");
445 fTreeD = fTOFLoader->TreeD();
446 }
447
d08a92dd 448 TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
449 Int_t bufsize = 32000;
450 fTreeD->Branch("TOF", &tofDigits, bufsize);
451
452
453 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
454
455 fRunLoader->GetEvent(iEvent);
456
457 AliDebug(2,Form(" Event number %2i ", iEvent));
458
d0eb8f39 459 Int_t ii = 0;
d08a92dd 460 Int_t indexDDL = 0;
461
462 Int_t detectorIndex[5];
d0eb8f39 463 Float_t digit[4];
d08a92dd 464
465 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
466
467 rawReader->Reset();
468 AliTOFRawStream tofInput(rawReader);
362c9d61 469 rawReader->Select("TOF", indexDDL, indexDDL);
d08a92dd 470
471 while(tofInput.Next()) {
472
d0eb8f39 473 for (ii=0; ii<5; ii++) detectorIndex[ii] = -1;
474
d08a92dd 475 detectorIndex[0] = tofInput.GetSector();
476 detectorIndex[1] = tofInput.GetPlate();
477 detectorIndex[2] = tofInput.GetStrip();
478 detectorIndex[3] = tofInput.GetPadX();
479 detectorIndex[4] = tofInput.GetPadZ();
480
d0eb8f39 481 //AliInfo(Form(" %2i %1i %2i %1i %2i ",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
482
483 if (detectorIndex[0]==-1 ||
484 detectorIndex[1]==-1 ||
485 detectorIndex[2]==-1 ||
486 detectorIndex[3]==-1 ||
487 detectorIndex[4]==-1) continue;
488
d08a92dd 489 digit[0] = (Float_t)tofInput.GetTofBin();
d0eb8f39 490 digit[1] = (Float_t)tofInput.GetToTbin();
491 digit[2] = (Float_t)tofInput.GetToTbin();
492 digit[3] = -1.;
d08a92dd 493
494 Int_t tracknum[3]={-1,-1,-1};
495
496 TClonesArray &aDigits = *tofDigits;
497 Int_t last=tofDigits->GetEntriesFast();
498 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
499
500 } // while loop
501
502 } // DDL Loop
503
504 fTreeD->Fill();
505
506 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
507 fTOFLoader->WriteDigits("OVERWRITE");
508
509}
510//______________________________________________________________________________
511
512Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
513 //---------------------------------------------------------------------------//
514 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
515 //---------------------------------------------------------------------------//
516 if (fNumberOfTofClusters==kTofMaxCluster) {
517 AliError("Too many clusters !");
518 return 1;
519 }
520
521 if (fNumberOfTofClusters==0) {
522 fTofClusters[fNumberOfTofClusters++] = tofCluster;
523 return 0;
524 }
525
526 Int_t ii = FindClusterIndex(tofCluster->GetZ());
527 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
528 fTofClusters[ii] = tofCluster;
529 fNumberOfTofClusters++;
530
531 return 0;
532
533}
534//_________________________________________________________________________
535
536Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
537 //--------------------------------------------------------------------
538 // This function returns the index of the nearest cluster
539 //--------------------------------------------------------------------
540 if (fNumberOfTofClusters==0) return 0;
541 if (z <= fTofClusters[0]->GetZ()) return 0;
542 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
543 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
544 for (; b<e; m=(b+e)/2) {
545 if (z > fTofClusters[m]->GetZ()) b=m+1;
546 else e=m;
547 }
548
549 return m;
550
551}
552//_________________________________________________________________________
553
554void AliTOFClusterFinder::FillRecPoint()
555{
556 //
557 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
558 // in Z) in the global TClonesArray of AliTOFcluster,
559 // i.e. fRecPoints.
560 //
561
562 Int_t ii, jj;
563
564 Int_t detectorIndex[5];
565 Double_t cylindricalPosition[5];
566 Int_t trackLabels[3];
567 Int_t digitIndex = -1;
58d8d9a3 568 Float_t tToT=0.;
569 Double_t tTdcND=0.;
570 Bool_t cStatus = kTRUE;
d08a92dd 571
572 TClonesArray &lRecPoints = *fRecPoints;
573
574 for (ii=0; ii<fNumberOfTofClusters; ii++) {
575
576 digitIndex = fTofClusters[ii]->GetIndex();
577 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
578 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
579 cylindricalPosition[0] = fTofClusters[ii]->GetR();
580 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
581 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
582 cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
583 cylindricalPosition[4] = fTofClusters[ii]->GetADC();
340693af 584 tToT = fTofClusters[ii]->GetToT();
585 tTdcND = fTofClusters[ii]->GetTDCND();
58d8d9a3 586 cStatus=fTofClusters[ii]->GetStatus();
d08a92dd 587
58d8d9a3 588 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND, cStatus);
d08a92dd 589
590 //AliInfo(Form("%3i %3i %f %f %f %f %f %2i %2i %2i %1i %2i",ii,digitIndex, cylindricalPosition[2],cylindricalPosition[0],cylindricalPosition[1],cylindricalPosition[3],cylindricalPosition[4],detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
591
592 } // loop on clusters
593
594}
37879eed 595
596//_________________________________________________________________________
597void AliTOFClusterFinder::CalibrateRecPoint()
598{
599 //
600 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
601 // in Z) in the global TClonesArray of AliTOFcluster,
602 // i.e. fRecPoints.
603 //
604
605 Int_t ii, jj;
606
607 Int_t detectorIndex[5];
608 Int_t digitIndex = -1;
340693af 609 Float_t tToT;
610 Float_t tdcCorr;
37879eed 611 AliInfo(" Calibrating TOF Clusters: ")
612 AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
3f2bddfa 613 // calib->ReadParFromCDB("TOF/Calib",0); // original
614 calib->ReadParFromCDB("TOF/Calib",-1); // Use AliCDBManager's run number
340693af 615 AliTOFCal *calTOFArray = calib->GetTOFCalArray();
37879eed 616
617 for (ii=0; ii<fNumberOfTofClusters; ii++) {
618 digitIndex = fTofClusters[ii]->GetIndex();
619 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
620
621 Int_t index = calib->GetIndex(detectorIndex);
622
340693af 623 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
58d8d9a3 624
625 // Get channel status
626 Bool_t status=calChannel->GetStatus();
627 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
628
629 // Get Rough channel online equalization
630 Float_t roughDelay=calChannel->GetDelay();
631
632 // Get Refined channel offline calibration parameters
37879eed 633 Float_t par[6];
634 for (Int_t j = 0; j<6; j++){
340693af 635 par[j]=calChannel->GetSlewPar(j);
37879eed 636 }
340693af 637 tToT = fTofClusters[ii]->GetToT();
58d8d9a3 638 Float_t timeCorr=par[0]+par[1]*tToT+par[2]*tToT*tToT+par[3]*tToT*tToT*tToT+par[4]*tToT*tToT*tToT*tToT+par[5]*tToT*tToT*tToT*tToT*tToT+roughDelay;
340693af 639 tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
640 tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
641 fTofClusters[ii]->SetTDC(tdcCorr);
37879eed 642
643 } // loop on clusters
644
645 delete calib;
646}
d08a92dd 647//______________________________________________________________________________
648
649void AliTOFClusterFinder::ResetRecpoint()
650{
651 //
652 // Clear the list of reconstructed points
653 //
654
655 fNumberOfTofClusters = 0;
656 if (fRecPoints) fRecPoints->Clear();
657
658}
659//______________________________________________________________________________
660
661void AliTOFClusterFinder::Load()
662{
663 //
664 // Load TOF.Digits.root and TOF.RecPoints.root files
665 //
666
667 fTOFLoader->LoadDigits("READ");
668 fTOFLoader->LoadRecPoints("recreate");
669
670}
671//______________________________________________________________________________
672
673void AliTOFClusterFinder::LoadClusters()
674{
675 //
676 // Load TOF.RecPoints.root file
677 //
678
679 fTOFLoader->LoadRecPoints("recreate");
680
681}
682//______________________________________________________________________________
683
684void AliTOFClusterFinder::UnLoad()
685{
686 //
687 // Unload TOF.Digits.root and TOF.RecPoints.root files
688 //
689
690 fTOFLoader->UnloadDigits();
691 fTOFLoader->UnloadRecPoints();
692
693}
694//______________________________________________________________________________
695
696void AliTOFClusterFinder::UnLoadClusters()
697{
698 //
699 // Unload TOF.RecPoints.root file
700 //
701
702 fTOFLoader->UnloadRecPoints();
703
704}
705//______________________________________________________________________________