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