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