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