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