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