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