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