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