]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOF.cxx
New TOF version with holes for PHOS/RICH
[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
16 /*
17 $Log$
18 Revision 1.11.2.1  2000/05/10 09:37:15  vicinanz
19 New version with Holes for PHOS/RICH
20
21 Revision 1.11  1999/11/05 22:39:06  fca
22 New hits structure
23
24 Revision 1.10  1999/11/01 20:41:57  fca
25 Added protections against using the wrong version of FRAME
26
27 Revision 1.9  1999/10/15 15:35:19  fca
28 New version for frame1099 with and without holes
29
30 Revision 1.9  1999/09/29 09:24:33  fca
31 Introduction of the Copyright and cvs Log
32
33 */
34 ///////////////////////////////////////////////////////////////////////////////
35 //                                                                           //
36 //  Time Of Flight                               FCA                         //
37 //  This class contains the basic functions for the Time Of Flight           //
38 //  detector. Functions specific to one particular geometry are              //
39 //  contained in the derived classes                                         //
40 //
41 //  VERSIONE WITH 5 SYMMETRIC MODULES ALONG Z AXIS
42 //  ==============================================
43 //  
44 //  VERSION WITH HOLES FOR PHOS AND TRD IN SPACEFRAME WITH HOLES
45 //
46 //  Volume sensibile : FPAD
47 //
48 //
49 //
50 // Begin_Html
51 /*
52 <img src="picts/AliTOFClass.gif">
53 */
54 //End_Html
55 //             
56 //
57 //                                                                           //
58 ///////////////////////////////////////////////////////////////////////////////
59
60 #include <iostream.h>
61
62 #include "AliTOF.h"
63 #include "AliTOFD.h"
64 #include "TBRIK.h"
65 #include "TNode.h"
66 #include "TObject.h"
67 #include "TRandom.h"
68 #include "AliRun.h"
69 #include "AliConst.h"
70
71  
72 ClassImp(AliTOF)
73  
74 //_____________________________________________________________________________
75 AliTOF::AliTOF()
76 {
77   //
78   // Default constructor
79   //
80   fIshunt   = 0;
81 }
82  
83 //_____________________________________________________________________________
84 AliTOF::AliTOF(const char *name, const char *title)
85        : AliDetector(name,title)
86 {
87   //
88   // AliTOF standard constructor
89   // 
90   // Here are fixed some important parameters
91   //
92
93   // Initialization of hits and digits array
94   //
95   fHits   = new TClonesArray("AliTOFhit",  405);
96   fIshunt  = 0;
97   fDigits = new TClonesArray("AliTOFdigit",405);
98   //
99   // Digitization parameters
100   //
101   // (Transfer Functions to be inserted here)
102   //
103   SetMarkerColor(7);
104   SetMarkerStyle(2);
105   SetMarkerSize(0.4);
106
107 // General Geometrical Parameters
108   fNTof    =  18;
109   fRmax    = 399.0;//cm
110   fRmin    = 370.0;//cm
111   fZlenC   = 177.5;//cm
112   fZlenB   = 141.0;//cm
113   fZlenA   = 106.0;//cm
114   fZtof    = 370.5;//cm
115
116 // Strip Parameters
117   fStripLn = 122.0;//cm  Strip Length
118   fSpace   =   5.5;//cm  Space Beetween the strip and the bottom of the plate 
119   fDeadBndZ=   1.5;//cm  Dead Boundaries of a Strip along Z direction (width)
120   fDeadBndX=   1.0;//cm  Dead Boundaries of a Strip along X direction (length)
121   fXpad    =   2.5;//cm  X size of a pad
122   fZpad    =   3.5;//cm  Z size of a pad
123   fGapA    =   4.; //cm  Gap beetween tilted strip in A-type plate
124   fGapB    =   6.; //cm  Gap beetween tilted strip in B-type plate
125   fOverSpc =  15.3;//cm Space available for sensitive layers in radial direction
126   fNpadX   =  48;  // Number of pads in a strip along the X direction
127   fNpadZ   =   2;  // Number of pads in a strip along the Z direction
128   fPadXStr = fNpadX*fNpadZ; //Number of pads per strip
129   fNStripA = 0;
130   fNStripB = 0;
131   fNStripC = 0;
132  
133 // Physical performances
134   fTimeRes = 100.;//ps
135   fChrgRes = 100.;//pC
136
137 // DAQ characteristics
138   fPadXSector = 1932;
139   fNRoc       = 14;
140   fNFec       = 32;
141   fNTdc       = 32;
142   fNPadXRoc   = (Int_t)fPadXSector/fNRoc;
143 }
144
145 //_____________________________________________________________________________
146 void AliTOF::AddHit(Int_t track, Int_t *vol, Float_t *hits)
147 {
148   //
149   // Add a TOF hit
150   //
151   TClonesArray &lhits = *fHits;
152   new(lhits[fNhits++]) AliTOFhit(fIshunt, track, vol, hits);
153 }
154  
155 //_____________________________________________________________________________
156 void AliTOF::AddDigit(Int_t *tracks, Int_t *vol, Float_t *digits)
157 {
158   //
159   // Add a TOF digit
160   //
161   TClonesArray &ldigits = *fDigits;
162   new (ldigits[fNdigits++]) AliTOFdigit(tracks, vol, digits);
163 }
164
165
166 //_____________________________________________________________________________
167 void AliTOF::CreateGeometry()
168 {
169   //
170   // Common geometry code 
171   //
172   //Begin_Html
173   /*
174     <img src="picts/AliTOFv23.gif">
175   */
176   //End_Html
177   //
178   const Double_t kPi=TMath::Pi();
179   const Double_t kDegrad=kPi/180.;
180   //
181   Float_t xTof, yTof, Wall;
182
183   // frame inbetween TOF modules
184   Wall = 4.;//cm
185
186   // Sizes of TOF module with its support etc..
187   xTof = 2.*(fRmin*TMath::Tan(10*kDegrad)-Wall/2-.5);
188   yTof = fRmax-fRmin;
189
190 //  TOF module internal definitions 
191   TOFpc(xTof, yTof, fZlenC, fZlenB, fZlenA, fZtof);
192 }
193
194 //_____________________________________________________________________________
195 void AliTOF::DrawModule()
196 {
197   //
198   // Draw a shaded view of the common part of the TOF geometry
199   //
200
201    cout << " Drawing of AliTOF"<< endl; 
202   // Set everything unseen
203   gMC->Gsatt("*", "seen", -1);
204   // 
205   // Set ALIC mother transparent
206   gMC->Gsatt("ALIC","SEEN",0);
207   //
208   // Set the volumes visible
209   gMC->Gsatt("FTOA","SEEN",1);
210   gMC->Gsatt("FTOB","SEEN",1);
211   gMC->Gsatt("FTOC","SEEN",1);
212   gMC->Gsatt("FLTA","SEEN",1);
213   gMC->Gsatt("FLTB","SEEN",1);
214   gMC->Gsatt("FLTC","SEEN",1);
215   gMC->Gsatt("FSTR","SEEN",1);
216   //
217   gMC->Gdopt("hide", "on");
218   gMC->Gdopt("shad", "on");
219   gMC->Gsatt("*", "fill", 7);
220   gMC->SetClipBox(".");
221   gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
222   gMC->DefaultRange();
223   gMC->Gdraw("alic", 40, 30, 0, 12, 9.5, .02, .02);
224   gMC->Gdhead(1111, "Time Of Flight");
225   gMC->Gdman(18, 4, "MAN");
226   gMC->Gdopt("hide","off");
227 }
228
229 //_____________________________________________________________________________
230 void AliTOF::CreateMaterials()
231 {
232   //
233   // Defines TOF materials for all versions
234   // Authors :   Maxim Martemianov, Boris Zagreev (ITEP) 
235   //            18/09/98 
236   //
237   Int_t   ISXFLD = gAlice->Field()->Integ();
238   Float_t SXMGMX = gAlice->Field()->Max();
239   //
240   //--- Quartz (SiO2) 
241   Float_t   aq[2] = { 28.0855,15.9994 };
242   Float_t   zq[2] = { 14.,8. };
243   Float_t   wq[2] = { 1.,2. };
244   Float_t   dq = 2.20;
245   Int_t nq = -2;
246   // --- Freon
247   Float_t afre[2]  = {12.011,18.9984032 };
248   Float_t zfre[2]  = { 6., 9.};
249   Float_t wfre[2]  = { 5.,12.};
250   Float_t densfre  = 1.5;
251   Int_t nfre = -2;
252   // --- CO2 
253   Float_t ac[2]   = {12.,16.};
254   Float_t zc[2]   = { 6., 8.};
255   Float_t wc[2]   = { 1., 2.};
256   Float_t dc = .001977;
257   Int_t nc = -2;
258    // For mylar (C5H4O2) 
259   Float_t amy[3] = { 12., 1., 16. };
260   Float_t zmy[3] = {  6., 1.,  8. };
261   Float_t wmy[3] = {  5., 4.,  2. };
262   Float_t dmy    = 1.39;
263   Int_t nmy = -3;
264  // For polyethilene (CH2) for honeycomb!!!!
265   Float_t ape[2] = { 12., 1. };
266   Float_t zpe[2] = {  6., 1. };
267   Float_t wpe[2] = {  1., 2. };
268   Float_t dpe    = 0.935*0.479; //To have 1%X0 for 1cm as for honeycomb
269   Int_t npe = -2;
270   // --- G10 
271   Float_t ag10[4] = { 12.,1.,16.,28. };
272   Float_t zg10[4] = {  6.,1., 8.,14. };
273   Float_t wmatg10[4] = { .259,.288,.248,.205 };
274   Float_t densg10  = 1.7;
275   Int_t nlmatg10 = -4;
276   // --- DME 
277   Float_t adme[5] = { 12.,1.,16.,19.,79. };
278   Float_t zdme[5] = {  6.,1., 8., 9.,35. };
279   Float_t wmatdme[5] = { .4056,.0961,.2562,.1014,.1407 };
280   Float_t densdme  = .00205;
281   Int_t nlmatdme = 5;
282   // ---- ALUMINA (AL203) 
283   Float_t aal[2] = { 27.,16.};
284   Float_t zal[2] = { 13., 8.};
285   Float_t wmatal[2] = { 2.,3. };
286   Float_t densal  = 2.3;
287   Int_t nlmatal = -2;
288   // -- Water
289   Float_t awa[2] = {  1., 16. };
290   Float_t zwa[2] = {  1.,  8. };
291   Float_t wwa[2] = {  2.,  1. };
292   Float_t dwa    = 1.0;
293   Int_t nwa = -2;
294   //
295   //AliMaterial(0, "Vacuum$", 1e-16, 1e-16, 1e-16, 1e16, 1e16);
296   AliMaterial( 1, "Air$",14.61,7.3,0.001205,30423.24,67500.);
297   AliMaterial( 2, "Cu $",  63.54, 29.0, 8.96, 1.43, 14.8);
298   AliMaterial( 3, "C  $",  12.01,  6.0, 2.265,18.8, 74.4);
299   AliMixture ( 4, "Polyethilene$", ape, zpe, dpe, npe, wpe);
300   AliMixture ( 5, "G10$", ag10, zg10, densg10, nlmatg10, wmatg10);
301   AliMixture ( 6, "DME ", adme, zdme, densdme, nlmatdme, wmatdme);
302   AliMixture ( 7, "CO2$", ac, zc, dc, nc, wc);
303   AliMixture ( 8, "ALUMINA$", aal, zal, densal, nlmatal, wmatal);
304   AliMaterial( 9, "Al $", 26.98, 13., 2.7, 8.9, 37.2);
305   AliMaterial(10, "C-TRD$", 12.01, 6., 2.265*18.8/69.282*15./100, 18.8, 74.4); // for 15%
306   AliMixture (11, "Mylar$",  amy, zmy, dmy, nmy, wmy);
307   AliMixture (12, "Freon$",  afre, zfre, densfre, nfre, wfre);
308   AliMixture (13, "Quartz$", aq, zq, dq, nq, wq);
309   AliMixture (14, "Water$",  awa, zwa, dwa, nwa, wwa);
310
311   Float_t epsil, stmin, deemax, stemax;
312  
313   //   Previous data
314   //       EPSIL  =  0.1   ! Tracking precision, 
315   //       STEMAX = 0.1      ! Maximum displacement for multiple scattering
316   //       DEEMAX = 0.1    ! Maximum fractional energy loss, DLS 
317   //       STMIN  = 0.1 
318   //
319   //   New data  
320   epsil  = .001;  // Tracking precision,
321   stemax = -1.;   // Maximum displacement for multiple scattering
322   deemax = -.3;   // Maximum fractional energy loss, DLS
323   stmin  = -.8;
324
325   AliMedium( 1, "Air$"  ,  1, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
326   AliMedium( 2, "Cu $"  ,  2, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
327   AliMedium( 3, "C  $"  ,  3, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
328   AliMedium( 4, "Pol$"  ,  4, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
329   AliMedium( 5, "G10$"  ,  5, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
330   AliMedium( 6, "DME$"  ,  6, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
331   AliMedium( 7, "CO2$"  ,  7, 0, ISXFLD, SXMGMX, 10., -.01, -.1, .01, -.01);
332   AliMedium( 8,"ALUMINA$", 8, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
333   AliMedium( 9,"Al Frame$",9, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
334   AliMedium(10, "DME-S$",  6, 1, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
335   AliMedium(11, "C-TRD$", 10, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
336   AliMedium(12, "Myl$"  , 11, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
337   AliMedium(13, "Fre$"  , 12, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
338   AliMedium(14, "Fre-S$", 12, 1, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
339   AliMedium(15, "Glass$", 13, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
340   AliMedium(16, "Water$", 14, 0, ISXFLD, SXMGMX, 10., stemax, deemax, epsil, stmin);
341 }
342
343 //_____________________________________________________________________________
344 Int_t AliTOF::DistancetoPrimitive(Int_t , Int_t )
345 {
346   //
347   // Returns distance from mouse pointer to detector, default version
348   //
349   return 9999;
350 }
351  
352 //_____________________________________________________________________________
353 void AliTOF::Init()
354 {
355   //
356   // Initialise TOF detector after it has been built
357   //
358   // Set id of TOF sensitive volume
359   if (IsVersion() !=0) fIdSens=gMC->VolId("FPAD");
360   //
361 }
362
363 //____________________________________________________________________________
364 void AliTOF::MakeBranch(Option_t* option)
365 //
366 // Initializes the Branches of the TOF inside the 
367 // trees written for each event. 
368 //  AliDetector::MakeBranch initializes just the 
369 // Branch inside TreeH. Here we add the branch in 
370 // TreeD.
371 //
372 {
373
374   AliDetector::MakeBranch(option);
375
376   Int_t buffersize = 4000;
377   Char_t branchname[10];
378   sprintf(branchname,"%s",GetName());
379   char *D = strstr(option,"D");
380
381   if (fDigits && gAlice->TreeD() && D){
382      gAlice->TreeD()->Branch(branchname,&fDigits,buffersize);
383      printf("Making Branch %s for digits \n",branchname);
384   }
385 }
386
387 //____________________________________________________________________________
388 void AliTOF::FinishEvent()
389 {
390 //  Hits2Digits();
391 }
392
393
394 //____________________________________________________________________________
395 void AliTOF::Hits2Digits(Int_t evNumber)
396 //
397 //  Starting from the Hits Tree (TreeH), this
398 // function writes the Digits Tree (TreeD) storing 
399 // the digits informations.
400 // Has to be called just at the end of an event or 
401 // at the end of a whole run.
402 //  It could  also be called by AliTOF::Finish Event()
403 // but it can be too heavy.
404 // Just for MC events. 
405 //
406 // Called by the macro H2D.C
407 //
408 {
409   AliTOFhit* currentHit;
410   TTree    *TD, *TH;
411   Int_t    tracks[3];
412   Int_t    vol[5];
413   Float_t  digit[2];
414   TClonesArray* TOFhits=this->Hits();
415
416   Int_t nparticles = gAlice->GetEvent(evNumber); 
417   if (nparticles <= 0) return;
418
419   TD = gAlice->TreeD();
420   TH = gAlice->TreeH();
421   Int_t    ntracks =(Int_t) TH->GetEntries();
422   Int_t    nbytes, nhits;
423   TRandom *rnd = new TRandom();
424
425   for (Int_t ntk=0; ntk<ntracks; ntk++){
426
427      nbytes = TH->GetEvent(ntk);
428      nhits  = TOFhits->GetEntriesFast();
429
430      for (Int_t hit=0; hit<nhits; hit++){
431
432         currentHit = (AliTOFhit*)(TOFhits->At(hit));
433
434         vol[0] = currentHit->GetSector();
435         vol[1] = currentHit->GetPlate();
436         vol[2] = currentHit->GetPad_x();
437         vol[3] = currentHit->GetPad_z();
438         vol[4] = currentHit->GetStrip();
439
440         Float_t IdealTime = currentHit->GetTof();
441         Float_t TDCTime   = rnd->Gaus(IdealTime, fTimeRes);     
442         digit[0] = TDCTime;
443
444         Float_t IdealCharge = currentHit->GetEdep();
445         Float_t ADCcharge = rnd->Gaus(IdealCharge, fChrgRes);
446         digit[1] = ADCcharge;
447
448         Int_t Track = currentHit -> GetTrack();
449         tracks[0] = Track;
450         tracks[1] = 0;
451         tracks[2] = 0;
452
453         Bool_t Overlap = CheckOverlap(vol, digit, Track);
454         if(!Overlap) AddDigit(tracks, vol, digit);
455      }
456   }
457   TD->Fill();
458   TD->Write();  
459 }
460
461 //___________________________________________________________________________
462 Bool_t AliTOF::CheckOverlap(Int_t* vol, Float_t* digit,Int_t Track)
463 //
464 // Checks if 2 or more hits belong to the same pad.
465 // In this case the data assigned to the digit object
466 // are the ones of the first hit in order of Time.
467 //
468 // Called only by Hits2Digits.
469 //
470 {
471         Bool_t Overlap = 0;
472         Int_t  vol2[5];
473
474         for (Int_t ndig=0; ndig<fNdigits; ndig++){
475            AliTOFdigit* currentDigit = (AliTOFdigit*)(fDigits->UncheckedAt(ndig));
476            currentDigit->GetLocation(vol2);
477            Bool_t Idem=1;
478            for (Int_t i=0;i<=4;i++){
479                if (vol[i]!=vol2[i]) Idem=0;}
480            if (Idem){
481               Float_t TDC2 = digit[0];
482               Float_t TDC1 = currentDigit->GetTdc();
483               if (TDC1>TDC2){
484                   currentDigit->SetTdc(TDC2); 
485                   currentDigit->SetAdc(digit[1]);
486               }
487               currentDigit->AddTrack(Track);
488               Overlap = 1;
489            }
490         }
491         return Overlap;
492 }
493
494
495 //____________________________________________________________________________
496 void AliTOF::Digits2Raw(Int_t evNumber)
497 //
498 // Starting from digits, writes the 
499 // Raw Data objects, i.e. a 
500 // TClonesArray of 18 AliTOFRawSector objects
501 //
502
503 {
504   TTree* TD;
505
506   Int_t nparticles = gAlice->GetEvent(evNumber); 
507   if (nparticles <= 0) return;
508
509   TD = gAlice->TreeD();
510   
511   TClonesArray* TOFdigits = this->Digits();
512   Int_t ndigits = TOFdigits->GetEntriesFast();
513
514   TClonesArray* Raw = new TClonesArray("AliTOFRawSector",fNTof+2); 
515
516   for (Int_t isect=1;isect<=fNTof;isect++){
517      AliTOFRawSector* currentSector = (AliTOFRawSector*)Raw->UncheckedAt(isect);
518      TClonesArray* RocData = (TClonesArray*)currentSector->GetRocData();
519
520      for (Int_t digit=0; digit<ndigits; digit++){
521         AliTOFdigit* currentDigit = (AliTOFdigit*)TOFdigits->UncheckedAt(digit);
522         Int_t sector = currentDigit->GetSector();
523         if (sector==isect){
524             Int_t   Pad    = currentDigit -> GetTotPad();
525             Int_t   Roc    = (Int_t)(Pad/fNPadXRoc)-1;
526             if (Roc>=fNRoc) printf("Wrong n. of ROC ! Roc = %i",Roc);
527             Int_t   PadRoc = (Int_t) Pad%fNPadXRoc;
528             Int_t   Fec    = (Int_t)(PadRoc/fNFec)-1;
529             Int_t   Tdc    = (Int_t)(PadRoc%fNFec)-1;
530             Float_t Time   = currentDigit->GetTdc();
531             Float_t Charge = currentDigit->GetAdc();
532             AliTOFRoc* currentROC = (AliTOFRoc*)RocData->UncheckedAt(Roc);
533             Int_t Error    = 0;
534             currentROC->AddItem(Fec, Tdc, Error, Charge, Time);
535         }
536      }
537      
538      UInt_t TotSize=16,RocSize=0;
539      UInt_t RocHead[14],RocChek[14];
540      UInt_t GlobalCheckSum=0;
541
542      for (UInt_t iRoc = 1; iRoc<(UInt_t)fNRoc; iRoc++){
543         AliTOFRoc* currentRoc = (AliTOFRoc*)RocData->UncheckedAt(iRoc); 
544         RocSize  = currentRoc->Items*2+1;
545         TotSize += RocSize*4;
546         if (RocSize>=pow(2,16)) RocSize=0;
547         RocHead[iRoc]   = iRoc<<28;
548         RocHead[iRoc]  += RocSize;
549         RocChek[iRoc]   = currentRoc->GetCheckSum();
550         Int_t HeadCheck = currentRoc->BitCount(RocHead[iRoc]);
551         GlobalCheckSum += HeadCheck;
552         GlobalCheckSum += RocChek[iRoc];
553      }
554      
555      AliTOFRoc* DummyRoc = new AliTOFRoc();
556      TotSize *= 4;
557      if (TotSize>=pow(2,24)) TotSize=0;
558      UInt_t Header = TotSize;
559      UInt_t SectId = ((UInt_t)isect)<<24;
560      Header += SectId;
561      GlobalCheckSum += DummyRoc->BitCount(Header);
562      currentSector->SetGlobalCS(GlobalCheckSum);
563      currentSector->SetHeader(Header);
564   }  
565 }
566  
567 //____________________________________________________________________________
568 void AliTOF::Raw2Digits(Int_t evNumber)
569 //
570 //  Converts Raw Data objects into digits objects.
571 //  We schematize the raw data with a 
572 //  TClonesArray of 18 AliTOFRawSector objects
573 //
574 {
575   TTree    *TD;
576   Int_t    vol[5];
577   Int_t    tracks[3];
578   Float_t  digit[2];
579  
580   tracks[0]=0;
581   tracks[1]=0;
582   tracks[2]=0;
583  
584   Int_t nparticles = gAlice->GetEvent(evNumber); 
585   if (nparticles <= 0) return;
586
587   TD = gAlice->TreeD();
588   
589   TClonesArray* Raw = new TClonesArray("AliTOFRawSector",fNTof+2);
590   
591   for(Int_t nSec=1; nSec<=fNTof; nSec++){
592      AliTOFRawSector* currentSector = (AliTOFRawSector*)Raw->UncheckedAt(nSec);
593      TClonesArray* RocData = (TClonesArray*)currentSector->GetRocData();
594      for(Int_t nRoc=1; nRoc<=14; nRoc++){
595         AliTOFRoc* currentRoc = (AliTOFRoc*)RocData->UncheckedAt(nRoc);
596         Int_t currentItems = currentRoc->GetItems();
597         for(Int_t item=1; item<currentItems; item++){ 
598            Int_t nPad = currentRoc->GetTotPad(item);        
599            vol[0] = nSec;
600            Int_t nStrip = (Int_t)(nPad/fPadXStr)+1;
601            Int_t nPlate = 5;
602            if (nStrip<=fNStripC+2*fNStripB+fNStripA) nPlate = 4;
603            if (nStrip<=fNStripC+fNStripB+fNStripA)   nPlate = 3;
604            if (nStrip<=fNStripC+fNStripB)            nPlate = 2;
605            if (nStrip<=fNStripC)                     nPlate=1;
606            vol[1] = nPlate;
607            switch (nPlate){
608            case 1: break;
609            case 2: nStrip -= (fNStripC);
610                    break;
611            case 3: nStrip -= (fNStripC+fNStripB);
612                    break;
613            case 4: nStrip -= (fNStripC+fNStripB+fNStripA);
614                    break;
615            case 5: nStrip -= (fNStripC+2*fNStripB+fNStripA);
616                    break;
617            }
618            vol[2] = nStrip;
619            Int_t Pad = nPad%fPadXStr;
620            if (Pad==0) Pad=fPadXStr;
621            Int_t nPadX=0, nPadZ=0;
622            (Pad>fNpadX)? nPadX -= fNpadX : nPadX = Pad ;
623            vol[3] = nPadX;
624            (Pad>fNpadX)? nPadZ = 2 : nPadZ = 1 ;
625            vol[4] = nPadZ;
626            UInt_t error=0;
627            Float_t TDC = currentRoc->GetTime(item,error);
628            if (!error) digit[0]=TDC;
629            digit[1] = currentRoc->GetCharge(item);
630            AddDigit(tracks,vol,digit);
631         }
632      }
633   }
634   TD->Fill();
635   TD->Write();
636
637
638
639 /******************************************************************************/
640
641 ClassImp(AliTOFhit)
642  
643 //______________________________________________________________________________
644 AliTOFhit::AliTOFhit(Int_t shunt, Int_t track, Int_t *vol,
645                      Float_t *hits)
646 :AliHit(shunt, track)
647 //
648 // Constructor of hit object
649 //
650 {
651   //
652   // Store a TOF hit
653   // _______________
654   //
655   // Hit Volume
656   // 
657   fSector= vol[0];
658   fPlate = vol[1];
659   fStrip = vol[2];
660   fPad_x = vol[3];
661   fPad_z = vol[4];
662   //
663   //Position
664   fX = hits[0];
665   fY = hits[1];
666   fZ = hits[2];
667   //
668   // Momentum
669   fPx  = hits[3];
670   fPy  = hits[4];
671   fPz  = hits[5];
672   fPmom= hits[6];
673   //
674   // Time Of Flight
675   fTof = hits[7];   //TOF[s]
676   //
677   // Other Data
678   fDx  = hits[8];   //Distance from the edge along x axis
679   fDy  = hits[9];   //Y cohordinate of the hit
680   fDz  = hits[10];  //Distance from the edge along z axis
681   fIncA= hits[11];  //Incidence angle
682   fEdep= hits[12];  //Energy loss in TOF pad
683 }
684
685 //******************************************************************************
686
687 ClassImp(AliTOFdigit)
688
689 //______________________________________________________________________________
690 AliTOFdigit::AliTOFdigit(Int_t *tracks, Int_t *vol,Float_t *digit)
691 :AliDigit(tracks)
692 //
693 // Constructor of digit object
694 //
695 {
696   fSector = vol[0];
697   fPlate  = vol[1];
698   fStrip  = vol[2];
699   fPad_x  = vol[3];
700   fPad_z  = vol[4];
701   fTdc    = digit[0];
702   fAdc    = digit[1];
703 }
704
705 //______________________________________________________________________________
706 void AliTOFdigit::GetLocation(Int_t *Loc)
707 //
708 // Get the cohordinates of the digit
709 // in terms of Sector - Plate - Strip - Pad
710 //
711 {
712    Loc[0]=fSector;
713    Loc[1]=fPlate;
714    Loc[2]=fStrip;
715    Loc[3]=fPad_x;
716    Loc[4]=fPad_z;
717 }
718
719 //______________________________________________________________________________
720 Int_t AliTOFdigit::GetTotPad()
721 //
722 // Get the "total" index of the pad inside a Sector
723 // starting from the digits data.
724 //
725 {
726   AliTOF* TOF;
727   
728   if(gAlice){
729      TOF =(AliTOF*) gAlice->GetDetector("TOF");
730   }else{
731      printf("AliTOFdigit::GetTotPad - No AliRun object present, exiting");
732      return 0;
733   }
734   
735   Int_t Pad = fPad_x+TOF->fNpadX*(fPad_z-1);
736   Int_t Before=0;
737
738   switch(fPlate){ 
739   case 1: Before = 0;
740           break;
741   case 2: Before = TOF->fNStripC;
742           break;
743   case 3: Before = TOF->fNStripB + TOF->fNStripC;
744           break;
745   case 4: Before = TOF->fNStripA + TOF->fNStripB + TOF->fNStripC;
746           break;
747   case 5: Before = TOF->fNStripA + 2*TOF->fNStripB + TOF->fNStripC;
748           break;
749   }
750   
751   Int_t Strip = fStrip+Before;
752   Int_t PadTot = TOF->fPadXStr*(Strip-1)+Pad;
753   return PadTot;
754 }
755
756 //______________________________________________________________________________
757 void AliTOFdigit::AddTrack(Int_t track)
758 //
759 // Add a track to the digit 
760 //
761 {
762   if (fTracks[1]==0){
763      fTracks[1] = track;
764   }else if (fTracks[2]==0){
765      fTracks[2] = track;
766   }else{
767      printf("AliTOFdigit::AddTrack ERROR: Too many Tracks (>3) \n");
768   }
769 }
770