]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSsimulationSDD.cxx
Correct comparison wiht AliGeant3::Class() introduced. Thanks to I.Hrivnacova
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationSDD.cxx
CommitLineData
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 48ClassImp(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
69Int_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
78void 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
138AliITSsimulationSDD::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 168AliITSsimulationSDD::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 177AliITSsimulationSDD& 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 187AliITSsimulationSDD::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
288AliITSsimulationSDD::~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
316void 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
633void 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
746void 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 821void 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//____________________________________________
880void 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
ece86d9a 898 TRandom random;
e8189707 899 Int_t i,k,kk;
ece86d9a 900
44a312c3 901 Float_t maxadc = fResponse->MaxAdc();
ece86d9a 902 if(!fDoFFT) {
903 for (i=0;i<fNofMaps;i++) {
904 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
905 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
906 fInZR[k] = fHitMap2->GetSignal(i,k);
44a312c3 907 contrib = (baseline + noise*random.Gaus());
ece86d9a 908 fInZR[k] += contrib;
e8189707 909 }
ece86d9a 910 for(k=0; k<fMaxNofSamples; k++) {
44a312c3 911 Double_t newcont = 0.;
912 Double_t maxcont = 0.;
ece86d9a 913 for(kk=0;kk<fScaleSize;kk++) {
44a312c3 914
ece86d9a 915 newcont = fInZR[fScaleSize*k+kk];
916 if(newcont > maxcont) maxcont = newcont;
44a312c3 917
918 //newcont += (fInZR[fScaleSize*k+kk]/fScaleSize);
ece86d9a 919 }
920 newcont = maxcont;
44a312c3 921 if (newcont >= maxadc) newcont = maxadc -1;
922 if(newcont >= baseline) cout << "newcont: " << newcont << endl;
ece86d9a 923 // back to analog: ?
44a312c3 924 fHitMap2->SetHit(i,k,newcont);
ece86d9a 925 }
926 } // loop over anodes
927 return;
928 } // end if DoFFT
929
930 for (i=0;i<fNofMaps;i++) {
931 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
932 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
933 fInZR[k] = fHitMap2->GetSignal(i,k);
44a312c3 934 contrib = (baseline + noise*random.Gaus());
ece86d9a 935 fInZR[k] += contrib;
936 fInZI[k] = 0.;
937 }
e8189707 938 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
ece86d9a 939 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
940 Double_t rw = fElectronics->GetTraFunReal(k);
e8189707 941 Double_t iw = fElectronics->GetTraFunImag(k);
ece86d9a 942 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
943 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
e8189707 944 }
ece86d9a 945 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
946 for(k=0; k<fMaxNofSamples; k++) {
44a312c3 947 Double_t newcont1 = 0.;
948 Double_t maxcont1 = 0.;
ece86d9a 949 for(kk=0;kk<fScaleSize;kk++) {
44a312c3 950
951 newcont1 = fOutZR[fScaleSize*k+kk];
952 if(newcont1 > maxcont1) maxcont1 = newcont1;
953
954 //newcont1 += (fInZR[fScaleSize*k+kk]/fScaleSize);
ece86d9a 955 }
44a312c3 956 newcont1 = maxcont1;
957 //cout << "newcont1: " << newcont1 << endl;
958 if (newcont1 >= maxadc) newcont1 = maxadc -1;
959 fHitMap2->SetHit(i,k,newcont1);
ece86d9a 960 }
b0f5e3fc 961 } // loop over anodes
ece86d9a 962 return;
963
b0f5e3fc 964}
965
966//____________________________________________
967void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
968 Float_t &noise){
969 // Returns the Baseline for a particular anode.
970 baseline=fBaseline[i];
971 noise=fNoise[i];
972
973}
974
975//____________________________________________
976void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
977 Int_t &th){
978 // Returns the compression alogirthm parameters
979 Int_t size = fD.GetSize();
980 if (size > 2 ) {
981 db=fD[i]; tl=fT1[i]; th=fT2[i];
982 } else {
983 if (size <= 2 && i>=fNofMaps/2) {
984 db=fD[1]; tl=fT1[1]; th=fT2[1];
985 } else {
986 db=fD[0]; tl=fT1[0]; th=fT2[0];
987 }
988 }
989}
990//____________________________________________
991void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
992 // returns the compression alogirthm parameters
993 Int_t size = fD.GetSize();
994 if (size > 2 ) {
995 db=fD[i]; tl=fT1[i];
996 } else {
997 if (size <= 2 && i>=fNofMaps/2) {
998 db=fD[1]; tl=fT1[1];
999 } else {
1000 db=fD[0]; tl=fT1[0];
1001 }
1002 }
1003
1004}
1005//____________________________________________
1006void AliITSsimulationSDD::SetCompressParam(){
1007 // Sets the compression alogirthm parameters
1008 Int_t cp[8],i;
1009
1010 fResponse->GiveCompressParam(cp);
e8189707 1011 for (i=0; i<2; i++) {
b0f5e3fc 1012 fD[i] =cp[i];
1013 fT1[i] =cp[i+2];
1014 fT2[i] =cp[i+4];
1015 fTol[i]=cp[i+6];
b0f5e3fc 1016 printf("\n i, fD, fT1, fT2, fTol %d %d %d %d %d\n",
1017 i,fD[i],fT1[i],fT2[i],fTol[i]);
b0f5e3fc 1018 }
1019}
1020
1021//____________________________________________
1022void AliITSsimulationSDD::ReadBaseline(){
1023 // read baseline and noise from file - either a .root file and in this
1024 // case data should be organised in a tree with one entry for each
1025 // module => reading should be done accordingly
1026 // or a classic file and do smth. like this:
1027 //
1028 // Read baselines and noise for SDD
1029 //
1030
1031
1032 Int_t na,pos;
1033 Float_t bl,n;
e8189707 1034 char input[100], base[100], param[100];
b0f5e3fc 1035 char *filtmp;
1036
e8189707 1037 fResponse->Filenames(input,base,param);
b0f5e3fc 1038 fFileName=base;
1039//
1040 filtmp = gSystem->ExpandPathName(fFileName.Data());
1041 FILE *bline = fopen(filtmp,"r");
1042 printf("filtmp %s\n",filtmp);
1043 na = 0;
1044
1045 if(bline) {
1046 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
b0f5e3fc 1047 if (pos != na+1) {
1048 Error("ReadBaseline","Anode number not in increasing order!",
1049 filtmp);
1050 exit(1);
1051 }
1052 fBaseline[na]=bl;
1053 fNoise[na]=n;
1054 na++;
1055 }
1056 } else {
1057 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",
1058 filtmp);
1059 exit(1);
1060 } // end if(bline)
1061
1062 fclose(bline);
1063 delete [] filtmp;
b0f5e3fc 1064}
1065
1066//____________________________________________
1067Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) {
1068 // To the 10 to 8 bit lossive compression.
1069 // code from Davide C. and Albert W.
1070
1071 if (signal < 128) return signal;
1072 if (signal < 256) return (128+((signal-128)>>1));
1073 if (signal < 512) return (192+((signal-256)>>3));
1074 if (signal < 1024) return (224+((signal-512)>>4));
1075 return 0;
1076
1077}
1078
1079//____________________________________________
1080Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) {
1081 // Undo the lossive 10 to 8 bit compression.
1082 // code from Davide C. and Albert W.
1083 if (signal < 0 || signal > 255) {
1084 printf("<Convert8to10> out of range %d \n",signal);
1085 return 0;
1086 }
1087
1088 if (signal < 128) return signal;
1089 if (signal < 192) {
1090 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1091 else return (128+((signal-128)<<1)+1);
1092 }
1093 if (signal < 224) {
1094 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1095 else return (256+((signal-192)<<3)+4);
1096 }
1097 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1098 else return (512+((signal-224)<<4)+7);
1099 return 0;
1100
1101}
1102
1103//____________________________________________
1104AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1105 //Return the correct map.
1106 return ((i==0)? fHitMap1 : fHitMap2);
1107}
1108
1109
1110//____________________________________________
e8189707 1111void AliITSsimulationSDD::ZeroSuppression(const char *option) {
b0f5e3fc 1112 // perform the zero suppresion
1113 if (strstr(option,"2D")) {
1114 //Init2D(); // activate if param change module by module
1115 Compress2D();
1116 } else if (strstr(option,"1D")) {
1117 //Init1D(); // activate if param change module by module
1118 Compress1D();
1119 } else StoreAllDigits();
1120
1121}
1122
1123//____________________________________________
1124void AliITSsimulationSDD::Init2D(){
1125 // read in and prepare arrays: fD, fT1, fT2
1126 // savemu[nanodes], savesigma[nanodes]
1127 // read baseline and noise from file - either a .root file and in this
1128 // case data should be organised in a tree with one entry for each
1129 // module => reading should be done accordingly
1130 // or a classic file and do smth. like this ( code from Davide C. and
1131 // Albert W.) :
1132 //
1133 // Read 2D zero-suppression parameters for SDD
1134 //
1135
1136 if (!strstr(fParam,"file")) return;
1137
1138 Int_t na,pos,tempTh;
1139 Float_t mu,sigma;
5d18fa90 1140 Float_t *savemu = new Float_t [fNofMaps];
e8189707 1141 Float_t *savesigma = new Float_t [fNofMaps];
1142 char input[100],basel[100],par[100];
b0f5e3fc 1143 char *filtmp;
1144
1145
1146 Int_t minval = fResponse->MinVal();
1147
e8189707 1148 fResponse->Filenames(input,basel,par);
1149 fFileName=par;
b0f5e3fc 1150
1151//
1152 filtmp = gSystem->ExpandPathName(fFileName.Data());
1153 FILE *param = fopen(filtmp,"r");
1154 na = 0;
1155
1156 if(param) {
1157 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1158 if (pos != na+1) {
1159 Error("Init2D ","Anode number not in increasing order!",
1160 filtmp);
1161 exit(1);
1162 }
1163 savemu[na]=mu;
1164 savesigma[na]=sigma;
1165 if ((2.*sigma) < mu) {
1166 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1167 mu = 2.0 * sigma;
1168 } else fD[na] = 0;
1169 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1170 if (tempTh < 0) tempTh=0;
1171 fT1[na] = tempTh;
1172 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1173 if (tempTh < 0) tempTh=0;
1174 fT2[na] = tempTh;
1175 na++;
1176 } // end while
1177
1178 } else {
1179 Error("Init2D "," THE FILE %s DOES NOT EXIST !",
1180 filtmp);
1181 exit(1);
1182 } // end if(param)
1183
1184 fclose(param);
1185 delete [] filtmp;
5d18fa90 1186 delete [] savemu;
e8189707 1187 delete [] savesigma;
b0f5e3fc 1188}
e8189707 1189
b0f5e3fc 1190//____________________________________________
1191void AliITSsimulationSDD::Compress2D(){
1192 //
1193 // simple ITS cluster finder -- online zero-suppression conditions
1194 //
1195 //
1196
b0f5e3fc 1197 Int_t db,tl,th;
1198 Int_t minval = fResponse->MinVal();
1199 Bool_t write=fResponse->OutputOption();
ece86d9a 1200 Bool_t do10to8=fResponse->Do10to8();
b0f5e3fc 1201
1202 Int_t nz, nl, nh, low, i, j;
1203
e8189707 1204 for (i=0; i<fNofMaps; i++) {
b0f5e3fc 1205 CompressionParam(i,db,tl,th);
1206 nz=0;
1207 nl=0;
1208 nh=0;
1209 low=0;
e8189707 1210 for (j=0; j<fMaxNofSamples; j++) {
b0f5e3fc 1211 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1212 signal -= db; // if baseline eq. is done here
1213 if (signal <= 0) {nz++; continue;}
1214 if ((signal - tl) < minval) low++;
1215 if ((signal - th) >= minval) {
1216 nh++;
1217 Bool_t cond=kTRUE;
b0f5e3fc 1218 FindCluster(i,j,signal,minval,cond);
44a312c3 1219 if (cond && j && ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)) {
ece86d9a 1220 if(do10to8) signal = Convert10to8(signal);
1221 AddDigit(i,j,signal);
1222 }
b0f5e3fc 1223 } else if ((signal - tl) >= minval) nl++;
1224 } // loop time samples
1225 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
b0f5e3fc 1226 } // loop anodes
1227
1228 char hname[30];
1229 if (write) {
1230 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1231 TreeB()->Write(hname);
1232 // reset tree
1233 TreeB()->Reset();
1234 }
1235
1236}
1237
1238//_____________________________________________________________________________
1239void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
ece86d9a 1240 Int_t minval,Bool_t &cond){
b0f5e3fc 1241//
1242// Find clusters according to the online 2D zero-suppression algorithm
1243//
1244
ece86d9a 1245 Bool_t do10to8=fResponse->Do10to8();
1246
b0f5e3fc 1247 Bool_t high=kFALSE;
1248
1249 fHitMap2->FlagHit(i,j);
1250//
1251// check the online zero-suppression conditions
1252//
e8189707 1253 const Int_t maxNeighbours = 4;
1254
b0f5e3fc 1255 Int_t nn;
1256 Int_t dbx,tlx,thx;
e8189707 1257 Int_t xList[maxNeighbours], yList[maxNeighbours];
1258 fSegmentation->Neighbours(i,j,&nn,xList,yList);
ece86d9a 1259 Int_t in,ix,iy,qns;
e8189707 1260 for (in=0; in<nn; in++) {
1261 ix=xList[in];
1262 iy=yList[in];
b0f5e3fc 1263 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1264 CompressionParam(ix,dbx,tlx,thx);
1265 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1266 qn -= dbx; // if baseline eq. is done here
1267 if ((qn-tlx) < minval) {
1268 fHitMap2->FlagHit(ix,iy);
1269 continue;
1270 } else {
1271 if ((qn - thx) >= minval) high=kTRUE;
1272 if (cond) {
ece86d9a 1273 if(do10to8) signal = Convert10to8(signal);
b0f5e3fc 1274 AddDigit(i,j,signal);
1275 }
ece86d9a 1276 if(do10to8) qns = Convert10to8(qn);
1277 else qns=qn;
b0f5e3fc 1278 if (!high) AddDigit(ix,iy,qns);
1279 cond=kFALSE;
1280 if(!high) fHitMap2->FlagHit(ix,iy);
1281 }
1282 } // TestHit
1283 } // loop over neighbours
1284
1285}
1286
1287//____________________________________________
1288void AliITSsimulationSDD::Init1D(){
1289 // this is just a copy-paste of input taken from 2D algo
1290 // Torino people should give input
1291 //
1292 // Read 1D zero-suppression parameters for SDD
1293 //
1294
1295 if (!strstr(fParam,"file")) return;
1296
1297 Int_t na,pos,tempTh;
1298 Float_t mu,sigma;
e8189707 1299 Float_t *savemu = new Float_t [fNofMaps];
1300 Float_t *savesigma = new Float_t [fNofMaps];
1301 char input[100],basel[100],par[100];
b0f5e3fc 1302 char *filtmp;
1303
1304
1305 Int_t minval = fResponse->MinVal();
e8189707 1306 fResponse->Filenames(input,basel,par);
1307 fFileName=par;
b0f5e3fc 1308
1309// set first the disable and tol param
1310 SetCompressParam();
1311//
1312 filtmp = gSystem->ExpandPathName(fFileName.Data());
1313 FILE *param = fopen(filtmp,"r");
1314 na = 0;
1315
1316 if (param) {
1317 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1318 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1319 if (pos != na+1) {
1320 Error("Init1D ","Anode number not in increasing order!",
1321 filtmp);
1322 exit(1);
1323 }
1324 savemu[na]=mu;
1325 savesigma[na]=sigma;
1326 if ((2.*sigma) < mu) {
1327 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1328 mu = 2.0 * sigma;
1329 } else fD[na] = 0;
1330 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1331 if (tempTh < 0) tempTh=0;
1332 fT1[na] = tempTh;
1333 na++;
1334 } // end while
1335 } else {
1336 Error("Init1D "," THE FILE %s DOES NOT EXIST !",
1337 filtmp);
1338 exit(1);
1339 } // end if(param)
1340
1341 fclose(param);
1342 delete [] filtmp;
749bd21a 1343 delete [] savemu;
e8189707 1344 delete [] savesigma;
1345
1346
b0f5e3fc 1347
1348}
1349
1350//____________________________________________
1351void AliITSsimulationSDD::Compress1D(){
1352 // 1D zero-suppression algorithm (from Gianluca A.)
1353
1354 Int_t dis,tol,thres,decr,diff;
b0f5e3fc 1355
1356 UChar_t *str=fStream->Stream();
1357 Int_t counter=0;
1358
ece86d9a 1359 Bool_t do10to8=fResponse->Do10to8();
1360
e8189707 1361 Int_t last=0;
1362 Int_t k,i,j;
ece86d9a 1363 for (k=0; k<2; k++) {
1364 tol = Tolerance(k);
1365 dis = Disable(k);
e8189707 1366 for (i=0; i<fNofMaps/2; i++) {
b0f5e3fc 1367 Bool_t firstSignal=kTRUE;
ece86d9a 1368 Int_t idx=i+k*fNofMaps/2;
1369 CompressionParam(idx,decr,thres);
e8189707 1370 for (j=0; j<fMaxNofSamples; j++) {
ece86d9a 1371 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
b0f5e3fc 1372 signal -= decr; // if baseline eq.
ece86d9a 1373 if(do10to8) signal = Convert10to8(signal);
1374 if (signal <= thres) {
b0f5e3fc 1375 signal=0;
1376 diff=128;
1377 last=0;
1378 // write diff in the buffer for HuffT
1379 str[counter]=(UChar_t)diff;
1380 counter++;
1381 continue;
1382 }
1383 diff=signal-last;
1384 if (diff > 127) diff=127;
1385 if (diff < -128) diff=-128;
1386
1387 if (signal < dis) {
ece86d9a 1388 // tol has changed to 8 possible cases ? - one can write
1389 // this if(TMath::Abs(diff)<tol) ... else ...
1390 if(TMath::Abs(diff)<tol) diff=0;
1391 // or keep it as it was before
1392 /*
b0f5e3fc 1393 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1394 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1395 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
ece86d9a 1396 */
1397 AddDigit(idx,j,last+diff);
b0f5e3fc 1398 } else {
ece86d9a 1399 AddDigit(idx,j,signal);
b0f5e3fc 1400 }
1401
1402 diff += 128;
1403 // write diff in the buffer used to compute Huffman tables
1404 if (firstSignal) str[counter]=(UChar_t)signal;
1405 else str[counter]=(UChar_t)diff;
1406 counter++;
1407
1408 last=signal;
1409 firstSignal=kFALSE;
1410 } // loop time samples
1411 } // loop anodes one half of detector
1412 }
1413
1414 // check
1415 fStream->CheckCount(counter);
1416
1417 // open file and write out the stream of diff's
1418
1419 static Bool_t open=kTRUE;
e8189707 1420 static TFile *outFile;
b0f5e3fc 1421 Bool_t write = fResponse->OutputOption();
1422
1423 if (write ) {
1424 if(open) {
1425 SetFileName("stream.root");
1426 cout<<"filename "<<fFileName<<endl;
e8189707 1427 outFile=new TFile(fFileName,"recreate");
b0f5e3fc 1428 cout<<"I have opened "<<fFileName<<" file "<<endl;
1429 }
1430 open=kFALSE;
e8189707 1431 outFile->cd();
b0f5e3fc 1432 fStream->Write();
1433 } // endif write
1434
1435 fStream->ClearStream();
1436
1437 // back to galice.root file
1438
1439 TTree *fAli=gAlice->TreeK();
1440 TFile *file = 0;
1441
1442 if (fAli) file =fAli->GetCurrentFile();
1443 file->cd();
1444
1445
1446}
1447//____________________________________________
1448void AliITSsimulationSDD::StoreAllDigits(){
ece86d9a 1449 // if non-zero-suppressed data
b0f5e3fc 1450
ece86d9a 1451 Bool_t do10to8=fResponse->Do10to8();
b0f5e3fc 1452
ece86d9a 1453 Int_t i, j, digits[3];
e8189707 1454 for (i=0; i<fNofMaps; i++) {
1455 for (j=0; j<fMaxNofSamples; j++) {
b0f5e3fc 1456 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
ece86d9a 1457 if(do10to8) signal = Convert10to8(signal);
1458 if(do10to8) signal = Convert8to10(signal);
b0f5e3fc 1459 digits[0]=i;
1460 digits[1]=j;
1461 digits[2]=signal;
ece86d9a 1462 fITS->AddRealDigit(1,digits);
b0f5e3fc 1463 }
1464 }
1465}
1466//____________________________________________
1467
ece86d9a 1468void AliITSsimulationSDD::CreateHistograms(Int_t scale){
b0f5e3fc 1469 // Creates histograms of maps for debugging
1470
1471 Int_t i;
ece86d9a 1472
1473 fHis=new TObjArray(fNofMaps);
e8189707 1474 for (i=0;i<fNofMaps;i++) {
a14dba92 1475 TString sddName("sdd_");
b0f5e3fc 1476 Char_t candNum[4];
1477 sprintf(candNum,"%d",i+1);
ece86d9a 1478 sddName.Append(candNum);
1479 (*fHis)[i] = new TH1F(sddName.Data(),"SDD maps",
1480 scale*fMaxNofSamples,0.,(Float_t) scale*fMaxNofSamples);
b0f5e3fc 1481 }
1482
1483}
ece86d9a 1484//____________________________________________
1485void AliITSsimulationSDD::FillHistograms(){
1486 // fill 1D histograms from map
1487 if (!fHis) return;
1488
1489 for( Int_t i=0; i<fNofMaps; i++) {
1490 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1491 Int_t nsamples = hist->GetNbinsX();
1492 for( Int_t j=0; j<nsamples; j++) {
1493 Double_t signal=fHitMap2->GetSignal(i,j);
1494 hist->Fill((Float_t)j,signal);
1495 }
1496 }
1497}
1498
b0f5e3fc 1499//____________________________________________
1500
1501void AliITSsimulationSDD::ResetHistograms(){
1502 //
1503 // Reset histograms for this detector
1504 //
1505 Int_t i;
e8189707 1506 for (i=0;i<fNofMaps;i++ ) {
b0f5e3fc 1507 if ((*fHis)[i]) ((TH1F*)(*fHis)[i])->Reset();
1508 }
1509
1510}
1511
1512
1513//____________________________________________
1514
1515TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1516 // Fills a histogram from a give anode.
1517 if (!fHis) return 0;
1518
1519 if(wing <=0 || wing > 2) {
1520 cout << "Wrong wing number: " << wing << endl;
1521 return NULL;
1522 }
1523 if(anode <=0 || anode > fNofMaps/2) {
1524 cout << "Wrong anode number: " << anode << endl;
1525 return NULL;
1526 }
1527
1528 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1529 return (TH1F*)((*fHis)[index]);
1530}
1531
1532//____________________________________________
1533
1534void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1535 // Writes the histograms to a file
1536 if (!fHis) return;
1537
1538 hfile->cd();
1539 Int_t i;
1540 for(i=0; i<fNofMaps; i++) (*fHis)[i]->Write(); //fAdcs[i]->Write();
1541 return;
1542}
1543//____________________________________________
ece86d9a 1544Float_t AliITSsimulationSDD::GetNoise() {
b0f5e3fc 1545 // Returns the noise value
b0f5e3fc 1546
ece86d9a 1547 //Bool_t do10to8=fResponse->Do10to8();
1548 //noise will always be in the liniar part of the signal
1549
1550 Int_t decr;
1551 Int_t threshold=fT1[0];
1552
1553 char opt1[20], opt2[20];
1554 fResponse->ParamOptions(opt1,opt2);
1555 fParam=opt2;
1556 char *same = strstr(opt1,"same");
1557 Float_t noise,baseline;
1558 if (same) {
1559 fResponse->GetNoiseParam(noise,baseline);
1560 } else {
1561 static Bool_t readfile=kTRUE;
1562 //read baseline and noise from file
1563 if (readfile) ReadBaseline();
1564 readfile=kFALSE;
1565 }
1566
1567 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1568 if(c2) delete c2->GetPrimitive("noisehist");
1569 if(c2) delete c2->GetPrimitive("anode");
1570 else c2=new TCanvas("c2");
1571 c2->cd();
1572 c2->SetFillColor(0);
1573
1574 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1575 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,(float)fMaxNofSamples);
b0f5e3fc 1576 Int_t i,k;
e8189707 1577 for (i=0;i<fNofMaps;i++) {
ece86d9a 1578 CompressionParam(i,decr,threshold);
1579 if (!same) GetAnodeBaseline(i,baseline,noise);
1580 anode->Reset();
1581 for (k=0;k<fMaxNofSamples;k++) {
1582 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1583 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1584 if (signal <= (float)threshold) noisehist->Fill(signal);
1585 anode->Fill((float)k,signal);
b0f5e3fc 1586 }
ece86d9a 1587 anode->Draw();
1588 c2->Update();
b0f5e3fc 1589 }
ece86d9a 1590 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
b0f5e3fc 1591 noisehist->Fit("gnoise","RQ");
1592 noisehist->Draw();
ece86d9a 1593 c2->Update();
b0f5e3fc 1594 Float_t mnoise = gnoise->GetParameter(1);
1595 cout << "mnoise : " << mnoise << endl;
1596 Float_t rnoise = gnoise->GetParameter(2);
1597 cout << "rnoise : " << rnoise << endl;
1598 delete noisehist;
1599 return rnoise;
1600}
44a312c3 1601
1602//____________________________________________
1603
1604void AliITSsimulationSDD::Print() {
1605
1606 // Print SDD simulation Parameters
1607
1608 cout << "**************************************************" << endl;
1609 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1610 cout << "**************************************************" << endl;
1611 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1612 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1613 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1614 cout << "Number pf Anodes used: " << fNofMaps << endl;
1615 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1616 cout << "Scale size factor: " << fScaleSize << endl;
1617 cout << "**************************************************" << endl;
1618}