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