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 Revision 1.20 2007/03/09 09:57:23 arcelli
19 Remove a forgotten include of Riostrem
21 Revision 1.19 2007/03/08 15:41:20 arcelli
22 set uncorrected times when filling RecPoints
24 Revision 1.18 2007/03/06 16:31:20 arcelli
25 Add Uncorrected TOF Time signal
27 Revision 1.17 2007/02/28 18:09:11 arcelli
28 Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
30 Revision 1.16 2007/02/20 15:57:00 decaro
31 Raw data update: to read the TOF raw data defined in UNPACKED mode
34 Revision 0.03 2005/07/28 A. De Caro
35 Implement public method
36 Raw2Digits(Int_t, AliRawReader *)
37 to convert digits from raw data in MC digits
40 Revision 0.02 2005/07/27 A. De Caro
41 Implement public method
42 Digits2RecPoint(Int_t)
43 to convert digits in clusters
45 Revision 0.02 2005/07/26 A. De Caro
46 Implement private methods
47 InsertCluster(AliTOFcluster *)
48 FindClusterIndex(Double_t)
49 originally implemented in AliTOFtracker
50 by S. Arcelli and C. Zampolli
52 Revision 0.01 2005/07/25 A. De Caro
53 Implement public methods
54 Digits2RecPoint(AliRawReader *, TTree *)
55 Digits2RecPoint(Int_t, AliRawReader *)
56 to convert raw data in clusters
59 ////////////////////////////////////////////////////////////////
61 // Class for TOF cluster finder //
63 // Starting from Raw Data, create rec points, //
64 // fill TreeR for TOF, //
65 // write TOF.RecPoints.root file //
67 ////////////////////////////////////////////////////////////////
70 #include "TClonesArray.h"
75 #include "AliLoader.h"
77 #include "AliRawReader.h"
78 #include "AliRunLoader.h"
80 #include "AliTOFCal.h"
81 #include "AliTOFcalib.h"
82 #include "AliTOFChannel.h"
83 #include "AliTOFClusterFinder.h"
84 #include "AliTOFcluster.h"
85 #include "AliTOFdigit.h"
86 #include "AliTOFGeometry.h"
87 #include "AliTOFGeometryV5.h"
88 #include "AliTOFrawData.h"
89 #include "AliTOFRawStream.h"
90 #include "Riostream.h"
92 //extern TFile *gFile;
94 ClassImp(AliTOFClusterFinder)
96 AliTOFClusterFinder::AliTOFClusterFinder():
101 fTOFGeometry(new AliTOFGeometryV5()),
102 fDigits(new TClonesArray("AliTOFdigit", 4000)),
103 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
104 fNumberOfTofClusters(0),
111 AliInfo("V5 TOF Geometry is taken as the default");
114 //______________________________________________________________________________
116 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
117 fRunLoader(runLoader),
118 fTOFLoader(runLoader->GetLoader("TOFLoader")),
121 fTOFGeometry(new AliTOFGeometryV5()),
122 fDigits(new TClonesArray("AliTOFdigit", 4000)),
123 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
124 fNumberOfTofClusters(0),
131 // runLoader->CdGAFile();
132 // TFile *in=(TFile*)gFile;
134 // fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
138 //------------------------------------------------------------------------
139 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
145 fTOFGeometry(new AliTOFGeometryV5()),
146 fDigits(new TClonesArray("AliTOFdigit", 4000)),
147 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
148 fNumberOfTofClusters(0),
152 this->fDigits=source.fDigits;
153 this->fRecPoints=source.fRecPoints;
154 this->fTOFGeometry=source.fTOFGeometry;
158 //------------------------------------------------------------------------
159 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
162 this->fDigits=source.fDigits;
163 this->fRecPoints=source.fRecPoints;
164 this->fTOFGeometry=source.fTOFGeometry;
165 this->fVerbose=source.fVerbose;
169 //______________________________________________________________________________
171 AliTOFClusterFinder::~AliTOFClusterFinder()
186 fRecPoints->Delete();
194 //______________________________________________________________________________
196 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
199 // Converts digits to recpoints for TOF
202 fRunLoader->GetEvent(iEvent);
204 fTreeD = fTOFLoader->TreeD();
207 AliFatal("AliTOFClusterFinder: Can not get TreeD");
210 TBranch *branch = fTreeD->GetBranch("TOF");
212 AliError("can't get the branch with the TOF digits !");
216 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
217 branch->SetAddress(&digits);
221 fTreeR = fTOFLoader->TreeR();
224 fTOFLoader->MakeTree("R");
225 fTreeR = fTOFLoader->TreeR();
228 Int_t bufsize = 32000;
229 fTreeR->Branch("TOF", &fRecPoints, bufsize);
232 Int_t nDigits = digits->GetEntriesFast();
233 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
241 for (ii=0; ii<nDigits; ii++) {
242 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
243 dig[0]=d->GetSector();
244 dig[1]=d->GetPlate();
245 dig[2]=d->GetStrip();
249 //AliInfo(Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
251 for (jj=0; jj<3; jj++) g[jj] = 0.;
252 fTOFGeometry->GetPos(dig,g);
254 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
255 h[1] = TMath::ATan2(g[1],g[0]);
260 tTdcND = d->GetTdcND();
262 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
263 tofCluster->SetTDCRAW(d->GetTdc());
264 InsertCluster(tofCluster);
268 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
276 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
277 fTOFLoader->WriteRecPoints("OVERWRITE");
280 //______________________________________________________________________________
282 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
286 // Converts RAW data to recpoints for TOF
289 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
290 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
294 Int_t bufsize = 32000;
295 clustersTree->Branch("TOF", &fRecPoints, bufsize);
297 TClonesArray * clonesRawData;
302 Int_t detectorIndex[5];
304 Double_t cylindricalPosition[5];
309 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
312 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
315 AliTOFRawStream tofInput(rawReader);
316 tofInput.LoadRawData(indexDDL);
318 clonesRawData = (TClonesArray*)tofInput.GetRawData();
320 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
322 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
324 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
327 if (indexDDL<10) ftxt << " " << indexDDL;
328 else ftxt << " " << indexDDL;
329 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
330 else ftxt << " " << tofRawDatum->GetTRM();
331 ftxt << " " << tofRawDatum->GetTRMchain();
332 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
333 else ftxt << " " << tofRawDatum->GetTDC();
334 ftxt << " " << tofRawDatum->GetTDCchannel();
337 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
338 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
339 dummy = detectorIndex[3];
340 detectorIndex[3] = detectorIndex[4];
341 detectorIndex[4] = dummy;
344 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
345 else ftxt << " -> " << detectorIndex[0];
346 ftxt << " " << detectorIndex[1];
347 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
348 else ftxt << " " << detectorIndex[2];
349 ftxt << " " << detectorIndex[3];
350 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
351 else ftxt << " " << detectorIndex[4];
354 for (ii=0; ii<3; ii++) position[ii] = 0.;
355 fTOFGeometry->GetPos(detectorIndex, position);
357 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
358 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
359 cylindricalPosition[2] = position[2];
360 cylindricalPosition[3] = tofRawDatum->GetTOF();
361 cylindricalPosition[4] = tofRawDatum->GetTOT();
362 tToT = tofRawDatum->GetTOT();
364 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
365 tofCluster->SetToT(tToT);
366 tofCluster->SetTDCND(tTdcND);
367 tofCluster->SetTDCRAW(tofRawDatum->GetTOF());
368 InsertCluster(tofCluster);
371 if (cylindricalPosition[4]<10) ftxt << " " << cylindricalPosition[4];
372 else if (cylindricalPosition[4]>=10 && cylindricalPosition[4]<100) ftxt << " " << cylindricalPosition[4];
373 else ftxt << " " << cylindricalPosition[4];
374 if (cylindricalPosition[3]<10) ftxt << " " << cylindricalPosition[3] << endl;
375 else if (cylindricalPosition[3]>=10 && cylindricalPosition[3]<100) ftxt << " " << cylindricalPosition[3] << endl;
376 else if (cylindricalPosition[3]>=100 && cylindricalPosition[3]<1000) ftxt << " " << cylindricalPosition[3] << endl;
377 else ftxt << " " << cylindricalPosition[3] << endl;
380 } // closed loop on TOF raw data per current DDL file
382 clonesRawData->Clear();
384 } // closed loop on DDL index
388 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
391 AliTOFRawStream tofInput(rawReader);
392 rawReader->Select("TOF", indexDDL, indexDDL);
394 while(tofInput.Next()) {
396 for (ii=0; ii<5; ii++) detectorIndex[ii] = -1;
398 detectorIndex[0] = tofInput.GetSector();
399 detectorIndex[1] = tofInput.GetPlate();
400 detectorIndex[2] = tofInput.GetStrip();
401 detectorIndex[3] = tofInput.GetPadZ();
402 detectorIndex[4] = tofInput.GetPadX();
404 //AliInfo(Form(" %2i %1i %2i %1i %2i ",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
406 if (detectorIndex[0]==-1 ||
407 detectorIndex[1]==-1 ||
408 detectorIndex[2]==-1 ||
409 detectorIndex[3]==-1 ||
410 detectorIndex[4]==-1) continue;
412 for (ii=0; ii<3; ii++) position[ii] = 0.;
414 fTOFGeometry->GetPos(detectorIndex, position);
416 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
417 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
418 cylindricalPosition[2] = position[2];
419 cylindricalPosition[3] = tofInput.GetTofBin();
420 cylindricalPosition[4] = tofInput.GetToTbin();
421 tToT = tofInput.GetToTbin();
423 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
424 tofCluster->SetToT(tToT);
425 tofCluster->SetTDCND(tTdcND);
426 InsertCluster(tofCluster);
430 } // loop on DDL files
433 if (fVerbose==2) ftxt.close();
435 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
440 clustersTree->Fill();
445 //______________________________________________________________________________
447 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
450 // Converts RAW data to recpoints for TOF
453 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
454 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
456 fRunLoader->GetEvent(iEvent);
458 AliDebug(2,Form(" Event number %2i ", iEvent));
460 fTreeR = fTOFLoader->TreeR();
463 fTOFLoader->MakeTree("R");
464 fTreeR = fTOFLoader->TreeR();
467 Int_t bufsize = 32000;
468 fTreeR->Branch("TOF", &fRecPoints, bufsize);
470 TClonesArray * clonesRawData;
475 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
477 Double_t cylindricalPosition[5];
482 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
485 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
488 AliTOFRawStream tofInput(rawReader);
489 tofInput.LoadRawData(indexDDL);
491 clonesRawData = (TClonesArray*)tofInput.GetRawData();
493 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
495 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
497 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
500 if (indexDDL<10) ftxt << " " << indexDDL;
501 else ftxt << " " << indexDDL;
502 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
503 else ftxt << " " << tofRawDatum->GetTRM();
504 ftxt << " " << tofRawDatum->GetTRMchain();
505 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
506 else ftxt << " " << tofRawDatum->GetTDC();
507 ftxt << " " << tofRawDatum->GetTDCchannel();
510 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
511 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
512 dummy = detectorIndex[3];
513 detectorIndex[3] = detectorIndex[4];
514 detectorIndex[4] = dummy;
517 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
518 else ftxt << " -> " << detectorIndex[0];
519 ftxt << " " << detectorIndex[1];
520 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
521 else ftxt << " " << detectorIndex[2];
522 ftxt << " " << detectorIndex[3];
523 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
524 else ftxt << " " << detectorIndex[4];
527 for (ii=0; ii<3; ii++) position[ii] = 0.;
528 fTOFGeometry->GetPos(detectorIndex, position);
530 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
531 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
532 cylindricalPosition[2] = position[2];
533 cylindricalPosition[3] = tofRawDatum->GetTOF();
534 cylindricalPosition[4] = tofRawDatum->GetTOT();
535 tToT = tofRawDatum->GetTOT();
537 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
538 tofCluster->SetToT(tToT);
539 tofCluster->SetTDCND(tTdcND);
540 tofCluster->SetTDCRAW(tofRawDatum->GetTOF());
541 InsertCluster(tofCluster);
544 if (cylindricalPosition[4]<10) ftxt << " " << cylindricalPosition[4];
545 else if (cylindricalPosition[4]>=10 && cylindricalPosition[4]<100) ftxt << " " << cylindricalPosition[4];
546 else ftxt << " " << cylindricalPosition[4];
547 if (cylindricalPosition[3]<10) ftxt << " " << cylindricalPosition[3] << endl;
548 else if (cylindricalPosition[3]>=10 && cylindricalPosition[3]<100) ftxt << " " << cylindricalPosition[3] << endl;
549 else if (cylindricalPosition[3]>=100 && cylindricalPosition[3]<1000) ftxt << " " << cylindricalPosition[3] << endl;
550 else ftxt << " " << cylindricalPosition[3] << endl;
553 } // closed loop on TOF raw data per current DDL file
555 clonesRawData->Clear();
557 } // closed loop on DDL index
559 if (fVerbose==2) ftxt.close();
561 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
569 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
570 fTOFLoader->WriteRecPoints("OVERWRITE");
573 //______________________________________________________________________________
575 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
578 // Converts RAW data to MC digits for TOF
580 // (temporary solution)
583 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
584 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
586 fRunLoader->GetEvent(iEvent);
588 fTreeD = fTOFLoader->TreeD();
591 AliInfo("TreeD re-creation");
593 fTOFLoader->MakeTree("D");
594 fTreeD = fTOFLoader->TreeD();
597 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
598 Int_t bufsize = 32000;
599 fTreeD->Branch("TOF", &tofDigits, bufsize);
601 fRunLoader->GetEvent(iEvent);
603 AliDebug(2,Form(" Event number %2i ", iEvent));
605 TClonesArray * clonesRawData;
609 Int_t detectorIndex[5];
613 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
616 AliTOFRawStream tofInput(rawReader);
617 tofInput.LoadRawData(indexDDL);
619 clonesRawData = (TClonesArray*)tofInput.GetRawData();
621 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
623 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
625 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
627 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
628 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
629 dummy = detectorIndex[3];
630 detectorIndex[3] = detectorIndex[4];
631 detectorIndex[4] = dummy;
633 digit[0] = (Float_t)tofInput.GetTofBin();
634 digit[1] = (Float_t)tofInput.GetToTbin();
635 digit[2] = (Float_t)tofInput.GetToTbin();
638 Int_t tracknum[3]={-1,-1,-1};
640 TClonesArray &aDigits = *tofDigits;
641 Int_t last=tofDigits->GetEntriesFast();
642 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
646 clonesRawData->Clear();
652 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
653 fTOFLoader->WriteDigits("OVERWRITE");
656 //______________________________________________________________________________
658 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
659 //---------------------------------------------------------------------------//
660 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
661 //---------------------------------------------------------------------------//
662 if (fNumberOfTofClusters==kTofMaxCluster) {
663 AliError("Too many clusters !");
667 if (fNumberOfTofClusters==0) {
668 fTofClusters[fNumberOfTofClusters++] = tofCluster;
672 Int_t ii = FindClusterIndex(tofCluster->GetZ());
673 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
674 fTofClusters[ii] = tofCluster;
675 fNumberOfTofClusters++;
680 //_________________________________________________________________________
682 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
683 //--------------------------------------------------------------------
684 // This function returns the index of the nearest cluster
685 //--------------------------------------------------------------------
686 if (fNumberOfTofClusters==0) return 0;
687 if (z <= fTofClusters[0]->GetZ()) return 0;
688 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
689 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
690 for (; b<e; m=(b+e)/2) {
691 if (z > fTofClusters[m]->GetZ()) b=m+1;
698 //_________________________________________________________________________
700 void AliTOFClusterFinder::FillRecPoint()
703 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
704 // in Z) in the global TClonesArray of AliTOFcluster,
710 Int_t detectorIndex[5];
711 Double_t cylindricalPosition[5];
712 Int_t trackLabels[3];
713 Int_t digitIndex = -1;
717 Bool_t cStatus = kTRUE;
719 TClonesArray &lRecPoints = *fRecPoints;
721 for (ii=0; ii<fNumberOfTofClusters; ii++) {
723 digitIndex = fTofClusters[ii]->GetIndex();
724 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
725 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
726 cylindricalPosition[0] = fTofClusters[ii]->GetR();
727 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
728 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
729 cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
730 cylindricalPosition[4] = fTofClusters[ii]->GetADC();
731 tToT = fTofClusters[ii]->GetToT();
732 tTdcND = fTofClusters[ii]->GetTDCND();
733 cStatus=fTofClusters[ii]->GetStatus();
734 tTdcRAW=fTofClusters[ii]->GetTDCRAW();
735 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND, tTdcRAW,cStatus);
737 //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]));
739 } // loop on clusters
743 //_________________________________________________________________________
744 void AliTOFClusterFinder::CalibrateRecPoint()
747 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
748 // in Z) in the global TClonesArray of AliTOFcluster,
754 Int_t detectorIndex[5];
755 Int_t digitIndex = -1;
758 AliInfo(" Calibrating TOF Clusters: ")
759 AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
760 // calib->ReadParFromCDB("TOF/Calib",0); // original
761 // Use AliCDBManager's run number
762 if(!calib->ReadParFromCDB("TOF/Calib",-1)) {AliFatal("Exiting, no CDB object found!!!");exit(0);}
764 AliTOFCal *calTOFArray = calib->GetTOFCalArray();
766 for (ii=0; ii<fNumberOfTofClusters; ii++) {
767 digitIndex = fTofClusters[ii]->GetIndex();
768 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
770 Int_t index = calib->GetIndex(detectorIndex);
772 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
774 // Get channel status
775 Bool_t status=calChannel->GetStatus();
776 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
778 // Get Rough channel online equalization
779 Float_t roughDelay=calChannel->GetDelay();
781 // Get Refined channel offline calibration parameters
783 for (Int_t j = 0; j<6; j++){
784 par[j]=calChannel->GetSlewPar(j);
786 tToT = fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()*1.E-3;
787 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;
788 tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
789 tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
790 fTofClusters[ii]->SetTDC(tdcCorr);
792 } // loop on clusters
796 //______________________________________________________________________________
798 void AliTOFClusterFinder::ResetRecpoint()
801 // Clear the list of reconstructed points
804 fNumberOfTofClusters = 0;
805 if (fRecPoints) fRecPoints->Clear();
808 //______________________________________________________________________________
810 void AliTOFClusterFinder::Load()
813 // Load TOF.Digits.root and TOF.RecPoints.root files
816 fTOFLoader->LoadDigits("READ");
817 fTOFLoader->LoadRecPoints("recreate");
820 //______________________________________________________________________________
822 void AliTOFClusterFinder::LoadClusters()
825 // Load TOF.RecPoints.root file
828 fTOFLoader->LoadRecPoints("recreate");
831 //______________________________________________________________________________
833 void AliTOFClusterFinder::UnLoad()
836 // Unload TOF.Digits.root and TOF.RecPoints.root files
839 fTOFLoader->UnloadDigits();
840 fTOFLoader->UnloadRecPoints();
843 //______________________________________________________________________________
845 void AliTOFClusterFinder::UnLoadClusters()
848 // Unload TOF.RecPoints.root file
851 fTOFLoader->UnloadRecPoints();
854 //______________________________________________________________________________