]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFClusterFinder.cxx
Change TOF digit Time, Tot etc to int type
[u/mrichter/AliRoot.git] / TOF / AliTOFClusterFinder.cxx
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 $Log$
18 Revision 1.24  2007/04/27 11:19:31  arcelli
19 updates for the new decoder
20
21 Revision 1.23  2007/04/23 16:51:39  decaro
22 Digits-to-raw_data conversion: correction for a more real description (A.De Caro, R.Preghenella)
23
24 Revision 1.22  2007/04/19 17:26:32  arcelli
25 Fix a bug (add some debug printout
26
27 Revision 1.21  2007/04/18 17:28:12  arcelli
28 Set the ToT bin width to the one actually used...
29
30 Revision 1.20  2007/03/09 09:57:23  arcelli
31  Remove a forgotten include of Riostrem
32
33 Revision 1.19  2007/03/08 15:41:20  arcelli
34 set uncorrected times when filling RecPoints
35
36 Revision 1.18  2007/03/06 16:31:20  arcelli
37 Add Uncorrected TOF Time signal
38
39 Revision 1.17  2007/02/28 18:09:11  arcelli
40 Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
41
42 Revision 1.16  2007/02/20 15:57:00  decaro
43 Raw data update: to read the TOF raw data defined in UNPACKED mode
44
45
46 Revision 0.03  2005/07/28 A. De Caro
47          Implement public method
48          Raw2Digits(Int_t, AliRawReader *)
49          to convert digits from raw data in MC digits
50          (temporary solution)
51
52 Revision 0.02  2005/07/27 A. De Caro
53          Implement public method
54          Digits2RecPoint(Int_t)
55          to convert digits in clusters
56
57 Revision 0.02  2005/07/26 A. De Caro
58          Implement private methods
59          InsertCluster(AliTOFcluster *)
60          FindClusterIndex(Double_t)
61          originally implemented in AliTOFtracker
62          by S. Arcelli and C. Zampolli
63
64 Revision 0.01  2005/07/25 A. De Caro
65          Implement public methods
66          Digits2RecPoint(AliRawReader *, TTree *)
67          Digits2RecPoint(Int_t, AliRawReader *)
68          to convert raw data in clusters
69  */
70
71 ////////////////////////////////////////////////////////////////
72 //                                                            //
73 //         Class for TOF cluster finder                       //
74 //                                                            //
75 // Starting from Raw Data, create rec points,                 //
76 //                         fill TreeR for TOF,                //
77 //                         write TOF.RecPoints.root file      //
78 //                                                            //
79 ////////////////////////////////////////////////////////////////
80
81
82 #include "TClonesArray.h"
83 //#include "TFile.h"
84 #include "TStopwatch.h"
85 #include "TTree.h"
86
87 #include "AliDAQ.h"
88 #include "AliLoader.h"
89 #include "AliLog.h"
90 #include "AliRawReader.h"
91 #include "AliRunLoader.h"
92
93 #include "AliTOFCal.h"
94 #include "AliTOFcalib.h"
95 #include "AliTOFChannel.h"
96 #include "AliTOFClusterFinder.h"
97 #include "AliTOFcluster.h"
98 #include "AliTOFdigit.h"
99 #include "AliTOFGeometry.h"
100 #include "AliTOFGeometryV5.h"
101 #include "AliTOFrawData.h"
102 #include "AliTOFRawStream.h"
103 #include "Riostream.h"
104
105 //extern TFile *gFile;
106
107 ClassImp(AliTOFClusterFinder)
108
109 AliTOFClusterFinder::AliTOFClusterFinder():
110   fRunLoader(0),
111   fTOFLoader(0),
112   fTreeD(0),
113   fTreeR(0),
114   fTOFGeometry(new AliTOFGeometryV5()),
115   fDigits(new TClonesArray("AliTOFdigit", 4000)),
116   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
117   fNumberOfTofClusters(0),
118   fVerbose(0),
119   fDecoderVersion(0)
120 {
121 //
122 // Constructor
123 //
124
125   AliInfo("V5 TOF Geometry is taken as the default");
126
127 }
128 //______________________________________________________________________________
129
130 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
131   fRunLoader(runLoader),
132   fTOFLoader(runLoader->GetLoader("TOFLoader")),
133   fTreeD(0),
134   fTreeR(0),
135   fTOFGeometry(new AliTOFGeometryV5()),
136   fDigits(new TClonesArray("AliTOFdigit", 4000)),
137   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
138   fNumberOfTofClusters(0),
139   fVerbose(0),
140   fDecoderVersion(0)
141 {
142 //
143 // Constructor
144 //
145
146 //  runLoader->CdGAFile();
147 //  TFile *in=(TFile*)gFile;
148 //  in->cd();
149 //  fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
150
151 }
152
153 //------------------------------------------------------------------------
154 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
155   :TObject(),
156   fRunLoader(0),
157   fTOFLoader(0),
158   fTreeD(0),
159   fTreeR(0),
160   fTOFGeometry(new AliTOFGeometryV5()),
161   fDigits(new TClonesArray("AliTOFdigit", 4000)),
162   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
163   fNumberOfTofClusters(0),
164    fVerbose(0),
165   fDecoderVersion(0)
166 {
167   // copy constructor
168   this->fDigits=source.fDigits;
169   this->fRecPoints=source.fRecPoints;
170   this->fTOFGeometry=source.fTOFGeometry;
171   this->fDecoderVersion=source.fDecoderVersion;
172
173 }
174
175 //------------------------------------------------------------------------
176 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
177 {
178   // ass. op.
179   this->fDigits=source.fDigits;
180   this->fRecPoints=source.fRecPoints;
181   this->fTOFGeometry=source.fTOFGeometry;
182   this->fVerbose=source.fVerbose;
183   this->fDecoderVersion=source.fDecoderVersion;
184   return *this;
185
186 }
187 //______________________________________________________________________________
188
189 AliTOFClusterFinder::~AliTOFClusterFinder()
190 {
191
192   //
193   // Destructor
194   //
195
196   if (fDigits)
197     {
198       fDigits->Delete();
199       delete fDigits;
200       fDigits=0;
201     }
202   if (fRecPoints)
203     {
204       fRecPoints->Delete();
205       delete fRecPoints;
206       fRecPoints=0;
207     }
208
209   delete fTOFGeometry;
210
211 }
212 //______________________________________________________________________________
213
214 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
215 {
216   //
217   // Converts digits to recpoints for TOF
218   //
219
220   TStopwatch stopwatch;
221   stopwatch.Start();
222
223   fRunLoader->GetEvent(iEvent);
224
225   fTreeD = fTOFLoader->TreeD();
226   if (fTreeD == 0x0)
227     {
228       AliFatal("AliTOFClusterFinder: Can not get TreeD");
229     }
230
231   TBranch *branch = fTreeD->GetBranch("TOF");
232   if (!branch) { 
233     AliError("can't get the branch with the TOF digits !");
234     return;
235   }
236
237   TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
238   branch->SetAddress(&digits);
239
240   ResetRecpoint();
241
242   fTreeR = fTOFLoader->TreeR();
243   if (fTreeR == 0x0)
244     {
245       fTOFLoader->MakeTree("R");
246       fTreeR = fTOFLoader->TreeR();
247     }
248
249   Int_t bufsize = 32000;
250   fTreeR->Branch("TOF", &fRecPoints, bufsize);
251
252   fTreeD->GetEvent(0);
253   Int_t nDigits = digits->GetEntriesFast();
254   AliDebug(2,Form("Number of TOF digits: %d",nDigits));
255
256   Int_t ii, jj;
257   Int_t dig[5]; //cluster detector indeces
258   Float_t g[3]; //cluster cartesian coord
259   Double_t h[3]; // the cluster spatial cyl. coordinates
260   Int_t  parTOF[5]; //The TOF signal parameters
261   Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
262   for (ii=0; ii<nDigits; ii++) {
263     AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
264     dig[0]=d->GetSector();
265     dig[1]=d->GetPlate();
266     dig[2]=d->GetStrip();
267     dig[3]=d->GetPadz();
268     dig[4]=d->GetPadx();
269
270     //    AliDebug(2,Form(" %2i  %1i  %2i  %1i  %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
271
272     for (jj=0; jj<3; jj++) g[jj] = 0.;
273     fTOFGeometry->GetPos(dig,g);
274
275     h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
276     h[1] = TMath::ATan2(g[1],g[0]);
277     h[2] = g[2];
278
279     parTOF[0] = d->GetTdc(); //the TDC signal
280     parTOF[1] = d->GetToT(); //the ToT signal
281     parTOF[2] = d->GetAdc(); // the adc charge
282     parTOF[3] = d->GetTdcND(); // non decalibrated sim time
283     parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
284     AliTOFcluster *tofCluster = new AliTOFcluster(h,dig,parTOF,status,d->GetTracks(),ii);
285     InsertCluster(tofCluster);
286
287   }
288
289   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
290
291   CalibrateRecPoint();
292   FillRecPoint();
293
294   fTreeR->Fill();
295   ResetRecpoint();
296
297   fTOFLoader = fRunLoader->GetLoader("TOFLoader");  
298   fTOFLoader->WriteRecPoints("OVERWRITE");
299
300   AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
301                stopwatch.RealTime(),stopwatch.CpuTime()));
302
303 }
304 //______________________________________________________________________________
305
306 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
307                                            TTree *clustersTree)
308 {
309   //
310   // Converts RAW data to recpoints for TOF
311   //
312
313   TStopwatch stopwatch;
314   stopwatch.Start();
315
316   //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
317   const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
318
319   ResetRecpoint();
320
321   Int_t bufsize = 32000;
322   clustersTree->Branch("TOF", &fRecPoints, bufsize);
323
324   TClonesArray * clonesRawData;
325
326   Int_t ii = 0;
327   Int_t dummy = -1;
328
329   Int_t detectorIndex[5];
330   Float_t position[3];
331   Double_t cylindricalPosition[3];
332   Int_t parTOF[5];
333
334   ofstream ftxt;
335   if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
336
337   Int_t indexDDL = 0;
338   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
339
340     rawReader->Reset();
341     AliTOFRawStream tofInput(rawReader);
342     if (fDecoderVersion) {
343       AliInfo("Using New Decoder \n"); 
344       tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
345     }
346     else tofInput.LoadRawData(indexDDL);
347
348     clonesRawData = (TClonesArray*)tofInput.GetRawData();
349
350     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
351
352       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
353
354       if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
355
356       if (fVerbose==2) {
357         if (indexDDL<10) ftxt << "  " << indexDDL;
358         else         ftxt << " " << indexDDL;
359         if (tofRawDatum->GetTRM()<10) ftxt << "  " << tofRawDatum->GetTRM();
360         else         ftxt << " " << tofRawDatum->GetTRM();
361         ftxt << "  " << tofRawDatum->GetTRMchain();
362         if (tofRawDatum->GetTDC()<10) ftxt << "  " << tofRawDatum->GetTDC();
363         else         ftxt << " " << tofRawDatum->GetTDC();
364         ftxt << "  " << tofRawDatum->GetTDCchannel();
365       }
366
367       tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
368                                     tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
369       dummy = detectorIndex[3];
370       detectorIndex[3] = detectorIndex[4];
371       detectorIndex[4] = dummy;
372
373       if (fVerbose==2) {
374         if (detectorIndex[0]<10) ftxt  << "  ->  " << detectorIndex[0];
375         else              ftxt  << "  -> " << detectorIndex[0];
376         ftxt << "  " << detectorIndex[1];
377         if (detectorIndex[2]<10) ftxt << "  " << detectorIndex[2];
378         else              ftxt << " " << detectorIndex[2];
379         ftxt << "  " << detectorIndex[3];
380         if (detectorIndex[4]<10) ftxt << "  " << detectorIndex[4];
381         else              ftxt << " " << detectorIndex[4];
382       }
383
384       for (ii=0; ii<3; ii++) position[ii] =  0.;
385       fTOFGeometry->GetPos(detectorIndex, position);
386
387       cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
388       cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
389       cylindricalPosition[2] = position[2];
390
391       parTOF[0] = tofRawDatum->GetTOF(); //TDC
392       parTOF[1] = tofRawDatum->GetTOT(); // TOT
393       parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
394       parTOF[3] = -1;//raw data: no track of undecalib sim time
395       parTOF[4] = tofRawDatum->GetTOF(); // RAW time
396       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex, parTOF);
397       InsertCluster(tofCluster);
398
399       if (fVerbose==2) {
400         if (parTOF[1]<10)ftxt << "        " << parTOF[1];
401         else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << "      " << parTOF[1];
402         else ftxt << "      " << parTOF[1];
403         if (parTOF[0]<10) ftxt << "      " << parTOF[0] << endl;
404         else if (parTOF[0]>=10 && parTOF[0]<100)   ftxt << "    " << parTOF[0] << endl;
405         else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << "    " << parTOF[0] << endl;
406         else ftxt << "   " << parTOF[3] << endl;
407       }
408
409     } // closed loop on TOF raw data per current DDL file
410
411     clonesRawData->Clear();
412
413   } // closed loop on DDL index
414
415   if (fVerbose==2) ftxt.close();
416
417   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
418
419   CalibrateRecPoint();
420   FillRecPoint();
421
422   clustersTree->Fill();
423
424   ResetRecpoint();
425
426   AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
427                    stopwatch.RealTime(),stopwatch.CpuTime()));
428
429 }
430 //______________________________________________________________________________
431
432 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
433 {
434   //
435   // Converts RAW data to recpoints for TOF
436   //
437
438   TStopwatch stopwatch;
439   stopwatch.Start();
440
441   //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
442   const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
443
444   fRunLoader->GetEvent(iEvent);
445
446   AliDebug(2,Form(" Event number %2i ", iEvent));
447
448   fTreeR = fTOFLoader->TreeR();
449
450   if (fTreeR == 0x0){
451     fTOFLoader->MakeTree("R");
452     fTreeR = fTOFLoader->TreeR();
453   }
454
455   Int_t bufsize = 32000;
456   fTreeR->Branch("TOF", &fRecPoints, bufsize);
457
458   TClonesArray * clonesRawData;
459
460   Int_t ii = 0;
461   Int_t dummy = -1;
462
463   Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
464   Float_t position[3];
465   Double_t cylindricalPosition[5];
466   Int_t parTOF[5];
467   ofstream ftxt;
468   if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
469
470   Int_t indexDDL = 0;
471   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
472
473     rawReader->Reset();
474     AliTOFRawStream tofInput(rawReader);
475     if (fDecoderVersion) {
476       AliInfo("Using New Decoder \n"); 
477       tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
478     }
479     else tofInput.LoadRawData(indexDDL);
480
481     clonesRawData = (TClonesArray*)tofInput.GetRawData();
482
483     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
484
485       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
486
487       if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
488
489       if (fVerbose==2) {
490         if (indexDDL<10) ftxt << "  " << indexDDL;
491         else         ftxt << " " << indexDDL;
492         if (tofRawDatum->GetTRM()<10) ftxt << "  " << tofRawDatum->GetTRM();
493         else         ftxt << " " << tofRawDatum->GetTRM();
494         ftxt << "  " << tofRawDatum->GetTRMchain();
495         if (tofRawDatum->GetTDC()<10) ftxt << "  " << tofRawDatum->GetTDC();
496         else         ftxt << " " << tofRawDatum->GetTDC();
497         ftxt << "  " << tofRawDatum->GetTDCchannel();
498       }
499
500       tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
501                                     tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
502       dummy = detectorIndex[3];
503       detectorIndex[3] = detectorIndex[4];
504       detectorIndex[4] = dummy;
505
506       if (fVerbose==2) {
507         if (detectorIndex[0]<10) ftxt  << "  ->  " << detectorIndex[0];
508         else              ftxt  << "  -> " << detectorIndex[0];
509         ftxt << "  " << detectorIndex[1];
510         if (detectorIndex[2]<10) ftxt << "  " << detectorIndex[2];
511         else              ftxt << " " << detectorIndex[2];
512         ftxt << "  " << detectorIndex[3];
513         if (detectorIndex[4]<10) ftxt << "  " << detectorIndex[4];
514         else              ftxt << " " << detectorIndex[4];
515       }
516
517       for (ii=0; ii<3; ii++) position[ii] =  0.;
518       fTOFGeometry->GetPos(detectorIndex, position);
519
520       cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
521       cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
522       cylindricalPosition[2] = position[2];
523       parTOF[0] = tofRawDatum->GetTOF(); // TDC
524       parTOF[1] = tofRawDatum->GetTOT(); // TOT
525       parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
526       parTOF[3] = -1; //raw data: no track of the undecalib sim time
527       parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
528       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex, parTOF);
529       InsertCluster(tofCluster);
530
531       if (fVerbose==2) {
532         if (parTOF[1]<10)ftxt << "        " << parTOF[1];
533         else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << "      " << parTOF[1];
534         else ftxt << "      " << parTOF[1];
535         if (parTOF[0]<10) ftxt << "      " << parTOF[0] << endl;
536         else if (parTOF[0]>=10 && parTOF[0]<100)   ftxt << "    " << parTOF[0] << endl;
537         else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << "    " << parTOF[0] << endl;
538         else ftxt << "   " << parTOF[3] << endl;
539       }
540
541     } // closed loop on TOF raw data per current DDL file
542
543     clonesRawData->Clear();
544
545   } // closed loop on DDL index
546
547   if (fVerbose==2) ftxt.close();
548
549   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
550
551   CalibrateRecPoint();
552   FillRecPoint();
553
554   fTreeR->Fill();
555   ResetRecpoint();
556
557   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
558   fTOFLoader->WriteRecPoints("OVERWRITE");
559   
560   AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
561                stopwatch.RealTime(),stopwatch.CpuTime()));
562
563 }
564 //______________________________________________________________________________
565
566 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
567 {
568   //
569   // Converts RAW data to MC digits for TOF
570   //
571   //             (temporary solution)
572   //
573
574   TStopwatch stopwatch;
575   stopwatch.Start();
576
577   //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
578   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
579
580   fRunLoader->GetEvent(iEvent);
581
582   fTreeD = fTOFLoader->TreeD();
583   if (fTreeD)
584     {
585     AliInfo("TreeD re-creation");
586     fTreeD = 0x0;
587     fTOFLoader->MakeTree("D");
588     fTreeD = fTOFLoader->TreeD();
589     }
590
591   TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
592   Int_t bufsize = 32000;
593   fTreeD->Branch("TOF", &tofDigits, bufsize);
594
595   fRunLoader->GetEvent(iEvent);
596
597   AliDebug(2,Form(" Event number %2i ", iEvent));
598
599   TClonesArray * clonesRawData;
600
601   Int_t dummy = -1;
602
603   Int_t detectorIndex[5];
604   Int_t digit[4];
605
606   Int_t indexDDL = 0;
607   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
608
609     rawReader->Reset();
610     AliTOFRawStream tofInput(rawReader);
611     if (fDecoderVersion) {
612       AliInfo("Using New Decoder \n"); 
613       tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
614     }
615     else tofInput.LoadRawData(indexDDL);
616
617     clonesRawData = (TClonesArray*)tofInput.GetRawData();
618
619     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
620
621       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
622
623       if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
624
625       tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
626                                     tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
627       dummy = detectorIndex[3];
628       detectorIndex[3] = detectorIndex[4];
629       detectorIndex[4] = dummy;
630
631       digit[0] = tofInput.GetTofBin();
632       digit[1] = tofInput.GetToTbin();
633       digit[2] = tofInput.GetToTbin();
634       digit[3] = -1;
635
636       Int_t tracknum[3]={-1,-1,-1};
637
638       TClonesArray &aDigits = *tofDigits;
639       Int_t last=tofDigits->GetEntriesFast();
640       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
641
642     } // while loop
643
644     clonesRawData->Clear();
645
646   } // DDL Loop
647
648   fTreeD->Fill();
649
650   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
651   fTOFLoader->WriteDigits("OVERWRITE");
652   
653   AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
654                    stopwatch.RealTime(),stopwatch.CpuTime()));
655
656 }
657 //______________________________________________________________________________
658
659 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
660   //---------------------------------------------------------------------------//
661   // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
662   //---------------------------------------------------------------------------//
663   if (fNumberOfTofClusters==kTofMaxCluster) {
664     AliError("Too many clusters !");
665     return 1;
666   }
667
668   if (fNumberOfTofClusters==0) {
669     fTofClusters[fNumberOfTofClusters++] = tofCluster;
670     return 0;
671   }
672
673   Int_t ii = FindClusterIndex(tofCluster->GetZ());
674   memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
675   fTofClusters[ii] = tofCluster;
676   fNumberOfTofClusters++;
677   
678   return 0;
679
680 }
681 //_________________________________________________________________________
682
683 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
684   //--------------------------------------------------------------------
685   // This function returns the index of the nearest cluster 
686   //--------------------------------------------------------------------
687   if (fNumberOfTofClusters==0) return 0;
688   if (z <= fTofClusters[0]->GetZ()) return 0;
689   if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
690   Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
691   for (; b<e; m=(b+e)/2) {
692     if (z > fTofClusters[m]->GetZ()) b=m+1;
693     else e=m;
694   }
695
696   return m;
697
698 }
699 //_________________________________________________________________________
700
701 void AliTOFClusterFinder::FillRecPoint()
702 {
703   //
704   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
705   // in Z) in the global TClonesArray of AliTOFcluster,
706   // i.e. fRecPoints.
707   //
708
709   Int_t ii, jj;
710
711   Int_t detectorIndex[5];
712   Double_t cylindricalPosition[3];
713   Int_t parTOF[5];
714   Int_t trackLabels[3];
715   Int_t digitIndex = -1;
716   Bool_t status=kTRUE;
717
718   TClonesArray &lRecPoints = *fRecPoints;
719   
720   for (ii=0; ii<fNumberOfTofClusters; ii++) {
721
722     digitIndex = fTofClusters[ii]->GetIndex();
723     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
724     for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
725     cylindricalPosition[0] = fTofClusters[ii]->GetR();
726     cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
727     cylindricalPosition[2] = fTofClusters[ii]->GetZ();
728     parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
729     parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
730     parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
731     parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
732     parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
733     status=fTofClusters[ii]->GetStatus();
734     new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, detectorIndex, parTOF,status,trackLabels,digitIndex);
735
736   } // loop on clusters
737
738 }
739
740 //_________________________________________________________________________
741 void AliTOFClusterFinder::CalibrateRecPoint()
742 {
743   //
744   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
745   // in Z) in the global TClonesArray of AliTOFcluster,
746   // i.e. fRecPoints.
747   //
748
749   Int_t ii, jj;
750
751   Int_t detectorIndex[5];
752   Int_t digitIndex = -1;
753   Float_t tToT;
754   Float_t timeCorr;
755   Int_t   tdcCorr;
756   AliInfo(" Calibrating TOF Clusters: ")
757   AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
758   if(!calib->ReadParFromCDB("TOF/Calib",-1)) {AliFatal("Exiting, no CDB object found!!!");exit(0);}  
759   
760   AliTOFCal *calTOFArray = calib->GetTOFCalArray();  
761
762   for (ii=0; ii<fNumberOfTofClusters; ii++) {
763     digitIndex = fTofClusters[ii]->GetIndex();
764     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
765
766     Int_t index = calib->GetIndex(detectorIndex);
767      
768     AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
769
770     // Get channel status 
771     Bool_t status=calChannel->GetStatus();
772     if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
773     // Get Rough channel online equalization 
774     Float_t roughDelay=calChannel->GetDelay();
775     AliDebug(2,Form(" channel delay = %f", roughDelay));
776     // Get Refined channel offline calibration parameters
777     Float_t par[6];
778     for (Int_t j = 0; j<6; j++){
779       par[j]=calChannel->GetSlewPar(j);
780       AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
781      }
782     AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %i , %i", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
783     tToT = (Float_t)(fTofClusters[ii]->GetToT())*AliTOFGeometry::ToTBinWidth()*1.E-3; //ToT in ns
784     AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %f, %f",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
785     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
786     AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
787     timeCorr=((Float_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;//redefine the time
788     AliDebug(2,Form(" The channel time, corr (ns)= %f",timeCorr ));
789     tdcCorr=(Int_t)((timeCorr*1E3-32)/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
790     fTofClusters[ii]->SetTDC(tdcCorr);
791     AliDebug(2,Form(" The channel time, corr (counts) counts= %i",fTofClusters[ii]->GetTDC()));
792
793   } // loop on clusters
794
795   delete calib;
796 }
797 //______________________________________________________________________________
798
799 void AliTOFClusterFinder::ResetRecpoint()
800 {
801   //
802   // Clear the list of reconstructed points
803   //
804
805   fNumberOfTofClusters = 0;
806   if (fRecPoints) fRecPoints->Clear();
807
808 }
809 //______________________________________________________________________________
810
811 void AliTOFClusterFinder::Load()
812 {
813   //
814   // Load TOF.Digits.root and TOF.RecPoints.root files
815   //
816
817   fTOFLoader->LoadDigits("READ");
818   fTOFLoader->LoadRecPoints("recreate");
819
820 }
821 //______________________________________________________________________________
822
823 void AliTOFClusterFinder::LoadClusters()
824 {
825   //
826   // Load TOF.RecPoints.root file
827   //
828
829   fTOFLoader->LoadRecPoints("recreate");
830
831 }
832 //______________________________________________________________________________
833
834 void AliTOFClusterFinder::UnLoad()
835 {
836   //
837   // Unload TOF.Digits.root and TOF.RecPoints.root files
838   //
839
840   fTOFLoader->UnloadDigits();
841   fTOFLoader->UnloadRecPoints();
842
843 }
844 //______________________________________________________________________________
845
846 void AliTOFClusterFinder::UnLoadClusters()
847 {
848   //
849   // Unload TOF.RecPoints.root file
850   //
851
852   fTOFLoader->UnloadRecPoints();
853
854 }
855 //______________________________________________________________________________