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