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