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