]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFClusterFinder.cxx
Remove old raw reader class
[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: AliTOFClusterFinder.cxx,v $
18 Revision 1.31  2007/11/24 14:53:19  zampolli
19 Status flag implemented as UChar_t
20
21 Revision 1.30  2007/10/04 13:08:52  arcelli
22 updates to comply with AliTOFGeometryV5 becoming AliTOFGeometry
23
24 Revision 1.29  2007/10/03 10:42:33  arcelli
25 updates to handle new AliTOFcluster, inheriting form AliCluster3D
26
27 Revision 1.28  2007/05/31 16:06:05  arcelli
28 move instance of AliRawStream outside loop on DDL
29
30 Revision 1.27  2007/05/02 16:31:49  arcelli
31 Add methods to handle single event reconstruction.  retrieval of Calib
32 info moved to AliTOFReconstructor ctor, and passed via a pointer to
33 AliTOFcalib
34
35 Revision 1.26  2007/04/30 19:02:24  arcelli
36 hopefully the last refinements for correct type conversion in calibration
37
38 Revision 1.25  2007/04/30 15:22:17  arcelli
39 Change TOF digit Time, Tot etc to int type
40
41 Revision 1.24  2007/04/27 11:19:31  arcelli
42 updates for the new decoder
43
44 Revision 1.23  2007/04/23 16:51:39  decaro
45 Digits-to-raw_data conversion: correction for a more real description
46 (A.De Caro, R.Preghenella)
47
48 Revision 1.22  2007/04/19 17:26:32  arcelli
49 Fix a bug (add some debug printout
50
51 Revision 1.21  2007/04/18 17:28:12  arcelli
52 Set the ToT bin width to the one actually used...
53
54 Revision 1.20  2007/03/09 09:57:23  arcelli
55  Remove a forgotten include of Riostrem
56
57 Revision 1.19  2007/03/08 15:41:20  arcelli
58 set uncorrected times when filling RecPoints
59
60 Revision 1.18  2007/03/06 16:31:20  arcelli
61 Add Uncorrected TOF Time signal
62
63 Revision 1.17  2007/02/28 18:09:11  arcelli
64 Add protection against failed retrieval of the CDB cal object,
65 now Reconstruction exits with AliFatal
66
67 Revision 1.16  2007/02/20 15:57:00  decaro
68 Raw data update: to read the TOF raw data defined in UNPACKED mode
69
70
71 Revision 0.03  2005/07/28 A. De Caro
72          Implement public method
73          Raw2Digits(Int_t, AliRawReader *)
74          to convert digits from raw data in MC digits
75          (temporary solution)
76
77 Revision 0.02  2005/07/27 A. De Caro
78          Implement public method
79          Digits2RecPoint(Int_t)
80          to convert digits in clusters
81
82 Revision 0.02  2005/07/26 A. De Caro
83          Implement private methods
84          InsertCluster(AliTOFcluster *)
85          FindClusterIndex(Double_t)
86          originally implemented in AliTOFtracker
87          by S. Arcelli and C. Zampolli
88
89 Revision 0.01  2005/07/25 A. De Caro
90          Implement public methods
91          Digits2RecPoint(AliRawReader *, TTree *)
92          Digits2RecPoint(Int_t, AliRawReader *)
93          to convert raw data in clusters
94  */
95
96 ////////////////////////////////////////////////////////////////
97 //                                                            //
98 //         Class for TOF cluster finder                       //
99 //                                                            //
100 // Starting from Raw Data, create rec points,                 //
101 //                         fill TreeR for TOF,                //
102 //                         write TOF.RecPoints.root file      //
103 //                                                            //
104 ////////////////////////////////////////////////////////////////
105
106 #include "Riostream.h"
107
108 #include "TClonesArray.h"
109 #include "TStopwatch.h"
110 #include "TTree.h"
111 //#include <TGeoManager.h>
112 #include <TGeoMatrix.h>
113 //#include <TGeoPhysicalNode.h>
114
115 #include "AliDAQ.h"
116 #include "AliLoader.h"
117 #include "AliLog.h"
118 #include "AliRawReader.h"
119 #include "AliRunLoader.h"
120 //#include "AliAlignObj.h"
121 #include <AliGeomManager.h>
122
123 #include "AliTOFcalib.h"
124 #include "AliTOFChannelOnlineArray.h"
125 #include "AliTOFChannelOnlineStatusArray.h"
126 #include "AliTOFChannelOffline.h"
127 #include "AliTOFClusterFinder.h"
128 #include "AliTOFcluster.h"
129 #include "AliTOFdigit.h"
130 #include "AliTOFGeometry.h"
131 #include "AliTOFrawData.h"
132
133 #include "AliTOFDeltaBCOffset.h"
134 #include "AliTOFCTPLatency.h"
135 #include "AliTOFRunParams.h"
136
137 //extern TFile *gFile;
138
139 ClassImp(AliTOFClusterFinder)
140
141 AliTOFClusterFinder::AliTOFClusterFinder(AliTOFcalib *calib):
142   TTask("AliTOFClusterFinder",""),
143   fRunLoader(0),
144   fTOFLoader(0),
145   fTreeD(0),
146   fTreeR(0),
147   fDigits(new TClonesArray("AliTOFdigit", 4000)),
148   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
149   fNumberOfTofClusters(0),
150   fVerbose(0),
151   fDecoderVersion(0),
152   fTOFcalib(calib),
153   fTOFRawStream(AliTOFRawStream())
154 {
155 //
156 // Constructor
157 //
158
159   TString validity = (TString)fTOFcalib->GetOfflineValidity();
160   if (validity.CompareTo("valid")==0) {
161     AliInfo(Form(" validity = %s - Using offline calibration parameters",validity.Data()));
162   } else {
163     AliInfo(Form(" validity = %s - Using online calibration parameters",validity.Data()));
164   }
165
166 }
167
168 //______________________________________________________________________________
169
170 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader, AliTOFcalib *calib):
171   TTask("AliTOFClusterFinder",""),
172   fRunLoader(runLoader),
173   fTOFLoader(runLoader->GetLoader("TOFLoader")),
174   fTreeD(0),
175   fTreeR(0),
176   fDigits(new TClonesArray("AliTOFdigit", 4000)),
177   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
178   fNumberOfTofClusters(0),
179   fVerbose(0),
180   fDecoderVersion(0),
181   fTOFcalib(calib),
182   fTOFRawStream(AliTOFRawStream())
183 {
184 //
185 // Constructor
186 //
187
188   TString validity = (TString)fTOFcalib->GetOfflineValidity();
189   if (validity.CompareTo("valid")==0) {
190     AliInfo(Form(" validity = %s - Using offline calibration parameters",validity.Data()));
191   } else {
192     AliInfo(Form(" validity = %s - Using online calibration parameters",validity.Data()));
193   }
194
195 }
196
197 //------------------------------------------------------------------------
198 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source) :
199   TTask(source),
200   fRunLoader(0),
201   fTOFLoader(0),
202   fTreeD(0),
203   fTreeR(0),
204   fDigits(source.fDigits),
205   fRecPoints(source.fRecPoints),
206   fNumberOfTofClusters(0),
207   fVerbose(0),
208   fDecoderVersion(source.fDecoderVersion),
209   fTOFcalib(source.fTOFcalib),
210   fTOFRawStream(source.fTOFRawStream)
211 {
212   // copy constructor
213 }
214
215 //------------------------------------------------------------------------
216 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
217 {
218   // ass. op.
219
220   if (this == &source)
221     return *this;
222
223   TTask::operator=(source);  
224   fDigits=source.fDigits;
225   fRecPoints=source.fRecPoints;
226   fVerbose=source.fVerbose;
227   fDecoderVersion=source.fDecoderVersion;
228   fTOFcalib=source.fTOFcalib;
229   fTOFRawStream=source.fTOFRawStream;
230   return *this;
231
232 }
233 //______________________________________________________________________________
234
235 AliTOFClusterFinder::~AliTOFClusterFinder()
236 {
237
238   //
239   // Destructor
240   //
241
242   if (fDigits)
243     {
244       fDigits->Delete();
245       delete fDigits;
246       fDigits=0;
247     }
248   if (fRecPoints)
249     {
250       fRecPoints->Delete();
251       delete fRecPoints;
252       fRecPoints=0;
253     }
254
255   if (fTofClusters || fNumberOfTofClusters) {
256     for (Int_t ii=0; ii<fNumberOfTofClusters; ii++)
257       if (fTofClusters[ii]) fTofClusters[ii]->Delete();
258     fNumberOfTofClusters=0;
259    }
260
261 }
262 //______________________________________________________________________________
263
264 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
265 {
266   //
267   // Converts digits to recpoints for TOF
268   //
269
270   TStopwatch stopwatch;
271   stopwatch.Start();
272
273   Int_t inholes = 0;
274
275   fRunLoader->GetEvent(iEvent);
276
277   fTreeD = fTOFLoader->TreeD();
278   if (fTreeD == 0x0)
279     {
280       AliFatal("AliTOFClusterFinder: Can not get TreeD");
281     }
282
283   fDigits->Clear();
284   fTreeD->GetBranch("TOF")->SetAutoDelete(kFALSE);
285   fTreeD->SetBranchAddress("TOF",&fDigits);
286
287   ResetRecpoint();
288
289   fTreeR = fTOFLoader->TreeR();
290   if (fTreeR == 0x0)
291     {
292       fTOFLoader->MakeTree("R");
293       fTreeR = fTOFLoader->TreeR();
294     }
295
296   Int_t bufsize = 32000;
297   fTreeR->Branch("TOF", &fRecPoints, bufsize);
298
299   fTreeD->GetEvent(0);
300   Int_t nDigits = fDigits->GetEntriesFast();
301   AliDebug(2,Form("Number of TOF digits: %d",nDigits));
302
303   Int_t ii;
304   Int_t dig[5]={-1,-1,-1,-1,-1}; //cluster detector indeces
305   Int_t parTOF[7]={0,0,0,0,0,0,0}; //The TOF signal parameters
306   Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
307   for (ii=0; ii<nDigits; ii++) {
308     AliTOFdigit *d = (AliTOFdigit*)fDigits->UncheckedAt(ii);
309     dig[0]=d->GetSector();
310     dig[1]=d->GetPlate();
311     dig[2]=d->GetStrip();
312     dig[3]=d->GetPadz();
313     dig[4]=d->GetPadx();
314
315     /* check valid index */
316     if (dig[0]==-1||dig[1]==-1||dig[2]==-1||dig[3]==-1||dig[4]==-1) continue;
317
318     // Do not reconstruct anything in the holes
319     if (dig[0]==13 || dig[0]==14 || dig[0]==15 ) { // sectors with holes
320       if (dig[1]==2) { // plate with holes
321         inholes++;
322         continue;
323       }
324     }
325
326     AliDebug(2,Form(" %2d  %1d  %2d  %1d  %2d ",dig[0],dig[1],dig[2],dig[3],dig[4]));
327
328     parTOF[0] = d->GetTdc(); //the TDC signal
329     parTOF[1] = d->GetToT(); //the ToT signal
330     parTOF[2] = d->GetAdc(); // the adc charge
331     parTOF[3] = d->GetTdcND(); // non decalibrated sim time
332     parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
333     parTOF[5] = 0; // deltaBC
334     parTOF[6] = 0; // L0-L1 latency
335     Double_t posClus[3];
336     Double_t covClus[6];
337     UShort_t volIdClus=GetClusterVolIndex(dig);
338     GetClusterPars(dig, posClus,covClus);
339     AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],d->GetTracks(),dig,parTOF,status,ii);
340     InsertCluster(tofCluster);
341
342   }
343
344   AliDebug(1,Form("Number of found clusters: %d for event: %d", fNumberOfTofClusters, iEvent));
345
346   CalibrateRecPoint();
347   FillRecPoint();
348
349   fTreeR->Fill();
350   ResetRecpoint();
351
352   fTOFLoader = fRunLoader->GetLoader("TOFLoader");  
353   fTOFLoader->WriteRecPoints("OVERWRITE");
354
355   AliDebug(1,Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
356                   stopwatch.RealTime(),stopwatch.CpuTime()));
357   if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
358
359 }
360
361 //______________________________________________________________________________
362
363 void AliTOFClusterFinder::Digits2RecPoints(TTree* digitsTree, TTree* clusterTree)
364 {
365   //
366   // Converts digits to recpoints for TOF
367   //
368
369   TStopwatch stopwatch;
370   stopwatch.Start();
371
372   Int_t inholes = 0;
373
374   if (digitsTree == 0x0)
375     {
376       AliFatal("AliTOFClusterFinder: Can not get TreeD");
377     }
378
379   fDigits->Clear();
380   digitsTree->GetBranch("TOF")->SetAutoDelete(kFALSE);
381   digitsTree->SetBranchAddress("TOF",&fDigits);
382
383   ResetRecpoint();
384   Int_t bufsize = 32000;
385   clusterTree->Branch("TOF", &fRecPoints, bufsize);
386
387   digitsTree->GetEvent(0);
388   Int_t nDigits = fDigits->GetEntriesFast();
389   AliDebug(2,Form("Number of TOF digits: %d",nDigits));
390
391   Int_t ii;
392   Int_t dig[5]={-1,-1,-1,-1,-1}; //cluster detector indeces
393   Int_t parTOF[7]={0,0,0,0,0,0,0}; //The TOF signal parameters
394   Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
395   for (ii=0; ii<nDigits; ii++) {
396     AliTOFdigit *d = (AliTOFdigit*)fDigits->UncheckedAt(ii);
397     dig[0]=d->GetSector();
398     dig[1]=d->GetPlate();
399     dig[2]=d->GetStrip();
400     dig[3]=d->GetPadz();
401     dig[4]=d->GetPadx();
402
403     /* check valid index */
404     if (dig[0]==-1||dig[1]==-1||dig[2]==-1||dig[3]==-1||dig[4]==-1) continue;
405
406     // Do not reconstruct anything in the holes
407     if (dig[0]==13 || dig[0]==14 || dig[0]==15 ) { // sectors with holes
408       if (dig[1]==2) { // plate with holes
409         inholes++;
410         continue;
411       }
412     }
413
414     //    AliDebug(2,Form(" %2d  %1d  %2d  %1d  %2d ",dig[0],dig[1],dig[2],dig[3],dig[4]));
415
416     parTOF[0] = d->GetTdc(); //the TDC signal
417     parTOF[1] = d->GetToT(); //the ToT signal
418     parTOF[2] = d->GetAdc(); // the adc charge
419     parTOF[3] = d->GetTdcND(); // non decalibrated sim time
420     parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
421     parTOF[5] = 0; // deltaBC
422     parTOF[6] = 0; // L0-L1 latency
423     
424     Double_t posClus[3];
425     Double_t covClus[6];
426     UShort_t volIdClus=GetClusterVolIndex(dig);
427     GetClusterPars(dig,posClus,covClus);
428     AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],d->GetTracks(),dig,parTOF,status,ii);
429     InsertCluster(tofCluster);
430
431   }
432
433   AliDebug(1,Form("Number of found clusters: %d", fNumberOfTofClusters));
434
435   CalibrateRecPoint();
436   FillRecPoint();
437
438   clusterTree->Fill();
439   ResetRecpoint();
440
441   AliDebug(1,Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
442                   stopwatch.RealTime(),stopwatch.CpuTime()));
443   if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
444
445 }
446 //______________________________________________________________________________
447
448 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
449                                            TTree *clustersTree)
450 {
451   //
452   // Converts RAW data to recpoints for TOF
453   //
454
455   TStopwatch stopwatch;
456   stopwatch.Start();
457
458   Int_t inholes = 0;
459
460   const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
461
462   ResetRecpoint();
463
464   Int_t bufsize = 32000;
465   clustersTree->Branch("TOF", &fRecPoints, bufsize);
466
467   TClonesArray * clonesRawData;
468   Int_t dummy = -1;
469
470   Int_t detectorIndex[5];
471   Int_t parTOF[7];
472
473   ofstream ftxt;
474   if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
475
476   fTOFRawStream.Clear();
477   fTOFRawStream.SetRawReader(rawReader);
478
479   if (fDecoderVersion == 1) {
480     AliInfo("Using New Decoder");
481   }
482   else if (fDecoderVersion == 2) {
483     AliInfo("Using Enhanced Decoder");
484   }
485   else {
486     AliInfo("Using Old Decoder");
487   }
488
489   Int_t indexDDL = 0;
490   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
491
492     rawReader->Reset();
493     if (fDecoderVersion == 1) {
494       fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
495     }
496     else if (fDecoderVersion == 2) {
497       fTOFRawStream.LoadRawDataBuffersV2(indexDDL,fVerbose);
498     }
499     else  {
500       fTOFRawStream.LoadRawData(indexDDL);
501     }
502     
503     clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
504
505     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
506
507       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
508
509       //if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
510       if (tofRawDatum->GetTOF()==-1) continue;
511
512       if (fVerbose==2) {
513         if (indexDDL<10) ftxt << "  " << indexDDL;
514         else         ftxt << " " << indexDDL;
515         if (tofRawDatum->GetTRM()<10) ftxt << "  " << tofRawDatum->GetTRM();
516         else         ftxt << " " << tofRawDatum->GetTRM();
517         ftxt << "  " << tofRawDatum->GetTRMchain();
518         if (tofRawDatum->GetTDC()<10) ftxt << "  " << tofRawDatum->GetTDC();
519         else         ftxt << " " << tofRawDatum->GetTDC();
520         ftxt << "  " << tofRawDatum->GetTDCchannel();
521       }
522
523       fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
524                                          tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
525       dummy = detectorIndex[3];
526       detectorIndex[3] = detectorIndex[4];
527       detectorIndex[4] = dummy;
528
529       if (fVerbose==2) {
530         if (detectorIndex[0]<10) ftxt  << "  ->  " << detectorIndex[0];
531         else              ftxt  << "  -> " << detectorIndex[0];
532         ftxt << "  " << detectorIndex[1];
533         if (detectorIndex[2]<10) ftxt << "  " << detectorIndex[2];
534         else              ftxt << " " << detectorIndex[2];
535         ftxt << "  " << detectorIndex[3];
536         if (detectorIndex[4]<10) ftxt << "  " << detectorIndex[4];
537         else              ftxt << " " << detectorIndex[4];
538       }
539
540     /* check valid index */
541     if (detectorIndex[0]==-1||detectorIndex[1]==-1||detectorIndex[2]==-1||detectorIndex[3]==-1||detectorIndex[4]==-1) continue;
542
543       // Do not reconstruct anything in the holes
544       if (detectorIndex[0]==13 || detectorIndex[0]==14 || detectorIndex[0]==15 ) { // sectors with holes
545         if (detectorIndex[1]==2) { // plate with holes
546           inholes++;
547           continue;
548         }
549       }
550
551       parTOF[0] = tofRawDatum->GetTOF(); //TDC
552       parTOF[1] = tofRawDatum->GetTOT(); // TOT
553       parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
554       parTOF[3] = 0;//raw data: no track of undecalib sim time
555       parTOF[4] = tofRawDatum->GetTOF(); // RAW time
556       parTOF[5] = tofRawDatum->GetDeltaBC(); // deltaBC
557       parTOF[6] = tofRawDatum->GetL0L1Latency(); // L0-L1 latency
558       Double_t posClus[3];
559       Double_t covClus[6];
560       UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
561       Int_t lab[3]={-1,-1,-1};
562       Bool_t status=kTRUE;
563       GetClusterPars(detectorIndex,posClus,covClus);
564       AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],lab,detectorIndex,parTOF,status,-1);
565       InsertCluster(tofCluster);
566
567       if (fVerbose==2) {
568         if (parTOF[1]<10)ftxt << "        " << parTOF[1];
569         else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << "      " << parTOF[1];
570         else ftxt << "      " << parTOF[1];
571         if (parTOF[0]<10) ftxt << "      " << parTOF[0] << endl;
572         else if (parTOF[0]>=10 && parTOF[0]<100)   ftxt << "    " << parTOF[0] << endl;
573         else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << "    " << parTOF[0] << endl;
574         else ftxt << "   " << parTOF[3] << endl;
575       }
576
577     } // closed loop on TOF raw data per current DDL file
578
579     clonesRawData->Clear("C");
580
581   } // closed loop on DDL index
582
583   if (fVerbose==2) ftxt.close();
584
585   AliDebug(1,Form("Number of found clusters: %d", fNumberOfTofClusters));
586
587   CalibrateRecPoint(rawReader->GetTimestamp());
588   FillRecPoint();
589
590   clustersTree->Fill();
591
592   ResetRecpoint();
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   if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
597
598 }
599 //______________________________________________________________________________
600
601 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
602 {
603   //
604   // Converts RAW data to recpoints for TOF
605   //
606
607   TStopwatch stopwatch;
608   stopwatch.Start();
609
610   Int_t inholes = 0;
611
612   const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
613
614   fRunLoader->GetEvent(iEvent);
615
616   AliDebug(2,Form(" Event number %2d ", iEvent));
617
618   fTreeR = fTOFLoader->TreeR();
619
620   if (fTreeR == 0x0){
621     fTOFLoader->MakeTree("R");
622     fTreeR = fTOFLoader->TreeR();
623   }
624
625   Int_t bufsize = 32000;
626   fTreeR->Branch("TOF", &fRecPoints, bufsize);
627
628   TClonesArray * clonesRawData;
629
630   Int_t dummy = -1;
631
632   Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
633   Int_t parTOF[7];
634   ofstream ftxt;
635   if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
636
637   fTOFRawStream.Clear();
638   fTOFRawStream.SetRawReader(rawReader);
639
640   if (fDecoderVersion == 1) {
641     AliInfo("Using New Decoder");
642   }
643   else if (fDecoderVersion == 2) {
644     AliInfo("Using Enhanced Decoder");
645   }
646   else {
647     AliInfo("Using Old Decoder");
648   }
649
650   Int_t indexDDL = 0;
651   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
652
653     rawReader->Reset();
654     if (fDecoderVersion == 1) {
655       fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
656     }
657     else if (fDecoderVersion == 2) {
658       fTOFRawStream.LoadRawDataBuffersV2(indexDDL,fVerbose);
659     }
660     else {
661       fTOFRawStream.LoadRawData(indexDDL);
662     }
663
664     clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
665
666     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
667
668       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
669
670       //if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
671       if (tofRawDatum->GetTOF()==-1) continue;
672
673       if (fVerbose==2) {
674         if (indexDDL<10) ftxt << "  " << indexDDL;
675         else         ftxt << " " << indexDDL;
676         if (tofRawDatum->GetTRM()<10) ftxt << "  " << tofRawDatum->GetTRM();
677         else         ftxt << " " << tofRawDatum->GetTRM();
678         ftxt << "  " << tofRawDatum->GetTRMchain();
679         if (tofRawDatum->GetTDC()<10) ftxt << "  " << tofRawDatum->GetTDC();
680         else         ftxt << " " << tofRawDatum->GetTDC();
681         ftxt << "  " << tofRawDatum->GetTDCchannel();
682       }
683
684       fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
685                                          tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
686       dummy = detectorIndex[3];
687       detectorIndex[3] = detectorIndex[4];
688       detectorIndex[4] = dummy;
689
690       if (fVerbose==2) {
691         if (detectorIndex[0]<10) ftxt  << "  ->  " << detectorIndex[0];
692         else              ftxt  << "  -> " << detectorIndex[0];
693         ftxt << "  " << detectorIndex[1];
694         if (detectorIndex[2]<10) ftxt << "  " << detectorIndex[2];
695         else              ftxt << " " << detectorIndex[2];
696         ftxt << "  " << detectorIndex[3];
697         if (detectorIndex[4]<10) ftxt << "  " << detectorIndex[4];
698         else              ftxt << " " << detectorIndex[4];
699       }
700
701     /* check valid index */
702     if (detectorIndex[0]==-1||detectorIndex[1]==-1||detectorIndex[2]==-1||detectorIndex[3]==-1||detectorIndex[4]==-1) continue;
703
704       // Do not reconstruct anything in the holes
705       if (detectorIndex[0]==13 || detectorIndex[0]==14 || detectorIndex[0]==15 ) { // sectors with holes
706         if (detectorIndex[1]==2) { // plate with holes
707           inholes++;
708           continue;
709         }
710       }
711
712       parTOF[0] = tofRawDatum->GetTOF(); // TDC
713       parTOF[1] = tofRawDatum->GetTOT(); // TOT
714       parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
715       parTOF[3] = 0; //raw data: no track of the undecalib sim time
716       parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
717       parTOF[5] = tofRawDatum->GetDeltaBC(); // deltaBC
718       parTOF[6] = tofRawDatum->GetL0L1Latency(); // L0-L1 latency
719       Double_t posClus[3];
720       Double_t covClus[6];
721       UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
722       Int_t lab[3]={-1,-1,-1};
723       Bool_t status=kTRUE;
724       GetClusterPars(detectorIndex,posClus,covClus);
725       AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],lab,detectorIndex,parTOF,status,-1);
726       InsertCluster(tofCluster);
727
728       if (fVerbose==2) {
729         if (parTOF[1]<10)ftxt << "        " << parTOF[1];
730         else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << "      " << parTOF[1];
731         else ftxt << "      " << parTOF[1];
732         if (parTOF[0]<10) ftxt << "      " << parTOF[0] << endl;
733         else if (parTOF[0]>=10 && parTOF[0]<100)   ftxt << "    " << parTOF[0] << endl;
734         else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << "    " << parTOF[0] << endl;
735         else ftxt << "   " << parTOF[3] << endl;
736       }
737
738     } // closed loop on TOF raw data per current DDL file
739
740     clonesRawData->Clear("C");
741
742   } // closed loop on DDL index
743
744   if (fVerbose==2) ftxt.close();
745
746   AliDebug(1,Form("Number of found clusters: %d for event: %d", fNumberOfTofClusters, iEvent));
747
748   CalibrateRecPoint(rawReader->GetTimestamp());
749   FillRecPoint();
750
751   fTreeR->Fill();
752   ResetRecpoint();
753
754   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
755   fTOFLoader->WriteRecPoints("OVERWRITE");
756   
757   AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
758                stopwatch.RealTime(),stopwatch.CpuTime()));
759   if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
760
761 }
762 //______________________________________________________________________________
763
764 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
765 {
766   //
767   // Converts RAW data to MC digits for TOF
768   //
769   //             (temporary solution)
770   //
771
772   TStopwatch stopwatch;
773   stopwatch.Start();
774
775   const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
776
777   fRunLoader->GetEvent(iEvent);
778
779   fTreeD = fTOFLoader->TreeD();
780   if (fTreeD)
781     {
782     AliInfo("TreeD re-creation");
783     fTreeD = 0x0;
784     fTOFLoader->MakeTree("D");
785     fTreeD = fTOFLoader->TreeD();
786     }
787
788   Int_t bufsize = 32000;
789   fDigits->Clear();
790   fTreeD->Branch("TOF", &fDigits, bufsize);
791
792   fRunLoader->GetEvent(iEvent);
793
794   AliDebug(2,Form(" Event number %2d ", iEvent));
795
796   TClonesArray * clonesRawData;
797
798   Int_t dummy = -1;
799
800   Int_t detectorIndex[5];
801   Int_t digit[4];
802
803   fTOFRawStream.Clear();
804   fTOFRawStream.SetRawReader(rawReader);
805
806   if (fDecoderVersion == 1) {
807     AliInfo("Using New Decoder");
808   }
809   else if (fDecoderVersion == 2) {
810     AliInfo("Using Enhanced Decoder");
811   }
812   else {
813     AliInfo("Using Old Decoder");
814   }
815
816   Int_t indexDDL = 0;
817   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
818
819     rawReader->Reset();
820     if (fDecoderVersion == 1) {
821       fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
822     }
823     else if (fDecoderVersion == 2) {
824       fTOFRawStream.LoadRawDataBuffersV2(indexDDL,fVerbose);
825     }
826     else {
827       fTOFRawStream.LoadRawData(indexDDL);
828     }
829
830     clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
831
832     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
833
834       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
835
836       //if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
837       if (tofRawDatum->GetTOF()==-1) continue;
838
839       fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
840                                          tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
841       dummy = detectorIndex[3];
842       detectorIndex[3] = detectorIndex[4];
843       detectorIndex[4] = dummy;
844
845       digit[0] = fTOFRawStream.GetTofBin();
846       digit[1] = fTOFRawStream.GetToTbin();
847       digit[2] = fTOFRawStream.GetToTbin();
848       digit[3] = 0;
849
850       Int_t tracknum[3]={-1,-1,-1};
851
852       TClonesArray &aDigits = *fDigits;
853       Int_t last=fDigits->GetEntriesFast();
854       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
855
856     } // while loop
857
858     clonesRawData->Clear("C");
859
860   } // DDL Loop
861
862   fTreeD->Fill();
863
864   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
865   fTOFLoader->WriteDigits("OVERWRITE");
866
867   AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
868                    stopwatch.RealTime(),stopwatch.CpuTime()));
869
870 }
871
872 //______________________________________________________________________________
873
874 void AliTOFClusterFinder::Raw2Digits(AliRawReader *rawReader, TTree* digitsTree)
875 {
876   //
877   // Converts RAW data to MC digits for TOF for the current event
878   //
879   //
880
881   TStopwatch stopwatch;
882   stopwatch.Start();
883
884   const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
885
886   if (!digitsTree)
887     {
888     AliError("No input digits Tree");
889     return;
890     }
891
892   Int_t bufsize = 32000;
893   digitsTree->Branch("TOF", &fDigits, bufsize);
894
895   TClonesArray * clonesRawData;
896   Int_t dummy = -1;
897
898   Int_t detectorIndex[5];
899   Int_t digit[4];
900
901   fTOFRawStream.Clear();
902   fTOFRawStream.SetRawReader(rawReader);
903
904   if (fDecoderVersion == 1) {
905     AliInfo("Using New Decoder");
906   }
907   else if (fDecoderVersion == 2) {
908     AliInfo("Using Enhanced Decoder");
909   }
910   else {
911     AliInfo("Using Old Decoder");
912   }
913
914   Int_t indexDDL = 0;
915   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
916
917     rawReader->Reset();
918     if (fDecoderVersion == 1) {
919       fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
920     }
921     else if (fDecoderVersion == 2) {
922       fTOFRawStream.LoadRawDataBuffersV2(indexDDL,fVerbose);
923     }
924     else {
925       fTOFRawStream.LoadRawData(indexDDL);
926     }
927
928     clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
929
930     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
931
932       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
933
934       //if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
935       if (tofRawDatum->GetTOF()==-1) continue;
936
937       fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
938                                          tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
939       dummy = detectorIndex[3];
940       detectorIndex[3] = detectorIndex[4];
941       detectorIndex[4] = dummy;
942
943       digit[0] = fTOFRawStream.GetTofBin();
944       digit[1] = fTOFRawStream.GetToTbin();
945       digit[2] = fTOFRawStream.GetToTbin();
946       digit[3] = 0;
947
948       Int_t tracknum[3]={-1,-1,-1};
949
950       TClonesArray &aDigits = *fDigits;
951       Int_t last=fDigits->GetEntriesFast();
952       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
953
954     } // while loop
955
956     clonesRawData->Clear("C");
957
958   } // DDL Loop
959
960   digitsTree->Fill();
961
962   AliDebug(1, Form("Got %d digits: ", fDigits->GetEntries()));
963   AliDebug(1, Form("Execution time to read TOF raw data and fill TOF digit tree : R:%.2fs C:%.2fs",
964                    stopwatch.RealTime(),stopwatch.CpuTime()));
965
966 }
967 //______________________________________________________________________________
968
969 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
970   //---------------------------------------------------------------------------//
971   // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
972   //---------------------------------------------------------------------------//
973   if (fNumberOfTofClusters==kTofMaxCluster) {
974     AliError("Too many clusters !");
975     return 1;
976   }
977
978   if (fNumberOfTofClusters==0) {
979     fTofClusters[fNumberOfTofClusters++] = tofCluster;
980     return 0;
981   }
982
983   Int_t ii = FindClusterIndex(tofCluster->GetZ());
984   memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
985   fTofClusters[ii] = tofCluster;
986   fNumberOfTofClusters++;
987   
988   return 0;
989
990 }
991 //_________________________________________________________________________
992
993 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
994   //--------------------------------------------------------------------
995   // This function returns the index of the nearest cluster in z
996   //--------------------------------------------------------------------
997   if (fNumberOfTofClusters==0) return 0;
998   if (z <= fTofClusters[0]->GetZ()) return 0;
999   if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
1000   Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
1001   for (; b<e; m=(b+e)/2) {
1002     if (z > fTofClusters[m]->GetZ()) b=m+1;
1003     else e=m;
1004   }
1005
1006   return m;
1007
1008 }
1009 //_________________________________________________________________________
1010
1011 void AliTOFClusterFinder::FillRecPoint()
1012 {
1013   //
1014   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
1015   // in Z) in the global TClonesArray of AliTOFcluster,
1016   // i.e. fRecPoints.
1017   //
1018
1019   Int_t ii, jj;
1020
1021   Int_t detectorIndex[5];
1022   Int_t parTOF[7];
1023   Int_t trackLabels[3];
1024   Int_t digitIndex = -1;
1025   Bool_t status=kTRUE;
1026
1027   TClonesArray &lRecPoints = *fRecPoints;
1028   
1029   for (ii=0; ii<fNumberOfTofClusters; ii++) {
1030
1031     digitIndex = fTofClusters[ii]->GetIndex();
1032     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
1033     for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
1034     parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
1035     parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
1036     parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
1037     parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
1038     parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
1039     parTOF[5] = fTofClusters[ii]->GetDeltaBC();//deltaBC
1040     parTOF[6] = fTofClusters[ii]->GetL0L1Latency();//L0-L1 latency
1041     status=fTofClusters[ii]->GetStatus();
1042     Double_t posClus[3];
1043     Double_t covClus[6];
1044     UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
1045     GetClusterPars(detectorIndex,posClus,covClus);
1046     new(lRecPoints[ii]) AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],trackLabels,detectorIndex, parTOF,status,digitIndex);
1047
1048     AliDebug(2, Form(" %4d  %4d  %f %f %f  %f %f %f %f %f %f  %3d %3d %3d  %2d %1d %2d %1d %2d  %4d %3d %3d %4d %4d  %1d  %4d", 
1049                      ii, volIdClus, posClus[0], posClus[1], posClus[2],
1050                      fTofClusters[ii]->GetSigmaX2(),
1051                      fTofClusters[ii]->GetSigmaXY(),
1052                      fTofClusters[ii]->GetSigmaXZ(),
1053                      fTofClusters[ii]->GetSigmaY2(),
1054                      fTofClusters[ii]->GetSigmaYZ(),
1055                      fTofClusters[ii]->GetSigmaZ2(),
1056                      trackLabels[0], trackLabels[1], trackLabels[2],
1057                      detectorIndex[0], detectorIndex[1], detectorIndex[2], detectorIndex[3], detectorIndex[4],
1058                      parTOF[0], parTOF[1], parTOF[2], parTOF[3], parTOF[4],
1059                      status, digitIndex));
1060
1061   } // loop on clusters
1062
1063 }
1064
1065 //_________________________________________________________________________
1066
1067 /*
1068  * OLD CALIBRATE REC POINTS FUNCTION
1069  */
1070
1071 #if 0
1072 void AliTOFClusterFinder::CalibrateRecPoint(UInt_t timestamp)
1073 {
1074   //
1075   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
1076   // in Z) in the global TClonesArray of AliTOFcluster,
1077   // i.e. fRecPoints.
1078   //
1079
1080   Int_t ii, jj;
1081
1082   Int_t detectorIndex[5];
1083   Int_t digitIndex = -1;
1084   Double_t tToT;
1085   Double_t timeCorr;
1086   Int_t   tdcCorr;
1087   Float_t tdcLatencyWindow;
1088   AliDebug(1," Calibrating TOF Clusters");
1089
1090   AliTOFChannelOnlineArray *calDelay = fTOFcalib->GetTOFOnlineDelay();  
1091   AliTOFChannelOnlineStatusArray *calStatus = fTOFcalib->GetTOFOnlineStatus();  
1092   TObjArray *calTOFArrayOffline = fTOFcalib->GetTOFCalArrayOffline();
1093   
1094   AliTOFDeltaBCOffset *deltaBCOffsetObj = fTOFcalib->GetDeltaBCOffset();
1095   Int_t deltaBCOffset = deltaBCOffsetObj->GetDeltaBCOffset();
1096   AliTOFCTPLatency *ctpLatencyObj = fTOFcalib->GetCTPLatency();
1097   Float_t ctpLatency = ctpLatencyObj->GetCTPLatency();
1098   AliTOFRunParams *runParamsObj = fTOFcalib->GetRunParams();
1099   Float_t t0 = runParamsObj->EvalT0(timestamp);
1100
1101   TString validity = (TString)fTOFcalib->GetOfflineValidity();
1102   Int_t calibration = -1;
1103   if (validity.CompareTo("valid")==0) {
1104     //AliInfo(Form(" validity = %s - Using offline calibration parameters",validity.Data()));
1105     calibration = 1;
1106   } else {
1107     //AliInfo(Form(" validity = %s - Using online calibration parameters",validity.Data()));
1108     calibration = 0 ;
1109   }
1110
1111   for (ii=0; ii<fNumberOfTofClusters; ii++) {
1112     digitIndex = fTofClusters[ii]->GetIndex();
1113     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
1114
1115     Int_t index = AliTOFGeometry::GetIndex(detectorIndex);
1116      
1117     UChar_t statusPulser=calStatus->GetPulserStatus(index);
1118     UChar_t statusNoise=calStatus->GetNoiseStatus(index);
1119     UChar_t statusHW=calStatus->GetHWStatus(index);
1120     UChar_t status=calStatus->GetStatus(index);
1121     tdcLatencyWindow = calStatus->GetLatencyWindow(index) * 1.e3; /* ns -> ps */
1122     
1123     //check the status, also unknown is fine!!!!!!!
1124
1125     AliDebug(2, Form(" Status for channel %d = %d",index, (Int_t)status));
1126     if((statusPulser & AliTOFChannelOnlineStatusArray::kTOFPulserBad)==(AliTOFChannelOnlineStatusArray::kTOFPulserBad)||(statusNoise & AliTOFChannelOnlineStatusArray::kTOFNoiseBad)==(AliTOFChannelOnlineStatusArray::kTOFNoiseBad)||(statusHW & AliTOFChannelOnlineStatusArray::kTOFHWBad)==(AliTOFChannelOnlineStatusArray::kTOFHWBad)){
1127       AliDebug(2, Form(" Bad Status for channel %d",index));
1128       fTofClusters[ii]->SetStatus(kFALSE); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
1129     }
1130     else {
1131       AliDebug(2, Form(" Good Status for channel %d",index));
1132     }
1133     // Get Rough channel online equalization 
1134     Double_t roughDelay=(Double_t)calDelay->GetDelay(index);  // in ns
1135     AliDebug(2,Form(" channel delay (ns) = %f", roughDelay));
1136     // Get Refined channel offline calibration parameters
1137     if (calibration ==1){
1138       AliTOFChannelOffline * calChannelOffline = (AliTOFChannelOffline*)calTOFArrayOffline->At(index);
1139       Double_t par[6];
1140       for (Int_t j = 0; j<6; j++){
1141         par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
1142      } 
1143       AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
1144       AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %d , %d", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
1145       tToT = (Double_t)(fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth());
1146       tToT*=1.E-3; //ToT in ns
1147
1148       /* check TOT limits and set new TOT in case */
1149       if (tToT < AliTOFGeometry::SlewTOTMin()) tToT = AliTOFGeometry::SlewTOTMin();
1150       if (tToT > AliTOFGeometry::SlewTOTMax()) tToT = AliTOFGeometry::SlewTOTMax();
1151
1152       AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
1153       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; // the time correction (ns)
1154     }
1155     else {
1156       timeCorr = roughDelay; // correction in ns
1157     }
1158     AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()));
1159     AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
1160     timeCorr=(Double_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()*1.E-3-timeCorr;//redefine the time
1161     timeCorr*=1.E3;
1162     AliDebug(2,Form(" The channel time, corr (ps)= %e",timeCorr ));
1163
1164     /* here timeCorr should be already corrected for calibration. 
1165      * we now go into further corrections keeping in mind that timeCorr
1166      * is in ps.
1167      *
1168      * the following corrections are performed in this way:
1169      *
1170      *    time = timeRaw - deltaBC + L0L1Latency + CTPLatency - TDCLatencyWindow - T0
1171      *
1172      */
1173
1174     AliDebug(2, Form("applying further corrections (DeltaBC): DeltaBC=%d (BC bins), DeltaBCoffset=%d (BC bins)", fTofClusters[ii]->GetDeltaBC(), deltaBCOffset));
1175     AliDebug(2, Form("applying further corrections (L0L1Latency): L0L1Latency=%d (BC bins)", fTofClusters[ii]->GetL0L1Latency()));
1176     AliDebug(2, Form("applying further corrections (CTPLatency): CTPLatency=%f (ps)", ctpLatency));
1177     AliDebug(2, Form("applying further corrections (TDCLatencyWindow): TDCLatencyWindow=%f (ps)", tdcLatencyWindow));
1178     AliDebug(2, Form("applying further corrections (T0): T0=%f (ps)", t0));
1179
1180     /* deltaBC correction (inhibited for the time being) */
1181     //    timeCorr -= (fTofClusters[ii]->GetDeltaBC() - deltaBCOffset) * AliTOFGeometry::BunchCrossingBinWidth();
1182     /* L0L1-latency correction */
1183     timeCorr += fTofClusters[ii]->GetL0L1Latency() * AliTOFGeometry::BunchCrossingBinWidth();
1184     /* CTP-latency correction (from OCDB) */
1185     timeCorr += ctpLatency;
1186     /* TDC latency-window correction (from OCDB) */
1187     timeCorr -= tdcLatencyWindow;
1188     /* T0 correction (from OCDB) */
1189     timeCorr -= t0;
1190
1191     /*
1192      * end of furhter corrections
1193      */
1194
1195     tdcCorr=(Int_t)(timeCorr/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
1196     fTofClusters[ii]->SetTDC(tdcCorr);
1197   } // loop on clusters
1198
1199 }
1200 #endif
1201
1202 //_________________________________________________________________________
1203
1204 void AliTOFClusterFinder::CalibrateRecPoint(UInt_t timestamp)
1205 {
1206   //
1207   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
1208   // in Z) in the global TClonesArray of AliTOFcluster,
1209   // i.e. fRecPoints.
1210   //
1211
1212   Int_t detectorIndex[5];
1213   Double_t time, tot, corr;
1214   Int_t deltaBC, l0l1, tdcBin;
1215   for (Int_t ii = 0; ii < fNumberOfTofClusters; ii++) {
1216     for(Int_t jj = 0; jj < 5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
1217
1218     Int_t index = AliTOFGeometry::GetIndex(detectorIndex);
1219
1220     /* check channel enabled */
1221     if (!fTOFcalib->IsChannelEnabled(index)) fTofClusters[ii]->SetStatus(kFALSE);
1222     
1223     /* get cluster info */
1224     time = fTofClusters[ii]->GetTDC() * AliTOFGeometry::TdcBinWidth(); /* ps */
1225     tot = fTofClusters[ii]->GetToT() * AliTOFGeometry::ToTBinWidth() * 1.e-3; /* ns */
1226     deltaBC = fTofClusters[ii]->GetDeltaBC();
1227     l0l1 = fTofClusters[ii]->GetL0L1Latency();
1228
1229     /* get correction */
1230     corr = fTOFcalib->GetTimeCorrection(index, tot, deltaBC, l0l1, timestamp); /* ps */
1231     AliDebug(2, Form("calibrate index %d: time=%f (ps) tot=%f (ns) deltaBC=%d l0l1=%d timestamp=%d corr=%f (ps)", index, time, tot, deltaBC, l0l1, timestamp, corr));
1232
1233     /* apply time correction */
1234     time -= corr;
1235
1236     /* convert in TDC bins and set cluster */
1237     //tdcBin = (Int_t)(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
1238     tdcBin = TMath::Nint(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
1239     fTofClusters[ii]->SetTDC(tdcBin);
1240
1241   } // loop on clusters
1242
1243 }
1244
1245 //______________________________________________________________________________
1246
1247 void AliTOFClusterFinder::ResetRecpoint()
1248 {
1249   //
1250   // Clear the list of reconstructed points
1251   //
1252
1253   fNumberOfTofClusters = 0;
1254   if (fRecPoints) fRecPoints->Clear();
1255
1256 }
1257 //______________________________________________________________________________
1258
1259 void AliTOFClusterFinder::Load()
1260 {
1261   //
1262   // Load TOF.Digits.root and TOF.RecPoints.root files
1263   //
1264
1265   fTOFLoader->LoadDigits("READ");
1266   fTOFLoader->LoadRecPoints("recreate");
1267
1268 }
1269 //______________________________________________________________________________
1270
1271 void AliTOFClusterFinder::LoadClusters()
1272 {
1273   //
1274   // Load TOF.RecPoints.root file
1275   //
1276
1277   fTOFLoader->LoadRecPoints("recreate");
1278
1279 }
1280 //______________________________________________________________________________
1281
1282 void AliTOFClusterFinder::UnLoad()
1283 {
1284   //
1285   // Unload TOF.Digits.root and TOF.RecPoints.root files
1286   //
1287
1288   fTOFLoader->UnloadDigits();
1289   fTOFLoader->UnloadRecPoints();
1290
1291 }
1292 //______________________________________________________________________________
1293
1294 void AliTOFClusterFinder::UnLoadClusters()
1295 {
1296   //
1297   // Unload TOF.RecPoints.root file
1298   //
1299
1300   fTOFLoader->UnloadRecPoints();
1301
1302 }
1303 //-------------------------------------------------------------------------
1304 UShort_t AliTOFClusterFinder::GetClusterVolIndex(const Int_t * const ind) const {
1305
1306   //First of all get the volume ID to retrieve the l2t transformation...
1307   //
1308   // Detector numbering scheme
1309   Int_t nSector = 18;
1310   Int_t nPlate  = 5;
1311   Int_t nStripA = 15;
1312   Int_t nStripB = 19;
1313   Int_t nStripC = 19;
1314
1315   Int_t isector =ind[0];
1316   if (isector >= nSector)
1317     AliError(Form("Wrong sector number in TOF (%d) !",isector));
1318   Int_t iplate = ind[1];
1319   if (iplate >= nPlate)
1320     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1321   Int_t istrip = ind[2];
1322
1323   Int_t stripOffset = 0;
1324   switch (iplate) {
1325   case 0:
1326     stripOffset = 0;
1327     break;
1328   case 1:
1329     stripOffset = nStripC;
1330     break;
1331   case 2:
1332     stripOffset = nStripC+nStripB;
1333     break;
1334   case 3:
1335     stripOffset = nStripC+nStripB+nStripA;
1336     break;
1337   case 4:
1338     stripOffset = nStripC+nStripB+nStripA+nStripB;
1339     break;
1340   default:
1341     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1342     break;
1343   };
1344
1345   Int_t index= (2*(nStripC+nStripB)+nStripA)*isector +
1346                stripOffset +
1347                istrip;
1348
1349   UShort_t volIndex = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,index);
1350   return volIndex;
1351 }
1352 //
1353 //-------------------------------------------------------------------------
1354 void AliTOFClusterFinder::GetClusterPars(Int_t *ind, Double_t* pos,Double_t* cov) const {
1355
1356   //First of all get the volume ID to retrieve the l2t transformation...
1357   //
1358   UShort_t volIndex = GetClusterVolIndex(ind);
1359   //
1360   //
1361   //we now go in the system of the strip: determine the local coordinates
1362   //
1363   //
1364   // 47---------------------------------------------------0  ^ z
1365   // | | | | | | | | | | | | | | | | | | | | | | | | | | | 1 |
1366   // -----------------------------------------------------   | y going outwards
1367   // | | | | | | | | | | | | | | | | | | | | | | | | | | | 0 |  par[0]=0;
1368
1369   // -----------------------------------------------------   |
1370   // x <-----------------------------------------------------
1371
1372   /*
1373   Float_t localX = (ind[4]-23.5)*AliTOFGeometry::XPad();
1374   Float_t localY = 0;
1375   Float_t localZ = (ind[3]- 0.5)*AliTOFGeometry::ZPad();
1376   */
1377   Float_t localX = (ind[4]-AliTOFGeometry::NpadX()/2)*AliTOFGeometry::XPad()+AliTOFGeometry::XPad()/2.;
1378   Float_t localY = 0;
1379   Float_t localZ = (ind[3]-AliTOFGeometry::NpadZ()/2)*AliTOFGeometry::ZPad()+AliTOFGeometry::ZPad()/2.;
1380   //move to the tracking ref system
1381
1382   Double_t lpos[3];
1383   lpos[0] = localX;
1384   lpos[1] = localY;
1385   lpos[2] = localZ; 
1386
1387   const TGeoHMatrix *l2t= AliGeomManager::GetTracking2LocalMatrix(volIndex);
1388   // Get The position in the track ref system
1389   Double_t tpos[3];
1390   l2t->MasterToLocal(lpos,tpos);
1391   pos[0] = tpos[0];
1392   pos[1] = tpos[1];
1393   pos[2] = tpos[2];
1394
1395   //Get The cluster covariance in the track ref system
1396   Double_t lcov[9];
1397
1398   //cluster covariance in the local system:
1399   // sx2   0   0
1400   // 0     0   0
1401   // 0     0   sz2
1402
1403   lcov[0] = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
1404   lcov[1] = 0;
1405   lcov[2] = 0;
1406   lcov[3] = 0;
1407   lcov[4] = 0;
1408   lcov[5] = 0;
1409   lcov[6] = 0;
1410   lcov[7] = 0;
1411   lcov[8] = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
1412
1413   //cluster covariance in the tracking system:
1414   TGeoHMatrix m;
1415   m.SetRotation(lcov);
1416   m.Multiply(l2t);
1417   m.MultiplyLeft(&l2t->Inverse());
1418   Double_t *tcov = m.GetRotationMatrix();
1419   cov[0] = tcov[0]; cov[1] = tcov[1]; cov[2] = tcov[2];
1420   cov[3] = tcov[4]; cov[4] = tcov[5];
1421   cov[5] = tcov[8];
1422
1423   return;
1424
1425 }