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