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