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