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