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