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