]>
Commit | Line | Data |
---|---|---|
b0f5e3fc | 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 | ||
b0f5e3fc | 16 | #include <iostream.h> |
17 | #include <stdlib.h> | |
18 | #include <stdio.h> | |
1ca7869b | 19 | #include <string.h> |
20 | ||
94de3818 | 21 | #include <TSystem.h> |
22 | #include <TROOT.h> | |
608f25d8 | 23 | #include <TStopwatch.h> |
ece86d9a | 24 | #include <TCanvas.h> |
25 | #include <TF1.h> | |
26 | #include <TRandom.h> | |
1ca7869b | 27 | #include <TH1.h> |
28 | #include <TFile.h> | |
29 | #include <TVector.h> | |
30 | #include <TArrayI.h> | |
31 | #include <TArrayF.h> | |
ece86d9a | 32 | |
b0f5e3fc | 33 | #include "AliRun.h" |
e8189707 | 34 | #include "AliITS.h" |
ece86d9a | 35 | #include "AliITShit.h" |
36 | #include "AliITSdigit.h" | |
37 | #include "AliITSmodule.h" | |
c7a4dac0 | 38 | #include "AliITSpList.h" |
e8189707 | 39 | #include "AliITSMapA1.h" |
40 | #include "AliITSMapA2.h" | |
e8189707 | 41 | #include "AliITSetfSDD.h" |
42 | #include "AliITSRawData.h" | |
b0f5e3fc | 43 | #include "AliITSHuffman.h" |
ece86d9a | 44 | #include "AliITSsegmentation.h" |
45 | #include "AliITSresponse.h" | |
c7a4dac0 | 46 | #include "AliITSsegmentationSDD.h" |
47 | #include "AliITSresponseSDD.h" | |
1ca7869b | 48 | #include "AliITSsimulationSDD.h" |
b0f5e3fc | 49 | |
b0f5e3fc | 50 | ClassImp(AliITSsimulationSDD) |
51 | //////////////////////////////////////////////////////////////////////// | |
52 | // Version: 0 | |
53 | // Written by Piergiorgio Cerello | |
54 | // November 23 1999 | |
55 | // | |
56 | // AliITSsimulationSDD is the simulation of SDDs. | |
57 | // | |
58 | //Begin_Html | |
59 | /* | |
60 | <img src="picts/ITS/AliITShit_Class_Diagram.gif"> | |
61 | </pre> | |
62 | <br clear=left> | |
63 | <font size=+2 color=red> | |
64 | <p>This show the relasionships between the ITS hit class and the rest of Aliroot. | |
65 | </font> | |
66 | <pre> | |
67 | */ | |
68 | //End_Html | |
8a33ae9e | 69 | //______________________________________________________________________ |
b0f5e3fc | 70 | Int_t power(Int_t b, Int_t e) { |
8a33ae9e | 71 | // compute b to the e power, where both b and e are Int_ts. |
72 | Int_t power = 1,i; | |
b0f5e3fc | 73 | |
8a33ae9e | 74 | for(i=0; i<e; i++) power *= b; |
75 | return power; | |
76 | } | |
77 | //______________________________________________________________________ | |
b0f5e3fc | 78 | void FastFourierTransform(AliITSetfSDD *alisddetf,Double_t *real, |
79 | Double_t *imag,Int_t direction) { | |
8a33ae9e | 80 | // Do a Fast Fourier Transform |
8a33ae9e | 81 | |
82 | Int_t samples = alisddetf->GetSamples(); | |
83 | Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5); | |
84 | Int_t m1 = samples; | |
85 | Int_t m = samples/2; | |
86 | Int_t m2 = samples/m1; | |
87 | Int_t i,j,k; | |
88 | for(i=1; i<=l; i++) { | |
89 | for(j=0; j<samples; j += m1) { | |
90 | Int_t p = 0; | |
91 | for(k=j; k<= j+m-1; k++) { | |
92 | Double_t wsr = alisddetf->GetWeightReal(p); | |
93 | Double_t wsi = alisddetf->GetWeightImag(p); | |
94 | if(direction == -1) wsi = -wsi; | |
95 | Double_t xr = *(real+k+m); | |
96 | Double_t xi = *(imag+k+m); | |
97 | *(real+k+m) = wsr*(*(real+k)-xr) - wsi*(*(imag+k)-xi); | |
98 | *(imag+k+m) = wsr*(*(imag+k)-xi) + wsi*(*(real+k)-xr); | |
99 | *(real+k) += xr; | |
100 | *(imag+k) += xi; | |
101 | p += m2; | |
102 | } // end for k | |
103 | } // end for j | |
104 | m1 = m; | |
105 | m /= 2; | |
106 | m2 += m2; | |
107 | } // end for i | |
b0f5e3fc | 108 | |
8a33ae9e | 109 | for(j=0; j<samples; j++) { |
110 | Int_t j1 = j; | |
111 | Int_t p = 0; | |
112 | Int_t i1; | |
113 | for(i1=1; i1<=l; i1++) { | |
114 | Int_t j2 = j1; | |
115 | j1 /= 2; | |
116 | p = p + p + j2 - j1 - j1; | |
117 | } // end for i1 | |
118 | if(p >= j) { | |
119 | Double_t xr = *(real+j); | |
120 | Double_t xi = *(imag+j); | |
121 | *(real+j) = *(real+p); | |
122 | *(imag+j) = *(imag+p); | |
123 | *(real+p) = xr; | |
124 | *(imag+p) = xi; | |
125 | } // end if p>=j | |
126 | } // end for j | |
127 | if(direction == -1) { | |
128 | for(i=0; i<samples; i++) { | |
129 | *(real+i) /= samples; | |
130 | *(imag+i) /= samples; | |
131 | } // end for i | |
132 | } // end if direction == -1 | |
b0f5e3fc | 133 | return; |
134 | } | |
8a33ae9e | 135 | //______________________________________________________________________ |
b0f5e3fc | 136 | AliITSsimulationSDD::AliITSsimulationSDD(){ |
8a33ae9e | 137 | // Default constructor |
138 | ||
139 | fResponse = 0; | |
140 | fSegmentation = 0; | |
141 | fHis = 0; | |
142 | fHitMap1 = 0; | |
143 | fHitMap2 = 0; | |
144 | fElectronics = 0; | |
145 | fStream = 0; | |
146 | fInZR = 0; | |
147 | fInZI = 0; | |
148 | fOutZR = 0; | |
149 | fOutZI = 0; | |
150 | fNofMaps = 0; | |
151 | fMaxNofSamples = 0; | |
152 | fITS = 0; | |
153 | fTreeB = 0; | |
154 | fD.Set(0); | |
155 | fT1.Set(0); | |
156 | fT2.Set(0); | |
157 | fTol.Set(0); | |
158 | fNoise.Set(0); | |
159 | fBaseline.Set(0); | |
160 | SetScaleFourier(); | |
161 | SetPerpendTracksFlag(); | |
162 | SetDoFFT(); | |
163 | SetCheckNoise(); | |
b0f5e3fc | 164 | } |
8a33ae9e | 165 | //______________________________________________________________________ |
166 | AliITSsimulationSDD::AliITSsimulationSDD(AliITSsimulationSDD &source){ | |
167 | // Copy constructor to satify Coding roules only. | |
168 | ||
169 | if(this==&source) return; | |
c7a4dac0 | 170 | Error("AliITSsimulationSSD","Not allowed to make a copy of " |
171 | "AliITSsimulationSDD Using default creater instead"); | |
8a33ae9e | 172 | AliITSsimulationSDD(); |
b0f5e3fc | 173 | } |
8a33ae9e | 174 | //______________________________________________________________________ |
c7a4dac0 | 175 | AliITSsimulationSDD& AliITSsimulationSDD::operator=(AliITSsimulationSDD &src){ |
8a33ae9e | 176 | // Assignment operator to satify Coding roules only. |
177 | ||
c7a4dac0 | 178 | if(this==&src) return *this; |
179 | Error("AliITSsimulationSSD","Not allowed to make a = with " | |
180 | "AliITSsimulationSDD Using default creater instead"); | |
8a33ae9e | 181 | return *this ; |
b0f5e3fc | 182 | } |
8a33ae9e | 183 | //______________________________________________________________________ |
c7a4dac0 | 184 | AliITSsimulationSDD::AliITSsimulationSDD(AliITSsegmentation *seg, |
185 | AliITSresponse *resp){ | |
186 | // Standard Constructor | |
187 | ||
188 | ||
189 | fResponse = 0; | |
190 | fSegmentation = 0; | |
191 | fHis = 0; | |
192 | fHitMap1 = 0; | |
193 | fHitMap2 = 0; | |
194 | fElectronics = 0; | |
195 | fStream = 0; | |
196 | fInZR = 0; | |
197 | fInZI = 0; | |
198 | fOutZR = 0; | |
199 | fOutZI = 0; | |
200 | fNofMaps = 0; | |
201 | fMaxNofSamples = 0; | |
202 | fITS = 0; | |
203 | fTreeB = 0; | |
204 | ||
205 | Init((AliITSsegmentationSDD*)seg,(AliITSresponseSDD*)resp); | |
206 | } | |
207 | //______________________________________________________________________ | |
208 | void AliITSsimulationSDD::Init(AliITSsegmentationSDD *seg, | |
209 | AliITSresponseSDD *resp){ | |
38867c90 | 210 | // Standard Constructor |
ece86d9a | 211 | |
38867c90 | 212 | fResponse = resp; |
213 | fSegmentation = seg; | |
214 | SetScaleFourier(); | |
215 | SetPerpendTracksFlag(); | |
216 | SetDoFFT(); | |
217 | SetCheckNoise(); | |
b0f5e3fc | 218 | |
38867c90 | 219 | fHitMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1); |
220 | fHitMap1 = new AliITSMapA1(fSegmentation); | |
b0f5e3fc | 221 | |
8a33ae9e | 222 | fNofMaps = fSegmentation->Npz(); |
223 | fMaxNofSamples = fSegmentation->Npx(); | |
224 | ||
38867c90 | 225 | Float_t sddLength = fSegmentation->Dx(); |
8a33ae9e | 226 | Float_t sddWidth = fSegmentation->Dz(); |
b0f5e3fc | 227 | |
8a33ae9e | 228 | Int_t dummy = 0; |
38867c90 | 229 | Float_t anodePitch = fSegmentation->Dpz(dummy); |
8a33ae9e | 230 | Double_t timeStep = (Double_t)fSegmentation->Dpx(dummy); |
231 | Float_t driftSpeed = fResponse->DriftSpeed(); | |
b0f5e3fc | 232 | |
38867c90 | 233 | if(anodePitch*(fNofMaps/2) > sddWidth) { |
234 | Warning("AliITSsimulationSDD", | |
235 | "Too many anodes %d or too big pitch %f \n", | |
236 | fNofMaps/2,anodePitch); | |
237 | } // end if | |
b0f5e3fc | 238 | |
38867c90 | 239 | if(timeStep*fMaxNofSamples < sddLength/driftSpeed) { |
240 | Error("AliITSsimulationSDD", | |
241 | "Time Interval > Allowed Time Interval: exit\n"); | |
242 | return; | |
243 | } // end if | |
b0f5e3fc | 244 | |
38867c90 | 245 | fElectronics = new AliITSetfSDD(timeStep/fScaleSize, |
246 | fResponse->Electronics()); | |
b0f5e3fc | 247 | |
38867c90 | 248 | char opt1[20], opt2[20]; |
249 | fResponse->ParamOptions(opt1,opt2); | |
8a33ae9e | 250 | fParam = opt2; |
38867c90 | 251 | char *same = strstr(opt1,"same"); |
252 | if (same) { | |
253 | fNoise.Set(0); | |
254 | fBaseline.Set(0); | |
255 | } else { | |
256 | fNoise.Set(fNofMaps); | |
257 | fBaseline.Set(fNofMaps); | |
258 | } // end if | |
b0f5e3fc | 259 | |
38867c90 | 260 | const char *kopt=fResponse->ZeroSuppOption(); |
261 | if (strstr(fParam,"file") ) { | |
262 | fD.Set(fNofMaps); | |
263 | fT1.Set(fNofMaps); | |
264 | if (strstr(kopt,"2D")) { | |
b0f5e3fc | 265 | fT2.Set(fNofMaps); |
266 | fTol.Set(0); | |
267 | Init2D(); // desactivate if param change module by module | |
38867c90 | 268 | } else if(strstr(kopt,"1D")) { |
b0f5e3fc | 269 | fT2.Set(2); |
270 | fTol.Set(2); | |
271 | Init1D(); // desactivate if param change module by module | |
38867c90 | 272 | } // end if strstr |
273 | } else { | |
274 | fD.Set(2); | |
275 | fTol.Set(2); | |
276 | fT1.Set(2); | |
277 | fT2.Set(2); | |
278 | SetCompressParam(); | |
279 | } // end if else strstr | |
b0f5e3fc | 280 | |
8a33ae9e | 281 | Bool_t write = fResponse->OutputOption(); |
38867c90 | 282 | if(write && strstr(kopt,"2D")) MakeTreeB(); |
b0f5e3fc | 283 | |
38867c90 | 284 | // call here if baseline does not change by module |
285 | // ReadBaseline(); | |
b0f5e3fc | 286 | |
8a33ae9e | 287 | fITS = (AliITS*)gAlice->GetModule("ITS"); |
288 | Int_t size = fNofMaps*fMaxNofSamples; | |
289 | fStream = new AliITSInStream(size); | |
b0f5e3fc | 290 | |
38867c90 | 291 | fInZR = new Double_t [fScaleSize*fMaxNofSamples]; |
292 | fInZI = new Double_t [fScaleSize*fMaxNofSamples]; | |
293 | fOutZR = new Double_t [fScaleSize*fMaxNofSamples]; | |
294 | fOutZI = new Double_t [fScaleSize*fMaxNofSamples]; | |
5d18fa90 | 295 | |
b0f5e3fc | 296 | } |
8a33ae9e | 297 | //______________________________________________________________________ |
b0f5e3fc | 298 | AliITSsimulationSDD::~AliITSsimulationSDD() { |
8a33ae9e | 299 | // destructor |
300 | ||
301 | delete fHitMap1; | |
302 | delete fHitMap2; | |
303 | delete fStream; | |
304 | delete fElectronics; | |
305 | ||
306 | fD.Set(0); | |
307 | fT1.Set(0); | |
308 | fT2.Set(0); | |
309 | fTol.Set(0); | |
310 | fNoise.Set(0); | |
311 | fBaseline.Set(0); | |
312 | fITS = 0; | |
313 | ||
314 | if (fHis) { | |
315 | fHis->Delete(); | |
316 | delete fHis; | |
317 | } // end if fHis | |
318 | if(fTreeB) delete fTreeB; | |
319 | if(fInZR) delete [] fInZR; | |
320 | if(fInZI) delete [] fInZI; | |
321 | if(fOutZR) delete [] fOutZR; | |
322 | if(fOutZI) delete [] fOutZI; | |
b0f5e3fc | 323 | } |
8a33ae9e | 324 | //______________________________________________________________________ |
c7a4dac0 | 325 | void AliITSsimulationSDD::SDigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){ |
326 | // create maps to build the lists of tracks for each summable digit | |
327 | ||
328 | TObjArray *fHits = mod->GetHits(); | |
329 | Int_t nhits = fHits->GetEntriesFast(); | |
330 | fModule = md; | |
331 | fEvent = ev; | |
332 | ||
333 | if(!nhits) return; | |
334 | ||
335 | AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(), | |
336 | fScaleSize*fSegmentation->Npx()); | |
337 | ||
338 | // inputs to ListOfFiredCells. | |
339 | TObjArray *alist = new TObjArray(); | |
340 | fHitMap1->SetArray(alist); | |
341 | static TClonesArray *padr = 0; | |
342 | if(!padr) padr = new TClonesArray("TVector",1000); | |
343 | ||
344 | HitsToAnalogDigits(mod,alist,padr,pList); | |
345 | ||
346 | WriteSDigits(pList); | |
347 | ||
348 | // clean memory | |
2fa318ee | 349 | delete pList; |
c7a4dac0 | 350 | alist->Delete(); |
351 | delete alist; | |
352 | padr->Delete(); | |
353 | fHitMap1->ClearMap(); | |
354 | fHitMap2->ClearMap(); | |
355 | } | |
356 | //______________________________________________________________________ | |
b0f5e3fc | 357 | void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){ |
8a33ae9e | 358 | // create maps to build the lists of tracks for each digit |
b0f5e3fc | 359 | |
360 | TObjArray *fHits = mod->GetHits(); | |
8a33ae9e | 361 | Int_t nhits = fHits->GetEntriesFast(); |
362 | fModule = md; | |
363 | fEvent = ev; | |
364 | ||
ece86d9a | 365 | if (!nhits && fCheckNoise) { |
366 | ChargeToSignal(); | |
367 | GetNoise(); | |
368 | fHitMap2->ClearMap(); | |
369 | return; | |
370 | } else if (!nhits) return; | |
371 | ||
c7a4dac0 | 372 | AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(), |
373 | fScaleSize*fSegmentation->Npx()); | |
374 | ||
375 | // inputs to ListOfFiredCells. | |
376 | TObjArray *alist = new TObjArray(); | |
377 | fHitMap1->SetArray(alist); | |
378 | static TClonesArray *padr = 0; | |
379 | if(!padr) padr = new TClonesArray("TVector",1000); | |
380 | ||
381 | HitsToAnalogDigits(mod,alist,padr,pList); | |
382 | ||
383 | FinishDigits(alist); | |
384 | ||
385 | // clean memory | |
2fa318ee | 386 | delete pList; |
c7a4dac0 | 387 | alist->Delete(); |
388 | delete alist; | |
389 | padr->Delete(); | |
390 | fHitMap1->ClearMap(); | |
391 | fHitMap2->ClearMap(); | |
392 | } | |
393 | //______________________________________________________________________ | |
394 | void AliITSsimulationSDD::SDigitsToDigits(AliITSpList *pList){ | |
395 | // Take Summable digits and create Digits. | |
396 | ||
8a33ae9e | 397 | // inputs to ListOfFiredCells. |
c7a4dac0 | 398 | TObjArray *alist = new TObjArray(); |
8a33ae9e | 399 | fHitMap1->SetArray(alist); |
400 | static TClonesArray *padr = 0; | |
401 | if(!padr) padr = new TClonesArray("TVector",1000); | |
c7a4dac0 | 402 | Int_t arg[6] = {0,0,0,0,0,0}; |
403 | Double_t timeAmplitude; | |
404 | Int_t i,j; | |
405 | ||
406 | // Fill maps from pList. | |
407 | for(i=0;i<pList->GetMaxIndex();i++){ | |
408 | pList->GetMapIndex(i,arg[0],arg[1]); | |
409 | for(j=0;j<pList->GetNEnteries();j++){ | |
410 | timeAmplitude = pList->GetTSignal(arg[0],arg[1],j); | |
411 | if(timeAmplitude>0.0) continue; | |
412 | arg[2] = pList->GetTrack(arg[0],arg[1],j); | |
413 | arg[3] = pList->GetHit(arg[0],arg[1],j); | |
414 | ListOfFiredCells(arg,timeAmplitude,alist,padr); | |
415 | } // end for j | |
416 | // Make sure map has full signal in it. | |
417 | fHitMap2->SetHit(arg[0],arg[1],pList->GetSignal(arg[0],arg[1])); | |
418 | } // end for i | |
419 | ||
420 | FinishDigits(alist); | |
421 | ||
422 | // clean memory | |
423 | alist->Delete(); | |
424 | delete alist; | |
425 | padr->Delete(); | |
426 | fHitMap1->ClearMap(); | |
427 | fHitMap2->ClearMap(); | |
428 | } | |
429 | //______________________________________________________________________ | |
430 | void AliITSsimulationSDD::FinishDigits(TObjArray *alist){ | |
431 | // introduce the electronics effects and do zero-suppression if required | |
432 | Int_t nentries=alist->GetEntriesFast(); | |
8a33ae9e | 433 | |
c7a4dac0 | 434 | if(!nentries) return; |
435 | ChargeToSignal(); | |
436 | const char *kopt=fResponse->ZeroSuppOption(); | |
437 | ZeroSuppression(kopt); | |
438 | } | |
439 | //______________________________________________________________________ | |
440 | void AliITSsimulationSDD::HitsToAnalogDigits(AliITSmodule *mod,TObjArray *alst, | |
441 | TClonesArray *padr, | |
442 | AliITSpList *pList){ | |
443 | // create maps to build the lists of tracks for each digit | |
444 | ||
445 | TObjArray *fHits = mod->GetHits(); | |
446 | Int_t nhits = fHits->GetEntriesFast(); | |
447 | Int_t arg[6] = {0,0,0,0,0,0}; | |
8a33ae9e | 448 | Int_t dummy = 0; |
449 | Int_t nofAnodes = fNofMaps/2; | |
450 | Float_t sddLength = fSegmentation->Dx(); | |
451 | Float_t sddWidth = fSegmentation->Dz(); | |
452 | Float_t anodePitch = fSegmentation->Dpz(dummy); | |
453 | Float_t timeStep = fSegmentation->Dpx(dummy); | |
454 | Float_t driftSpeed = fResponse->DriftSpeed(); | |
455 | Float_t maxadc = fResponse->MaxAdc(); | |
456 | Float_t topValue = fResponse->DynamicRange(); | |
457 | Float_t cHloss = fResponse->ChargeLoss(); | |
458 | Float_t norm = maxadc/topValue; | |
459 | Float_t dfCoeff, s1; fResponse->DiffCoeff(dfCoeff,s1); // Signal 2d Shape | |
460 | Double_t eVpairs = 3.6; // electron pair energy eV. | |
461 | Float_t nsigma = fResponse->NSigmaIntegration(); // | |
462 | Int_t nlookups = fResponse->GausNLookUp(); // | |
44a312c3 | 463 | |
b0f5e3fc | 464 | // Piergiorgio's part (apart for few variables which I made float |
465 | // when i thought that can be done | |
b0f5e3fc | 466 | // Fill detector maps with GEANT hits |
467 | // loop over hits in the module | |
468 | ||
8a33ae9e | 469 | const Float_t kconv = 1.0e+6; // GeV->KeV |
470 | Int_t itrack = 0; | |
471 | Int_t hitDetector; // detector number (lay,lad,hitDetector) | |
472 | Int_t iWing; // which detector wing/side. | |
473 | Int_t detector; // 2*(detector-1)+iWing | |
474 | Int_t ii,kk,ka,kt; // loop indexs | |
475 | Int_t ia,it,index; // sub-pixel integration indexies | |
476 | Int_t iAnode; // anode number. | |
477 | Int_t timeSample; // time buckett. | |
478 | Int_t anodeWindow; // anode direction charge integration width | |
479 | Int_t timeWindow; // time direction charge integration width | |
480 | Int_t jamin,jamax; // anode charge integration window | |
481 | Int_t jtmin,jtmax; // time charge integration window | |
482 | Int_t ndiv; // Anode window division factor. | |
483 | Int_t nsplit; // the number of splits in anode and time windows==1. | |
484 | Int_t nOfSplits; // number of times track length is split into | |
485 | Float_t nOfSplitsF; // Floating point version of nOfSplits. | |
486 | Float_t kkF; // Floating point version of loop index kk. | |
487 | Float_t pathInSDD; // Track length in SDD. | |
488 | Float_t drPath; // average position of track in detector. in microns | |
489 | Float_t drTime; // Drift time | |
490 | Float_t nmul; // drift time window multiplication factor. | |
491 | Float_t avDrft; // x position of path length segment in cm. | |
492 | Float_t avAnode; // Anode for path length segment in Anode number (float) | |
493 | Float_t xAnode; // Floating point anode number. | |
494 | Float_t driftPath; // avDrft in microns. | |
495 | Float_t width; // width of signal at anodes. | |
496 | Double_t depEnergy; // Energy deposited in this GEANT step. | |
497 | Double_t xL[3],dxL[3]; // local hit coordinates and diff. | |
498 | Double_t sigA; // sigma of signal at anode. | |
499 | Double_t sigT; // sigma in time/drift direction for track segment | |
500 | Double_t aStep,aConst; // sub-pixel size and offset anode | |
501 | Double_t tStep,tConst; // sub-pixel size and offset time | |
502 | Double_t amplitude; // signal amplitude for track segment in nanoAmpere | |
503 | Double_t chargeloss; // charge loss for track segment. | |
504 | Double_t anodeAmplitude; // signal amplitude in anode direction | |
505 | Double_t aExpo; // exponent of Gaussian anode direction | |
506 | Double_t timeAmplitude; // signal amplitude in time direction | |
507 | Double_t tExpo; // exponent of Gaussian time direction | |
508 | // Double_t tof; // Time of flight in ns of this step. | |
44a312c3 | 509 | |
b0f5e3fc | 510 | for(ii=0; ii<nhits; ii++) { |
8a33ae9e | 511 | if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2], |
512 | depEnergy,itrack)) continue; | |
513 | depEnergy *= kconv; | |
514 | hitDetector = mod->GetDet(); | |
515 | //tof = 1.E+09*(mod->GetHit(ii)->GetTOF()); // tof in ns. | |
516 | //if(tof>sddLength/driftSpeed) continue; // hit happed too late. | |
517 | ||
518 | // scale path to simulate a perpendicular track | |
519 | // continue if the particle did not lose energy | |
520 | // passing through detector | |
521 | if (!depEnergy) { | |
c7a4dac0 | 522 | Warning("HitsToAnalogDigits", |
523 | "This particle has passed without losing energy!"); | |
8a33ae9e | 524 | continue; |
525 | } // end if !depEnergy | |
526 | ||
527 | pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]); | |
528 | ||
529 | if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); } | |
530 | drPath = 10000.*(dxL[0]+2.*xL[0])*0.5; | |
531 | if(drPath < 0) drPath = -drPath; | |
532 | drPath = sddLength-drPath; | |
533 | if(drPath < 0) { | |
c7a4dac0 | 534 | Warning("HitsToAnalogDigits","negative drift path %e",drPath); |
8a33ae9e | 535 | continue; |
536 | } // end if drPath < 0 | |
537 | ||
538 | // Compute number of segments to brake step path into | |
539 | drTime = drPath/driftSpeed; // Drift Time | |
540 | sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes | |
541 | // calcuate the number of time the path length should be split into. | |
542 | nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA); | |
543 | if(fFlag) nOfSplits = 1; | |
544 | ||
545 | // loop over path segments, init. some variables. | |
546 | depEnergy /= nOfSplits; | |
547 | nOfSplitsF = (Float_t) nOfSplits; | |
548 | for(kk=0;kk<nOfSplits;kk++) { // loop over path segments | |
549 | kkF = (Float_t) kk + 0.5; | |
550 | avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF; | |
551 | avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF; | |
552 | driftPath = 10000.*avDrft; | |
553 | ||
554 | iWing = 2; // Assume wing is 2 | |
555 | if(driftPath < 0) { // if wing is not 2 it is 1. | |
556 | iWing = 1; | |
557 | driftPath = -driftPath; | |
558 | } // end if driftPath < 0 | |
559 | driftPath = sddLength-driftPath; | |
560 | detector = 2*(hitDetector-1) + iWing; | |
561 | if(driftPath < 0) { | |
c7a4dac0 | 562 | Warning("HitsToAnalogDigits","Warning: negative drift path %e", |
563 | driftPath); | |
8a33ae9e | 564 | continue; |
565 | } // end if driftPath < 0 | |
566 | ||
567 | // Drift Time | |
568 | drTime = driftPath/driftSpeed; // drift time for segment. | |
569 | timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1); | |
570 | // compute time Sample including tof information. The tof only | |
571 | // effects the time of the signal is recoreded and not the | |
572 | // the defusion. | |
573 | // timeSample = (Int_t) (fScaleSize*(drTime+tof)/timeStep + 1); | |
574 | if(timeSample > fScaleSize*fMaxNofSamples) { | |
c7a4dac0 | 575 | Warning("HItsToAnalogDigits","Wrong Time Sample: %e", |
576 | timeSample); | |
8a33ae9e | 577 | continue; |
578 | } // end if timeSample > fScaleSize*fMaxNoofSamples | |
579 | ||
580 | // Anode | |
581 | xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1? | |
582 | if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.) | |
c7a4dac0 | 583 | Warning("HitsToAnalogDigits","Z = %e", |
584 | xAnode*anodePitch); | |
8a33ae9e | 585 | iAnode = (Int_t) (1.+xAnode); // xAnode? |
586 | if(iAnode < 1 || iAnode > nofAnodes) { | |
c7a4dac0 | 587 | Warning("HitToAnalogDigits","Wrong iAnode: %d",iAnode); |
8a33ae9e | 588 | continue; |
589 | } // end if iAnode < 1 || iAnode > nofAnodes | |
590 | ||
591 | // store straight away the particle position in the array | |
592 | // of particles and take idhit=ii only when part is entering (this | |
593 | // requires FillModules() in the macro for analysis) : | |
b0f5e3fc | 594 | |
8a33ae9e | 595 | // Sigma along the anodes for track segment. |
596 | sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1); | |
597 | sigT = sigA/driftSpeed; | |
598 | // Peak amplitude in nanoAmpere | |
599 | amplitude = fScaleSize*160.*depEnergy/ | |
600 | (timeStep*eVpairs*2.*acos(-1.)*sigT*sigA); | |
601 | amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to | |
602 | // account for clock variations | |
603 | // (reference value: 40 MHz) | |
604 | chargeloss = 1.-cHloss*driftPath/1000; | |
605 | amplitude *= chargeloss; | |
606 | width = 2.*nsigma/(nlookups-1); | |
607 | // Spread the charge | |
608 | // Pixel index | |
609 | ndiv = 2; | |
610 | nmul = 3.; | |
611 | if(drTime > 1200.) { | |
612 | ndiv = 4; | |
613 | nmul = 1.5; | |
614 | } // end if drTime > 1200. | |
615 | // Sub-pixel index | |
616 | nsplit = 4; // hard-wired //nsplit=4;nsplit = (nsplit+1)/2*2; | |
617 | // Sub-pixel size see computation of aExpo and tExpo. | |
618 | aStep = anodePitch/(nsplit*fScaleSize*sigA); | |
619 | aConst = xAnode*anodePitch/sigA; | |
620 | tStep = timeStep/(nsplit*fScaleSize*sigT); | |
621 | tConst = drTime/sigT; | |
622 | // Define SDD window corresponding to the hit | |
623 | anodeWindow = (Int_t)(fScaleSize*nsigma*sigA/anodePitch+1); | |
624 | timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.); | |
625 | jamin = (iAnode - anodeWindow/ndiv - 1)*fScaleSize*nsplit +1; | |
626 | jamax = (iAnode + anodeWindow/ndiv)*fScaleSize*nsplit; | |
627 | if(jamin <= 0) jamin = 1; | |
628 | if(jamax > fScaleSize*nofAnodes*nsplit) | |
629 | jamax = fScaleSize*nofAnodes*nsplit; | |
630 | // jtmin and jtmax are Hard-wired | |
631 | jtmin = (Int_t)(timeSample-timeWindow*nmul-1)*nsplit+1; | |
632 | jtmax = (Int_t)(timeSample+timeWindow*nmul)*nsplit; | |
633 | if(jtmin <= 0) jtmin = 1; | |
634 | if(jtmax > fScaleSize*fMaxNofSamples*nsplit) | |
635 | jtmax = fScaleSize*fMaxNofSamples*nsplit; | |
636 | // Spread the charge in the anode-time window | |
637 | for(ka=jamin; ka <=jamax; ka++) { | |
638 | ia = (ka-1)/(fScaleSize*nsplit) + 1; | |
c7a4dac0 | 639 | if(ia <= 0) { |
640 | Warning("HitsToAnalogDigits","ia < 1: "); | |
641 | continue; | |
642 | } // end if | |
8a33ae9e | 643 | if(ia > nofAnodes) ia = nofAnodes; |
644 | aExpo = (aStep*(ka-0.5)-aConst); | |
645 | if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.; | |
646 | else { | |
647 | dummy = (Int_t) ((aExpo+nsigma)/width); | |
648 | anodeAmplitude = amplitude*fResponse->GausLookUp(dummy); | |
649 | } // end if TMath::Abs(aEspo) > nsigma | |
650 | // index starts from 0 | |
651 | index = ((detector+1)%2)*nofAnodes+ia-1; | |
652 | if(anodeAmplitude) for(kt=jtmin; kt<=jtmax; kt++) { | |
653 | it = (kt-1)/nsplit+1; // it starts from 1 | |
c7a4dac0 | 654 | if(it<=0){ |
655 | Warning("HitsToAnalogDigits","it < 1:"); | |
656 | continue; | |
657 | } // end if | |
8a33ae9e | 658 | if(it>fScaleSize*fMaxNofSamples) |
659 | it = fScaleSize*fMaxNofSamples; | |
660 | tExpo = (tStep*(kt-0.5)-tConst); | |
661 | if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.; | |
662 | else { | |
663 | dummy = (Int_t) ((tExpo+nsigma)/width); | |
664 | timeAmplitude = anodeAmplitude* | |
665 | fResponse->GausLookUp(dummy); | |
666 | } // end if TMath::Abs(tExpo) > nsigma | |
667 | // build the list of digits for this module | |
668 | arg[0] = index; | |
669 | arg[1] = it; | |
670 | arg[2] = itrack; // track number | |
671 | arg[3] = ii-1; // hit number. | |
672 | timeAmplitude *= norm; | |
673 | timeAmplitude *= 10; | |
c7a4dac0 | 674 | ListOfFiredCells(arg,timeAmplitude,alst,padr); |
675 | pList->AddSignal(index,it,itrack,ii-1, | |
676 | mod->GetIndex(),timeAmplitude); | |
8a33ae9e | 677 | } // end if anodeAmplitude and loop over time in window |
678 | } // loop over anodes in window | |
679 | } // end loop over "sub-hits" | |
44a312c3 | 680 | } // end loop over hits |
b0f5e3fc | 681 | } |
8a33ae9e | 682 | //______________________________________________________________________ |
b0f5e3fc | 683 | void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude, |
8a33ae9e | 684 | TObjArray *alist,TClonesArray *padr){ |
685 | // Returns the list of "fired" cells. | |
686 | ||
687 | Int_t index = arg[0]; | |
688 | Int_t ik = arg[1]; | |
689 | Int_t idtrack = arg[2]; | |
690 | Int_t idhit = arg[3]; | |
691 | Int_t counter = arg[4]; | |
692 | Int_t countadr = arg[5]; | |
693 | Double_t charge = timeAmplitude; | |
694 | charge += fHitMap2->GetSignal(index,ik-1); | |
695 | fHitMap2->SetHit(index, ik-1, charge); | |
696 | ||
697 | Int_t digits[3]; | |
698 | Int_t it = (Int_t)((ik-1)/fScaleSize); | |
699 | digits[0] = index; | |
700 | digits[1] = it; | |
701 | digits[2] = (Int_t)timeAmplitude; | |
702 | Float_t phys; | |
703 | if (idtrack >= 0) phys = (Float_t)timeAmplitude; | |
704 | else phys = 0; | |
705 | ||
706 | Double_t cellcharge = 0.; | |
707 | AliITSTransientDigit* pdigit; | |
708 | // build the list of fired cells and update the info | |
709 | if (!fHitMap1->TestHit(index, it)) { | |
710 | new((*padr)[countadr++]) TVector(3); | |
711 | TVector &trinfo=*((TVector*) (*padr)[countadr-1]); | |
712 | trinfo(0) = (Float_t)idtrack; | |
713 | trinfo(1) = (Float_t)idhit; | |
714 | trinfo(2) = (Float_t)timeAmplitude; | |
715 | ||
716 | alist->AddAtAndExpand(new AliITSTransientDigit(phys,digits),counter); | |
717 | fHitMap1->SetHit(index, it, counter); | |
718 | counter++; | |
719 | pdigit=(AliITSTransientDigit*)alist->At(alist->GetLast()); | |
720 | // list of tracks | |
721 | TObjArray *trlist=(TObjArray*)pdigit->TrackList(); | |
722 | trlist->Add(&trinfo); | |
723 | } else { | |
724 | pdigit = (AliITSTransientDigit*) fHitMap1->GetHit(index, it); | |
725 | for(Int_t kk=0;kk<fScaleSize;kk++) { | |
726 | cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk); | |
727 | } // end for kk | |
728 | // update charge | |
729 | (*pdigit).fSignal = (Int_t)cellcharge; | |
730 | (*pdigit).fPhysics += phys; | |
731 | // update list of tracks | |
732 | TObjArray* trlist = (TObjArray*)pdigit->TrackList(); | |
733 | Int_t lastentry = trlist->GetLast(); | |
734 | TVector *ptrkp = (TVector*)trlist->At(lastentry); | |
735 | TVector &trinfo = *ptrkp; | |
736 | Int_t lasttrack = Int_t(trinfo(0)); | |
737 | Float_t lastcharge=(trinfo(2)); | |
738 | if (lasttrack==idtrack ) { | |
739 | lastcharge += (Float_t)timeAmplitude; | |
740 | trlist->RemoveAt(lastentry); | |
741 | trinfo(0) = lasttrack; | |
742 | trinfo(1) = idhit; | |
743 | trinfo(2) = lastcharge; | |
744 | trlist->AddAt(&trinfo,lastentry); | |
745 | } else { | |
746 | new((*padr)[countadr++]) TVector(3); | |
747 | TVector &trinfo=*((TVector*) (*padr)[countadr-1]); | |
748 | trinfo(0) = (Float_t)idtrack; | |
749 | trinfo(1) = (Float_t)idhit; | |
750 | trinfo(2) = (Float_t)timeAmplitude; | |
751 | trlist->Add(&trinfo); | |
752 | } // end if lasttrack==idtrack | |
b0f5e3fc | 753 | |
754 | #ifdef print | |
8a33ae9e | 755 | // check the track list - debugging |
756 | Int_t trk[20], htrk[20]; | |
757 | Float_t chtrk[20]; | |
758 | Int_t nptracks = trlist->GetEntriesFast(); | |
759 | if (nptracks > 2) { | |
760 | Int_t tr; | |
761 | for (tr=0;tr<nptracks;tr++) { | |
762 | TVector *pptrkp = (TVector*)trlist->At(tr); | |
763 | TVector &pptrk = *pptrkp; | |
764 | trk[tr] = Int_t(pptrk(0)); | |
765 | htrk[tr] = Int_t(pptrk(1)); | |
766 | chtrk[tr] = (pptrk(2)); | |
c7a4dac0 | 767 | cout << "nptracks "<<nptracks << endl; |
8a33ae9e | 768 | } // end for tr |
769 | } // end if nptracks | |
b0f5e3fc | 770 | #endif |
8a33ae9e | 771 | } // end if pdigit |
b0f5e3fc | 772 | |
8a33ae9e | 773 | // update counter and countadr for next call. |
774 | arg[4] = counter; | |
775 | arg[5] = countadr; | |
b0f5e3fc | 776 | } |
b0f5e3fc | 777 | //____________________________________________ |
778 | ||
779 | void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){ | |
8a33ae9e | 780 | // Adds a Digit. |
b0f5e3fc | 781 | // tag with -1 signals coming from background tracks |
782 | // tag with -2 signals coming from pure electronic noise | |
783 | ||
e8189707 | 784 | Int_t digits[3], tracks[3], hits[3]; |
b0f5e3fc | 785 | Float_t phys, charges[3]; |
786 | ||
e8189707 | 787 | Int_t trk[20], htrk[20]; |
b0f5e3fc | 788 | Float_t chtrk[20]; |
789 | ||
ece86d9a | 790 | Bool_t do10to8=fResponse->Do10to8(); |
791 | ||
792 | if(do10to8) signal=Convert8to10(signal); | |
b0f5e3fc | 793 | AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j); |
8a33ae9e | 794 | digits[0] = i; |
795 | digits[1] = j; | |
796 | digits[2] = signal; | |
b0f5e3fc | 797 | if (!obj) { |
798 | phys=0; | |
799 | Int_t k; | |
e8189707 | 800 | for (k=0;k<3;k++) { |
8a33ae9e | 801 | tracks[k]=-2; |
802 | charges[k]=0; | |
803 | hits[k]=-1; | |
804 | } // end for k | |
e8189707 | 805 | fITS->AddSimDigit(1,phys,digits,tracks,hits,charges); |
b0f5e3fc | 806 | } else { |
8a33ae9e | 807 | phys=obj->fPhysics; |
808 | TObjArray* trlist=(TObjArray*)obj->TrackList(); | |
809 | Int_t nptracks=trlist->GetEntriesFast(); | |
810 | if (nptracks > 20) { | |
c7a4dac0 | 811 | Warning("AddDigit","nptracks=%d > 20 nptracks set to 20",nptracks); |
8a33ae9e | 812 | nptracks=20; |
813 | } // end if nptracks > 20 | |
814 | Int_t tr; | |
815 | for (tr=0;tr<nptracks;tr++) { | |
816 | TVector &pp =*((TVector*)trlist->At(tr)); | |
817 | trk[tr]=Int_t(pp(0)); | |
818 | htrk[tr]=Int_t(pp(1)); | |
819 | chtrk[tr]=(pp(2)); | |
820 | } // end for tr | |
821 | if (nptracks > 1) { | |
822 | SortTracks(trk,chtrk,htrk,nptracks); | |
823 | } // end if nptracks > 1 | |
824 | Int_t i; | |
825 | if (nptracks < 3 ) { | |
826 | for (i=0; i<nptracks; i++) { | |
827 | tracks[i]=trk[i]; | |
828 | charges[i]=chtrk[i]; | |
829 | hits[i]=htrk[i]; | |
830 | } // end for i | |
831 | for (i=nptracks; i<3; i++) { | |
832 | tracks[i]=-3; | |
833 | hits[i]=-1; | |
834 | charges[i]=0; | |
835 | } // end for i | |
836 | } else { | |
837 | for (i=0; i<3; i++) { | |
838 | tracks[i]=trk[i]; | |
839 | charges[i]=chtrk[i]; | |
840 | hits[i]=htrk[i]; | |
841 | } // end for i | |
842 | } // end if/else nptracks < 3 | |
843 | ||
844 | fITS->AddSimDigit(1,phys,digits,tracks,hits,charges); | |
b0f5e3fc | 845 | |
8a33ae9e | 846 | } // end if/else !obj |
b0f5e3fc | 847 | } |
8a33ae9e | 848 | //______________________________________________________________________ |
849 | void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges, | |
850 | Int_t *hits,Int_t ntr){ | |
851 | // Sort the list of tracks contributing to a given digit | |
852 | // Only the 3 most significant tracks are acctually sorted | |
853 | // Loop over signals, only 3 times | |
854 | ||
855 | Float_t qmax; | |
856 | Int_t jmax; | |
857 | Int_t idx[3] = {-3,-3,-3}; | |
858 | Float_t jch[3] = {-3,-3,-3}; | |
859 | Int_t jtr[3] = {-3,-3,-3}; | |
860 | Int_t jhit[3] = {-3,-3,-3}; | |
861 | Int_t i,j,imax; | |
862 | ||
863 | if (ntr<3) imax = ntr; | |
864 | else imax = 3; | |
865 | for(i=0;i<imax;i++){ | |
866 | qmax = 0; | |
867 | jmax = 0; | |
868 | for(j=0;j<ntr;j++){ | |
869 | if((i == 1 && j == idx[i-1] ) | |
870 | ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue; | |
871 | if(charges[j] > qmax) { | |
872 | qmax = charges[j]; | |
873 | jmax=j; | |
874 | } // end if charges[j]>qmax | |
875 | } // end for j | |
876 | if(qmax > 0) { | |
877 | idx[i] = jmax; | |
878 | jch[i] = charges[jmax]; | |
879 | jtr[i] = tracks[jmax]; | |
880 | jhit[i] = hits[jmax]; | |
881 | } // end if qmax > 0 | |
882 | } // end for i | |
883 | ||
884 | for(i=0;i<3;i++){ | |
885 | if (jtr[i] == -3) { | |
886 | charges[i] = 0; | |
887 | tracks[i] = -3; | |
888 | hits[i] = -1; | |
889 | } else { | |
890 | charges[i] = jch[i]; | |
891 | tracks[i] = jtr[i]; | |
892 | hits[i] = jhit[i]; | |
893 | } // end if jtr[i] == -3 | |
894 | } // end for i | |
b0f5e3fc | 895 | } |
8a33ae9e | 896 | //______________________________________________________________________ |
b0f5e3fc | 897 | void AliITSsimulationSDD::ChargeToSignal() { |
8a33ae9e | 898 | // add baseline, noise, electronics and ADC saturation effects |
b0f5e3fc | 899 | |
8a33ae9e | 900 | char opt1[20], opt2[20]; |
901 | fResponse->ParamOptions(opt1,opt2); | |
902 | char *read = strstr(opt1,"file"); | |
903 | Float_t baseline, noise; | |
904 | ||
905 | if (read) { | |
906 | static Bool_t readfile=kTRUE; | |
907 | //read baseline and noise from file | |
908 | if (readfile) ReadBaseline(); | |
909 | readfile=kFALSE; | |
910 | } else fResponse->GetNoiseParam(noise,baseline); | |
911 | ||
912 | Float_t contrib=0; | |
913 | Int_t i,k,kk; | |
914 | Float_t maxadc = fResponse->MaxAdc(); | |
915 | if(!fDoFFT) { | |
916 | for (i=0;i<fNofMaps;i++) { | |
917 | if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise); | |
918 | for(k=0; k<fScaleSize*fMaxNofSamples; k++) { | |
919 | fInZR[k] = fHitMap2->GetSignal(i,k); | |
920 | contrib = (baseline + noise*gRandom->Gaus()); | |
921 | fInZR[k] += contrib; | |
922 | } // end for k | |
923 | for(k=0; k<fMaxNofSamples; k++) { | |
924 | Double_t newcont = 0.; | |
925 | Double_t maxcont = 0.; | |
926 | for(kk=0;kk<fScaleSize;kk++) { | |
927 | newcont = fInZR[fScaleSize*k+kk]; | |
928 | if(newcont > maxcont) maxcont = newcont; | |
929 | } // end for kk | |
930 | newcont = maxcont; | |
931 | if (newcont >= maxadc) newcont = maxadc -1; | |
c7a4dac0 | 932 | if(newcont >= baseline){ |
933 | Warning("","newcont=%d>=baseline=%d",newcont,baseline); | |
934 | } // end if | |
8a33ae9e | 935 | // back to analog: ? |
936 | fHitMap2->SetHit(i,k,newcont); | |
937 | } // end for k | |
938 | } // end for i loop over anodes | |
939 | return; | |
940 | } // end if DoFFT | |
ece86d9a | 941 | |
ece86d9a | 942 | for (i=0;i<fNofMaps;i++) { |
8a33ae9e | 943 | if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise); |
944 | for(k=0; k<fScaleSize*fMaxNofSamples; k++) { | |
945 | fInZR[k] = fHitMap2->GetSignal(i,k); | |
946 | contrib = (baseline + noise*gRandom->Gaus()); | |
947 | fInZR[k] += contrib; | |
948 | fInZI[k] = 0.; | |
949 | } // end for k | |
950 | FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1); | |
ece86d9a | 951 | for(k=0; k<fScaleSize*fMaxNofSamples; k++) { |
8a33ae9e | 952 | Double_t rw = fElectronics->GetTraFunReal(k); |
953 | Double_t iw = fElectronics->GetTraFunImag(k); | |
954 | fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw; | |
955 | fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw; | |
956 | } // end for k | |
957 | FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1); | |
ece86d9a | 958 | for(k=0; k<fMaxNofSamples; k++) { |
8a33ae9e | 959 | Double_t newcont1 = 0.; |
960 | Double_t maxcont1 = 0.; | |
961 | for(kk=0;kk<fScaleSize;kk++) { | |
962 | newcont1 = fOutZR[fScaleSize*k+kk]; | |
963 | if(newcont1 > maxcont1) maxcont1 = newcont1; | |
964 | } // end for kk | |
965 | newcont1 = maxcont1; | |
966 | if (newcont1 >= maxadc) newcont1 = maxadc -1; | |
967 | fHitMap2->SetHit(i,k,newcont1); | |
968 | } // end for k | |
969 | } // end for i loop over anodes | |
ece86d9a | 970 | return; |
b0f5e3fc | 971 | } |
8a33ae9e | 972 | //______________________________________________________________________ |
b0f5e3fc | 973 | void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline, |
974 | Float_t &noise){ | |
8a33ae9e | 975 | // Returns the Baseline for a particular anode. |
976 | baseline = fBaseline[i]; | |
977 | noise = fNoise[i]; | |
b0f5e3fc | 978 | } |
8a33ae9e | 979 | //______________________________________________________________________ |
b0f5e3fc | 980 | void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl, |
981 | Int_t &th){ | |
8a33ae9e | 982 | // Returns the compression alogirthm parameters |
983 | Int_t size = fD.GetSize(); | |
984 | if (size > 2 ) { | |
985 | db=fD[i]; tl=fT1[i]; th=fT2[i]; | |
986 | } else { | |
987 | if (size <= 2 && i>=fNofMaps/2) { | |
988 | db=fD[1]; tl=fT1[1]; th=fT2[1]; | |
989 | } else { | |
990 | db=fD[0]; tl=fT1[0]; th=fT2[0]; | |
991 | } // end if size <=2 && i>=fNofMaps/2 | |
992 | } // end if size >2 | |
b0f5e3fc | 993 | } |
8a33ae9e | 994 | //______________________________________________________________________ |
b0f5e3fc | 995 | void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){ |
8a33ae9e | 996 | // returns the compression alogirthm parameters |
997 | Int_t size = fD.GetSize(); | |
b0f5e3fc | 998 | |
8a33ae9e | 999 | if (size > 2 ) { |
1000 | db=fD[i]; tl=fT1[i]; | |
1001 | } else { | |
1002 | if (size <= 2 && i>=fNofMaps/2) { | |
1003 | db=fD[1]; tl=fT1[1]; | |
1004 | } else { | |
1005 | db=fD[0]; tl=fT1[0]; | |
1006 | } // end if size <=2 && i>=fNofMaps/2 | |
1007 | } // end if size > 2 | |
b0f5e3fc | 1008 | } |
8a33ae9e | 1009 | //______________________________________________________________________ |
b0f5e3fc | 1010 | void AliITSsimulationSDD::SetCompressParam(){ |
8a33ae9e | 1011 | // Sets the compression alogirthm parameters |
1012 | Int_t cp[8],i; | |
1013 | ||
1014 | fResponse->GiveCompressParam(cp); | |
1015 | for (i=0; i<2; i++) { | |
1016 | fD[i] = cp[i]; | |
1017 | fT1[i] = cp[i+2]; | |
1018 | fT2[i] = cp[i+4]; | |
1019 | fTol[i] = cp[i+6]; | |
8a33ae9e | 1020 | } // end for i |
b0f5e3fc | 1021 | } |
8a33ae9e | 1022 | //______________________________________________________________________ |
b0f5e3fc | 1023 | void AliITSsimulationSDD::ReadBaseline(){ |
8a33ae9e | 1024 | // read baseline and noise from file - either a .root file and in this |
1025 | // case data should be organised in a tree with one entry for each | |
1026 | // module => reading should be done accordingly | |
1027 | // or a classic file and do smth. like this: | |
1028 | // Read baselines and noise for SDD | |
b0f5e3fc | 1029 | |
1030 | Int_t na,pos; | |
1031 | Float_t bl,n; | |
e8189707 | 1032 | char input[100], base[100], param[100]; |
b0f5e3fc | 1033 | char *filtmp; |
1034 | ||
e8189707 | 1035 | fResponse->Filenames(input,base,param); |
b0f5e3fc | 1036 | fFileName=base; |
1037 | // | |
1038 | filtmp = gSystem->ExpandPathName(fFileName.Data()); | |
1039 | FILE *bline = fopen(filtmp,"r"); | |
b0f5e3fc | 1040 | na = 0; |
1041 | ||
1042 | if(bline) { | |
8a33ae9e | 1043 | while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) { |
1044 | if (pos != na+1) { | |
1045 | Error("ReadBaseline","Anode number not in increasing order!", | |
1046 | filtmp); | |
1047 | exit(1); | |
1048 | } // end if pos != na+1 | |
1049 | fBaseline[na]=bl; | |
1050 | fNoise[na]=n; | |
1051 | na++; | |
1052 | } // end while | |
b0f5e3fc | 1053 | } else { |
8a33ae9e | 1054 | Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",filtmp); |
1055 | exit(1); | |
b0f5e3fc | 1056 | } // end if(bline) |
8a33ae9e | 1057 | |
b0f5e3fc | 1058 | fclose(bline); |
1059 | delete [] filtmp; | |
b0f5e3fc | 1060 | } |
8a33ae9e | 1061 | //______________________________________________________________________ |
1062 | Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const { | |
1063 | // To the 10 to 8 bit lossive compression. | |
1064 | // code from Davide C. and Albert W. | |
1065 | ||
1066 | if (signal < 128) return signal; | |
1067 | if (signal < 256) return (128+((signal-128)>>1)); | |
1068 | if (signal < 512) return (192+((signal-256)>>3)); | |
1069 | if (signal < 1024) return (224+((signal-512)>>4)); | |
1070 | return 0; | |
b0f5e3fc | 1071 | } |
8a33ae9e | 1072 | //______________________________________________________________________ |
1073 | Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const { | |
1074 | // Undo the lossive 10 to 8 bit compression. | |
1075 | // code from Davide C. and Albert W. | |
1076 | if (signal < 0 || signal > 255) { | |
c7a4dac0 | 1077 | Warning("Convert8to10","out of range signal=%d",signal); |
8a33ae9e | 1078 | return 0; |
1079 | } // end if signal <0 || signal >255 | |
1080 | ||
1081 | if (signal < 128) return signal; | |
1082 | if (signal < 192) { | |
1083 | if (TMath::Odd(signal)) return (128+((signal-128)<<1)); | |
1084 | else return (128+((signal-128)<<1)+1); | |
1085 | } // end if signal < 192 | |
1086 | if (signal < 224) { | |
1087 | if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3); | |
1088 | else return (256+((signal-192)<<3)+4); | |
1089 | } // end if signal < 224 | |
1090 | if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7); | |
9226cc6a | 1091 | return (512+((signal-224)<<4)+7); |
8a33ae9e | 1092 | } |
1093 | //______________________________________________________________________ | |
b0f5e3fc | 1094 | AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){ |
8a33ae9e | 1095 | //Return the correct map. |
1096 | ||
b0f5e3fc | 1097 | return ((i==0)? fHitMap1 : fHitMap2); |
1098 | } | |
8a33ae9e | 1099 | //______________________________________________________________________ |
e8189707 | 1100 | void AliITSsimulationSDD::ZeroSuppression(const char *option) { |
8a33ae9e | 1101 | // perform the zero suppresion |
1102 | ||
1103 | if (strstr(option,"2D")) { | |
1104 | //Init2D(); // activate if param change module by module | |
1105 | Compress2D(); | |
1106 | } else if (strstr(option,"1D")) { | |
1107 | //Init1D(); // activate if param change module by module | |
1108 | Compress1D(); | |
1109 | } else StoreAllDigits(); | |
b0f5e3fc | 1110 | } |
8a33ae9e | 1111 | //______________________________________________________________________ |
b0f5e3fc | 1112 | void AliITSsimulationSDD::Init2D(){ |
8a33ae9e | 1113 | // read in and prepare arrays: fD, fT1, fT2 |
1114 | // savemu[nanodes], savesigma[nanodes] | |
1115 | // read baseline and noise from file - either a .root file and in this | |
1116 | // case data should be organised in a tree with one entry for each | |
1117 | // module => reading should be done accordingly | |
1118 | // or a classic file and do smth. like this ( code from Davide C. and | |
1119 | // Albert W.) : | |
1120 | // Read 2D zero-suppression parameters for SDD | |
b0f5e3fc | 1121 | |
1122 | if (!strstr(fParam,"file")) return; | |
1123 | ||
1124 | Int_t na,pos,tempTh; | |
1125 | Float_t mu,sigma; | |
8a33ae9e | 1126 | Float_t *savemu = new Float_t [fNofMaps]; |
e8189707 | 1127 | Float_t *savesigma = new Float_t [fNofMaps]; |
1128 | char input[100],basel[100],par[100]; | |
b0f5e3fc | 1129 | char *filtmp; |
b0f5e3fc | 1130 | Int_t minval = fResponse->MinVal(); |
1131 | ||
e8189707 | 1132 | fResponse->Filenames(input,basel,par); |
8a33ae9e | 1133 | fFileName = par; |
b0f5e3fc | 1134 | // |
1135 | filtmp = gSystem->ExpandPathName(fFileName.Data()); | |
1136 | FILE *param = fopen(filtmp,"r"); | |
1137 | na = 0; | |
1138 | ||
1139 | if(param) { | |
8a33ae9e | 1140 | while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) { |
1141 | if (pos != na+1) { | |
c7a4dac0 | 1142 | Error("Init2D","Anode number not in increasing order!",filtmp); |
8a33ae9e | 1143 | exit(1); |
1144 | } // end if pos != na+1 | |
1145 | savemu[na] = mu; | |
1146 | savesigma[na] = sigma; | |
b0f5e3fc | 1147 | if ((2.*sigma) < mu) { |
1148 | fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5); | |
1149 | mu = 2.0 * sigma; | |
1150 | } else fD[na] = 0; | |
1151 | tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval; | |
1152 | if (tempTh < 0) tempTh=0; | |
1153 | fT1[na] = tempTh; | |
1154 | tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval; | |
1155 | if (tempTh < 0) tempTh=0; | |
1156 | fT2[na] = tempTh; | |
1157 | na++; | |
8a33ae9e | 1158 | } // end while |
b0f5e3fc | 1159 | } else { |
c7a4dac0 | 1160 | Error("Init2D","THE FILE %s DOES NOT EXIST !",filtmp); |
8a33ae9e | 1161 | exit(1); |
b0f5e3fc | 1162 | } // end if(param) |
8a33ae9e | 1163 | |
b0f5e3fc | 1164 | fclose(param); |
1165 | delete [] filtmp; | |
5d18fa90 | 1166 | delete [] savemu; |
e8189707 | 1167 | delete [] savesigma; |
8a33ae9e | 1168 | } |
1169 | //______________________________________________________________________ | |
b0f5e3fc | 1170 | void AliITSsimulationSDD::Compress2D(){ |
8a33ae9e | 1171 | // simple ITS cluster finder -- online zero-suppression conditions |
b0f5e3fc | 1172 | |
b0f5e3fc | 1173 | Int_t db,tl,th; |
8a33ae9e | 1174 | Int_t minval = fResponse->MinVal(); |
1175 | Bool_t write = fResponse->OutputOption(); | |
1176 | Bool_t do10to8 = fResponse->Do10to8(); | |
b0f5e3fc | 1177 | Int_t nz, nl, nh, low, i, j; |
1178 | ||
e8189707 | 1179 | for (i=0; i<fNofMaps; i++) { |
b0f5e3fc | 1180 | CompressionParam(i,db,tl,th); |
8a33ae9e | 1181 | nz = 0; |
1182 | nl = 0; | |
1183 | nh = 0; | |
1184 | low = 0; | |
e8189707 | 1185 | for (j=0; j<fMaxNofSamples; j++) { |
b0f5e3fc | 1186 | Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j)); |
1187 | signal -= db; // if baseline eq. is done here | |
1188 | if (signal <= 0) {nz++; continue;} | |
1189 | if ((signal - tl) < minval) low++; | |
1190 | if ((signal - th) >= minval) { | |
1191 | nh++; | |
1192 | Bool_t cond=kTRUE; | |
b0f5e3fc | 1193 | FindCluster(i,j,signal,minval,cond); |
c7a4dac0 | 1194 | if(cond && j && |
1195 | ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)){ | |
8a33ae9e | 1196 | if(do10to8) signal = Convert10to8(signal); |
1197 | AddDigit(i,j,signal); | |
1198 | } // end if cond&&j&&() | |
b0f5e3fc | 1199 | } else if ((signal - tl) >= minval) nl++; |
8a33ae9e | 1200 | } // end for j loop time samples |
1201 | if (write) TreeB()->Fill(nz,nl,nh,low,i+1); | |
1202 | } //end for i loop anodes | |
b0f5e3fc | 1203 | |
8a33ae9e | 1204 | char hname[30]; |
1205 | if (write) { | |
b0f5e3fc | 1206 | sprintf(hname,"TNtuple%d_%d",fModule,fEvent); |
1207 | TreeB()->Write(hname); | |
1208 | // reset tree | |
1209 | TreeB()->Reset(); | |
8a33ae9e | 1210 | } // end if write |
1211 | } | |
1212 | //______________________________________________________________________ | |
b0f5e3fc | 1213 | void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal, |
ece86d9a | 1214 | Int_t minval,Bool_t &cond){ |
8a33ae9e | 1215 | // Find clusters according to the online 2D zero-suppression algorithm |
1216 | Bool_t do10to8 = fResponse->Do10to8(); | |
1217 | Bool_t high = kFALSE; | |
b0f5e3fc | 1218 | |
1219 | fHitMap2->FlagHit(i,j); | |
1220 | // | |
1221 | // check the online zero-suppression conditions | |
1222 | // | |
8a33ae9e | 1223 | const Int_t kMaxNeighbours = 4; |
b0f5e3fc | 1224 | Int_t nn; |
1225 | Int_t dbx,tlx,thx; | |
8a33ae9e | 1226 | Int_t xList[kMaxNeighbours], yList[kMaxNeighbours]; |
e8189707 | 1227 | fSegmentation->Neighbours(i,j,&nn,xList,yList); |
ece86d9a | 1228 | Int_t in,ix,iy,qns; |
e8189707 | 1229 | for (in=0; in<nn; in++) { |
1230 | ix=xList[in]; | |
1231 | iy=yList[in]; | |
b0f5e3fc | 1232 | if (fHitMap2->TestHit(ix,iy)==kUnused) { |
8a33ae9e | 1233 | CompressionParam(ix,dbx,tlx,thx); |
1234 | Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy)); | |
1235 | qn -= dbx; // if baseline eq. is done here | |
1236 | if ((qn-tlx) < minval) { | |
1237 | fHitMap2->FlagHit(ix,iy); | |
1238 | continue; | |
1239 | } else { | |
1240 | if ((qn - thx) >= minval) high=kTRUE; | |
1241 | if (cond) { | |
1242 | if(do10to8) signal = Convert10to8(signal); | |
1243 | AddDigit(i,j,signal); | |
1244 | } // end if cond | |
1245 | if(do10to8) qns = Convert10to8(qn); | |
1246 | else qns=qn; | |
1247 | if (!high) AddDigit(ix,iy,qns); | |
1248 | cond=kFALSE; | |
1249 | if(!high) fHitMap2->FlagHit(ix,iy); | |
1250 | } // end if qn-tlx < minval | |
1251 | } // end if TestHit | |
1252 | } // end for in loop over neighbours | |
b0f5e3fc | 1253 | } |
8a33ae9e | 1254 | //______________________________________________________________________ |
b0f5e3fc | 1255 | void AliITSsimulationSDD::Init1D(){ |
8a33ae9e | 1256 | // this is just a copy-paste of input taken from 2D algo |
1257 | // Torino people should give input | |
1258 | // Read 1D zero-suppression parameters for SDD | |
b0f5e3fc | 1259 | |
1260 | if (!strstr(fParam,"file")) return; | |
1261 | ||
1262 | Int_t na,pos,tempTh; | |
1263 | Float_t mu,sigma; | |
8a33ae9e | 1264 | Float_t *savemu = new Float_t [fNofMaps]; |
e8189707 | 1265 | Float_t *savesigma = new Float_t [fNofMaps]; |
1266 | char input[100],basel[100],par[100]; | |
b0f5e3fc | 1267 | char *filtmp; |
b0f5e3fc | 1268 | Int_t minval = fResponse->MinVal(); |
8a33ae9e | 1269 | |
e8189707 | 1270 | fResponse->Filenames(input,basel,par); |
1271 | fFileName=par; | |
b0f5e3fc | 1272 | |
1273 | // set first the disable and tol param | |
1274 | SetCompressParam(); | |
1275 | // | |
1276 | filtmp = gSystem->ExpandPathName(fFileName.Data()); | |
1277 | FILE *param = fopen(filtmp,"r"); | |
1278 | na = 0; | |
1279 | ||
1280 | if (param) { | |
8a33ae9e | 1281 | fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]); |
1282 | while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) { | |
1283 | if (pos != na+1) { | |
c7a4dac0 | 1284 | Error("Init1D","Anode number not in increasing order!",filtmp); |
8a33ae9e | 1285 | exit(1); |
1286 | } // end if pos != na+1 | |
1287 | savemu[na]=mu; | |
1288 | savesigma[na]=sigma; | |
1289 | if ((2.*sigma) < mu) { | |
1290 | fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5); | |
1291 | mu = 2.0 * sigma; | |
1292 | } else fD[na] = 0; | |
1293 | tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval; | |
1294 | if (tempTh < 0) tempTh=0; | |
1295 | fT1[na] = tempTh; | |
1296 | na++; | |
1297 | } // end while | |
b0f5e3fc | 1298 | } else { |
c7a4dac0 | 1299 | Error("Init1D","THE FILE %s DOES NOT EXIST !",filtmp); |
8a33ae9e | 1300 | exit(1); |
b0f5e3fc | 1301 | } // end if(param) |
8a33ae9e | 1302 | |
b0f5e3fc | 1303 | fclose(param); |
1304 | delete [] filtmp; | |
749bd21a | 1305 | delete [] savemu; |
e8189707 | 1306 | delete [] savesigma; |
8a33ae9e | 1307 | } |
1308 | //______________________________________________________________________ | |
b0f5e3fc | 1309 | void AliITSsimulationSDD::Compress1D(){ |
1310 | // 1D zero-suppression algorithm (from Gianluca A.) | |
8a33ae9e | 1311 | Int_t dis,tol,thres,decr,diff; |
b0f5e3fc | 1312 | UChar_t *str=fStream->Stream(); |
8a33ae9e | 1313 | Int_t counter=0; |
1314 | Bool_t do10to8=fResponse->Do10to8(); | |
1315 | Int_t last=0; | |
1316 | Int_t k,i,j; | |
ece86d9a | 1317 | |
ece86d9a | 1318 | for (k=0; k<2; k++) { |
8a33ae9e | 1319 | tol = Tolerance(k); |
1320 | dis = Disable(k); | |
1321 | for (i=0; i<fNofMaps/2; i++) { | |
1322 | Bool_t firstSignal=kTRUE; | |
1323 | Int_t idx=i+k*fNofMaps/2; | |
1324 | CompressionParam(idx,decr,thres); | |
1325 | for (j=0; j<fMaxNofSamples; j++) { | |
1326 | Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j)); | |
1327 | signal -= decr; // if baseline eq. | |
1328 | if(do10to8) signal = Convert10to8(signal); | |
1329 | if (signal <= thres) { | |
1330 | signal=0; | |
1331 | diff=128; | |
1332 | last=0; | |
1333 | // write diff in the buffer for HuffT | |
1334 | str[counter]=(UChar_t)diff; | |
1335 | counter++; | |
1336 | continue; | |
1337 | } // end if signal <= thres | |
1338 | diff=signal-last; | |
1339 | if (diff > 127) diff=127; | |
1340 | if (diff < -128) diff=-128; | |
1341 | if (signal < dis) { | |
1342 | // tol has changed to 8 possible cases ? - one can write | |
1343 | // this if(TMath::Abs(diff)<tol) ... else ... | |
ece86d9a | 1344 | if(TMath::Abs(diff)<tol) diff=0; |
1345 | // or keep it as it was before | |
1346 | /* | |
b0f5e3fc | 1347 | if (tol==1 && (diff >= -2 && diff <= 1)) diff=0; |
1348 | if (tol==2 && (diff >= -4 && diff <= 3)) diff=0; | |
1349 | if (tol==3 && (diff >= -16 && diff <= 15)) diff=0; | |
ece86d9a | 1350 | */ |
1351 | AddDigit(idx,j,last+diff); | |
8a33ae9e | 1352 | } else { |
1353 | AddDigit(idx,j,signal); | |
1354 | } // end if singal < dis | |
1355 | diff += 128; | |
1356 | // write diff in the buffer used to compute Huffman tables | |
1357 | if (firstSignal) str[counter]=(UChar_t)signal; | |
1358 | else str[counter]=(UChar_t)diff; | |
1359 | counter++; | |
1360 | last=signal; | |
1361 | firstSignal=kFALSE; | |
1362 | } // end for j loop time samples | |
1363 | } // end for i loop anodes one half of detector | |
1364 | } // end for k | |
b0f5e3fc | 1365 | |
1366 | // check | |
1367 | fStream->CheckCount(counter); | |
1368 | ||
1369 | // open file and write out the stream of diff's | |
b0f5e3fc | 1370 | static Bool_t open=kTRUE; |
e8189707 | 1371 | static TFile *outFile; |
b0f5e3fc | 1372 | Bool_t write = fResponse->OutputOption(); |
1373 | ||
1374 | if (write ) { | |
1375 | if(open) { | |
1376 | SetFileName("stream.root"); | |
1377 | cout<<"filename "<<fFileName<<endl; | |
e8189707 | 1378 | outFile=new TFile(fFileName,"recreate"); |
b0f5e3fc | 1379 | cout<<"I have opened "<<fFileName<<" file "<<endl; |
8a33ae9e | 1380 | } // end if open |
1381 | open = kFALSE; | |
e8189707 | 1382 | outFile->cd(); |
b0f5e3fc | 1383 | fStream->Write(); |
1384 | } // endif write | |
1385 | ||
8a33ae9e | 1386 | fStream->ClearStream(); |
b0f5e3fc | 1387 | |
8a33ae9e | 1388 | // back to galice.root file |
b0f5e3fc | 1389 | |
8a33ae9e | 1390 | TTree *fAli=gAlice->TreeK(); |
1391 | TFile *file = 0; | |
b0f5e3fc | 1392 | |
8a33ae9e | 1393 | if (fAli) file =fAli->GetCurrentFile(); |
1394 | file->cd(); | |
1395 | } | |
1396 | //______________________________________________________________________ | |
b0f5e3fc | 1397 | void AliITSsimulationSDD::StoreAllDigits(){ |
8a33ae9e | 1398 | // if non-zero-suppressed data |
1399 | Bool_t do10to8 = fResponse->Do10to8(); | |
ece86d9a | 1400 | Int_t i, j, digits[3]; |
8a33ae9e | 1401 | |
e8189707 | 1402 | for (i=0; i<fNofMaps; i++) { |
1403 | for (j=0; j<fMaxNofSamples; j++) { | |
8a33ae9e | 1404 | Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j)); |
1405 | if(do10to8) signal = Convert10to8(signal); | |
1406 | if(do10to8) signal = Convert8to10(signal); | |
1407 | digits[0] = i; | |
1408 | digits[1] = j; | |
1409 | digits[2] = signal; | |
1410 | fITS->AddRealDigit(1,digits); | |
1411 | } // end for j | |
1412 | } // end for i | |
b0f5e3fc | 1413 | } |
8a33ae9e | 1414 | //______________________________________________________________________ |
ece86d9a | 1415 | void AliITSsimulationSDD::CreateHistograms(Int_t scale){ |
8a33ae9e | 1416 | // Creates histograms of maps for debugging |
1417 | Int_t i; | |
ece86d9a | 1418 | |
1419 | fHis=new TObjArray(fNofMaps); | |
e8189707 | 1420 | for (i=0;i<fNofMaps;i++) { |
a14dba92 | 1421 | TString sddName("sdd_"); |
b0f5e3fc | 1422 | Char_t candNum[4]; |
1423 | sprintf(candNum,"%d",i+1); | |
ece86d9a | 1424 | sddName.Append(candNum); |
8a33ae9e | 1425 | fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples, |
1426 | 0.,(Float_t) scale*fMaxNofSamples), i); | |
1427 | } // end for i | |
b0f5e3fc | 1428 | } |
8a33ae9e | 1429 | //______________________________________________________________________ |
ece86d9a | 1430 | void AliITSsimulationSDD::FillHistograms(){ |
8a33ae9e | 1431 | // fill 1D histograms from map |
1432 | ||
1433 | if (!fHis) return; | |
1434 | ||
1435 | for( Int_t i=0; i<fNofMaps; i++) { | |
1436 | TH1F *hist =(TH1F *)fHis->UncheckedAt(i); | |
1437 | Int_t nsamples = hist->GetNbinsX(); | |
1438 | for( Int_t j=0; j<nsamples; j++) { | |
1439 | Double_t signal=fHitMap2->GetSignal(i,j); | |
1440 | hist->Fill((Float_t)j,signal); | |
1441 | } // end for j | |
1442 | } // end for i | |
ece86d9a | 1443 | } |
8a33ae9e | 1444 | //______________________________________________________________________ |
b0f5e3fc | 1445 | void AliITSsimulationSDD::ResetHistograms(){ |
b0f5e3fc | 1446 | // Reset histograms for this detector |
b0f5e3fc | 1447 | Int_t i; |
8a33ae9e | 1448 | |
e8189707 | 1449 | for (i=0;i<fNofMaps;i++ ) { |
2682e810 | 1450 | if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset(); |
8a33ae9e | 1451 | } // end for i |
b0f5e3fc | 1452 | } |
8a33ae9e | 1453 | //______________________________________________________________________ |
b0f5e3fc | 1454 | TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) { |
8a33ae9e | 1455 | // Fills a histogram from a give anode. |
1456 | ||
1457 | if (!fHis) return 0; | |
1458 | ||
1459 | if(wing <=0 || wing > 2) { | |
c7a4dac0 | 1460 | Warning("GetAnode","Wrong wing number: %d",wing); |
8a33ae9e | 1461 | return NULL; |
1462 | } // end if wing <=0 || wing >2 | |
1463 | if(anode <=0 || anode > fNofMaps/2) { | |
c7a4dac0 | 1464 | Warning("GetAnode","Wrong anode number: %d",anode); |
8a33ae9e | 1465 | return NULL; |
1466 | } // end if ampde <=0 || andoe > fNofMaps/2 | |
1467 | ||
1468 | Int_t index = (wing-1)*fNofMaps/2 + anode-1; | |
8a33ae9e | 1469 | return (TH1F*)(fHis->At(index)); |
b0f5e3fc | 1470 | } |
8a33ae9e | 1471 | //______________________________________________________________________ |
b0f5e3fc | 1472 | void AliITSsimulationSDD::WriteToFile(TFile *hfile) { |
8a33ae9e | 1473 | // Writes the histograms to a file |
b0f5e3fc | 1474 | |
8a33ae9e | 1475 | if (!fHis) return; |
1476 | ||
1477 | hfile->cd(); | |
1478 | Int_t i; | |
8a33ae9e | 1479 | for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write(); |
1480 | return; | |
b0f5e3fc | 1481 | } |
8a33ae9e | 1482 | //______________________________________________________________________ |
ece86d9a | 1483 | Float_t AliITSsimulationSDD::GetNoise() { |
8a33ae9e | 1484 | // Returns the noise value |
1485 | //Bool_t do10to8=fResponse->Do10to8(); | |
1486 | //noise will always be in the liniar part of the signal | |
1487 | Int_t decr; | |
1488 | Int_t threshold = fT1[0]; | |
1489 | char opt1[20], opt2[20]; | |
1490 | ||
1491 | fResponse->ParamOptions(opt1,opt2); | |
1492 | fParam=opt2; | |
1493 | char *same = strstr(opt1,"same"); | |
1494 | Float_t noise,baseline; | |
1495 | if (same) { | |
1496 | fResponse->GetNoiseParam(noise,baseline); | |
1497 | } else { | |
1498 | static Bool_t readfile=kTRUE; | |
1499 | //read baseline and noise from file | |
1500 | if (readfile) ReadBaseline(); | |
1501 | readfile=kFALSE; | |
1502 | } // end if same | |
1503 | ||
1504 | TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2"); | |
1505 | if(c2) delete c2->GetPrimitive("noisehist"); | |
1506 | if(c2) delete c2->GetPrimitive("anode"); | |
1507 | else c2=new TCanvas("c2"); | |
1508 | c2->cd(); | |
1509 | c2->SetFillColor(0); | |
1510 | ||
1511 | TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold); | |
1512 | TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0., | |
1513 | (float)fMaxNofSamples); | |
1514 | Int_t i,k; | |
1515 | for (i=0;i<fNofMaps;i++) { | |
1516 | CompressionParam(i,decr,threshold); | |
1517 | if (!same) GetAnodeBaseline(i,baseline,noise); | |
1518 | anode->Reset(); | |
1519 | for (k=0;k<fMaxNofSamples;k++) { | |
1520 | Float_t signal=(Float_t)fHitMap2->GetSignal(i,k); | |
1521 | //if (signal <= (float)threshold) noisehist->Fill(signal-baseline); | |
1522 | if (signal <= (float)threshold) noisehist->Fill(signal); | |
1523 | anode->Fill((float)k,signal); | |
1524 | } // end for k | |
1525 | anode->Draw(); | |
1526 | c2->Update(); | |
1527 | } // end for i | |
1528 | TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold); | |
1529 | noisehist->Fit("gnoise","RQ"); | |
1530 | noisehist->Draw(); | |
ece86d9a | 1531 | c2->Update(); |
8a33ae9e | 1532 | Float_t mnoise = gnoise->GetParameter(1); |
1533 | cout << "mnoise : " << mnoise << endl; | |
1534 | Float_t rnoise = gnoise->GetParameter(2); | |
1535 | cout << "rnoise : " << rnoise << endl; | |
1536 | delete noisehist; | |
1537 | return rnoise; | |
c7a4dac0 | 1538 | }//______________________________________________________________________ |
1539 | void AliITSsimulationSDD::WriteSDigits(AliITSpList *pList){ | |
1540 | // Fills the Summable digits Tree | |
1541 | Int_t i,ni,j,nj; | |
1542 | static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS"); | |
1543 | ||
1544 | pList->GetMaxMapIndex(ni,nj); | |
1545 | for(i=0;i<ni;i++)for(j=0;j<nj;j++){ | |
1546 | if(pList->GetSignalOnly(i,j)>0.5*fT1[0]){ // above small threshold. | |
1547 | aliITS->AddSumDigit(*(pList->GetpListItem(i,j))); | |
1548 | // cout << "pListSDD: " << *(pList->GetpListItem(i,j)) << endl; | |
1549 | } // end if | |
1550 | } // end for i,j | |
1551 | return; | |
b0f5e3fc | 1552 | } |
8a33ae9e | 1553 | //______________________________________________________________________ |
44a312c3 | 1554 | void AliITSsimulationSDD::Print() { |
8a33ae9e | 1555 | // Print SDD simulation Parameters |
1556 | ||
1557 | cout << "**************************************************" << endl; | |
1558 | cout << " Silicon Drift Detector Simulation Parameters " << endl; | |
1559 | cout << "**************************************************" << endl; | |
1560 | cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl; | |
1561 | cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl; | |
1562 | cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl; | |
1563 | cout << "Number pf Anodes used: " << fNofMaps << endl; | |
1564 | cout << "Number of Time Samples: " << fMaxNofSamples << endl; | |
1565 | cout << "Scale size factor: " << fScaleSize << endl; | |
1566 | cout << "**************************************************" << endl; | |
44a312c3 | 1567 | } |