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