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