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.23 2007/04/23 16:51:39 decaro
19 Digits-to-raw_data conversion: correction for a more real description (A.De Caro, R.Preghenella)
21 Revision 1.22 2007/04/19 17:26:32 arcelli
22 Fix a bug (add some debug printout
24 Revision 1.21 2007/04/18 17:28:12 arcelli
25 Set the ToT bin width to the one actually used...
27 Revision 1.20 2007/03/09 09:57:23 arcelli
28 Remove a forgotten include of Riostrem
30 Revision 1.19 2007/03/08 15:41:20 arcelli
31 set uncorrected times when filling RecPoints
33 Revision 1.18 2007/03/06 16:31:20 arcelli
34 Add Uncorrected TOF Time signal
36 Revision 1.17 2007/02/28 18:09:11 arcelli
37 Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
39 Revision 1.16 2007/02/20 15:57:00 decaro
40 Raw data update: to read the TOF raw data defined in UNPACKED mode
43 Revision 0.03 2005/07/28 A. De Caro
44 Implement public method
45 Raw2Digits(Int_t, AliRawReader *)
46 to convert digits from raw data in MC digits
49 Revision 0.02 2005/07/27 A. De Caro
50 Implement public method
51 Digits2RecPoint(Int_t)
52 to convert digits in clusters
54 Revision 0.02 2005/07/26 A. De Caro
55 Implement private methods
56 InsertCluster(AliTOFcluster *)
57 FindClusterIndex(Double_t)
58 originally implemented in AliTOFtracker
59 by S. Arcelli and C. Zampolli
61 Revision 0.01 2005/07/25 A. De Caro
62 Implement public methods
63 Digits2RecPoint(AliRawReader *, TTree *)
64 Digits2RecPoint(Int_t, AliRawReader *)
65 to convert raw data in clusters
68 ////////////////////////////////////////////////////////////////
70 // Class for TOF cluster finder //
72 // Starting from Raw Data, create rec points, //
73 // fill TreeR for TOF, //
74 // write TOF.RecPoints.root file //
76 ////////////////////////////////////////////////////////////////
79 #include "TClonesArray.h"
81 #include "TStopwatch.h"
85 #include "AliLoader.h"
87 #include "AliRawReader.h"
88 #include "AliRunLoader.h"
90 #include "AliTOFCal.h"
91 #include "AliTOFcalib.h"
92 #include "AliTOFChannel.h"
93 #include "AliTOFClusterFinder.h"
94 #include "AliTOFcluster.h"
95 #include "AliTOFdigit.h"
96 #include "AliTOFGeometry.h"
97 #include "AliTOFGeometryV5.h"
98 #include "AliTOFrawData.h"
99 #include "AliTOFRawStream.h"
100 #include "Riostream.h"
102 //extern TFile *gFile;
104 ClassImp(AliTOFClusterFinder)
106 AliTOFClusterFinder::AliTOFClusterFinder():
111 fTOFGeometry(new AliTOFGeometryV5()),
112 fDigits(new TClonesArray("AliTOFdigit", 4000)),
113 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
114 fNumberOfTofClusters(0),
122 AliInfo("V5 TOF Geometry is taken as the default");
125 //______________________________________________________________________________
127 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
128 fRunLoader(runLoader),
129 fTOFLoader(runLoader->GetLoader("TOFLoader")),
132 fTOFGeometry(new AliTOFGeometryV5()),
133 fDigits(new TClonesArray("AliTOFdigit", 4000)),
134 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
135 fNumberOfTofClusters(0),
143 // runLoader->CdGAFile();
144 // TFile *in=(TFile*)gFile;
146 // fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
150 //------------------------------------------------------------------------
151 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
157 fTOFGeometry(new AliTOFGeometryV5()),
158 fDigits(new TClonesArray("AliTOFdigit", 4000)),
159 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
160 fNumberOfTofClusters(0),
165 this->fDigits=source.fDigits;
166 this->fRecPoints=source.fRecPoints;
167 this->fTOFGeometry=source.fTOFGeometry;
168 this->fDecoderVersion=source.fDecoderVersion;
172 //------------------------------------------------------------------------
173 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
176 this->fDigits=source.fDigits;
177 this->fRecPoints=source.fRecPoints;
178 this->fTOFGeometry=source.fTOFGeometry;
179 this->fVerbose=source.fVerbose;
180 this->fDecoderVersion=source.fDecoderVersion;
184 //______________________________________________________________________________
186 AliTOFClusterFinder::~AliTOFClusterFinder()
201 fRecPoints->Delete();
209 //______________________________________________________________________________
211 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
214 // Converts digits to recpoints for TOF
217 TStopwatch stopwatch;
220 fRunLoader->GetEvent(iEvent);
222 fTreeD = fTOFLoader->TreeD();
225 AliFatal("AliTOFClusterFinder: Can not get TreeD");
228 TBranch *branch = fTreeD->GetBranch("TOF");
230 AliError("can't get the branch with the TOF digits !");
234 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
235 branch->SetAddress(&digits);
239 fTreeR = fTOFLoader->TreeR();
242 fTOFLoader->MakeTree("R");
243 fTreeR = fTOFLoader->TreeR();
246 Int_t bufsize = 32000;
247 fTreeR->Branch("TOF", &fRecPoints, bufsize);
250 Int_t nDigits = digits->GetEntriesFast();
251 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
259 for (ii=0; ii<nDigits; ii++) {
260 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
261 dig[0]=d->GetSector();
262 dig[1]=d->GetPlate();
263 dig[2]=d->GetStrip();
267 //AliInfo(Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
269 for (jj=0; jj<3; jj++) g[jj] = 0.;
270 fTOFGeometry->GetPos(dig,g);
272 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
273 h[1] = TMath::ATan2(g[1],g[0]);
278 tTdcND = d->GetTdcND();
280 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
281 tofCluster->SetTDCRAW(d->GetTdc());
282 InsertCluster(tofCluster);
286 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
294 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
295 fTOFLoader->WriteRecPoints("OVERWRITE");
297 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
298 stopwatch.RealTime(),stopwatch.CpuTime()));
301 //______________________________________________________________________________
303 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
307 // Converts RAW data to recpoints for TOF
310 TStopwatch stopwatch;
313 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
314 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
318 Int_t bufsize = 32000;
319 clustersTree->Branch("TOF", &fRecPoints, bufsize);
321 TClonesArray * clonesRawData;
326 Int_t detectorIndex[5];
328 Double_t cylindricalPosition[5];
333 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
336 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
339 AliTOFRawStream tofInput(rawReader);
340 if (fDecoderVersion) {
341 AliInfo("Using New Decoder \n");
342 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
344 else tofInput.LoadRawData(indexDDL);
346 clonesRawData = (TClonesArray*)tofInput.GetRawData();
348 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
350 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
352 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
355 if (indexDDL<10) ftxt << " " << indexDDL;
356 else ftxt << " " << indexDDL;
357 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
358 else ftxt << " " << tofRawDatum->GetTRM();
359 ftxt << " " << tofRawDatum->GetTRMchain();
360 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
361 else ftxt << " " << tofRawDatum->GetTDC();
362 ftxt << " " << tofRawDatum->GetTDCchannel();
365 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
366 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
367 dummy = detectorIndex[3];
368 detectorIndex[3] = detectorIndex[4];
369 detectorIndex[4] = dummy;
372 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
373 else ftxt << " -> " << detectorIndex[0];
374 ftxt << " " << detectorIndex[1];
375 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
376 else ftxt << " " << detectorIndex[2];
377 ftxt << " " << detectorIndex[3];
378 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
379 else ftxt << " " << detectorIndex[4];
382 for (ii=0; ii<3; ii++) position[ii] = 0.;
383 fTOFGeometry->GetPos(detectorIndex, position);
385 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
386 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
387 cylindricalPosition[2] = position[2];
388 cylindricalPosition[3] = tofRawDatum->GetTOF();
389 cylindricalPosition[4] = tofRawDatum->GetTOT();
390 tToT = tofRawDatum->GetTOT();
392 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
393 tofCluster->SetToT(tToT);
394 tofCluster->SetTDCND(tTdcND);
395 tofCluster->SetTDCRAW(tofRawDatum->GetTOF());
396 InsertCluster(tofCluster);
399 if (cylindricalPosition[4]<10) ftxt << " " << cylindricalPosition[4];
400 else if (cylindricalPosition[4]>=10 && cylindricalPosition[4]<100) ftxt << " " << cylindricalPosition[4];
401 else ftxt << " " << cylindricalPosition[4];
402 if (cylindricalPosition[3]<10) ftxt << " " << cylindricalPosition[3] << endl;
403 else if (cylindricalPosition[3]>=10 && cylindricalPosition[3]<100) ftxt << " " << cylindricalPosition[3] << endl;
404 else if (cylindricalPosition[3]>=100 && cylindricalPosition[3]<1000) ftxt << " " << cylindricalPosition[3] << endl;
405 else ftxt << " " << cylindricalPosition[3] << endl;
408 } // closed loop on TOF raw data per current DDL file
410 clonesRawData->Clear();
412 } // closed loop on DDL index
416 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
419 AliTOFRawStream tofInput(rawReader);
420 rawReader->Select("TOF", indexDDL, indexDDL);
422 while(tofInput.Next()) {
424 for (ii=0; ii<5; ii++) detectorIndex[ii] = -1;
426 detectorIndex[0] = tofInput.GetSector();
427 detectorIndex[1] = tofInput.GetPlate();
428 detectorIndex[2] = tofInput.GetStrip();
429 detectorIndex[3] = tofInput.GetPadZ();
430 detectorIndex[4] = tofInput.GetPadX();
432 //AliInfo(Form(" %2i %1i %2i %1i %2i ",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
434 if (detectorIndex[0]==-1 ||
435 detectorIndex[1]==-1 ||
436 detectorIndex[2]==-1 ||
437 detectorIndex[3]==-1 ||
438 detectorIndex[4]==-1) continue;
440 for (ii=0; ii<3; ii++) position[ii] = 0.;
442 fTOFGeometry->GetPos(detectorIndex, position);
444 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
445 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
446 cylindricalPosition[2] = position[2];
447 cylindricalPosition[3] = tofInput.GetTofBin();
448 cylindricalPosition[4] = tofInput.GetToTbin();
449 tToT = tofInput.GetToTbin();
451 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
452 tofCluster->SetToT(tToT);
453 tofCluster->SetTDCND(tTdcND);
454 InsertCluster(tofCluster);
458 } // loop on DDL files
461 if (fVerbose==2) ftxt.close();
463 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
468 clustersTree->Fill();
472 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
473 stopwatch.RealTime(),stopwatch.CpuTime()));
476 //______________________________________________________________________________
478 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
481 // Converts RAW data to recpoints for TOF
484 TStopwatch stopwatch;
487 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
488 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
490 fRunLoader->GetEvent(iEvent);
492 AliDebug(2,Form(" Event number %2i ", iEvent));
494 fTreeR = fTOFLoader->TreeR();
497 fTOFLoader->MakeTree("R");
498 fTreeR = fTOFLoader->TreeR();
501 Int_t bufsize = 32000;
502 fTreeR->Branch("TOF", &fRecPoints, bufsize);
504 TClonesArray * clonesRawData;
509 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
511 Double_t cylindricalPosition[5];
516 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
519 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
522 AliTOFRawStream tofInput(rawReader);
523 if (fDecoderVersion) {
524 AliInfo("Using New Decoder \n");
525 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
527 else tofInput.LoadRawData(indexDDL);
529 clonesRawData = (TClonesArray*)tofInput.GetRawData();
531 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
533 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
535 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
538 if (indexDDL<10) ftxt << " " << indexDDL;
539 else ftxt << " " << indexDDL;
540 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
541 else ftxt << " " << tofRawDatum->GetTRM();
542 ftxt << " " << tofRawDatum->GetTRMchain();
543 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
544 else ftxt << " " << tofRawDatum->GetTDC();
545 ftxt << " " << tofRawDatum->GetTDCchannel();
548 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
549 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
550 dummy = detectorIndex[3];
551 detectorIndex[3] = detectorIndex[4];
552 detectorIndex[4] = dummy;
555 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
556 else ftxt << " -> " << detectorIndex[0];
557 ftxt << " " << detectorIndex[1];
558 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
559 else ftxt << " " << detectorIndex[2];
560 ftxt << " " << detectorIndex[3];
561 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
562 else ftxt << " " << detectorIndex[4];
565 for (ii=0; ii<3; ii++) position[ii] = 0.;
566 fTOFGeometry->GetPos(detectorIndex, position);
568 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
569 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
570 cylindricalPosition[2] = position[2];
571 cylindricalPosition[3] = tofRawDatum->GetTOF();
572 cylindricalPosition[4] = tofRawDatum->GetTOT();
573 tToT = tofRawDatum->GetTOT();
575 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
576 tofCluster->SetToT(tToT);
577 tofCluster->SetTDCND(tTdcND);
578 tofCluster->SetTDCRAW(tofRawDatum->GetTOF());
579 InsertCluster(tofCluster);
582 if (cylindricalPosition[4]<10) ftxt << " " << cylindricalPosition[4];
583 else if (cylindricalPosition[4]>=10 && cylindricalPosition[4]<100) ftxt << " " << cylindricalPosition[4];
584 else ftxt << " " << cylindricalPosition[4];
585 if (cylindricalPosition[3]<10) ftxt << " " << cylindricalPosition[3] << endl;
586 else if (cylindricalPosition[3]>=10 && cylindricalPosition[3]<100) ftxt << " " << cylindricalPosition[3] << endl;
587 else if (cylindricalPosition[3]>=100 && cylindricalPosition[3]<1000) ftxt << " " << cylindricalPosition[3] << endl;
588 else ftxt << " " << cylindricalPosition[3] << endl;
591 } // closed loop on TOF raw data per current DDL file
593 clonesRawData->Clear();
595 } // closed loop on DDL index
597 if (fVerbose==2) ftxt.close();
599 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
607 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
608 fTOFLoader->WriteRecPoints("OVERWRITE");
610 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
611 stopwatch.RealTime(),stopwatch.CpuTime()));
614 //______________________________________________________________________________
616 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
619 // Converts RAW data to MC digits for TOF
621 // (temporary solution)
624 TStopwatch stopwatch;
627 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
628 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
630 fRunLoader->GetEvent(iEvent);
632 fTreeD = fTOFLoader->TreeD();
635 AliInfo("TreeD re-creation");
637 fTOFLoader->MakeTree("D");
638 fTreeD = fTOFLoader->TreeD();
641 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
642 Int_t bufsize = 32000;
643 fTreeD->Branch("TOF", &tofDigits, bufsize);
645 fRunLoader->GetEvent(iEvent);
647 AliDebug(2,Form(" Event number %2i ", iEvent));
649 TClonesArray * clonesRawData;
653 Int_t detectorIndex[5];
657 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
660 AliTOFRawStream tofInput(rawReader);
661 if (fDecoderVersion) {
662 AliInfo("Using New Decoder \n");
663 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
665 else tofInput.LoadRawData(indexDDL);
667 clonesRawData = (TClonesArray*)tofInput.GetRawData();
669 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
671 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
673 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
675 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
676 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
677 dummy = detectorIndex[3];
678 detectorIndex[3] = detectorIndex[4];
679 detectorIndex[4] = dummy;
681 digit[0] = (Float_t)tofInput.GetTofBin();
682 digit[1] = (Float_t)tofInput.GetToTbin();
683 digit[2] = (Float_t)tofInput.GetToTbin();
686 Int_t tracknum[3]={-1,-1,-1};
688 TClonesArray &aDigits = *tofDigits;
689 Int_t last=tofDigits->GetEntriesFast();
690 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
694 clonesRawData->Clear();
700 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
701 fTOFLoader->WriteDigits("OVERWRITE");
703 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
704 stopwatch.RealTime(),stopwatch.CpuTime()));
707 //______________________________________________________________________________
709 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
710 //---------------------------------------------------------------------------//
711 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
712 //---------------------------------------------------------------------------//
713 if (fNumberOfTofClusters==kTofMaxCluster) {
714 AliError("Too many clusters !");
718 if (fNumberOfTofClusters==0) {
719 fTofClusters[fNumberOfTofClusters++] = tofCluster;
723 Int_t ii = FindClusterIndex(tofCluster->GetZ());
724 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
725 fTofClusters[ii] = tofCluster;
726 fNumberOfTofClusters++;
731 //_________________________________________________________________________
733 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
734 //--------------------------------------------------------------------
735 // This function returns the index of the nearest cluster
736 //--------------------------------------------------------------------
737 if (fNumberOfTofClusters==0) return 0;
738 if (z <= fTofClusters[0]->GetZ()) return 0;
739 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
740 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
741 for (; b<e; m=(b+e)/2) {
742 if (z > fTofClusters[m]->GetZ()) b=m+1;
749 //_________________________________________________________________________
751 void AliTOFClusterFinder::FillRecPoint()
754 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
755 // in Z) in the global TClonesArray of AliTOFcluster,
761 Int_t detectorIndex[5];
762 Double_t cylindricalPosition[5];
763 Int_t trackLabels[3];
764 Int_t digitIndex = -1;
768 Bool_t cStatus = kTRUE;
770 TClonesArray &lRecPoints = *fRecPoints;
772 for (ii=0; ii<fNumberOfTofClusters; ii++) {
774 digitIndex = fTofClusters[ii]->GetIndex();
775 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
776 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
777 cylindricalPosition[0] = fTofClusters[ii]->GetR();
778 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
779 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
780 cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
781 cylindricalPosition[4] = fTofClusters[ii]->GetADC();
782 tToT = fTofClusters[ii]->GetToT();
783 tTdcND = fTofClusters[ii]->GetTDCND();
784 cStatus=fTofClusters[ii]->GetStatus();
785 tTdcRAW=fTofClusters[ii]->GetTDCRAW();
786 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND, tTdcRAW,cStatus);
788 //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]));
790 } // loop on clusters
794 //_________________________________________________________________________
795 void AliTOFClusterFinder::CalibrateRecPoint()
798 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
799 // in Z) in the global TClonesArray of AliTOFcluster,
805 Int_t detectorIndex[5];
806 Int_t digitIndex = -1;
809 AliInfo(" Calibrating TOF Clusters: ")
810 AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
811 if(!calib->ReadParFromCDB("TOF/Calib",-1)) {AliFatal("Exiting, no CDB object found!!!");exit(0);}
813 AliTOFCal *calTOFArray = calib->GetTOFCalArray();
815 for (ii=0; ii<fNumberOfTofClusters; ii++) {
816 digitIndex = fTofClusters[ii]->GetIndex();
817 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
819 Int_t index = calib->GetIndex(detectorIndex);
821 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
823 // Get channel status
824 Bool_t status=calChannel->GetStatus();
825 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
827 // Get Rough channel online equalization
828 Float_t roughDelay=calChannel->GetDelay();
829 AliDebug(2,Form(" channel delay = %f", roughDelay));
830 // Get Refined channel offline calibration parameters
832 for (Int_t j = 0; j<6; j++){
833 par[j]=calChannel->GetSlewPar(j);
835 tToT = fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()*1.E-3;
836 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;
837 AliDebug(2,Form(" time correction (ns) = %f", timeCorr));
838 AliDebug(2,Form(" channel time, uncorr (ns)= %f",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3 ));
839 tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
840 tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
841 fTofClusters[ii]->SetTDC(tdcCorr);
842 AliDebug(2,Form(" channel time, corr (ns)= %f",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3 ));
844 } // loop on clusters
848 //______________________________________________________________________________
850 void AliTOFClusterFinder::ResetRecpoint()
853 // Clear the list of reconstructed points
856 fNumberOfTofClusters = 0;
857 if (fRecPoints) fRecPoints->Clear();
860 //______________________________________________________________________________
862 void AliTOFClusterFinder::Load()
865 // Load TOF.Digits.root and TOF.RecPoints.root files
868 fTOFLoader->LoadDigits("READ");
869 fTOFLoader->LoadRecPoints("recreate");
872 //______________________________________________________________________________
874 void AliTOFClusterFinder::LoadClusters()
877 // Load TOF.RecPoints.root file
880 fTOFLoader->LoadRecPoints("recreate");
883 //______________________________________________________________________________
885 void AliTOFClusterFinder::UnLoad()
888 // Unload TOF.Digits.root and TOF.RecPoints.root files
891 fTOFLoader->UnloadDigits();
892 fTOFLoader->UnloadRecPoints();
895 //______________________________________________________________________________
897 void AliTOFClusterFinder::UnLoadClusters()
900 // Unload TOF.RecPoints.root file
903 fTOFLoader->UnloadRecPoints();
906 //______________________________________________________________________________