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