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