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.25 2007/04/30 15:22:17 arcelli
19 Change TOF digit Time, Tot etc to int type
21 Revision 1.24 2007/04/27 11:19:31 arcelli
22 updates for the new decoder
24 Revision 1.23 2007/04/23 16:51:39 decaro
25 Digits-to-raw_data conversion: correction for a more real description (A.De Caro, R.Preghenella)
27 Revision 1.22 2007/04/19 17:26:32 arcelli
28 Fix a bug (add some debug printout
30 Revision 1.21 2007/04/18 17:28:12 arcelli
31 Set the ToT bin width to the one actually used...
33 Revision 1.20 2007/03/09 09:57:23 arcelli
34 Remove a forgotten include of Riostrem
36 Revision 1.19 2007/03/08 15:41:20 arcelli
37 set uncorrected times when filling RecPoints
39 Revision 1.18 2007/03/06 16:31:20 arcelli
40 Add Uncorrected TOF Time signal
42 Revision 1.17 2007/02/28 18:09:11 arcelli
43 Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
45 Revision 1.16 2007/02/20 15:57:00 decaro
46 Raw data update: to read the TOF raw data defined in UNPACKED mode
49 Revision 0.03 2005/07/28 A. De Caro
50 Implement public method
51 Raw2Digits(Int_t, AliRawReader *)
52 to convert digits from raw data in MC digits
55 Revision 0.02 2005/07/27 A. De Caro
56 Implement public method
57 Digits2RecPoint(Int_t)
58 to convert digits in clusters
60 Revision 0.02 2005/07/26 A. De Caro
61 Implement private methods
62 InsertCluster(AliTOFcluster *)
63 FindClusterIndex(Double_t)
64 originally implemented in AliTOFtracker
65 by S. Arcelli and C. Zampolli
67 Revision 0.01 2005/07/25 A. De Caro
68 Implement public methods
69 Digits2RecPoint(AliRawReader *, TTree *)
70 Digits2RecPoint(Int_t, AliRawReader *)
71 to convert raw data in clusters
74 ////////////////////////////////////////////////////////////////
76 // Class for TOF cluster finder //
78 // Starting from Raw Data, create rec points, //
79 // fill TreeR for TOF, //
80 // write TOF.RecPoints.root file //
82 ////////////////////////////////////////////////////////////////
85 #include "TClonesArray.h"
87 #include "TStopwatch.h"
91 #include "AliLoader.h"
93 #include "AliRawReader.h"
94 #include "AliRunLoader.h"
96 #include "AliTOFCal.h"
97 #include "AliTOFcalib.h"
98 #include "AliTOFChannel.h"
99 #include "AliTOFClusterFinder.h"
100 #include "AliTOFcluster.h"
101 #include "AliTOFdigit.h"
102 #include "AliTOFGeometry.h"
103 #include "AliTOFGeometryV5.h"
104 #include "AliTOFrawData.h"
105 #include "AliTOFRawStream.h"
106 #include "Riostream.h"
108 //extern TFile *gFile;
110 ClassImp(AliTOFClusterFinder)
112 AliTOFClusterFinder::AliTOFClusterFinder():
117 fTOFGeometry(new AliTOFGeometryV5()),
118 fDigits(new TClonesArray("AliTOFdigit", 4000)),
119 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
120 fNumberOfTofClusters(0),
128 AliInfo("V5 TOF Geometry is taken as the default");
131 //______________________________________________________________________________
133 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
134 fRunLoader(runLoader),
135 fTOFLoader(runLoader->GetLoader("TOFLoader")),
138 fTOFGeometry(new AliTOFGeometryV5()),
139 fDigits(new TClonesArray("AliTOFdigit", 4000)),
140 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
141 fNumberOfTofClusters(0),
149 // runLoader->CdGAFile();
150 // TFile *in=(TFile*)gFile;
152 // fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
156 //------------------------------------------------------------------------
157 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
163 fTOFGeometry(new AliTOFGeometryV5()),
164 fDigits(new TClonesArray("AliTOFdigit", 4000)),
165 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
166 fNumberOfTofClusters(0),
171 this->fDigits=source.fDigits;
172 this->fRecPoints=source.fRecPoints;
173 this->fTOFGeometry=source.fTOFGeometry;
174 this->fDecoderVersion=source.fDecoderVersion;
178 //------------------------------------------------------------------------
179 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
182 this->fDigits=source.fDigits;
183 this->fRecPoints=source.fRecPoints;
184 this->fTOFGeometry=source.fTOFGeometry;
185 this->fVerbose=source.fVerbose;
186 this->fDecoderVersion=source.fDecoderVersion;
190 //______________________________________________________________________________
192 AliTOFClusterFinder::~AliTOFClusterFinder()
207 fRecPoints->Delete();
215 //______________________________________________________________________________
217 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
220 // Converts digits to recpoints for TOF
223 TStopwatch stopwatch;
226 fRunLoader->GetEvent(iEvent);
228 fTreeD = fTOFLoader->TreeD();
231 AliFatal("AliTOFClusterFinder: Can not get TreeD");
234 TBranch *branch = fTreeD->GetBranch("TOF");
236 AliError("can't get the branch with the TOF digits !");
240 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
241 branch->SetAddress(&digits);
245 fTreeR = fTOFLoader->TreeR();
248 fTOFLoader->MakeTree("R");
249 fTreeR = fTOFLoader->TreeR();
252 Int_t bufsize = 32000;
253 fTreeR->Branch("TOF", &fRecPoints, bufsize);
256 Int_t nDigits = digits->GetEntriesFast();
257 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
260 Int_t dig[5]; //cluster detector indeces
261 Float_t g[3]; //cluster cartesian coord
262 Double_t h[3]; // the cluster spatial cyl. coordinates
263 Int_t parTOF[5]; //The TOF signal parameters
264 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
265 for (ii=0; ii<nDigits; ii++) {
266 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
267 dig[0]=d->GetSector();
268 dig[1]=d->GetPlate();
269 dig[2]=d->GetStrip();
273 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
275 for (jj=0; jj<3; jj++) g[jj] = 0.;
276 fTOFGeometry->GetPos(dig,g);
278 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
279 h[1] = TMath::ATan2(g[1],g[0]);
282 parTOF[0] = d->GetTdc(); //the TDC signal
283 parTOF[1] = d->GetToT(); //the ToT signal
284 parTOF[2] = d->GetAdc(); // the adc charge
285 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
286 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
287 AliTOFcluster *tofCluster = new AliTOFcluster(h,dig,parTOF,status,d->GetTracks(),ii);
288 InsertCluster(tofCluster);
292 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
300 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
301 fTOFLoader->WriteRecPoints("OVERWRITE");
303 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
304 stopwatch.RealTime(),stopwatch.CpuTime()));
307 //______________________________________________________________________________
309 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
313 // Converts RAW data to recpoints for TOF
316 TStopwatch stopwatch;
319 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
320 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
324 Int_t bufsize = 32000;
325 clustersTree->Branch("TOF", &fRecPoints, bufsize);
327 TClonesArray * clonesRawData;
332 Int_t detectorIndex[5];
334 Double_t cylindricalPosition[3];
338 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
341 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
344 AliTOFRawStream tofInput(rawReader);
345 if (fDecoderVersion) {
346 AliInfo("Using New Decoder \n");
347 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
349 else tofInput.LoadRawData(indexDDL);
351 clonesRawData = (TClonesArray*)tofInput.GetRawData();
353 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
355 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
357 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
360 if (indexDDL<10) ftxt << " " << indexDDL;
361 else ftxt << " " << indexDDL;
362 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
363 else ftxt << " " << tofRawDatum->GetTRM();
364 ftxt << " " << tofRawDatum->GetTRMchain();
365 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
366 else ftxt << " " << tofRawDatum->GetTDC();
367 ftxt << " " << tofRawDatum->GetTDCchannel();
370 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
371 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
372 dummy = detectorIndex[3];
373 detectorIndex[3] = detectorIndex[4];
374 detectorIndex[4] = dummy;
377 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
378 else ftxt << " -> " << detectorIndex[0];
379 ftxt << " " << detectorIndex[1];
380 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
381 else ftxt << " " << detectorIndex[2];
382 ftxt << " " << detectorIndex[3];
383 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
384 else ftxt << " " << detectorIndex[4];
387 for (ii=0; ii<3; ii++) position[ii] = 0.;
388 fTOFGeometry->GetPos(detectorIndex, position);
390 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
391 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
392 cylindricalPosition[2] = position[2];
394 parTOF[0] = tofRawDatum->GetTOF(); //TDC
395 parTOF[1] = tofRawDatum->GetTOT(); // TOT
396 parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
397 parTOF[3] = -1;//raw data: no track of undecalib sim time
398 parTOF[4] = tofRawDatum->GetTOF(); // RAW time
399 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex, parTOF);
400 InsertCluster(tofCluster);
403 if (parTOF[1]<10)ftxt << " " << parTOF[1];
404 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
405 else ftxt << " " << parTOF[1];
406 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
407 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
408 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
409 else ftxt << " " << parTOF[3] << endl;
412 } // closed loop on TOF raw data per current DDL file
414 clonesRawData->Clear();
416 } // closed loop on DDL index
418 if (fVerbose==2) ftxt.close();
420 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
425 clustersTree->Fill();
429 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
430 stopwatch.RealTime(),stopwatch.CpuTime()));
433 //______________________________________________________________________________
435 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
438 // Converts RAW data to recpoints for TOF
441 TStopwatch stopwatch;
444 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
445 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
447 fRunLoader->GetEvent(iEvent);
449 AliDebug(2,Form(" Event number %2i ", iEvent));
451 fTreeR = fTOFLoader->TreeR();
454 fTOFLoader->MakeTree("R");
455 fTreeR = fTOFLoader->TreeR();
458 Int_t bufsize = 32000;
459 fTreeR->Branch("TOF", &fRecPoints, bufsize);
461 TClonesArray * clonesRawData;
466 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
468 Double_t cylindricalPosition[5];
471 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
474 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
477 AliTOFRawStream tofInput(rawReader);
478 if (fDecoderVersion) {
479 AliInfo("Using New Decoder \n");
480 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
482 else tofInput.LoadRawData(indexDDL);
484 clonesRawData = (TClonesArray*)tofInput.GetRawData();
486 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
488 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
490 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
493 if (indexDDL<10) ftxt << " " << indexDDL;
494 else ftxt << " " << indexDDL;
495 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
496 else ftxt << " " << tofRawDatum->GetTRM();
497 ftxt << " " << tofRawDatum->GetTRMchain();
498 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
499 else ftxt << " " << tofRawDatum->GetTDC();
500 ftxt << " " << tofRawDatum->GetTDCchannel();
503 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
504 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
505 dummy = detectorIndex[3];
506 detectorIndex[3] = detectorIndex[4];
507 detectorIndex[4] = dummy;
510 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
511 else ftxt << " -> " << detectorIndex[0];
512 ftxt << " " << detectorIndex[1];
513 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
514 else ftxt << " " << detectorIndex[2];
515 ftxt << " " << detectorIndex[3];
516 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
517 else ftxt << " " << detectorIndex[4];
520 for (ii=0; ii<3; ii++) position[ii] = 0.;
521 fTOFGeometry->GetPos(detectorIndex, position);
523 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
524 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
525 cylindricalPosition[2] = position[2];
526 parTOF[0] = tofRawDatum->GetTOF(); // TDC
527 parTOF[1] = tofRawDatum->GetTOT(); // TOT
528 parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
529 parTOF[3] = -1; //raw data: no track of the undecalib sim time
530 parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
531 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex, parTOF);
532 InsertCluster(tofCluster);
535 if (parTOF[1]<10)ftxt << " " << parTOF[1];
536 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
537 else ftxt << " " << parTOF[1];
538 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
539 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
540 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
541 else ftxt << " " << parTOF[3] << endl;
544 } // closed loop on TOF raw data per current DDL file
546 clonesRawData->Clear();
548 } // closed loop on DDL index
550 if (fVerbose==2) ftxt.close();
552 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
560 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
561 fTOFLoader->WriteRecPoints("OVERWRITE");
563 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
564 stopwatch.RealTime(),stopwatch.CpuTime()));
567 //______________________________________________________________________________
569 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
572 // Converts RAW data to MC digits for TOF
574 // (temporary solution)
577 TStopwatch stopwatch;
580 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
581 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
583 fRunLoader->GetEvent(iEvent);
585 fTreeD = fTOFLoader->TreeD();
588 AliInfo("TreeD re-creation");
590 fTOFLoader->MakeTree("D");
591 fTreeD = fTOFLoader->TreeD();
594 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
595 Int_t bufsize = 32000;
596 fTreeD->Branch("TOF", &tofDigits, bufsize);
598 fRunLoader->GetEvent(iEvent);
600 AliDebug(2,Form(" Event number %2i ", iEvent));
602 TClonesArray * clonesRawData;
606 Int_t detectorIndex[5];
610 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
613 AliTOFRawStream tofInput(rawReader);
614 if (fDecoderVersion) {
615 AliInfo("Using New Decoder \n");
616 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
618 else tofInput.LoadRawData(indexDDL);
620 clonesRawData = (TClonesArray*)tofInput.GetRawData();
622 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
624 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
626 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
628 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
629 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
630 dummy = detectorIndex[3];
631 detectorIndex[3] = detectorIndex[4];
632 detectorIndex[4] = dummy;
634 digit[0] = tofInput.GetTofBin();
635 digit[1] = tofInput.GetToTbin();
636 digit[2] = tofInput.GetToTbin();
639 Int_t tracknum[3]={-1,-1,-1};
641 TClonesArray &aDigits = *tofDigits;
642 Int_t last=tofDigits->GetEntriesFast();
643 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
647 clonesRawData->Clear();
653 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
654 fTOFLoader->WriteDigits("OVERWRITE");
656 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
657 stopwatch.RealTime(),stopwatch.CpuTime()));
660 //______________________________________________________________________________
662 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
663 //---------------------------------------------------------------------------//
664 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
665 //---------------------------------------------------------------------------//
666 if (fNumberOfTofClusters==kTofMaxCluster) {
667 AliError("Too many clusters !");
671 if (fNumberOfTofClusters==0) {
672 fTofClusters[fNumberOfTofClusters++] = tofCluster;
676 Int_t ii = FindClusterIndex(tofCluster->GetZ());
677 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
678 fTofClusters[ii] = tofCluster;
679 fNumberOfTofClusters++;
684 //_________________________________________________________________________
686 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
687 //--------------------------------------------------------------------
688 // This function returns the index of the nearest cluster
689 //--------------------------------------------------------------------
690 if (fNumberOfTofClusters==0) return 0;
691 if (z <= fTofClusters[0]->GetZ()) return 0;
692 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
693 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
694 for (; b<e; m=(b+e)/2) {
695 if (z > fTofClusters[m]->GetZ()) b=m+1;
702 //_________________________________________________________________________
704 void AliTOFClusterFinder::FillRecPoint()
707 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
708 // in Z) in the global TClonesArray of AliTOFcluster,
714 Int_t detectorIndex[5];
715 Double_t cylindricalPosition[3];
717 Int_t trackLabels[3];
718 Int_t digitIndex = -1;
721 TClonesArray &lRecPoints = *fRecPoints;
723 for (ii=0; ii<fNumberOfTofClusters; ii++) {
725 digitIndex = fTofClusters[ii]->GetIndex();
726 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
727 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
728 cylindricalPosition[0] = fTofClusters[ii]->GetR();
729 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
730 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
731 parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
732 parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
733 parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
734 parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
735 parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
736 status=fTofClusters[ii]->GetStatus();
737 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, detectorIndex, parTOF,status,trackLabels,digitIndex);
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;
759 AliInfo(" Calibrating TOF Clusters: ")
760 AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
761 if(!calib->ReadParFromCDB("TOF/Calib",-1)) {AliFatal("Exiting, no CDB object found!!!");exit(0);}
763 AliTOFCal *calTOFArray = calib->GetTOFCalArray();
765 for (ii=0; ii<fNumberOfTofClusters; ii++) {
766 digitIndex = fTofClusters[ii]->GetIndex();
767 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
769 Int_t index = calib->GetIndex(detectorIndex);
771 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
773 // Get channel status
774 Bool_t status=calChannel->GetStatus();
775 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
776 // Get Rough channel online equalization
777 Float_t roughDelay=calChannel->GetDelay();
778 AliDebug(2,Form(" channel delay = %f", roughDelay));
779 // Get Refined channel offline calibration parameters
781 for (Int_t j = 0; j<6; j++){
782 par[j]=calChannel->GetSlewPar(j);
783 AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
785 AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %i , %i", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
786 tToT = (Double_t)(fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()); tToT*=1.E-3; //ToT in ns
787 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
788 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; // the time correction
789 AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
790 timeCorr=(Double_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()*1.E-3-timeCorr;//redefine the time
792 AliDebug(2,Form(" The channel time, corr (ps)= %e",timeCorr ));
793 tdcCorr=(Int_t)(timeCorr/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
794 fTofClusters[ii]->SetTDC(tdcCorr);
795 AliDebug(2,Form(" The channel time, corr (counts) = %i",fTofClusters[ii]->GetTDC()));
797 } // loop on clusters
801 //______________________________________________________________________________
803 void AliTOFClusterFinder::ResetRecpoint()
806 // Clear the list of reconstructed points
809 fNumberOfTofClusters = 0;
810 if (fRecPoints) fRecPoints->Clear();
813 //______________________________________________________________________________
815 void AliTOFClusterFinder::Load()
818 // Load TOF.Digits.root and TOF.RecPoints.root files
821 fTOFLoader->LoadDigits("READ");
822 fTOFLoader->LoadRecPoints("recreate");
825 //______________________________________________________________________________
827 void AliTOFClusterFinder::LoadClusters()
830 // Load TOF.RecPoints.root file
833 fTOFLoader->LoadRecPoints("recreate");
836 //______________________________________________________________________________
838 void AliTOFClusterFinder::UnLoad()
841 // Unload TOF.Digits.root and TOF.RecPoints.root files
844 fTOFLoader->UnloadDigits();
845 fTOFLoader->UnloadRecPoints();
848 //______________________________________________________________________________
850 void AliTOFClusterFinder::UnLoadClusters()
853 // Unload TOF.RecPoints.root file
856 fTOFLoader->UnloadRecPoints();
859 //______________________________________________________________________________