]>
Commit | Line | Data |
---|---|---|
517b7f8f | 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 "AliTOFdigit.h" | |
31 | #include "AliTOFhit.h" | |
32 | #include "AliTOF.h" | |
33 | #include "AliTOFv1.h" | |
34 | #include "AliTOFv2.h" | |
35 | #include "AliTOFv3.h" | |
36 | #include "AliTOFv4.h" | |
37 | #include "AliTOFSDigitizer.h" | |
38 | #include "AliRun.h" | |
39 | #include "AliDetector.h" | |
40 | #include "AliMC.h" | |
41 | ||
42 | #include "TFile.h" | |
43 | #include "TTask.h" | |
44 | #include "TTree.h" | |
45 | #include "TSystem.h" | |
46 | #include "TROOT.h" | |
47 | #include "TFolder.h" | |
48 | #include <stdlib.h> | |
49 | #include <iostream.h> | |
50 | #include <fstream.h> | |
51 | ||
52 | ClassImp(AliTOFSDigitizer) | |
53 | ||
54 | //____________________________________________________________________________ | |
55 | AliTOFSDigitizer::AliTOFSDigitizer():TTask("AliTOFSDigitizer","") | |
56 | { | |
57 | // ctor | |
58 | fNevents = 0 ; | |
59 | fSDigits = 0 ; | |
60 | fHits = 0 ; | |
61 | ||
62 | } | |
63 | ||
64 | //____________________________________________________________________________ | |
65 | AliTOFSDigitizer::AliTOFSDigitizer(char* HeaderFile,char *SdigitsFile ):TTask("AliTOFSDigitizer","") | |
66 | { | |
67 | fNevents = 0 ; // Number of events to digitize, 0 means all evens in current file | |
68 | // add Task to //root/Tasks folder | |
69 | TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; | |
70 | roottasks->Add(this) ; | |
71 | } | |
72 | ||
73 | //____________________________________________________________________________ | |
74 | AliTOFSDigitizer::~AliTOFSDigitizer() | |
75 | { | |
76 | // dtor | |
77 | } | |
78 | ||
79 | //____________________________________________________________________________ | |
80 | void AliTOFSDigitizer::Exec(Option_t *option) { | |
81 | ||
82 | ||
83 | // Initialise Hit array | |
84 | fHits = new TClonesArray ("AliTOFhit", 1000); | |
85 | fSDigits = new TClonesArray ("AliTOFdigit", 1000); | |
86 | ||
87 | AliTOF *TOF = (AliTOF *) gAlice->GetDetector ("TOF"); | |
88 | ||
89 | if (fNevents == 0) | |
90 | fNevents = (Int_t) gAlice->TreeE ()->GetEntries (); | |
91 | ||
92 | for (Int_t ievent = 0; ievent < fNevents; ievent++) | |
93 | { | |
94 | gAlice->GetEvent (ievent); | |
95 | if (gAlice->TreeH () == 0) | |
96 | return; | |
97 | if (gAlice->TreeS () == 0) | |
98 | gAlice->MakeTree ("S"); | |
99 | ||
100 | ||
101 | //Make branches | |
102 | char branchname[20]; | |
103 | sprintf (branchname, "%s", TOF->GetName ()); | |
104 | //Make branch for digits | |
105 | TOF->MakeBranch ("S"); | |
106 | ||
107 | //Now made SDigits from hits | |
108 | ||
109 | Int_t tracks[3]; // track info | |
110 | Int_t vol[5]; // location for a digit | |
111 | Float_t digit[2]; // TOF digit variables | |
112 | Int_t hit, nbytes; | |
113 | TParticle *particle; | |
114 | AliTOFhit *tofHit; | |
115 | TClonesArray *TOFhits = TOF->Hits(); | |
116 | ||
117 | ||
118 | // Event ------------------------- LOOP | |
119 | ||
120 | ||
121 | if (TOF) | |
122 | { | |
123 | TOFhits = TOF->Hits (); | |
124 | TTree *TH = gAlice->TreeH (); | |
125 | Stat_t ntracks = TH->GetEntries (); | |
126 | for (Int_t track = 0; track < ntracks; track++) | |
127 | { | |
128 | gAlice->ResetHits (); | |
129 | nbytes += TH->GetEvent (track); | |
130 | particle = gAlice->Particle (track); | |
131 | Int_t nhits = TOFhits->GetEntriesFast (); | |
132 | ||
133 | for (hit = 0; hit < nhits; hit++) | |
134 | { | |
135 | tofHit = (AliTOFhit *) TOFhits->UncheckedAt(hit); | |
136 | vol[0] = tofHit->GetSector(); | |
137 | vol[1] = tofHit->GetPlate(); | |
138 | vol[2] = tofHit->GetPadx(); | |
139 | vol[3] = tofHit->GetPadz(); | |
140 | vol[4] = tofHit->GetStrip(); | |
141 | ||
142 | // 95% of efficiency to be inserted here | |
143 | // edge effect to be inserted here | |
144 | // cross talk to be inserted here | |
145 | ||
146 | Float_t idealtime = tofHit->GetTof(); // unit s | |
147 | idealtime *= 1.E+12; // conversion from s to ps | |
148 | // fTimeRes is given usually in ps | |
149 | // Float_t tdctime = gRandom->Gaus(idealtime, fTimeRes); | |
150 | Float_t tdctime = gRandom->Gaus(idealtime, TOF->GetTimeRes()); | |
151 | digit[0] = tdctime; | |
152 | ||
153 | // typical Landau Distribution to be inserted here | |
154 | // instead of Gaussian Distribution | |
155 | Float_t idealcharge = tofHit->GetEdep(); | |
156 | // Float_t adccharge = gRandom->Gaus(idealcharge, fChrgRes); | |
157 | Float_t adccharge = gRandom->Gaus(idealcharge, TOF->GetChrgRes()); | |
158 | digit[1] = adccharge; | |
159 | Int_t tracknum = tofHit->GetTrack(); | |
160 | tracks[0] = tracknum; | |
161 | tracks[1] = 0; | |
162 | tracks[2] = 0; | |
163 | ||
164 | // check if two digit are on the same pad; in that case we sum | |
165 | // the two or more digits | |
166 | // Bool_t overlap = AliTOF::CheckOverlap(vol, digit, tracknum); | |
167 | Bool_t overlap = TOF->CheckOverlap(vol, digit, tracknum); | |
168 | if(!overlap) | |
169 | // new ((*fSDigits)[nSdigits++]) AliTOFdigit(tracks, vol, digit); | |
170 | TOF->AddSDigit(tracks, vol, digit); | |
171 | cout << "nSdigits" << endl; | |
172 | } // end loop on hits for the current track | |
173 | ||
174 | } // end loop on ntracks | |
175 | ||
176 | } // close if TOF switched ON | |
177 | ||
178 | gAlice->TreeS()->Reset(); | |
179 | gAlice->TreeS()->Fill(); | |
180 | gAlice->TreeS()->Write(0,TObject::kOverwrite) ; | |
181 | } //event loop | |
182 | ||
183 | ||
184 | } | |
185 | ||
186 | //__________________________________________________________________ | |
187 | void AliTOFSDigitizer::SetSDigitsFile(char * file ){ | |
188 | if(!fSDigitsFile.IsNull()) | |
189 | cout << "Changing SDigits file from " <<(char *)fSDigitsFile.Data() << " to " << file << endl ; | |
190 | fSDigitsFile=file ; | |
191 | } | |
192 | //__________________________________________________________________ | |
193 | void AliTOFSDigitizer::Print(Option_t* option)const | |
194 | { | |
195 | cout << "------------------- "<< GetName() << " -------------" << endl ; | |
196 | if(fSDigitsFile.IsNull()) | |
197 | cout << " Writing SDigitis to file galice.root "<< endl ; | |
198 | else | |
199 | cout << " Writing SDigitis to file " << (char*) fSDigitsFile.Data() << endl ; | |
200 | ||
201 | } | |
202 | ||
203 | ||
204 | ||
205 | ||
206 | ||
207 | ||
208 | ||
209 | ||
210 | ||
211 | ||
212 | ||
213 | ||
214 | ||
215 | ||
216 | ||
217 | ||
218 | ||
219 | ||
220 | ||
221 | ||
222 | ||
223 | ||
224 | ||
225 | ||
226 | ||
227 | ||
228 | ||
229 | ||
230 | ||
231 | ||
232 | ||
233 | ||
234 | ||
235 | ||
236 | ||
237 | ||
238 | ||
239 | ||
240 | ||
241 | ||
242 | ||
243 | ||
244 | ||
245 | ||
246 | ||
247 | ||
248 | ||
249 | ||
250 | ||
251 | ||
252 | ||
253 |