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