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