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