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