]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOF.cxx
Added flag for MC qa
[u/mrichter/AliRoot.git] / TOF / AliTOF.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 /* $Id$ */
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  Time Of Flight                                                           //
20 //  This class contains the basic functions for the Time Of Flight           //
21 //  detector. Functions specific to one particular geometry are              //
22 //  contained in the derived classes                                         //
23 //                                                                           //
24 //  VERSIONE WITH 5 SYMMETRIC MODULES ALONG Z AXIS                           //
25 //  ============================================================             //
26 //                                                                           //
27 //  VERSION WITH HOLES FOR PHOS AND TRD IN SPACEFRAME WITH HOLES             //
28 //                                                                           //
29 //  Volume sensibile : FPAD                                                  //
30 //                                                                           //
31 //                                                                           //
32 //                                                                           //
33 //                                                                           //
34 ///////////////////////////////////////////////////////////////////////////////
35 // Begin_Html
36 /*
37 <img src="picts/AliTOFClass.gif">
38 */
39 //End_Html
40
41 #include <TClonesArray.h>
42 #include <TFile.h>
43 #include <TFolder.h>
44 #include <TROOT.h>
45 #include <TTree.h>
46 #include <TVirtualMC.h>
47 #include <TStopwatch.h>
48
49 #include "AliConst.h"
50 #include "AliLoader.h"
51 #include "AliLog.h"
52 #include "AliMC.h"
53 #include "AliRun.h"
54 #include "AliDAQ.h"
55 #include "AliRawReader.h"
56
57 //#include "AliTOFDDLRawData.h"
58 #include "AliTOFDigitizer.h"
59 #include "AliTOFdigit.h"
60 #include "AliTOFhitT0.h"
61 #include "AliTOFhit.h"
62 #include "AliTOFGeometry.h"
63 #include "AliTOFSDigitizer.h"
64 #include "AliTOFSDigit.h"
65 #include "AliTOF.h"
66 #include "AliTOFrawData.h"
67 #include "AliTOFRawStream.h"
68
69 class AliTOFcluster;
70
71 extern TROOT *gROOT;
72
73 extern AliRun *gAlice;
74
75 ClassImp(AliTOF)
76
77 //_____________________________________________________________________________
78 AliTOF::AliTOF():
79   fFGeom(0x0),
80   fSDigits(0x0),
81   fNSDigits(0),
82   fReconParticles(0x0),
83   fIdSens(-1),
84   fTZero(kFALSE),
85   fTOFHoles(kTRUE),
86   fTOFGeometry(0x0),
87   fTOFRawWriter(AliTOFDDLRawData())
88 {
89   //
90   // Default constructor
91   //
92
93   //by default all sectors switched on
94   for (Int_t ii=0; ii<18; ii++) fTOFSectors[ii]=0;
95
96   fDigits = 0;
97   fIshunt   = 0;
98   fName = "TOF";
99
100 }
101  
102 //_____________________________________________________________________________
103 AliTOF::AliTOF(const char *name, const char *title, Option_t *option)
104        : 
105   AliDetector(name,title),
106   fFGeom(0x0),
107   fSDigits(0x0),
108   fNSDigits(0),
109   fReconParticles(0x0),
110   fIdSens(-1),
111   fTZero(kFALSE),
112   fTOFHoles(kTRUE),
113   fTOFGeometry(0x0),
114   fTOFRawWriter(AliTOFDDLRawData())
115 {
116   //
117   // AliTOF standard constructor
118   // 
119   // Here are fixed some important parameters
120   //
121
122   // Initialization of hits, sdigits and digits array
123   // added option for time zero analysis
124   //skowron
125   fTOFGeometry = new AliTOFGeometry();
126
127   //by default all sectors switched on
128   for (Int_t ii=0; ii<18; ii++) fTOFSectors[ii]=0;
129
130   if (strstr(option,"tzero")){
131     fHits   = new TClonesArray("AliTOFhitT0",  1000);
132     fTZero = kTRUE;
133     //    AliWarning("tzero option requires AliTOFv4T0/AliTOFv5T0 as TOF version (check Your Config.C)");
134   }else{
135     fHits   = new TClonesArray("AliTOFhit",  1000);
136     fTZero = kFALSE;
137   }
138   if (gAlice==0) {
139      AliFatal("gAlice==0 !");
140   }
141
142   AliMC *mcApplication = (AliMC*)gAlice->GetMCApp();
143
144   if (mcApplication->GetHitLists())
145      mcApplication->AddHitList(fHits);
146   else AliError("gAlice->GetHitLists()==0");
147
148   fIshunt  = 0;
149   fSDigits = new TClonesArray("AliTOFSDigit", 1000);
150   fDigits  = new TClonesArray("AliTOFdigit",  1000);
151
152   //
153   // Digitization parameters
154   //
155   // (Transfer Functions to be inserted here)
156   //
157   //PH  SetMarkerColor(7);
158   //PH  SetMarkerStyle(2);
159   //PH  SetMarkerSize(0.4);
160
161 // Strip Parameters
162   //fGapA    =   4.; //cm  Gap beetween tilted strip in A-type plate
163   //fGapB    =   6.; //cm  Gap beetween tilted strip in B-type plate
164
165   // Physical performances
166   //fTimeRes = 100.;//ps
167   //fChrgRes = 100.;//pC
168
169 }
170
171 //____________________________________________________________________________
172 void AliTOF::SetTOFSectors(Int_t * const sectors)
173 {
174   // Setter for partial/full TOF configuration
175
176   for(Int_t isec=0;isec<18;isec++){
177     fTOFSectors[isec]=sectors[isec];
178   }
179 }
180 //____________________________________________________________________________
181 void AliTOF::GetTOFSectors(Int_t *sectors) const
182 {
183   // Getter for partial/full TOF configuration
184
185   for(Int_t isec=0;isec<18;isec++){
186     sectors[isec]=fTOFSectors[isec];
187   }
188 }
189
190 //_____________________________________________________________________________
191 void AliTOF::CreateTOFFolders()
192 {
193   // create the ALICE TFolder
194   // create the ALICE main TFolder
195   // to be done by AliRun
196
197   TFolder * alice = new TFolder();
198   alice->SetNameTitle("FPAlice", "Alice Folder") ;
199   gROOT->GetListOfBrowsables()->Add(alice) ;
200
201   TFolder * aliceF  = alice->AddFolder("folders", "Alice memory Folder") ;
202   //  make it the owner of the objects that it contains
203   aliceF->SetOwner() ;
204   // geometry folder
205   TFolder * geomF = aliceF->AddFolder("Geometry", "Geometry objects") ;
206  
207   // creates the TOF geometry  folder
208   geomF->AddFolder("TOF", "Geometry for TOF") ;
209 }
210
211 //_____________________________________________________________________________
212 AliTOF::~AliTOF()
213 {
214   // dtor:
215   // it remove also the alice folder 
216   /* PH Temporarily commented because of problems
217   TFolder * alice = (TFolder*)gROOT->GetListOfBrowsables()->FindObject("FPAlice") ;
218   delete alice;
219   alice = 0;
220   */
221   if (fHits)
222     {
223       fHits->Delete ();
224       delete fHits;
225       fHits = 0;
226     }
227   if (fDigits)
228     {
229       fDigits->Delete ();
230       delete fDigits;
231       fDigits = 0;
232     }
233   if (fSDigits)
234     {
235       fSDigits->Delete();
236       delete fSDigits;
237       fSDigits = 0;
238     }
239
240   if (fReconParticles)
241     {
242       fReconParticles->Delete ();
243       delete fReconParticles;
244       fReconParticles = 0;
245     }
246
247 }
248
249 //_____________________________________________________________________________
250 void AliTOF::AddHit(Int_t track, Int_t *vol, Float_t *hits)
251 {
252   //
253   // Add a TOF hit
254   // new with placement used
255   //
256   TClonesArray &lhits = *fHits;
257   new(lhits[fNhits++]) AliTOFhit(fIshunt, track, vol, hits);
258 }
259
260 //_____________________________________________________________________________
261 void AliTOF::AddT0Hit(Int_t track, Int_t *vol, Float_t *hits)
262 {
263   //
264   // Add a TOF hit
265   // new with placement used
266   //
267   TClonesArray &lhits = *fHits;
268   new(lhits[fNhits++]) AliTOFhitT0(fIshunt, track, vol, hits);
269 }
270
271 //_____________________________________________________________________________
272 void AliTOF::AddDigit(Int_t *tracks, Int_t *vol, Int_t *digits)
273 {
274   //
275   // Add a TOF digit
276   // new with placement used
277   //
278   TClonesArray &ldigits = *fDigits;
279   new (ldigits[fNdigits++]) AliTOFdigit(tracks, vol, digits);
280 }
281
282 //_____________________________________________________________________________
283 void AliTOF::AddSDigit(Int_t tracknum, Int_t *vol, Int_t *digits)
284 {
285      
286 //
287 // Add a TOF sdigit
288 //
289         
290   TClonesArray &lSDigits = *fSDigits;   
291   new(lSDigits[fNSDigits++]) AliTOFSDigit(tracknum, vol, digits);
292 }
293
294 //_____________________________________________________________________________
295 void AliTOF::SetTreeAddress ()
296 {
297   // Set branch address for the Hits and Digits Tree.
298   
299   if (fLoader->TreeH())
300    {
301      if (fHits == 0x0)
302       {
303         if (fTZero) fHits   = new TClonesArray("AliTOFhitT0", 1000);
304         else fHits   = new TClonesArray("AliTOFhit", 1000);
305       }
306    }
307   AliDetector::SetTreeAddress ();
308
309   TBranch *branch;
310
311   if (fLoader->TreeS () )
312     {
313       branch = fLoader->TreeS ()->GetBranch ("TOF");
314       if (branch) {
315         if (fSDigits == 0x0) fSDigits = new TClonesArray("AliTOFSDigit",  1000);
316         branch->SetAddress (&fSDigits);
317       }
318     }
319
320   if (fLoader->TreeR() ) 
321     {
322       branch = fLoader->TreeR()->GetBranch("TOF"); 
323       if (branch) 
324        {
325          if (fReconParticles == 0x0) fReconParticles = new TClonesArray("AliTOFcluster",  1000);
326          branch->SetAddress(&fReconParticles);
327        }
328     }
329
330   /*
331   if (fLoader->TreeR() && fReconParticles) //I do not know where this array is created - skowron
332     {
333       branch = fLoader->TreeR()->GetBranch("TOF"); 
334       if (branch) 
335        {
336          branch->SetAddress(&fReconParticles) ;
337        }
338     }
339   */
340 }
341
342 //_____________________________________________________________________________
343 void AliTOF::CreateGeometry()
344 {
345   //
346   // Common geometry code 
347   //
348   //Begin_Html
349   /*
350     <img src="picts/AliTOFv23.gif">
351   */
352   //End_Html
353   //
354
355   Float_t xTof, yTof;
356
357   if (IsVersion()==8) {
358
359     xTof = 124.5;//fTOFGeometry->StripLength()+2.*(0.3+0.03); // cm,  x-dimension of FTOA volume
360     yTof = fTOFGeometry->Rmax()-fTOFGeometry->Rmin(); // cm,  y-dimension of FTOA volume
361     Float_t zTof = fTOFGeometry->ZlenA();             // cm,  z-dimension of FTOA volume
362     
363     //  TOF module internal definitions
364     TOFpc(xTof, yTof, zTof);
365
366   } else if (IsVersion()==7) {
367
368     xTof = 124.5;//fTOFGeometry->StripLength()+2.*(0.3+0.03); // cm,  x-dimension of FTOA volume
369     yTof = fTOFGeometry->Rmax()-fTOFGeometry->Rmin(); // cm,  y-dimension of FTOA volume
370     Float_t zTof = fTOFGeometry->ZlenA();             // cm,  z-dimension of FTOA volume
371     
372     //  TOF module internal definitions
373     TOFpc(xTof, yTof, zTof, fTOFGeometry->ZlenB());
374
375   } else {
376
377     Float_t wall = 4.;//cm // frame inbetween TOF modules
378
379     // Sizes of TOF module with its support etc..
380     xTof = 2.*(fTOFGeometry->Rmin()*TMath::Tan(10.*kDegrad)-wall/2.-0.5);
381     yTof = fTOFGeometry->Rmax()-fTOFGeometry->Rmin();
382
383     //  TOF module internal definitions 
384     TOFpc(xTof, yTof, fTOFGeometry->ZlenC(), fTOFGeometry->ZlenB(), fTOFGeometry->ZlenA(), fTOFGeometry->MaxhZtof());
385   }
386
387 }
388
389
390 //___________________________________________
391 void AliTOF::ResetHits ()
392 {
393   // Reset number of clusters and the cluster array for this detector
394   AliDetector::ResetHits ();
395 }
396
397 //____________________________________________
398 void AliTOF::ResetDigits ()
399 {
400   //
401   // Reset number of digits and the digits array for this detector
402   AliDetector::ResetDigits ();
403   //
404
405 //____________________________________________
406 void AliTOF::ResetSDigits ()
407 {
408   //
409   // Reset number of sdigits and the sdigits array for this detector
410   fNSDigits = 0;
411   //fSDigits = 0x0;
412   //
413
414 //_____________________________________________________________________________
415 void AliTOF::Init()
416 {
417   //
418   // Initialise TOF detector after it has been built
419   //
420   // Set id of TOF sensitive volume
421   if (IsVersion() !=0) fIdSens=gMC->VolId("FPAD");
422
423   /*
424   // Save the geometry
425   TDirectory* saveDir = gDirectory;
426   AliRunLoader::Instance()->CdGAFile();
427   fTOFGeometry->Write("TOFGeometry");
428   saveDir->cd();
429   */
430 }
431
432 //____________________________________________________________________________
433 void AliTOF::MakeBranch(Option_t* option)
434 {
435  //
436  // Initializes the Branches of the TOF inside the 
437  // trees written for each event. 
438  // AliDetector::MakeBranch initializes just the 
439  // Branch inside TreeH. Here we add the branches in 
440  // TreeD, TreeS and TreeR.
441  //
442   const char *oH = strstr(option,"H");
443   if (fLoader->TreeH() && oH)
444    {
445      if (fHits == 0x0)
446       {
447         if (fTZero) fHits   = new TClonesArray("AliTOFhitT0", 1000);
448         else fHits   = new TClonesArray("AliTOFhit", 1000);
449       }
450    }
451   
452   AliDetector::MakeBranch(option);
453
454   Int_t buffersize = 4000;
455   const Int_t kSize=10;
456   Char_t branchname[kSize];
457   snprintf(branchname,kSize,"%s",GetName());
458   
459   const char *oD = strstr(option,"D");
460   const char *oS = strstr(option,"S");
461   const char *oR = strstr(option,"R");
462
463   if (fLoader->TreeD() && oD){
464     if (fDigits == 0x0) fDigits = new TClonesArray("AliTOFdigit",  1000); 
465     MakeBranchInTree(fLoader->TreeD(), branchname, &fDigits,buffersize, 0) ;
466   }
467
468   if (fLoader->TreeS() && oS){
469     if (fSDigits == 0x0) fSDigits = new TClonesArray("AliTOFSDigit",  1000);
470     MakeBranchInTree(fLoader->TreeS(), branchname, &fSDigits,buffersize, 0) ;
471   }
472
473   if (fLoader->TreeR() && oR){
474     if (fReconParticles == 0x0) fReconParticles = new TClonesArray("AliTOFcluster",  1000);
475     MakeBranchInTree(fLoader->TreeR(), branchname, &fReconParticles,buffersize, 0) ;
476   }
477
478   /*
479   if (fReconParticles && fLoader->TreeR() && oR){
480     MakeBranchInTree(fLoader->TreeR(), branchname, &fReconParticles,buffersize, 0) ;
481   }
482   */
483 }
484
485 //____________________________________________________________________________
486 void AliTOF::Makehits(Bool_t hits) 
487 {
488 // default argument used, see AliTOF.h
489 // Enable/Disable the writing of the TOF-hits branch 
490 // on TreeH
491 // by default :  enabled for TOFv1, v2, v3, v4, v5
492 //              disabled for TOFv0
493 // 
494    if (hits &&  (IsVersion()!=0))
495       fIdSens = gMC->VolId("FPAD");
496    else
497       AliInfo("Option for writing the TOF-hits branch on TreeH: disabled");
498 }
499
500 //____________________________________________________________________________
501 void AliTOF::FinishEvent()
502 {
503 // do nothing
504 }
505
506 //____________________________________________________________________________
507 void AliTOF::Hits2SDigits()
508 {
509 //
510 // Use the TOF SDigitizer to make TOF SDigits
511 //
512
513 //  AliInfo("Start...");
514   
515   AliRunLoader * rl = fLoader->GetRunLoader();
516   AliDebug(2,"Initialized runLoader");
517   AliTOFSDigitizer sd((rl->GetFileName()).Data());
518   AliDebug(2,"Initialized TOF sdigitizer");
519   //ToAliDebug(1, sd.Print(""));
520   //AliInfo("ToAliDebug");
521
522   //sd.Digitize("all") ;
523   sd.Digitize("partial") ;
524
525   AliDebug(2,"I am sorting from AliTOF class");
526
527 }
528
529 //____________________________________________________________________________
530 void AliTOF::Hits2SDigits(Int_t evNumber1, Int_t evNumber2)
531 {
532 //
533 // Use the TOF SDigitizer to make TOF SDigits
534 //
535
536   if ((evNumber2-evNumber1)==1) 
537     AliDebug(1, Form("I am making sdigits for the %dth event", evNumber1));
538   if ((evNumber2-evNumber1)>1)
539     AliDebug(1, Form("I am making sdigits for the events from the %dth to the %dth", evNumber1, evNumber2-1));
540  
541   AliRunLoader * rl = fLoader->GetRunLoader();
542   AliTOFSDigitizer sd((rl->GetFileName()).Data(),evNumber1,evNumber2) ;
543   ToAliDebug(1, sd.Print(""));
544
545   sd.Digitize("") ;
546
547 }
548
549 //___________________________________________________________________________
550 AliDigitizer* AliTOF::CreateDigitizer(AliDigitizationInput* digInput) const
551 {
552   AliDebug(2,"I am creating the TOF digitizer");
553   return new AliTOFDigitizer(digInput);
554 }
555
556 //___________________________________________________________________________
557 Bool_t AliTOF::CheckOverlap(const Int_t * const vol,
558                             Int_t* digit,Int_t Track)
559 {
560 //
561 // Checks if 2 or more hits belong to the same pad.
562 // In this case the data assigned to the digit object
563 // are the ones of the first hit in order of Time.
564 // 2 hits from the same track on the same pad are collected.
565 // Called only by Hits2SDigits.
566 // This procedure has to be optimized in the next TOF release.
567 //
568
569   Bool_t overlap = kFALSE;
570   Int_t  vol2[5];
571
572   for (Int_t ndig=0; ndig<fSDigits->GetEntries(); ndig++){
573     AliTOFdigit* currentDigit = (AliTOFdigit*)(fSDigits->UncheckedAt(ndig));
574     currentDigit->GetLocation(vol2);
575     Bool_t idem= kTRUE;
576     // check on digit volume
577     for (Int_t i=0;i<=4;i++){
578       if (!idem) break;
579       if (vol[i]!=vol2[i]) idem=kFALSE;}
580
581     if (idem){  // same pad fired
582       Int_t tdc2 = digit[0];
583       Int_t tdc1 = currentDigit->GetTdc();
584
585       // we separate two digits on the same pad if
586       // they are separated in time by at least 25 ns
587       // remember that tdc time is given in ps
588
589       if (TMath::Abs(tdc1-tdc2)<25000){
590         // in case of overlap we take the earliest
591         if (tdc1>tdc2){
592           currentDigit->SetTdc(tdc2); 
593           currentDigit->SetAdc(digit[1]);
594         }
595         else {
596           currentDigit->SetTdc(tdc1);
597           currentDigit->SetAdc(digit[1]);
598         }
599         currentDigit->AddTrack(Track); // add track number in the track array
600         overlap = kTRUE;
601         return overlap;
602       } else 
603         overlap= kFALSE;
604
605     } // close if (idem) -> two digits on the same TOF pad
606
607   } // end loop on existing sdigits
608
609   return overlap;
610 }
611 //____________________________________________________________________________
612 void AliTOF::Digits2Raw()
613 {
614 //
615 // Starting from the TOF digits, writes the Raw Data objects
616 //
617
618   TStopwatch stopwatch;
619   stopwatch.Start();
620
621   fLoader->LoadDigits();
622
623   TTree* digits = fLoader->TreeD();
624   if (!digits) {
625     AliError("no digits tree");
626     return;
627   }
628   
629   //AliTOFDDLRawData rawWriter;
630   fTOFRawWriter.Clear();
631   fTOFRawWriter.SetVerbose(0);
632   if (fTOFRawWriter.GetPackedAcquisitionMode()) {
633     if(fTOFRawWriter.GetMatchingWindow()>8192)
634       AliWarning(Form("You are running in packing mode and the matching window is %.2f ns, i.e. greater than 199.8848 ns",
635                       fTOFRawWriter.GetMatchingWindow()*AliTOFGeometry::TdcBinWidth()*1.e-03));
636   }
637   
638   AliDebug(1,"Formatting raw data for TOF");
639   digits->GetEvent(0);
640   fTOFRawWriter.RawDataTOF(digits->GetBranch("TOF"));  
641
642   fLoader->UnloadDigits();
643   
644   AliDebug(1, Form("Execution time to write TOF raw data : R:%.2fs C:%.2fs",
645                    stopwatch.RealTime(),stopwatch.CpuTime()));
646
647 }
648
649 //____________________________________________________________________________
650 void AliTOF::RecreateSDigitsArray() {
651 //
652 // delete TClonesArray fSDigits and create it again
653 //  needed for backward compatability with PPR test production
654 //
655   fSDigits->Clear();
656 }
657 //____________________________________________________________________________
658 void AliTOF::CreateSDigitsArray() {
659 //
660 // create TClonesArray fSDigits
661 //  needed for backward compatability with PPR test production
662 //
663   fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
664 }
665 //____________________________________________________________________________
666 Bool_t AliTOF::Raw2SDigits(AliRawReader* rawReader)
667 {
668   //
669   // Converts raw data to sdigits for TOF
670   //
671
672   TStopwatch stopwatch;
673   stopwatch.Start();
674
675   if(!GetLoader()->TreeS()) {MakeTree("S");  MakeBranch("S");}
676   //TClonesArray &aSDigits = *fSDigits;
677
678   AliTOFRawStream tofRawStream = AliTOFRawStream();
679   tofRawStream.Raw2SDigits(rawReader, fSDigits);
680
681   GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
682   Int_t nSDigits = fSDigits->GetEntries();
683
684   ResetSDigits();
685
686   AliDebug(1, Form("Got %d TOF sdigits", nSDigits));
687   AliDebug(1, Form("Execution time to read TOF raw data and fill TOF sdigit tree : R:%.2fs C:%.2fs",
688                    stopwatch.RealTime(),stopwatch.CpuTime()));
689
690   return kTRUE;
691
692 }
693
694 //____________________________________________________________________________
695 void AliTOF::Raw2Digits(AliRawReader* rawReader)
696 {
697   //
698   // Converts raw data to digits for TOF
699   //
700
701   TStopwatch stopwatch;
702   stopwatch.Start();
703
704   if(!GetLoader()->TreeD()) {MakeTree("D");  MakeBranch("D");}
705   //TClonesArray &aDigits = *fDigits;
706
707   AliTOFRawStream tofRawStream = AliTOFRawStream();
708   tofRawStream.Raw2Digits(rawReader, fDigits);
709
710   GetLoader()->TreeD()->Fill(); GetLoader()->WriteDigits("OVERWRITE");//write out digits
711   Int_t nDigits = fDigits->GetEntries();
712
713   ResetDigits();
714
715   AliDebug(1, Form("Got %d TOF digits", nDigits));
716   AliDebug(1, Form("Execution time to read TOF raw data and fill TOF digit tree : R:%.2fs C:%.2fs",
717                    stopwatch.RealTime(),stopwatch.CpuTime()));
718
719 }