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