]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOF.cxx
Track reconstruction code for TOF: updating
[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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  Time Of Flight                                                           //
21 //  This class contains the basic functions for the Time Of Flight           //
22 //  detector. Functions specific to one particular geometry are              //
23 //  contained in the derived classes                                         //
24 //
25 //  VERSIONE WITH 5 SYMMETRIC MODULES ALONG Z AXIS
26 //  ==============================================
27 //  
28 //  VERSION WITH HOLES FOR PHOS AND TRD IN SPACEFRAME WITH HOLES
29 //
30 //  Volume sensibile : FPAD
31 //
32 //
33 //
34 // Begin_Html
35 /*
36 <img src="picts/AliTOFClass.gif">
37 */
38 //End_Html
39 //             
40 //
41 //                                                                           //
42 ///////////////////////////////////////////////////////////////////////////////
43
44
45 #include <Riostream.h>
46 #include <Rstrstream.h>
47 #include <stdlib.h>
48
49 #include <TBRIK.h>
50 #include <TFile.h>
51 #include <TFolder.h>
52 #include <TNode.h>
53 #include <TObject.h>
54 #include <TROOT.h>
55 #include <TSystem.h>
56 #include <TTask.h>
57 #include <TTree.h>
58 #include <TVirtualMC.h>
59
60 #include "AliConst.h"
61 #include "AliLoader.h"
62 #include "AliMagF.h"
63 #include "AliRun.h"
64 #include "AliTOF.h"
65 #include "AliTOFRawDigit.h"
66 #include "AliTOFRawSector.h"
67 #include "AliTOFRoc.h"
68 #include "AliTOFSDigit.h"
69 #include "AliTOFdigit.h"
70 #include "AliTOFhit.h"
71 #include "AliTOFhitT0.h"
72 #include "AliMC.h"
73 #include "AliTOFDigitizer.h"
74 #include "AliTOFtracker.h"
75  
76 ClassImp(AliTOF)
77  
78 //_____________________________________________________________________________
79 AliTOF::AliTOF()
80 {
81   //
82   // Default constructor
83   //
84   fFGeom = 0x0;
85   fDTask = 0x0;
86   fReTask = 0x0;
87   fIshunt   = 0;
88   fSDigits       = 0 ;
89   fNSDigits = 0;
90   fDigits        = 0 ;
91   fReconParticles = 0x0;
92   fName="TOF";
93   fMerger = 0x0;
94   fTZero = kFALSE;
95   fTOFGeometry = 0;
96 }
97  
98 //_____________________________________________________________________________
99 AliTOF::AliTOF(const char *name, const char *title, Option_t *option)
100        : AliDetector(name,title)
101 {
102   //
103   // AliTOF standard constructor
104   // 
105   // Here are fixed some important parameters
106   //
107
108   // Initialization of hits, sdigits and digits array
109   // added option for time zero analysis
110   fFGeom = 0x0; //skowron
111   fDTask = 0x0;
112   fReTask= 0x0;
113   fReconParticles= 0x0;
114   fMerger = 0x0;
115   fTOFGeometry = 0;
116
117   if (strstr(option,"tzero")){
118     fHits   = new TClonesArray("AliTOFhitT0",  1000);
119     fTZero = kTRUE;
120     cout << "tzero option requires AliTOFv4T0 as TOF version (check Your Config.C)" << endl;
121   }else{
122     fHits   = new TClonesArray("AliTOFhit",  1000);
123     fTZero = kFALSE;
124   }
125   if (gAlice==0) {
126      Fatal("AliTOF","gAlice==0 !");
127   }
128   if (gAlice->GetMCApp()->GetHitLists())
129      gAlice->GetMCApp()->AddHitList(fHits);
130   else Error("AliTOF","gAlice->GetHitLists()==0");
131
132   fIshunt  = 0;
133   fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
134   fDigits        = new TClonesArray("AliTOFdigit",  1000);
135   fNSDigits = 0;
136
137   fFGeom = 0x0;
138   fDTask = 0x0;
139   fReTask = 0x0;
140   fReconParticles = 0x0;
141   fMerger = 0x0;
142
143   //
144   // Digitization parameters
145   //
146   // (Transfer Functions to be inserted here)
147   //
148   SetMarkerColor(7);
149   SetMarkerStyle(2);
150   SetMarkerSize(0.4);
151
152 // General Geometrical Parameters
153   fNTof    =  18;  // number of sectors
154   fRmax    = 399.0;//cm 
155   fRmin    = 370.0;//cm
156   fZlenC   = 177.5;//cm length of module C
157   fZlenB   = 141.0;//cm length of module B
158   fZlenA   = 106.0;//cm length of module A
159   fZtof    = 371.5;//cm total semi-length of TOF detector
160
161 // Strip Parameters
162   fStripLn = 122.0;//cm  Strip Length
163   fSpace   =   5.5;//cm  Space Beetween the strip and the bottom of the plate 
164   fDeadBndZ=   1.5;//cm  Dead Boundaries of a Strip along Z direction (width)
165   fDeadBndX=   1.0;//cm  Dead Boundaries of a Strip along X direction (length)
166   fXpad    =   2.5;//cm  X size of a pad
167   fZpad    =   3.5;//cm  Z size of a pad
168   fGapA    =   4.; //cm  Gap beetween tilted strip in A-type plate
169   fGapB    =   6.; //cm  Gap beetween tilted strip in B-type plate
170   fOverSpc =  15.3;//cm Space available for sensitive layers in radial direction
171   fNpadX   =  48;  // Number of pads in a strip along the X direction
172   fNpadZ   =   2;  // Number of pads in a strip along the Z direction
173   fPadXStr = fNpadX*fNpadZ; //Number of pads per strip
174   fNStripA = 15; // number of strips in A type module 
175   fNStripB = 19; // number of strips in B type module
176   fNStripC = 20; // number of strips in C type module
177  
178 // Physical performances
179   fTimeRes = 100.;//ps
180   fChrgRes = 100.;//pC
181
182 // DAQ characteristics
183 // cfr. TOF-TDR pag. 105 for Glossary
184 // TARODA : TOF-ALICE Read Out and Data Acquisition system
185   fPadXSector = 8928; // number of pad per sector -with no holes-
186                       // ((15+2*19+2*20)*(48*2))
187   fNRoc       = 14;   // number of Roc (Read Out Controller) (TARODA)
188   fNFec       = 32;   // number of Fec (Front-End electronic Card)
189                       // (TARODA)
190   fNTdc       = 32;   // number of Tdc (Time to Digital Converter)
191   fNPadXRoc   = (Int_t)fPadXSector/fNRoc; // number of pads for each ROC
192   
193 }
194
195 //_____________________________________________________________________________
196 void AliTOF::CreateTOFFolders()
197 {
198   // create the ALICE TFolder
199   // create the ALICE TTasks
200   // create the ALICE main TFolder
201   // to be done by AliRun
202
203   TFolder * alice = new TFolder();
204   alice->SetNameTitle("FPAlice", "Alice Folder") ;
205   gROOT->GetListOfBrowsables()->Add(alice) ;
206
207   TFolder * aliceF  = alice->AddFolder("folders", "Alice memory Folder") ;
208   //  make it the owner of the objects that it contains
209   aliceF->SetOwner() ;
210   // geometry folder
211   TFolder * geomF = aliceF->AddFolder("Geometry", "Geometry objects") ;
212   TFolder * aliceT  = alice->AddFolder("tasks", "Alice tasks Folder") ;   
213   //  make it the owner of the objects that it contains
214   aliceT->SetOwner() ;
215
216   TTask * aliceDi = new TTask("(S)Digitizer", "Alice SDigitizer & Digitizer") ;
217   aliceT->Add(aliceDi);
218
219   TTask * aliceRe = new TTask("Reconstructioner", "Alice Reconstructioner") ;
220   aliceT->Add(aliceRe);
221
222   char * tempo = new char[80] ;
223
224   // creates the TOF Digitizer and adds it to alice main (S)Digitizer task
225   sprintf(tempo, "%sDigitizers container",GetName() ) ;
226   fDTask = new TTask(GetName(), tempo);
227   aliceDi->Add(fDTask) ;
228
229   // creates the TOF reconstructioner and adds it to alice main Reconstructioner task
230   sprintf(tempo, "%sReconstructioner container",GetName() ) ;
231   fReTask = new TTask(GetName(), tempo);
232   aliceRe->Add(fReTask) ;
233
234   delete [] tempo ;
235  
236   // creates the TOF geometry  folder
237   geomF->AddFolder("TOF", "Geometry for TOF") ;
238 }
239
240 //_____________________________________________________________________________
241 AliTOF::~AliTOF()
242 {
243   // dtor:
244   // it remove also the alice folder 
245   // and task that TOF creates instead of AliRun
246   /* PH Temporarily commented because of problems
247   TFolder * alice = (TFolder*)gROOT->GetListOfBrowsables()->FindObject("FPAlice") ;
248   delete alice;
249   alice = 0;
250   */
251   if (fHits)
252     {
253       fHits->Delete ();
254       delete fHits;
255       fHits = 0;
256     }
257   if (fDigits)
258     {
259       fDigits->Delete ();
260       delete fDigits;
261       fDigits = 0;
262     }
263   if (fSDigits)
264     {
265       fSDigits->Delete();
266       delete fSDigits;
267       fSDigits = 0;
268     }
269   if (fReconParticles)
270     {
271       fReconParticles->Delete ();
272       delete fReconParticles;
273       fReconParticles = 0;
274     }
275
276 }
277
278 //_____________________________________________________________________________
279 void AliTOF::AddHit(Int_t track, Int_t *vol, Float_t *hits)
280 {
281   //
282   // Add a TOF hit
283   // new with placement used
284   //
285   TClonesArray &lhits = *fHits;
286   new(lhits[fNhits++]) AliTOFhit(fIshunt, track, vol, hits);
287 }
288
289 //_____________________________________________________________________________
290 void AliTOF::AddT0Hit(Int_t track, Int_t *vol, Float_t *hits)
291 {
292   //
293   // Add a TOF hit
294   // new with placement used
295   //
296   TClonesArray &lhits = *fHits;
297   new(lhits[fNhits++]) AliTOFhitT0(fIshunt, track, vol, hits);
298 }
299
300 //_____________________________________________________________________________
301 void AliTOF::AddDigit(Int_t *tracks, Int_t *vol, Float_t *digits)
302 {
303   //
304   // Add a TOF digit
305   // new with placement used
306   //
307   TClonesArray &ldigits = *fDigits;
308   new (ldigits[fNdigits++]) AliTOFdigit(tracks, vol, digits);
309 }
310
311 //___________________________________________
312 void AliTOF::AddSDigit(Int_t tracknum, Int_t *vol, Float_t *digits)
313 {
314      
315 //
316 // Add a TOF sdigit
317 //
318         
319   TClonesArray &lSDigits = *fSDigits;   
320   new(lSDigits[fNSDigits++]) AliTOFSDigit(tracknum, vol, digits);
321 }
322
323 //_____________________________________________________________________________
324 void AliTOF::SetTreeAddress ()
325 {
326   // Set branch address for the Hits and Digits Tree.
327   char branchname[30];
328   
329   if (fLoader->TreeH())
330    {
331      if (fHits == 0x0)
332       {
333         if (fTZero) fHits   = new TClonesArray("AliTOFhitT0", 1000);
334         else fHits   = new TClonesArray("AliTOFhit", 1000);
335       }
336    }
337   AliDetector::SetTreeAddress ();
338
339   TBranch *branch;
340   TTree *treeD = fLoader->TreeD();
341
342   if (treeD)
343     {
344       branch = treeD->GetBranch (branchname);
345       if (branch)
346        {
347          if (fDigits == 0x0) fDigits = new TClonesArray("AliTOFdigit",  1000); 
348          branch->SetAddress (&fDigits);
349        }
350     }
351
352   if (fLoader->TreeS () )
353     {
354       branch = fLoader->TreeS ()->GetBranch ("TOF");
355       if (branch) {
356         if (fSDigits == 0x0) fSDigits = new TClonesArray("AliTOFSDigit",  1000);
357         branch->SetAddress (&fSDigits);
358       }
359     }
360
361   if (fLoader->TreeR() && fReconParticles) //I do not know where this array is created - skowron
362     {
363       branch = fLoader->TreeR()->GetBranch("TOF"); 
364       if (branch) 
365        {
366          branch->SetAddress(&fReconParticles) ;
367        }
368     }   
369 }
370
371 //_____________________________________________________________________________
372 void AliTOF::CreateGeometry()
373 {
374   //
375   // Common geometry code 
376   //
377   //Begin_Html
378   /*
379     <img src="picts/AliTOFv23.gif">
380   */
381   //End_Html
382   //
383   const Double_t kPi=TMath::Pi();
384   const Double_t kDegrad=kPi/180.;
385   //
386   Float_t xTof, yTof, wall;
387
388   // frame inbetween TOF modules
389   wall = 4.;//cm
390
391   // Sizes of TOF module with its support etc..
392   xTof = 2.*(fRmin*TMath::Tan(10*kDegrad)-wall/2-.5);
393   yTof = fRmax-fRmin;
394
395 //  TOF module internal definitions 
396   TOFpc(xTof, yTof, fZlenC, fZlenB, fZlenA, fZtof);
397 }
398
399 //_____________________________________________________________________________
400 void AliTOF::DrawModule() const
401 {
402   //
403   // Draw a shaded view of the common part of the TOF geometry
404   //
405
406    cout << " Drawing of AliTOF"<< endl; 
407   // Set everything unseen
408   gMC->Gsatt("*", "seen", -1);
409   // 
410   // Set ALIC mother transparent
411   gMC->Gsatt("ALIC","SEEN",0);
412   //
413   // Set the volumes visible
414   gMC->Gsatt("FTOA","SEEN",1);
415   gMC->Gsatt("FTOB","SEEN",1);
416   gMC->Gsatt("FTOC","SEEN",1);
417   gMC->Gsatt("FLTA","SEEN",1);
418   gMC->Gsatt("FLTB","SEEN",1);
419   gMC->Gsatt("FLTC","SEEN",1);
420   gMC->Gsatt("FSTR","SEEN",1);
421   //
422   gMC->Gdopt("hide", "on");
423   gMC->Gdopt("shad", "on");
424   gMC->Gsatt("*", "fill", 7);
425   gMC->SetClipBox(".");
426   gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
427   gMC->DefaultRange();
428   gMC->Gdraw("alic", 40, 30, 0, 12, 9.5, .02, .02);
429   gMC->Gdhead(1111, "Time Of Flight");
430   gMC->Gdman(18, 4, "MAN");
431   gMC->Gdopt("hide","off");
432 }
433
434 //_____________________________________________________________________________
435 void AliTOF::CreateMaterials()
436 {
437   //
438   // Defines TOF materials for all versions
439   // Revision: F. Pierella 18-VI-2002
440   //
441
442   Int_t   isxfld = gAlice->Field()->Integ();
443   Float_t sxmgmx = gAlice->Field()->Max();
444   //
445   //--- Quartz (SiO2) to simulate float glass
446   //    density tuned to have correct float glass 
447   //    radiation length
448   Float_t   aq[2] = { 28.0855,15.9994 };
449   Float_t   zq[2] = { 14.,8. };
450   Float_t   wq[2] = { 1.,2. };
451   Float_t   dq = 2.55; // std value: 2.2
452   Int_t nq = -2;
453
454   // --- Freon C2F4H2 (TOF-TDR pagg.)
455   // Geant Manual CONS110-1, pag. 43 (Geant, Detector Description and Simulation Tool)
456   Float_t afre[3]  = {12.011,18.998,1.007};
457   Float_t zfre[3]  = { 6., 9., 1.}; 
458   Float_t wfre[3]  = { 2., 4., 2.};
459   Float_t densfre  = 0.00375;   
460 // http://www.fi.infn.it/sezione/prevprot/gas/freon.html
461   Int_t nfre = -3; 
462 /*
463   //-- Isobutane quencher C4H10 (5% in the sensitive mixture)
464   Float_t aiso[2]  = {12.011,1.007};
465   Float_t ziso[2]  = { 6.,  1.};
466   Float_t wiso[2]  = { 4., 10.};
467   Float_t densiso  = .......;  // (g/cm3) density
468   Int_t nfre = -2; // < 0 i.e. proportion by number of atoms of each kind
469   //-- SF6 (5% in the sensitive mixture)
470   Float_t asf[3]  = {32.066,18.998};
471   Float_t zsf[3]  = { 16., 9.};
472   Float_t wsf[3]  = {  1., 6.}; 
473   Float_t denssf  = .....;   // (g/cm3) density
474   Int_t nfre = -2; // < 0 i.e. proportion by number of atoms of each kind
475 */
476   // --- CO2 
477   Float_t ac[2]   = {12.,16.};
478   Float_t zc[2]   = { 6., 8.};
479   Float_t wc[2]   = { 1., 2.};
480   Float_t dc = .001977;
481   Int_t nc = -2;
482    // For mylar (C5H4O2) 
483   Float_t amy[3] = { 12., 1., 16. };
484   Float_t zmy[3] = {  6., 1.,  8. };
485   Float_t wmy[3] = {  5., 4.,  2. };
486   Float_t dmy    = 1.39;
487   Int_t nmy = -3;
488  // For polyethilene (CH2) - honeycomb -
489   Float_t ape[2] = { 12., 1. };
490   Float_t zpe[2] = {  6., 1. };
491   Float_t wpe[2] = {  1., 2. };
492   Float_t dpe    = 0.935*0.479; //To have 1%X0 for 1cm as for honeycomb
493   Int_t npe = -2;
494   // --- G10 
495   Float_t ag10[4] = { 12.,1.,16.,28. };
496   Float_t zg10[4] = {  6.,1., 8.,14. };
497   Float_t wmatg10[4] = { .259,.288,.248,.205 };
498   Float_t densg10  = 1.7;
499   Int_t nlmatg10 = -4;
500
501   // plexiglass CH2=C(CH3)CO2CH3
502   Float_t aplex[3] = { 12.,1.,16.};
503   Float_t zplex[3] = {  6.,1., 8.};
504   Float_t wmatplex[3] = {5.,8.,2.};
505   Float_t densplex  =1.16;
506   Int_t nplex = -3;
507
508   // ---- ALUMINA (AL203) 
509   Float_t aal[2] = { 27.,16.};
510   Float_t zal[2] = { 13., 8.};
511   Float_t wmatal[2] = { 2.,3. };
512   Float_t densal  = 2.3;
513   Int_t nlmatal = -2;
514   // -- Water
515   Float_t awa[2] = {  1., 16. };
516   Float_t zwa[2] = {  1.,  8. };
517   Float_t wwa[2] = {  2.,  1. };
518   Float_t dwa    = 1.0;
519   Int_t nwa = -2;
520
521 // stainless steel
522   Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
523   Float_t zsteel[4] = { 26.,24.,28.,14. };
524   Float_t wsteel[4] = { .715,.18,.1,.005 };
525
526   //
527   //AliMaterial(0, "Vacuum$", 1e-16, 1e-16, 1e-16, 1e16, 1e16);
528   AliMaterial( 1, "Air$",14.61,7.3,0.001205,30423.24,67500.);
529   AliMaterial( 2, "Cu $",  63.54, 29.0, 8.96, 1.43, 14.8);
530   AliMaterial( 3, "C  $",  12.01,  6.0, 2.265,18.8, 74.4);
531   AliMixture ( 4, "Polyethilene$", ape, zpe, dpe, npe, wpe);
532   AliMixture ( 5, "G10$", ag10, zg10, densg10, nlmatg10, wmatg10);
533   AliMixture ( 6, "PLE$", aplex, zplex, densplex, nplex, wmatplex);
534   AliMixture ( 7, "CO2$", ac, zc, dc, nc, wc);
535   AliMixture ( 8, "ALUMINA$", aal, zal, densal, nlmatal, wmatal);
536   AliMaterial( 9, "Al $", 26.98, 13., 2.7, 8.9, 37.2);
537   AliMaterial(10, "C-TRD$", 12.01, 6., 2.265*18.8/69.282*15./100, 18.8, 74.4); // for 15%
538   AliMixture (11, "Mylar$",  amy, zmy, dmy, nmy, wmy);
539   AliMixture (12, "Freon$",  afre, zfre, densfre, nfre, wfre);
540   AliMixture (13, "Glass$", aq, zq, dq, nq, wq);
541   AliMixture (14, "Water$",  awa, zwa, dwa, nwa, wwa);
542   AliMixture (15, "STAINLESS STEEL$", asteel, zsteel, 7.88, 4, wsteel);
543
544   Float_t epsil, stmin, deemax, stemax;
545  
546   //   Previous data
547   //       EPSIL  =  0.1   ! Tracking precision, 
548   //       STEMAX = 0.1      ! Maximum displacement for multiple scattering
549   //       DEEMAX = 0.1    ! Maximum fractional energy loss, DLS 
550   //       STMIN  = 0.1 
551   //
552   //   New data  
553   epsil  = .001;  // Tracking precision,
554   stemax = -1.;   // Maximum displacement for multiple scattering
555   deemax = -.3;   // Maximum fractional energy loss, DLS
556   stmin  = -.8;
557
558   AliMedium( 1, "Air$"  ,  1, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
559   AliMedium( 2, "Cu $"  ,  2, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
560   AliMedium( 3, "C  $"  ,  3, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
561   AliMedium( 4, "Pol$"  ,  4, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
562   AliMedium( 5, "G10$"  ,  5, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
563   AliMedium( 6, "PLE$"  ,  6, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
564   AliMedium( 7, "CO2$"  ,  7, 0, isxfld, sxmgmx, 10., -.01, -.1, .01, -.01);
565   AliMedium( 8,"ALUMINA$", 8, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
566   AliMedium( 9,"Al Frame$",9, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
567   AliMedium(10, "DME-S$",  6, 1, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
568   AliMedium(11, "C-TRD$", 10, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
569   AliMedium(12, "Myl$"  , 11, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
570   AliMedium(13, "Fre$"  , 12, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
571   AliMedium(14, "Fre-S$", 12, 1, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
572   AliMedium(15, "Glass$", 13, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
573   AliMedium(16, "Water$", 14, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
574   AliMedium(17, "STEEL$", 15, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
575 }
576
577 //_____________________________________________________________________________
578 Int_t AliTOF::DistancetoPrimitive(Int_t , Int_t ) const
579 {
580   //
581   // Returns distance from mouse pointer to detector, default version
582   //
583   return 9999;
584 }
585
586 //___________________________________________
587 void AliTOF::ResetHits ()
588 {
589   // Reset number of clusters and the cluster array for this detector
590   AliDetector::ResetHits ();
591 }
592
593 //____________________________________________
594 void AliTOF::ResetDigits ()
595 {
596   //
597   // Reset number of digits and the digits array for this detector
598   AliDetector::ResetDigits ();
599   //
600
601 //____________________________________________
602 void AliTOF::ResetSDigits ()
603 {
604   //
605   // Reset number of sdigits and the sdigits array for this detector
606   fNSDigits = 0;
607   //fSDigits = 0x0;
608   //
609
610 //_____________________________________________________________________________
611 void AliTOF::Init()
612 {
613   //
614   // Initialise TOF detector after it has been built
615   //
616   // Set id of TOF sensitive volume
617   if (IsVersion() !=0) fIdSens=gMC->VolId("FPAD");
618   //
619 }
620
621 //____________________________________________________________________________
622 void AliTOF::MakeBranch(Option_t* option)
623 {
624  //
625  // Initializes the Branches of the TOF inside the 
626  // trees written for each event. 
627  // AliDetector::MakeBranch initializes just the 
628  // Branch inside TreeH. Here we add the branches in 
629  // TreeD, TreeS and TreeR.
630  //
631   const char *oH = strstr(option,"H");
632   if (fLoader->TreeH() && oH)
633    {
634      if (fHits == 0x0)
635       {
636         if (fTZero) fHits   = new TClonesArray("AliTOFhitT0", 1000);
637         else fHits   = new TClonesArray("AliTOFhit", 1000);
638       }
639    }
640   
641   AliDetector::MakeBranch(option);
642
643   Int_t buffersize = 4000;
644   Char_t branchname[10];
645   sprintf(branchname,"%s",GetName());
646   
647   const char *oD = strstr(option,"D");
648   const char *oS = strstr(option,"S");
649   const char *oR = strstr(option,"R");
650
651   if (fLoader->TreeD() && oD){
652     if (fDigits == 0x0) fDigits = new TClonesArray("AliTOFdigit",  1000); 
653     MakeBranchInTree(fLoader->TreeD(), branchname, &fDigits,buffersize, 0) ;
654   }
655
656   if (fLoader->TreeS() && oS){
657     if (fSDigits == 0x0) fSDigits = new TClonesArray("AliTOFSDigit",  1000);
658     MakeBranchInTree(fLoader->TreeS(), branchname, &fSDigits,buffersize, 0) ;
659   }
660
661   if (fReconParticles && fLoader->TreeR() && oR){
662     MakeBranchInTree(fLoader->TreeR(), branchname, &fReconParticles,buffersize, 0) ;
663   }
664 }
665
666 //____________________________________________________________________________
667 void AliTOF::Makehits(Bool_t hits) 
668 {
669 // default argument used, see AliTOF.h
670 // Enable/Disable the writing of the TOF-hits branch 
671 // on TreeH
672 // by default : enabled for TOFv1, v2, v3, v4
673 //              disabled for TOFv0
674 // 
675    if (hits &&  (IsVersion()!=0))
676       fIdSens = gMC->VolId("FPAD");
677    else
678       cout << "Option for writing the TOF-hits branch on TreeH: disabled" << endl;
679 }
680
681 //____________________________________________________________________________
682 void AliTOF::FinishEvent()
683 {
684 // do nothing
685 }
686
687 //___________________________________________
688 void AliTOF::SDigits2Digits()
689 {
690 //
691 // Generate digits performing merging
692 //
693   /*
694     int nparticles = gAlice->GetNtrack();
695     cout << "Particles       :" <<nparticles<<endl;
696     if (nparticles > 0 ) {
697       
698       AliTOF::Hits2Digits();
699       
700     }
701   */
702   cout<<"AliTOF::SDigits2Digits"<<endl; 
703     if (fMerger) {
704       fMerger->Init();
705       cout<<"AliTOF::SDigits2Digits Init"<<endl; 
706       fMerger->Digitise();
707       cout<<"AliTOF::SDigits2Digits Digitise() "<<endl; 
708      }
709 }
710
711 //---------------------------------------------------------------------
712 void   AliTOF::SetMerger(AliTOFMerger* merger)
713 {
714 // Set pointer to merger
715     fMerger = merger;
716 }
717
718 //---------------------------------------------------------------------
719 AliTOFMerger*  AliTOF::Merger()
720 {
721 // Return pointer to merger
722     return fMerger;
723 }
724
725 //---------------------------------------------------------------------
726 void AliTOF::Hits2SDigits()
727 {
728 //
729 // Use the TOF SDigitizer to make TOF SDigits
730 //
731
732 //  cout<<"AliTOF::Hits2SDigits> start...\n";
733   
734   AliRunLoader * rl = fLoader->GetRunLoader();
735   AliTOFSDigitizer sd((rl->GetFileName()).Data());
736   if (GetDebug()) sd.Print("");
737
738   sd.Exec("") ;
739
740 }
741
742 //---------------------------------------------------------------------
743 void AliTOF::Hits2SDigits(Int_t evNumber1, Int_t evNumber2)
744 {
745 //
746 // Use the TOF SDigitizer to make TOF SDigits
747 //
748
749   if ((evNumber2-evNumber1)==1) 
750       cout << "<AliTOF::Hits2SDigits>: I am making sdigits for the " << evNumber1 << "th event \n";
751   else if ((evNumber2-evNumber1)>1)
752       cout << "<AliTOF::Hits2SDigits>: I am making sdigits for the events from the " 
753            << evNumber1 << "th to the " << evNumber2-1 << "th \n";
754  
755   AliRunLoader * rl = fLoader->GetRunLoader();
756   AliTOFSDigitizer sd((rl->GetFileName()).Data(),evNumber1,evNumber2) ;
757   if (GetDebug()) sd.Print("");
758
759   sd.Exec("") ;
760
761 }
762
763 //___________________________________________________________________________
764 AliDigitizer* AliTOF::CreateDigitizer(AliRunDigitizer* manager) const
765 {
766   return new AliTOFDigitizer(manager);
767 }
768
769 //___________________________________________________________________________
770 Bool_t AliTOF::CheckOverlap(Int_t* vol, Float_t* digit,Int_t Track)
771 {
772 //
773 // Checks if 2 or more hits belong to the same pad.
774 // In this case the data assigned to the digit object
775 // are the ones of the first hit in order of Time.
776 // 2 hits from the same track on the same pad are collected.
777 // Called only by Hits2SDigits.
778 // This procedure has to be optimized in the next TOF release.
779 //
780
781         Bool_t overlap = kFALSE;
782         Int_t  vol2[5];
783
784         for (Int_t ndig=0; ndig<fSDigits->GetEntries(); ndig++){
785            AliTOFdigit* currentDigit = (AliTOFdigit*)(fSDigits->UncheckedAt(ndig));
786            currentDigit->GetLocation(vol2);
787            Bool_t idem= kTRUE;
788            // check on digit volume
789            for (Int_t i=0;i<=4;i++){
790                if (!idem) break;
791                if (vol[i]!=vol2[i]) idem=kFALSE;}
792
793            if (idem){  // same pad fired
794               Float_t tdc2 = digit[0];
795               Float_t tdc1 = currentDigit->GetTdc();
796
797               // we separate two digits on the same pad if
798               // they are separated in time by at least 25 ns
799               // remember that tdc time is given in ps
800
801               if (TMath::Abs(tdc1-tdc2)<25000){
802                   // in case of overlap we take the earliest
803                   if (tdc1>tdc2){
804                    currentDigit->SetTdc(tdc2); 
805                    currentDigit->SetAdc(digit[1]);
806                   }
807                   else {
808                    currentDigit->SetTdc(tdc1);
809                    currentDigit->SetAdc(digit[1]);
810                   }
811                   currentDigit->AddTrack(Track); // add track number in the track array
812                   overlap = kTRUE;
813                   return overlap;
814               } else 
815                 overlap= kFALSE;
816
817            } // close if (idem) -> two digits on the same TOF pad
818
819         } // end loop on existing sdigits
820
821         return overlap;
822 }
823
824
825 //____________________________________________________________________________
826 void AliTOF::Digits2Raw(Int_t evNumber)
827 {
828 //
829 // Starting from digits, writes the 
830 // Raw Data objects, i.e. a 
831 // TClonesArray of 18 AliTOFRawSector objects
832 //
833
834   TTree* tD;
835
836   // do nothing if no particles
837   Int_t nparticles = gAlice->GetEvent(evNumber); 
838   if (nparticles <= 0) return;
839
840   tD = gAlice->TreeD();
841   
842   TClonesArray* tofdigits = this->Digits();
843   Int_t ndigits = tofdigits->GetEntriesFast();
844
845   TClonesArray* rawsectors = new TClonesArray("AliTOFRawSector",fNTof+2); 
846
847   for (Int_t isect=1;isect<=fNTof;isect++){
848      AliTOFRawSector* currentSector = (AliTOFRawSector*)rawsectors->UncheckedAt(isect);
849      TClonesArray* rocData = (TClonesArray*)currentSector->GetRocData();
850
851      for (Int_t digit=0; digit<ndigits; digit++){
852         AliTOFdigit* currentDigit = (AliTOFdigit*)tofdigits->UncheckedAt(digit);
853         Int_t sector = currentDigit->GetSector();
854         if (sector==isect){
855             Int_t   pad    = currentDigit -> GetTotPad();
856             Int_t   roc    = (Int_t)(pad/fNPadXRoc)-1;
857             if (roc>=fNRoc) printf("Wrong n. of ROC ! Roc = %i",roc);
858             Int_t   padRoc = (Int_t) pad%fNPadXRoc;
859             Int_t   fec    = (Int_t)(padRoc/fNFec)-1;
860             Int_t   tdc    = (Int_t)(padRoc%fNFec)-1;
861             Float_t time   = currentDigit->GetTdc();
862             Float_t charge = currentDigit->GetAdc();
863             AliTOFRoc* currentROC = (AliTOFRoc*)rocData->UncheckedAt(roc);
864             Int_t error    = 0;
865             currentROC->AddItem(fec, tdc, error, charge, time);
866         } // close if (sector==isect) i.e. end loop on digits for the current sector
867      } // end loop on TOF digits
868      
869      UInt_t totSize=16,rocSize=0;
870      UInt_t rocHead[14],rocChek[14];
871      UInt_t globalCheckSum=0;
872
873      for (UInt_t iRoc = 1; iRoc<(UInt_t)fNRoc; iRoc++){
874         AliTOFRoc* currentRoc = (AliTOFRoc*)rocData->UncheckedAt(iRoc); 
875         rocSize  = currentRoc->GetItems()*2+1;
876         totSize += rocSize*4;
877         if (rocSize>=TMath::Power(2,16)) rocSize=0;
878         rocHead[iRoc]   = iRoc<<28;
879         rocHead[iRoc]  += rocSize;
880         rocChek[iRoc]   = currentRoc->GetCheckSum();
881         Int_t headCheck = currentRoc->BitCount(rocHead[iRoc]);
882         globalCheckSum += headCheck;
883         globalCheckSum += rocChek[iRoc];
884      }
885      
886      AliTOFRoc* dummyRoc = new AliTOFRoc();
887      totSize *= 4;
888      if (totSize>=TMath::Power(2,24)) totSize=0;
889      UInt_t header = totSize;
890      UInt_t sectId = ((UInt_t)isect)<<24;
891      header += sectId;
892      globalCheckSum += dummyRoc->BitCount(header);
893      currentSector->SetGlobalCS(globalCheckSum);
894      currentSector->SetHeader(header);
895   }  
896 }
897  
898 //____________________________________________________________________________
899 void AliTOF::Raw2Digits(Int_t evNumber)
900 {
901 //
902 //  Converts Raw Data objects into digits objects.
903 //  We schematize the raw data with a 
904 //  TClonesArray of 18 AliTOFRawSector objects
905 //
906
907   TTree    *tD;
908   Int_t    vol[5];
909   Int_t    tracks[3];
910   Float_t  digit[2];
911  
912   tracks[0]=0;
913   tracks[1]=0;
914   tracks[2]=0;
915  
916   Int_t nparticles = gAlice->GetEvent(evNumber); 
917   if (nparticles <= 0) return;
918
919   tD = gAlice->TreeD();
920   
921   TClonesArray* rawsectors = new TClonesArray("AliTOFRawSector",fNTof+2);
922   
923   for(Int_t nSec=1; nSec<=fNTof; nSec++){
924      AliTOFRawSector* currentSector = (AliTOFRawSector*)rawsectors->UncheckedAt(nSec);
925      TClonesArray* rocData = (TClonesArray*)currentSector->GetRocData();
926      for(Int_t nRoc=1; nRoc<=14; nRoc++){
927         AliTOFRoc* currentRoc = (AliTOFRoc*)rocData->UncheckedAt(nRoc);
928         Int_t currentItems = currentRoc->GetItems();
929         for(Int_t item=1; item<currentItems; item++){ 
930            Int_t nPad = currentRoc->GetTotPad(item);        
931            vol[0] = nSec;
932            Int_t nStrip = (Int_t)(nPad/fPadXStr)+1;
933            Int_t nPlate = 5;
934            if (nStrip<=fNStripC+2*fNStripB+fNStripA) nPlate = 4;
935            if (nStrip<=fNStripC+fNStripB+fNStripA)   nPlate = 3;
936            if (nStrip<=fNStripC+fNStripB)            nPlate = 2;
937            if (nStrip<=fNStripC)                     nPlate = 1;
938            vol[1] = nPlate;
939            switch (nPlate){
940            case 1: break;
941            case 2: nStrip -= (fNStripC);
942                    break;
943            case 3: nStrip -= (fNStripC+fNStripB);
944                    break;
945            case 4: nStrip -= (fNStripC+fNStripB+fNStripA);
946                    break;
947            case 5: nStrip -= (fNStripC+2*fNStripB+fNStripA);
948                    break;
949            }
950            vol[2] = nStrip;
951            Int_t pad = nPad%fPadXStr;
952            if (pad==0) pad=fPadXStr;
953            Int_t nPadX=0, nPadZ=0;
954            (pad>fNpadX)? nPadX -= fNpadX : nPadX = pad ;
955            vol[3] = nPadX;
956            (pad>fNpadX)? nPadZ = 2 : nPadZ = 1 ;
957            vol[4] = nPadZ;
958            UInt_t error=0;
959            Float_t tdc = currentRoc->GetTime(item,error);
960            if (!error) digit[0]=tdc;
961            digit[1] = currentRoc->GetCharge(item);
962            AddDigit(tracks,vol,digit);
963         }
964      }
965   }
966   tD->Fill();
967   tD->Write(0,TObject::kOverwrite);
968
969
970 ////////////////////////////////////////////////////////////////////////
971 void AliTOF::RecreateSDigitsArray() {
972 //
973 // delete TClonesArray fSDigits and create it again
974 //  needed for backward compatability with PPR test production
975 //
976   delete fSDigits;
977   fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
978 }
979 ////////////////////////////////////////////////////////////////////////
980 void AliTOF::CreateSDigitsArray() {
981 //
982 // create TClonesArray fSDigits
983 //  needed for backward compatability with PPR test production
984 //
985   fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
986 }
987
988 ////////////////////////////////////////////////////////////////////////
989 AliTracker* AliTOF::CreateTracker() const
990 {
991 // create a TOF tracker
992
993   Double_t parPID[] = {130., 5.};
994   return new AliTOFtracker(GetGeometry(), parPID);
995 }
996