]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSsimulationSDD.cxx
o fix for coverity bug 19907
[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
18 #include <Riostream.h>
19 #include <stdlib.h>
20 #include <stdio.h>
21 #include <cstring>
22
23 #include <TCanvas.h>
24 #include <TF1.h>
25 #include <TH1.h>
26 #include <TFile.h>
27 #include <TRandom.h>
28 #include <TROOT.h>
29 #include "AliITS.h"
30 #include "AliITSMapA2.h"
31 #include "AliITSdigitSPD.h"
32 #include "AliITSetfSDD.h"
33 #include "AliITSmodule.h"
34 #include "AliITShit.h"
35 #include "AliITSpList.h"
36 #include "AliITSCalibrationSDD.h"
37 #include "AliITSresponseSDD.h"
38 #include "AliITSsimulationSDD.h"
39 #include "AliLog.h"
40 #include "AliRun.h"
41
42 ClassImp(AliITSsimulationSDD)
43 ////////////////////////////////////////////////////////////////////////
44 // Version: 0                                                         //
45 // Written by Piergiorgio Cerello                                     //
46 // November 23 1999                                                   //
47 //                                                                    //
48 // AliITSsimulationSDD is the simulation of SDDs.                     //
49 ////////////////////////////////////////////////////////////////////////
50
51 //______________________________________________________________________
52 AliITSsimulationSDD::AliITSsimulationSDD():
53 AliITSsimulation(),
54 fITS(0),
55 fHitMap2(0),
56 fHitSigMap2(0),
57 fHitNoiMap2(0),
58 fElectronics(0),
59 fInZR(0),
60 fInZI(0),
61 fOutZR(0),
62 fOutZI(0),
63 fAnodeFire(0),
64 fHis(0),
65 fFlag(kFALSE),
66 fCrosstalkFlag(kFALSE),
67 fDoFFT(1),
68 fNofMaps(0),
69 fMaxNofSamples(0),
70 fScaleSize(0){
71     // Default constructor
72     SetPerpendTracksFlag();
73     SetCrosstalkFlag();
74     SetDoFFT();
75 }
76
77 //______________________________________________________________________
78 AliITSsimulationSDD::AliITSsimulationSDD(AliITSDetTypeSim* dettyp):
79 AliITSsimulation(dettyp),
80 fITS(0),
81 fHitMap2(0),
82 fHitSigMap2(0),
83 fHitNoiMap2(0),
84 fElectronics(0),
85 fInZR(0),
86 fInZI(0),
87 fOutZR(0),
88 fOutZI(0),
89 fAnodeFire(0),
90 fHis(0),
91 fFlag(kFALSE),
92 fCrosstalkFlag(kFALSE),
93 fDoFFT(1),
94 fNofMaps(0),
95 fMaxNofSamples(0),
96 fScaleSize(0){
97     // Default Constructor
98   Init();
99 }
100 //______________________________________________________________________
101 void AliITSsimulationSDD::Init(){
102     // Standard Constructor
103
104     AliITSsegmentationSDD* seg = (AliITSsegmentationSDD*)GetSegmentationModel(1);
105     fScaleSize = ScaleFourier(seg);
106     SetPerpendTracksFlag();
107     SetCrosstalkFlag();
108     SetDoFFT();
109
110     AliITSSimuParam* simpar = fDetType->GetSimuParam();
111     fpList = new AliITSpList( seg->Npz(),
112                               fScaleSize*seg->Npx() );
113     fHitSigMap2 = new AliITSMapA2(seg,fScaleSize,1);
114     fHitNoiMap2 = new AliITSMapA2(seg,fScaleSize,1);
115     fHitMap2 = fHitSigMap2;
116
117     fNofMaps = seg->Npz();
118     fMaxNofSamples = seg->Npx();
119     fAnodeFire = new Bool_t [fNofMaps];
120     
121     Float_t sddWidth  = seg->Dz();
122     Float_t anodePitch = seg->Dpz(0);
123     Double_t timeStep  = (Double_t)seg->Dpx(0);
124
125     if(anodePitch*(fNofMaps/2) > sddWidth) {
126       AliWarning(Form("Too many anodes %d or too big pitch %f ",
127                 fNofMaps/2,anodePitch));
128     } // end if
129
130
131     fElectronics = new AliITSetfSDD(timeStep/fScaleSize,
132                                     simpar->GetSDDElectronics());
133
134
135     fITS       = (AliITS*)gAlice->GetModule("ITS");
136  
137     fInZR  = new Double_t [fScaleSize*fMaxNofSamples];
138     fInZI  = new Double_t [fScaleSize*fMaxNofSamples];
139     fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
140     fOutZI = new Double_t [fScaleSize*fMaxNofSamples];  
141 }
142 //______________________________________________________________________
143 AliITSsimulationSDD::~AliITSsimulationSDD() { 
144     // destructor
145
146     //    delete fpList;
147     delete fHitSigMap2;
148     delete fHitNoiMap2;
149     delete fElectronics;
150
151     fITS = 0;
152
153     if (fHis) {
154         fHis->Delete(); 
155         delete fHis;     
156     } // end if fHis
157     if(fInZR)  delete [] fInZR;
158     if(fInZI)  delete [] fInZI;        
159     if(fOutZR) delete [] fOutZR;
160     if(fOutZI) delete [] fOutZI;
161     if(fAnodeFire) delete [] fAnodeFire;
162 }
163 //______________________________________________________________________
164 void AliITSsimulationSDD::InitSimulationModule( Int_t module, Int_t event ) {
165     // create maps to build the lists of tracks for each summable digit
166     fModule = module;
167     fEvent  = event;
168     ClearMaps();
169     memset(fAnodeFire,0,sizeof(Bool_t)*fNofMaps);    
170 }
171 //______________________________________________________________________
172 void AliITSsimulationSDD::ClearMaps() {
173     // clear maps
174     fpList->ClearMap();
175     fHitSigMap2->ClearMap();
176     fHitNoiMap2->ClearMap();
177 }
178 //______________________________________________________________________
179 void AliITSsimulationSDD::FastFourierTransform(Double_t *real,
180                           Double_t *imag,Int_t direction) {
181     // Do a Fast Fourier Transform
182
183     Int_t samples = fElectronics->GetSamples();
184     Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5);
185     Int_t m1 = samples;
186     Int_t m  = samples/2;
187     Int_t m2 = samples/m1;
188     Int_t i,j,k;
189     for(i=1; i<=l; i++) {
190         for(j=0; j<samples; j += m1) {
191             Int_t p = 0;
192             for(k=j; k<= j+m-1; k++) {
193                 Double_t wsr = fElectronics->GetWeightReal(p);
194                 Double_t wsi = fElectronics->GetWeightImag(p);
195                 if(direction == -1) wsi = -wsi;
196                 Double_t xr = *(real+k+m);
197                 Double_t xi = *(imag+k+m);
198                 *(real+k+m) = wsr*(*(real+k)-xr) - wsi*(*(imag+k)-xi);
199                 *(imag+k+m) = wsr*(*(imag+k)-xi) + wsi*(*(real+k)-xr);
200                 *(real+k) += xr;
201                 *(imag+k) += xi;
202                 p += m2;
203             } // end for k
204         } // end for j
205         m1 = m;
206         m /= 2;
207         m2 += m2;
208     } // end for i
209     for(j=0; j<samples; j++) {
210         Int_t j1 = j;
211         Int_t p = 0;
212         Int_t i1;
213         for(i1=1; i1<=l; i1++) {
214             Int_t j2 = j1;
215             j1 /= 2;
216             p = p + p + j2 - j1 - j1;
217         } // end for i1
218         if(p >= j) {
219             Double_t xr = *(real+j);
220             Double_t xi = *(imag+j);
221             *(real+j) = *(real+p);
222             *(imag+j) = *(imag+p);
223             *(real+p) = xr;
224             *(imag+p) = xi;
225         } // end if p>=j
226     } // end for j
227     if(direction == -1) {
228         for(i=0; i<samples; i++) {
229             *(real+i) /= samples;
230             *(imag+i) /= samples;
231         } // end for i
232     } // end if direction == -1
233     return;
234 }
235
236 //______________________________________________________________________
237 void AliITSsimulationSDD::SDigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
238     // digitize module using the "slow" detector simulator creating
239     // summable digits.
240
241     TObjArray *fHits = mod->GetHits();
242     Int_t nhits      = fHits->GetEntriesFast();
243     if( !nhits ) return;
244
245     InitSimulationModule( md, ev );
246     HitsToAnalogDigits( mod );  // fills fHitMap2 which is = fHitSigmap2
247     ChargeToSignal( fModule,kFALSE,kTRUE ); // - Process signal adding gain without adding noise
248     fHitMap2 = fHitNoiMap2;   // - Swap to noise map
249     ChargeToSignal( fModule,kTRUE,kFALSE );  // - Process only noise
250     fHitMap2 = fHitSigMap2;   // - Return to signal map
251     WriteSDigits();
252     ClearMaps();
253 }
254 //______________________________________________________________________
255 Bool_t AliITSsimulationSDD::AddSDigitsToModule(TClonesArray *pItemArray,
256                                                Int_t mask ) {
257     // Add Summable digits to module maps.
258     AliITSSimuParam* simpar = fDetType->GetSimuParam();
259     Int_t    nItems = pItemArray->GetEntries();
260     Double_t maxadc = simpar->GetSDDMaxAdc();
261     Bool_t sig = kFALSE;
262     
263     // cout << "Adding "<< nItems <<" SDigits to module " << fModule << endl;
264     for( Int_t i=0; i<nItems; i++ ) {
265         AliITSpListItem * pItem = (AliITSpListItem *)(pItemArray->At( i ));
266         if( pItem->GetModule() != fModule ) {
267             Error( "AliITSsimulationSDD","Error reading, SDigits module "
268                    "%d != current module %d: exit",
269                    pItem->GetModule(), fModule );
270             return sig;
271         } // end if
272
273         if(pItem->GetSignal()>0.0 ) sig = kTRUE;
274         
275         fpList->AddItemTo( mask, pItem ); // Add SignalAfterElect + noise
276         AliITSpListItem * pItem2 = fpList->GetpListItem( pItem->GetIndex() );
277         Double_t sigAE = pItem2->GetSignalAfterElect();
278         if( sigAE >= maxadc ) sigAE = maxadc-1; // avoid overflow signal
279         Int_t ia;
280         Int_t it;
281         fpList->GetMapIndex( pItem->GetIndex(), ia, it );
282         fHitMap2->SetHit( ia, it, sigAE );
283         fAnodeFire[ia] = kTRUE;
284     }
285     return sig;
286 }
287 //______________________________________________________________________
288 void AliITSsimulationSDD::FinishSDigitiseModule() {
289     // digitize module using the "slow" detector simulator from
290     // the sum of summable digits.
291     FinishDigits() ;
292     ClearMaps();
293 }
294 //______________________________________________________________________
295 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
296     // create maps to build the lists of tracks for each digit
297
298     TObjArray *fHits = mod->GetHits();
299     Int_t nhits      = fHits->GetEntriesFast();
300
301     InitSimulationModule( md, ev );
302     if( !nhits ) return;
303         
304     HitsToAnalogDigits( mod );
305     ChargeToSignal( fModule,kTRUE,kTRUE );  // process signal + noise
306
307     for( Int_t i=0; i<fNofMaps; i++ ) {
308         for( Int_t j=0; j<fMaxNofSamples; j++ ) {
309             Int_t jdx = j*fScaleSize;
310             Int_t index = fpList->GetHitIndex( i, j );
311             AliITSpListItem pItemTmp2( fModule, index, 0. );
312             // put the fScaleSize analog digits in only one
313             for( Int_t ik=0; ik<fScaleSize; ik++ ) {
314                 AliITSpListItem *pItemTmp = fpList->GetpListItem( i, jdx+ik );
315                 if( pItemTmp == 0 ) continue;
316                 pItemTmp2.Add( pItemTmp );
317             }
318             fpList->DeleteHit( i, j );
319             fpList->AddItemTo( 0, &pItemTmp2 );
320         }
321     }
322     FinishDigits();
323     ClearMaps();
324 }
325 //______________________________________________________________________
326 void AliITSsimulationSDD::FinishDigits() {
327     // introduce the electronics effects and do zero-suppression if required
328
329     if( fCrosstalkFlag ) ApplyCrosstalk(fModule);
330
331     AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(fModule);
332     Bool_t isZeroSupp = res->GetZeroSupp();
333     if (isZeroSupp) Compress2D();
334     else StoreAllDigits();
335 }
336 //______________________________________________________________________
337 void AliITSsimulationSDD::HitsToAnalogDigits( AliITSmodule *mod ) {
338     // create maps to build the lists of tracks for each digit
339   AliITSsegmentationSDD* seg = (AliITSsegmentationSDD*)GetSegmentationModel(1);
340   AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(fModule);
341   AliITSSimuParam* simpar = fDetType->GetSimuParam();
342   TObjArray *hits     = mod->GetHits();
343   Int_t      nhits    = hits->GetEntriesFast();
344
345   //    Int_t      arg[6]   = {0,0,0,0,0,0};
346   Int_t     nofAnodes  = fNofMaps/2;
347   Double_t  sddLength  = seg->Dx();
348   Double_t  anodePitch = seg->Dpz(0);
349   Double_t  timeStep   = seg->Dpx(0);
350   Double_t  driftSpeed ;  // drift velocity (anode dependent)
351   Double_t  nanoampToADC       = simpar->GetSDDMaxAdc()/simpar->GetSDDDynamicRange(); //   maxadc/topValue;
352   Double_t  cHloss     = simpar->GetSDDChargeLoss();
353   Float_t   dfCoeff, s1; 
354   simpar->GetSDDDiffCoeff(dfCoeff,s1); // Signal 2d Shape
355   Double_t  eVpairs    = simpar->GetGeVToCharge()*1.0E9; // 3.6 eV by def.
356   Double_t  nsigma     = simpar->GetNSigmaIntegration(); //
357   Int_t     nlookups   = simpar->GetGausNLookUp();       //
358   Float_t   jitter     = simpar->GetSDDJitterError(); // 
359   Float_t   mapsmear   = simpar->GetSDDCorrMapPrecision(); // 
360   Float_t   trigDelay  = simpar->GetSDDTrigDelay(); // compensation for MC time zero
361   if(res->IsAMAt20MHz()) trigDelay+=12.5; // compensation for discretization step
362
363   Float_t   timeZero=fDetType->GetResponseSDD()->GetTimeZero(fModule);
364   Float_t   adcscale   = fDetType->GetResponseSDD()->GetADCtokeV(fModule);
365   adcscale/=simpar->GetSDDkeVtoADC();
366
367   // Piergiorgio's part (apart for few variables which I made float
368   // when i thought that can be done
369   // Fill detector maps with GEANT hits
370   // loop over hits in the module
371   
372   const Float_t kconv = 1.0e+6;  // GeV->KeV
373   Int_t     itrack      = 0;
374   Int_t     iWing;       // which detector wing/side.
375   Int_t     ii,kk,ka,kt; // loop indexs
376   Int_t     ia,it,index; // sub-pixel integration indexies
377   Int_t     iAnode;      // anode number.
378   Int_t     timeSample;  // time buckett.
379   Int_t     anodeWindow; // anode direction charge integration width
380   Int_t     timeWindow;  // time direction charge integration width
381   Int_t     jamin,jamax; // anode charge integration window
382   Int_t     jtmin,jtmax; // time charge integration window
383   Int_t     nsplitAn;    // the number of splits in anode and time windows
384   Int_t     nsplitTb;    // the number of splits in anode and time windows
385   Int_t     nOfSplits;   // number of times track length is split into
386   Float_t   nOfSplitsF;  // Floating point version of nOfSplits.
387   Float_t   kkF;         // Floating point version of loop index kk.
388   Double_t  pathInSDD; // Track length in SDD.
389   Double_t  drPath; // average position of track in detector. in microns
390   Double_t  drTime; // Drift time
391   Double_t  avDrft;  // x position of path length segment in cm.
392   Double_t  avAnode; // Anode for path length segment in Anode number (float)
393   Double_t  zAnode;  // Floating point anode number.
394   Double_t  driftPath; // avDrft in microns.
395   Double_t  width;     // width of signal at anodes.
396   Double_t  depEnergy; // Energy deposited in this GEANT step.
397   Double_t  xL[3],dxL[3]; // local hit coordinates and diff.
398   Double_t  sigA; // sigma of signal at anode.
399   Double_t  sigT; // sigma in time/drift direction for track segment
400   Double_t  aStep,aConst; // sub-pixel size and offset anode
401   Double_t  tStep,tConst; // sub-pixel size and offset time
402   Double_t  amplitude; // signal amplitude for track segment in nanoAmpere
403   Double_t  chargeloss; // charge loss for track segment.
404   Double_t  anodeAmplitude; // signal amplitude in anode direction
405   Double_t  aExpo;          // exponent of Gaussian anode direction
406   Double_t  timeAmplitude;  // signal amplitude in time direction
407   Double_t  tExpo;          // exponent of Gaussian time direction
408   Double_t  tof;            // Time of flight in ns of this step.    
409   
410   for(ii=0; ii<nhits; ii++) {
411     if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
412                           depEnergy,itrack)) continue;
413     Float_t xloc=xL[0];
414     if(xloc>0) iWing=0; // left side, carlos channel 0
415     else iWing=1; // right side
416     
417     Float_t zloc=xL[2]+0.5*dxL[2];
418     zAnode=seg->GetAnodeFromLocal(xloc,zloc); // anode number in the range 0.-511.
419     driftSpeed = res->GetDriftSpeedAtAnode(zAnode);
420     driftSpeed+= fDetType->GetResponseSDD()->GetDeltaVDrift(fModule,zAnode>255);
421
422     if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
423       AliWarning("Time Interval > Allowed Time Interval");
424     }
425     depEnergy  *= kconv;
426     if (!depEnergy) {
427       AliDebug(1,
428                Form("fTrack = %d hit=%d module=%d This particle has passed without losing energy!",
429                     itrack,ii,mod->GetIndex()));
430       continue;
431       // continue if the particle did not lose energy
432       // passing through detector
433     } // end if !depEnergy
434      
435     tof=0.;
436     AliITShit* h=(AliITShit*)hits->At(ii);
437     if(h){ 
438       tof=h->GetTOF()*1E9; 
439       AliDebug(1,Form("TOF for hit %d on mod %d (particle %d)=%g",ii,fModule,h->Track(),tof));
440     }
441
442     Float_t corrx=0, corrz=0;
443     res->GetShiftsForSimulation(xL[2],xL[0],corrz,corrx,seg);
444     xL[2]-=corrz;
445     xL[0]-=corrx;
446     xL[0] += 0.0001*gRandom->Gaus( 0, mapsmear); //
447     xL[0] += 0.0001*gRandom->Gaus( 0, jitter ); //
448
449     pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
450     
451     if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
452     drPath = TMath::Abs(10000.*(dxL[0]+2.*xL[0])*0.5);
453     drPath = sddLength-drPath;
454     if(drPath < 0) {
455       AliInfo( // this should be fixed at geometry level
456                Form("negative drift path drPath=%e sddLength=%e dxL[0]=%e xL[0]=%e",
457                     drPath,sddLength,dxL[0],xL[0]));
458       continue;
459     } // end if drPath < 0
460     
461     // Compute number of segments to brake step path into
462     drTime = drPath/driftSpeed;  //   Drift Time
463     sigA   = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
464     // calcuate the number of time the path length should be split into.
465     nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
466     if(fFlag) nOfSplits = 1;
467     
468     // loop over path segments, init. some variables.
469     depEnergy /= nOfSplits;
470     nOfSplitsF = (Float_t) nOfSplits;
471     Float_t theAverage=0.,theSteps=0.;
472     for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
473       kkF       = (Float_t) kk + 0.5;
474       avDrft    = xL[0]+dxL[0]*kkF/nOfSplitsF;
475       avAnode   = xL[2]+dxL[2]*kkF/nOfSplitsF;
476       theSteps+=1.;
477       theAverage+=avAnode;
478       zAnode = seg->GetAnodeFromLocal(avDrft,avAnode);
479       driftSpeed = res->GetDriftSpeedAtAnode(zAnode);   
480       driftSpeed+= fDetType->GetResponseSDD()->GetDeltaVDrift(fModule,zAnode>255);
481       driftPath = TMath::Abs(10000.*avDrft);
482       driftPath = sddLength-driftPath;
483       if(driftPath < 0) {
484         AliDebug(1, // this should be fixed at geometry level
485                  Form("negative drift path driftPath=%e sddLength=%e avDrft=%e dxL[0]=%e xL[0]=%e",
486                       driftPath,sddLength,avDrft,dxL[0],xL[0]));
487         continue;
488       } // end if driftPath < 0
489       drTime     = driftPath/driftSpeed; // drift time for segment.
490       // Sigma along the anodes for track segment.
491       sigA       = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
492       sigT       = sigA/driftSpeed;
493
494       drTime+=tof; // take into account Time Of Flight from production point
495       drTime-=trigDelay;
496       drTime+=timeZero;
497       timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1.001); // time bin in range 1-256 !!!
498       if(zAnode>nofAnodes) zAnode-=nofAnodes;  // to have the anode number between 0. and 256.
499       iAnode = (Int_t) (1.001+zAnode); // iAnode in range 1-256 !!!!
500       
501         // Peak amplitude in nanoAmpere
502       amplitude  = fScaleSize*160.*depEnergy/
503         (timeStep*eVpairs*2.*acos(-1.));
504       chargeloss = 1.-cHloss*driftPath/1000.;
505       amplitude *= chargeloss;
506       amplitude *= adcscale;
507       width  = 2.*nsigma/(nlookups-1);
508       // Spread the charge 
509       nsplitAn = 4; 
510       nsplitTb=4;
511       aStep  = anodePitch/(nsplitAn*sigA);
512       aConst = zAnode*anodePitch/sigA;
513       tStep  = timeStep/(nsplitTb*fScaleSize*sigT);
514       tConst = drTime/sigT;
515       // Define SDD window corresponding to the hit
516       anodeWindow = (Int_t)(nsigma*sigA/anodePitch+1);
517       timeWindow  = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
518       jamin = (iAnode - anodeWindow - 2)*nsplitAn+1;
519       if(jamin <= 0) jamin = 1;
520       if(jamin > nofAnodes*nsplitAn){ 
521         AliDebug(1,Form("Energy deposition completely outside anode acceptance: anode min=%d",jamin));
522         continue;
523       }
524       jamax = (iAnode + anodeWindow + 2)*nsplitAn;
525       if(jamax > nofAnodes*nsplitAn) jamax = nofAnodes*nsplitAn;
526       if(jamax <=0){ 
527         AliDebug(1,Form("Energy deposition completely outside anode acceptance: anode max=%d",jamax));
528         continue;
529       }
530       jtmin = (Int_t)(timeSample-timeWindow-2)*nsplitTb+1;
531       if(jtmin <= 0) jtmin = 1;
532       if(jtmin > fScaleSize*fMaxNofSamples*nsplitTb){ 
533         AliDebug(1,Form("Energy deposition completely outside time acceptance: time sample min=%d  tof=%f",jtmin,tof));
534         continue; 
535       }
536       jtmax = (Int_t)(timeSample+timeWindow+2)*nsplitTb;
537       if(jtmax > fScaleSize*fMaxNofSamples*nsplitTb) jtmax = fScaleSize*fMaxNofSamples*nsplitTb;
538       if(jtmax <= 0){
539         AliDebug(1,Form("Energy deposition completely outside time acceptance: time sample max=%d  tof=%f",jtmax,tof));
540         continue; 
541       }
542
543       // Spread the charge in the anode-time window
544       for(ka=jamin; ka <=jamax; ka++) {   
545         ia = (ka-1)/nsplitAn + 1;
546         if(ia <= 0) ia=1; 
547         if(ia > nofAnodes) ia = nofAnodes;
548         aExpo     = (aStep*(ka-0.5)-aConst);
549         if(TMath::Abs(aExpo) > nsigma)  anodeAmplitude = 0.;
550         else {
551           Int_t theBin = (Int_t) ((aExpo+nsigma)/width+0.5);
552           anodeAmplitude = amplitude*simpar->GetGausLookUp(theBin);
553         }
554         // index starts from 0
555         index = iWing*nofAnodes+ia-1;
556         if(anodeAmplitude){
557           for(kt=jtmin; kt<=jtmax; kt++) {
558             it = (kt-1)/nsplitTb+1;  // it starts from 1
559             if(it<=0) it=1;
560             if(it>fScaleSize*fMaxNofSamples)
561               it = fScaleSize*fMaxNofSamples;
562             tExpo    = (tStep*(kt-0.5)-tConst);
563             if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
564             else {
565               Int_t theBin = (Int_t) ((tExpo+nsigma)/width+0.5);
566               timeAmplitude = anodeAmplitude*simpar->GetGausLookUp(theBin)*aStep*tStep;
567             }
568             timeAmplitude *= nanoampToADC;
569             //         ListOfFiredCells(arg,timeAmplitude,alst,padr);
570             Double_t charge = timeAmplitude;
571             charge += fHitMap2->GetSignal(index,it-1);
572             fHitMap2->SetHit(index, it-1, charge);
573             fpList->AddSignal(index,it-1,itrack,ii-1,
574                               mod->GetIndex(),timeAmplitude);
575             fAnodeFire[index] = kTRUE;
576           }  // end loop over time in window               
577         } // end if anodeAmplitude 
578       } // loop over anodes in window
579     } // end loop over "sub-hits"
580   } // end loop over hits
581 }
582
583 //____________________________________________
584 void AliITSsimulationSDD::AddDigit( Int_t i, Int_t j, Int_t signalc, Int_t signale) {
585   // Adds a Digit.
586   Int_t size = AliITSdigit::GetNTracks();
587
588   Int_t digits[3];
589   Int_t * tracks = new Int_t[size];
590   Int_t * hits = new Int_t[size];
591   Float_t phys;
592   Float_t * charges = new Float_t[size];
593
594   digits[0] = i;
595   digits[1] = j;
596   digits[2] = signalc;
597
598   AliITSpListItem *pItem = fpList->GetpListItem( i, j );
599   if( pItem == 0 ) {
600     phys = 0.0;
601     for( Int_t l=0; l<size; l++ ) {
602       tracks[l]  = 0;
603       hits[l]    = 0;
604       charges[l] = 0.0;
605     }
606   } else {
607     Int_t idtrack =  pItem->GetTrack( 0 );
608     if( idtrack >= 0 ) phys = pItem->GetSignal();  
609     else phys = 0.0;
610
611     for( Int_t l=0; l<size; l++ ) if(l<pItem->GetMaxKept()) {
612       tracks[l]  = pItem->GetTrack( l );
613       hits[l]    = pItem->GetHit( l );
614       charges[l] = pItem->GetSignal( l );
615     }else{
616       tracks[l]  = -3;
617       hits[l]    = -1;
618       charges[l] = 0.0;
619     }// end for if
620   }
621
622   fITS->AddSimDigit( 1, phys, digits, tracks, hits, charges, signale ); 
623   delete [] tracks;
624   delete [] hits;
625   delete [] charges;
626 }
627 //______________________________________________________________________
628 void AliITSsimulationSDD::ChargeToSignal(Int_t mod,Bool_t bAddNoise, Bool_t bAddGain) {
629   // add baseline, noise, gain, electronics and ADC saturation effects
630   // apply dead channels
631
632   AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(mod);
633   Double_t baseline=0; 
634   Double_t noise=0; 
635   Double_t gain=0; 
636   Float_t contrib=0;
637   Int_t i,k,kk;
638   AliITSSimuParam* simpar = fDetType->GetSimuParam();
639   Float_t maxadc = simpar->GetSDDMaxAdc();    
640   Int_t nGroup=fScaleSize;
641   if(res->IsAMAt20MHz()){
642     nGroup=fScaleSize/2;
643   }
644
645   for (i=0;i<fNofMaps;i++) {
646     if( !fAnodeFire[i] ) continue;
647     baseline = res->GetBaseline(i);
648     noise = res->GetNoise(i);
649     gain = res->GetChannelGain(i)/fDetType->GetAverageGainSDD();
650     if(res->IsBad()) gain=0.;
651     if( res->IsChipBad(res->GetChip(i)) )gain=0.;
652     for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
653       fInZR[k]  = fHitMap2->GetSignal(i,k);
654       if(bAddGain) fInZR[k]*=gain;
655       if( bAddNoise ) {
656         contrib   = (baseline + noise*gRandom->Gaus());
657         fInZR[k] += contrib;
658       }
659       fInZI[k]  = 0.;
660     } // end for k
661     if(!fDoFFT) {      
662       for(k=0; k<fMaxNofSamples; k++) {
663         Double_t newcont = 0.;
664         Double_t maxcont = 0.;
665         for(kk=0;kk<fScaleSize;kk++) {
666           newcont = fInZR[fScaleSize*k+kk];
667           if(newcont > maxcont) maxcont = newcont;
668         } // end for kk
669         newcont = maxcont;
670         if (newcont >= maxadc) newcont = maxadc -1;
671         if(newcont >= baseline){
672           Warning("","newcont=%f>=baseline=%f",newcont,baseline);
673         } // end if
674           // back to analog: ?
675         fHitMap2->SetHit(i,k,newcont);
676       }  // end for k
677     }else{
678       FastFourierTransform(&fInZR[0],&fInZI[0],1);
679       for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
680         Double_t rw = fElectronics->GetTraFunReal(k);
681         Double_t iw = fElectronics->GetTraFunImag(k);
682         fOutZR[k]   = fInZR[k]*rw - fInZI[k]*iw;
683         fOutZI[k]   = fInZR[k]*iw + fInZI[k]*rw;
684       } // end for k
685       FastFourierTransform(&fOutZR[0],&fOutZI[0],-1);
686       for(k=0; k<fMaxNofSamples; k++) {
687         Double_t newcont1 = 0.;
688         Double_t maxcont1 = 0.;
689         for(kk=0;kk<nGroup;kk++) {
690           newcont1 = fOutZR[fScaleSize*k+kk];
691           if(newcont1 > maxcont1) maxcont1 = newcont1;
692         } // end for kk
693         newcont1 = maxcont1;
694         if (newcont1 >= maxadc) newcont1 = maxadc -1;
695         fHitMap2->SetHit(i,k,newcont1);
696       } // end for k
697     }
698   } // end for i loop over anodes
699   return;
700 }
701
702 //______________________________________________________________________
703 void AliITSsimulationSDD::ApplyCrosstalk(Int_t mod) {
704     // function add the crosstalk effect to signal
705     // temporal function, should be checked...!!!
706   
707     // create and inizialice crosstalk map
708     Float_t* ctk = new Float_t[fNofMaps*fMaxNofSamples+1];
709     memset( ctk, 0, sizeof(Float_t)*(fNofMaps*fMaxNofSamples+1) );
710     AliITSCalibrationSDD* calibr = (AliITSCalibrationSDD*)GetCalibrationModel(mod);
711     for( Int_t z=0; z<fNofMaps; z++ ) {
712       Double_t baseline = calibr->GetBaseline(z);
713         Bool_t on = kFALSE;
714         Int_t tstart = 0;
715         Int_t tstop = 0;
716         Int_t nTsteps = 0;
717         
718         for( Int_t l=0; l<fMaxNofSamples; l++ ) {
719             Float_t fadc = (Float_t)fHitMap2->GetSignal( z, l );
720             if( fadc > baseline ) {
721                 if( on == kFALSE && l<fMaxNofSamples-4 ) {
722                     Float_t fadc1 = (Float_t)fHitMap2->GetSignal( z, l+1 );
723                     if( fadc1 < fadc ) continue;
724                     on = kTRUE;
725                     nTsteps = 0;
726                     tstart = l;
727                 }
728                 nTsteps++;
729             }
730             else { // end fadc > baseline
731                 if( on == kTRUE ) {        
732                     if( nTsteps > 2 ) {
733                         tstop = l;
734                         // make smooth derivative
735                         Float_t* dev = new Float_t[fMaxNofSamples+1];
736                         memset( dev, 0, sizeof(Float_t)*(fMaxNofSamples+1) );
737                         for( Int_t i=tstart; i<tstop; i++ ) {   
738                             if( i > 2 && i < fMaxNofSamples-2 )
739                                 dev[i] = -0.2*fHitMap2->GetSignal( z,i-2 ) 
740                                     -0.1*fHitMap2->GetSignal( z,i-1 ) 
741                                     +0.1*fHitMap2->GetSignal( z,i+1 ) 
742                                     +0.2*fHitMap2->GetSignal( z,i+2 );
743                         }
744                         
745                         // add crosstalk contribution to neibourg anodes  
746                         for( Int_t i=tstart; i<tstop; i++ ) {
747                             Int_t anode = z - 1;
748                             Int_t i1 = (Int_t)((i-tstart)*.61+tstart+0.5); // 
749                             Float_t ctktmp =  -dev[i1] * 0.25;
750                             if( anode > 0 ) {
751                                 ctk[anode*fMaxNofSamples+i] += ctktmp;
752                             }
753                             anode = z + 1;
754                             if( anode < fNofMaps ) {
755                                 ctk[anode*fMaxNofSamples+i] += ctktmp;
756                             }
757                         }
758                         delete [] dev;
759                         
760                     } // if( nTsteps > 2 )
761                     on = kFALSE;
762                 }  // if( on == kTRUE )
763             }  // else
764         }
765     }
766     
767     for( Int_t a=0; a<fNofMaps; a++ )
768         for( Int_t t=0; t<fMaxNofSamples; t++ ) {     
769             Float_t signal = fHitMap2->GetSignal(a,t)+ctk[a*fMaxNofSamples+t];
770             fHitMap2->SetHit( a, t, signal );
771         }
772
773     delete [] ctk;
774 }
775
776 //______________________________________________________________________
777 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
778     // To the 10 to 8 bit lossive compression.
779     // code from Davide C. and Albert W.
780
781     if (signal < 128)  return signal;
782     if (signal < 256)  return (128+((signal-128)>>1));
783     if (signal < 512)  return (192+((signal-256)>>3));
784     if (signal < 1024) return (224+((signal-512)>>4));
785     return 0;
786 }
787 //______________________________________________________________________
788 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
789   // Decompression from 8 to 10 bit
790
791   if (signal < 0 || signal > 255) {
792     AliWarning(Form("Signal value %d out of range",signal));
793     return 0;
794   } // end if signal <0 || signal >255
795
796   if (signal < 128) return signal;
797   if (signal < 192) {
798     if (TMath::Odd(signal)) return (128+((signal-128)<<1));
799     else  return (128+((signal-128)<<1)+1);
800   } // end if signal < 192
801   if (signal < 224) {
802     if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
803     else  return (256+((signal-192)<<3)+4);
804   } // end if signal < 224
805   if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
806   return (512+((signal-224)<<4)+8);
807 }
808 //______________________________________________________________________
809 void AliITSsimulationSDD::Compress2D(){
810   // 2D zero-suppression algorithm as described in ALICE-INT-1999-28 V10
811   AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(fModule);  
812   for (Int_t iWing=0; iWing<2; iWing++) {
813     Int_t tL=res->GetZSLowThreshold(iWing);
814     Int_t tH=res->GetZSHighThreshold(iWing);
815     for (Int_t i=0; i<fNofMaps/2; i++) {  
816       Int_t ian=i+iWing*fNofMaps/2;
817       if( !fAnodeFire[ian] ) continue;
818       for (Int_t itb=0; itb<fMaxNofSamples; itb++) {
819         Int_t nLow=0, nHigh=0;      
820         Float_t cC=fHitMap2->GetSignal(ian,itb);
821         if(cC<=tL) continue;
822         nLow++; // cC is greater than tL
823         if(cC>tH) nHigh++;
824         //                     N
825         // Get "quintuple":   WCE
826         //                     S
827         Float_t wW=0.;
828         if(itb>0) wW=fHitMap2->GetSignal(ian,itb-1);
829         if(wW>tL) nLow++;
830         if(wW>tH) nHigh++;
831         Float_t eE=0.;
832         if(itb<fMaxNofSamples-1) eE=fHitMap2->GetSignal(ian,itb+1);
833         if(eE>tL) nLow++;
834         if(eE>tH) nHigh++;
835         Float_t nN=0.;
836         if(i<(fNofMaps/2-1)) nN=fHitMap2->GetSignal(ian+1,itb);
837         if(nN>tL) nLow++;
838         if(nN>tH) nHigh++;
839         Float_t sS=0.;
840         if(i>0) sS=fHitMap2->GetSignal(ian-1,itb);
841         if(sS>tL) nLow++;
842         if(sS>tH) nHigh++;
843         
844         if(nLow>=2 && nHigh>=1){
845           Int_t signal=(Int_t)cC;
846           Int_t signalc = Convert10to8(signal);
847           Int_t signale = Convert8to10(signalc);
848           signalc-=tL; // subtract low threshold after 10 to 8 bit compression
849           if(signalc>=4) AddDigit(ian,itb,signalc,signale);  // store C 
850         }
851       }
852     }
853   }
854 }
855
856
857 //______________________________________________________________________
858 void AliITSsimulationSDD::StoreAllDigits(){
859   // store digits for non-zero-suppressed data
860   for (Int_t ian=0; ian<fNofMaps; ian++) {
861     for (Int_t itb=0; itb<fMaxNofSamples; itb++){
862       Int_t signal=(Int_t)(fHitMap2->GetSignal(ian,itb));
863       Int_t signalc = Convert10to8(signal);
864       Int_t signale = Convert8to10(signalc);
865       AddDigit(ian,itb,signalc,signale);  
866     } 
867   }
868
869 //______________________________________________________________________
870 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
871   // Creates histograms of maps for debugging
872   Int_t i;
873   
874   fHis=new TObjArray(fNofMaps);
875   for (i=0;i<fNofMaps;i++) {
876     TString sddName;
877     sddName.Form("sdd_%d",i+1);
878     fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
879                          0.,(Float_t) scale*fMaxNofSamples), i);
880   } // end for i
881 }
882 //______________________________________________________________________
883 void AliITSsimulationSDD::FillHistograms(){
884     // fill 1D histograms from map
885
886     if (!fHis) return;
887
888     for( Int_t i=0; i<fNofMaps; i++) {
889         TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
890         Int_t nsamples = hist->GetNbinsX();
891         for( Int_t j=0; j<nsamples; j++) {
892             Double_t signal=fHitMap2->GetSignal(i,j);
893             hist->Fill((Float_t)j,signal);
894         } // end for j
895     } // end for i
896 }
897 //______________________________________________________________________
898 void AliITSsimulationSDD::ResetHistograms(){
899     // Reset histograms for this detector
900     Int_t i;
901
902     for (i=0;i<fNofMaps;i++ ) {
903         if (fHis->At(i))    ((TH1F*)fHis->At(i))->Reset();
904     } // end for i
905 }
906 //______________________________________________________________________
907 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) { 
908     // Fills a histogram from a give anode.  
909
910     if (!fHis) return 0;
911
912     if(wing <=0 || wing > 2) {
913         Warning("GetAnode","Wrong wing number: %d",wing);
914         return NULL;
915     } // end if wing <=0 || wing >2
916     if(anode <=0 || anode > fNofMaps/2) {
917         Warning("GetAnode","Wrong anode number: %d",anode);
918         return NULL;
919     } // end if ampde <=0 || andoe > fNofMaps/2
920
921     Int_t index = (wing-1)*fNofMaps/2 + anode-1;
922     return (TH1F*)(fHis->At(index));
923 }
924 //______________________________________________________________________
925 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
926     // Writes the histograms to a file
927
928     if (!fHis) return;
929
930     hfile->cd();
931     Int_t i;
932     for(i=0; i<fNofMaps; i++)  fHis->At(i)->Write(); //fAdcs[i]->Write();
933     return;
934 }
935 //______________________________________________________________________
936 void AliITSsimulationSDD::WriteSDigits(){
937     // Fills the Summable digits Tree
938     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
939
940     for( Int_t i=0; i<fNofMaps; i++ ) {
941         if( !fAnodeFire[i] ) continue;
942         for( Int_t j=0; j<fMaxNofSamples; j++ ) {
943             Double_t sig = fHitMap2->GetSignal( i, j );
944             if( sig > 0.2 ) {
945                 Int_t jdx = j*fScaleSize;
946                 Int_t index = fpList->GetHitIndex( i, j );
947                 AliITSpListItem pItemTmp2( fModule, index, 0. );
948                 // put the fScaleSize analog digits in only one
949                 for( Int_t ik=0; ik<fScaleSize; ik++ ) {
950                     AliITSpListItem *pItemTmp = fpList->GetpListItem(i,jdx+ik);
951                     if( pItemTmp == 0 ) continue;
952                     pItemTmp2.Add( pItemTmp );
953                 }
954                 pItemTmp2.AddSignalAfterElect( fModule, index, sig );
955                 pItemTmp2.AddNoise(fModule,index,fHitNoiMap2->GetSignal(i,j));
956                 aliITS->AddSumDigit( pItemTmp2 );
957             } // end if (sig > 0.2)
958         }
959     }
960     return;
961 }
962 //______________________________________________________________________
963 void AliITSsimulationSDD::PrintStatus() const {
964     // Print SDD simulation Parameters
965
966     cout << "**************************************************" << endl;
967     cout << "   Silicon Drift Detector Simulation Parameters   " << endl;
968     cout << "**************************************************" << endl;
969     cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
970     cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
971     cout << "Number of Anodes used: " << fNofMaps << endl;
972     cout << "Number of Time Samples: " << fMaxNofSamples << endl;
973     cout << "Scale size factor: " << fScaleSize << endl;
974     cout << "**************************************************" << endl;
975 }