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