]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFSDigitizer.cxx
Updated version of TOF digitization, N^2 problem solved (J.Chudoba)
[u/mrichter/AliRoot.git] / TOF / AliTOFSDigitizer.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 // This is a TTask that constructs SDigits out of Hits
18 // A Summable Digits is the sum of all hits in a pad
19 // 
20 //
21 //-- Author: F. Pierella
22 //////////////////////////////////////////////////////////////////////////////
23
24
25 #include "TTask.h"
26 #include "TTree.h"
27 #include "TSystem.h"
28 #include "TFile.h"
29
30 #include "AliTOFHitMap.h"
31 #include "AliTOFSDigit.h"
32 #include "AliTOFhit.h"
33 #include "AliTOF.h"
34 #include "AliTOFv1.h"
35 #include "AliTOFv2.h"
36 #include "AliTOFv3.h"
37 #include "AliTOFv4.h"
38 #include "AliTOFSDigitizer.h"
39 #include "AliRun.h"
40 #include "AliDetector.h"
41 #include "AliMC.h"
42
43 #include "TFile.h"
44 #include "TTask.h"
45 #include "TTree.h"
46 #include "TSystem.h"
47 #include "TROOT.h"
48 #include "TFolder.h"
49 #include <stdlib.h>
50 #include <iostream.h>
51 #include <fstream.h>
52
53 ClassImp(AliTOFSDigitizer)
54
55 //____________________________________________________________________________ 
56   AliTOFSDigitizer::AliTOFSDigitizer():TTask("AliTOFSDigitizer","") 
57 {
58   // ctor
59   fNevents = 0 ;     
60 //  fSDigits = 0 ;
61   fHits = 0 ;
62
63 }
64            
65 //____________________________________________________________________________ 
66   AliTOFSDigitizer::AliTOFSDigitizer(char* HeaderFile,char *SdigitsFile ):TTask("AliTOFSDigitizer","") 
67 {
68   fNevents = 0 ;     // Number of events to digitize, 0 means all evens in current file
69   // add Task to //root/Tasks folder
70   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
71   roottasks->Add(this) ; 
72 }
73
74 //____________________________________________________________________________ 
75   AliTOFSDigitizer::~AliTOFSDigitizer()
76 {
77   // dtor
78 }
79
80 //____________________________________________________________________________
81 void AliTOFSDigitizer::Exec(Option_t *option) { 
82
83
84   AliTOF *TOF = (AliTOF *) gAlice->GetDetector ("TOF");
85
86   if (!TOF) {
87     Error("AliTOFSDigitizer","TOF not found");
88     return;
89   }
90
91   if (fNevents == 0)
92     fNevents = (Int_t) gAlice->TreeE()->GetEntries();
93
94   for (Int_t ievent = 0; ievent < fNevents; ievent++) {
95     gAlice->GetEvent(ievent);
96     TTree *TH = gAlice->TreeH ();
97     if (!TH)
98       return;
99     if (gAlice->TreeS () == 0)
100       gAlice->MakeTree ("S");
101
102       
103     //Make branches
104     char branchname[20];
105     sprintf (branchname, "%s", TOF->GetName ());
106     //Make branch for digits
107     TOF->MakeBranch ("S");
108     
109     //Now made SDigits from hits
110
111     Int_t    vol[5];       // location for a digit
112     Float_t  digit[2];     // TOF digit variables
113     TParticle *particle;
114     AliTOFhit *tofHit;
115     TClonesArray *TOFhits = TOF->Hits();
116
117 // create hit map
118     AliTOFHitMap *hitMap = new AliTOFHitMap(TOF->SDigits());
119
120     Int_t ntracks = static_cast<Int_t>(TH->GetEntries());
121     for (Int_t track = 0; track < ntracks; track++)
122     {
123       gAlice->ResetHits();
124       TH->GetEvent(track);
125       particle = gAlice->Particle(track);
126       Int_t nhits = TOFhits->GetEntriesFast();
127
128       for (Int_t hit = 0; hit < nhits; hit++)
129       {
130         tofHit = (AliTOFhit *) TOFhits->UncheckedAt(hit);
131         vol[0] = tofHit->GetSector();
132         vol[1] = tofHit->GetPlate();
133         vol[2] = tofHit->GetStrip();
134         vol[3] = tofHit->GetPadx();
135         vol[4] = tofHit->GetPadz();
136
137         // 95% of efficiency to be inserted here
138         // edge effect to be inserted here
139         // cross talk  to be inserted here
140
141         Float_t idealtime = tofHit->GetTof(); // unit s
142         idealtime *= 1.E+12;  // conversion from s to ps
143         // fTimeRes is given usually in ps
144         Float_t tdctime   = gRandom->Gaus(idealtime, TOF->GetTimeRes());
145         digit[0] = tdctime;
146
147         // typical Landau Distribution to be inserted here
148         // instead of Gaussian Distribution
149         Float_t idealcharge = tofHit->GetEdep();
150         Float_t adccharge = gRandom->Gaus(idealcharge, TOF->GetChrgRes());
151         digit[1] = adccharge;
152         Int_t tracknum = tofHit->GetTrack();
153
154         // check if two digit are on the same pad; in that case we sum
155         // the two or more digits
156         if (hitMap->TestHit(vol) != kEmpty) {
157           AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(hitMap->GetHit(vol));
158           sdig->Update(tdctime,adccharge,tracknum);
159         } else {
160           TOF->AddSDigit(tracknum, vol, digit);
161           hitMap->SetHit(vol);
162         }
163       } // end loop on hits for the current track
164     } // end loop on ntracks
165
166     delete hitMap;
167       
168     gAlice->TreeS()->Reset();
169     gAlice->TreeS()->Fill();
170     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
171   }                             //event loop
172
173
174 }
175  
176 //__________________________________________________________________
177 void AliTOFSDigitizer::SetSDigitsFile(char * file ){
178   if(!fSDigitsFile.IsNull())
179     cout << "Changing SDigits file from " <<(char *)fSDigitsFile.Data() << " to " << file << endl ;
180   fSDigitsFile=file ;
181 }
182 //__________________________________________________________________
183 void AliTOFSDigitizer::Print(Option_t* option)const
184 {
185   cout << "------------------- "<< GetName() << " -------------" << endl ;
186   if(fSDigitsFile.IsNull())
187     cout << " Writing SDigitis to file galice.root "<< endl ;
188   else
189     cout << "    Writing SDigitis to file  " << (char*) fSDigitsFile.Data() << endl ;
190
191 }