]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSsimulationSPDdubna.cxx
370c96e14d4334be0355d248b267ca776c099793
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationSPDdubna.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 $Log$
18 Revision 1.6  2002/08/21 22:09:58  nilsen
19 Updated SPD simulation with difusion effects. ReWritten Hit to SDigits
20 code.
21
22 */
23 #include <iostream.h>
24 #include <TRandom.h>
25 #include <TH1.h>
26 #include <TMath.h>
27 #include <TString.h>
28 #include <TParticle.h>
29
30 #include "AliRun.h"
31 #include "AliITS.h"
32 #include "AliITShit.h"
33 #include "AliITSdigit.h"
34 #include "AliITSmodule.h"
35 #include "AliITSMapA2.h" 
36 #include "AliITSpList.h"
37 #include "AliITSsimulationSPDdubna.h"
38 #include "AliITSsegmentationSPD.h"
39 #include "AliITSresponseSPDdubna.h"
40
41 //#define DEBUG
42
43 ClassImp(AliITSsimulationSPDdubna)
44 ////////////////////////////////////////////////////////////////////////
45 // Version: 0
46 // Written by Boris Batyunya
47 // December 20 1999
48 //
49 // AliITSsimulationSPDdubna is the simulation of SPDs
50 //______________________________________________________________________
51
52
53 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(){
54     // constructor
55
56     fResponse = 0;
57     fSegmentation = 0;
58     fMapA2 = 0;
59     fpList = 0;
60     fModule = 0;
61     fEvent = 0;
62     fHis = 0;
63     fNoise = 0.;
64     fBaseline = 0.;
65     fNPixelsZ = 0;
66     fNPixelsX = 0;
67 }
68 //______________________________________________________________________
69 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(AliITSsegmentation *seg,
70                                                    AliITSresponse *resp){
71     // standard constructor
72     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
73
74     fHis = 0;
75     fResponse = resp;
76     fSegmentation = seg;
77     fModule = 0;
78     fEvent = 0;
79
80     fNPixelsZ=fSegmentation->Npz();
81     fNPixelsX=fSegmentation->Npx();
82
83     fResponse->GetNoiseParam(fNoise,fBaseline);
84     fResponse->SetDistanceOverVoltage(kmictocm*fSegmentation->Dy(),50.0);
85
86 //    fMapA2 = new AliITSMapA2(fSegmentation);
87     fMapA2 = 0;
88
89     fpList = new AliITSpList(fNPixelsZ+1,fNPixelsX+1);
90
91 }
92 //______________________________________________________________________
93 AliITSsimulationSPDdubna::~AliITSsimulationSPDdubna(){
94     // destructor
95
96     if(fMapA2) delete fMapA2;
97
98     if (fHis) {
99         fHis->Delete(); 
100         delete fHis;     
101     } // end if fHis
102 }
103 //______________________________________________________________________
104 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(const 
105                                                    AliITSsimulationSPDdubna 
106                                                    &source){
107     //     Copy Constructor 
108     if(&source == this) return;
109     this->fMapA2 = source.fMapA2;
110     this->fNoise = source.fNoise;
111     this->fBaseline = source.fBaseline;
112     this->fNPixelsX = source.fNPixelsX;
113     this->fNPixelsZ = source.fNPixelsZ;
114     this->fHis = source.fHis;
115     return;
116 }
117 //______________________________________________________________________
118 AliITSsimulationSPDdubna&  AliITSsimulationSPDdubna::operator=(const 
119                                            AliITSsimulationSPDdubna &source){
120     //    Assignment operator
121     if(&source == this) return *this;
122     this->fMapA2 = source.fMapA2;
123     this->fNoise = source.fNoise;
124     this->fBaseline = source.fBaseline;
125     this->fNPixelsX = source.fNPixelsX;
126     this->fNPixelsZ = source.fNPixelsZ;
127     this->fHis = source.fHis;
128     return *this;
129 }
130 //______________________________________________________________________
131 void AliITSsimulationSPDdubna::InitSimulationModule(Int_t module, Int_t event){
132     //  This function creates maps to build the list of tracks for each
133     //  summable digit.
134     //
135     //  Inputs:
136     //    Int_t module   // Module number to be simulated
137     //    Int_t event    // Event number to be simulated
138     //
139     //  Outputs:
140     //    none
141     //
142     //  Returns:
143     //    none
144
145     fModule = module;
146     fEvent  = event;
147 //    fMapA2->ClearMap();
148     fpList->ClearMap();
149 }
150 //_____________________________________________________________________
151 void AliITSsimulationSPDdubna::SDigitiseModule(AliITSmodule *mod, Int_t mask,
152                                                Int_t event){
153     //  This function begins the work of creating S-Digits
154     //
155     //  Inputs:
156     //    AliITSmodule *mod  //  module
157     //    Int_t mask         //  mask to be applied to the module
158     //
159     //  Outputs:
160     //    none
161     //
162     //  Return:
163     //    test              //  test returns kTRUE if the module contained hits
164     //                      //  test returns kFALSE if it did not contain hits
165
166     Int_t module = 0;
167
168     if(!(mod->GetNhits())) return;// if module has no hits don't create Sdigits
169     fModule = mod->GetIndex();
170     HitToSDigit(mod, module, mask, fpList);
171     WriteSDigits(fpList);
172 //    fMapA2->ClearMap();
173     fpList->ClearMap();
174 }
175 //______________________________________________________________________
176 void AliITSsimulationSPDdubna::WriteSDigits(AliITSpList *pList){
177     //  This function adds each S-Digit to pList
178     //
179     //  Inputs:
180     //    AliITSpList *pList
181     //
182     //  Outputs:
183     //    none
184     //
185     //  Return:
186     //    none
187     Int_t ix, nix, iz, niz;
188     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
189
190     pList->GetMaxMapIndex(niz, nix);
191     for(iz=0; iz<niz; iz++)for(ix=0; ix<nix; ix++){
192         if(pList->GetSignalOnly(iz,ix)>0.0){
193             aliITS->AddSumDigit(*(pList->GetpListItem(iz,ix)));
194 #ifdef DEBUG
195             cout <<"SDigits " << iz << "," << ix << "," << 
196                 *(pList->GetpListItem(iz,ix)) << endl;
197 #endif
198         } // end if pList
199     } // end for iz,ix
200     return; 
201 }
202 //______________________________________________________________________
203 void AliITSsimulationSPDdubna::FinishSDigitiseModule(){
204     //  This function calls SDigitsToDigits which creates Digits from SDigits
205     //
206     //  Inputs:
207     //    none
208     //
209     //  Outputs:
210     //    none
211     //  Return
212     //    none
213
214     SDigitsToDigits(fModule, fpList);
215     return;
216 }
217 //______________________________________________________________________
218 void AliITSsimulationSPDdubna::SDigitsToDigits(Int_t module,
219                                                AliITSpList *pList){
220     //  This function adds electronic noise to the S-Digits and then adds them
221     // to  a new pList
222     //
223     //  Inputs:
224     //    Int_t       module  // module number
225     //    AliITSpList *pList  // pList
226     //
227     //  Outputs:
228     //    pList is passed along to the functions ChargeToSignal and GetList
229     //
230     //  Return:
231     //    none
232
233     fModule = module;
234     ChargeToSignal(pList); // Charge To Signal both adds noise and
235 //    fMapA2->ClearMap();
236     pList->ClearMap();
237 }
238 //______________________________________________________________________
239 void AliITSsimulationSPDdubna::DigitiseModule(AliITSmodule *mod, Int_t module,
240                                               Int_t dummy){
241     //  This function creates Digits straight from the hits and then adds
242     //  electronic noise to the digits before adding them to pList
243     //
244     //  Inputs:
245     //    AliITSmodule *mod    // module
246     //    Int_t        module  // module number  Dummy.
247     //    Int_t        dummy
248     //
249     //  Outputs:
250     //    Each of the input variables is passed along to HitToSDigit
251     //
252     //  Return:
253     //    none
254
255     fModule = mod->GetIndex();  //This calls the module for HitToSDigit
256     HitToSDigit(mod,fModule, dummy, fpList);
257     ChargeToSignal(fpList);
258 //    fMapA2->ClearMap();
259     fpList->ClearMap();
260 }
261 //______________________________________________________________________
262 void AliITSsimulationSPDdubna::UpdateMapSignal(Int_t iz, Int_t ix, Int_t trk,
263                                                Int_t ht, Int_t module,
264                                                Double_t signal,
265                                                AliITSpList *pList){
266     //  This function adds a signal to the pList from the pList class
267     //
268     //  Inputs:
269     //    Int_t       iz     // row number
270     //    Int_t       ix     // column number
271     //    Int_t       trk    // track number
272     //    Int_t       ht     // hit number
273     //    Double_t    signal // signal strength
274     //    AliITSpList *pList // pList
275     //
276     //  Outputs:
277     //    All of the inputs are passed to AliITSpList::AddSignal
278     //    Int_t    ix  // row number
279     //    Int_t    iz  // column number
280     //    Double_t sig // signal strength
281     //          // These three variables are defined to preserve the
282     //          // assignments used in the function AliITSMapA2::AddSignal
283     //
284     //  Return:
285     //    none
286
287 //    fMapA2->AddSignal(iz, ix, signal);
288     pList->AddSignal(iz,ix, trk, ht, fModule, signal);
289 }
290 //______________________________________________________________________
291 void AliITSsimulationSPDdubna::UpdateMapNoise(Int_t iz,
292                                               Int_t ix, Int_t fModule,
293                                               Double_t sig, Float_t noise,
294                                               AliITSpList *pList){
295     //  This function adds noise to data in the MapA2 as well as the pList
296     //
297     //  Inputs:
298     //    Int_t       iz       // row number
299     //    Int_t       ix       // column number
300     //    Int_t       mod     // module number
301     //    Double_t    sig     // signal strength
302     //    Double_t    noise   // electronic noise generated by ChargeToSignal
303     //    AliITSpList *pList  // pList
304     //
305     //  Outputs:
306     //    All of the inputs are passed to AliITSMapA2::AddSignal or
307     //    AliITSpList::AddNoise
308     //
309     //  Return:
310     //    none
311
312 //    fMapA2->AddSignal(iz, ix, noise);
313     pList->AddNoise(iz,ix, fModule, noise);
314 }
315 //______________________________________________________________________
316 void AliITSsimulationSPDdubna::HitToDigit(AliITSmodule *mod, Int_t module,
317                                           Int_t dummy){
318     DigitiseModule(mod, module, dummy);
319 }
320 //______________________________________________________________________
321 void AliITSsimulationSPDdubna::HitToSDigit(AliITSmodule *mod, Int_t module,
322                                             Int_t dummy,AliITSpList *pList){
323     // Does the charge distributions using Gaussian diffusion charge charing.
324     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
325     TObjArray *hits = mod->GetHits();
326     Int_t nhits = hits->GetEntriesFast();
327     Int_t h,ix,iz;
328     Int_t idtrack;
329     Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
330     Double_t x,y,z,t,tp,st,dt=0.2,el,sig;
331     Double_t thick = kmictocm*GetSeg()->Dy();
332
333     if(nhits<=0) return;
334     for(h=0;h<nhits;h++){
335 #ifdef DEBUG
336         cout << "Hits=" << h << "," << *(mod->GetHit(h)) << endl;
337 #endif
338         if(mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)){
339         st =TMath::Sqrt(x1*x1+y1*y1+z1*z1);
340         if(st>0.0) for(t=0;t<1.0;t+=dt){ // Integrate over t
341             tp = t+0.5*dt;
342             el = GetResp()->GeVToCharge((Float_t)(dt*de));
343 #ifdef DEBUG
344             if(el<=0.0) cout << "el="<<el<<" dt="<<dt<<" de="<<de<<endl;
345 #endif
346             x = x0+x1*tp;
347             y = y0+y1*tp;
348             z = z0+z1*tp;
349             GetSeg()->LocalToDet(x,z,ix,iz);
350             sig = GetResp()->SigmaDiffusion1D(thick + y);
351             SpreadCharge(x,y,z,ix,iz,el,sig,
352                          idtrack,mod->GetHitTrackIndex(h),h,mod->GetIndex());
353         } else { // st == 0.0 deposit it at this point
354             el = GetResp()->GeVToCharge((Float_t)de);
355             x = x0;
356             y = y0;
357             z = z0;
358             GetSeg()->LocalToDet(x,z,ix,iz);
359             sig = GetResp()->SigmaDiffusion1D(thick + y);
360             SpreadCharge(x,y,z,ix,iz,el,sig,
361                          idtrack,mod->GetHitTrackIndex(h),h,mod->GetIndex());
362         } // end if st>0.0
363     }} // Loop over all hits h
364 }
365 //______________________________________________________________________
366 void AliITSsimulationSPDdubna::SpreadCharge(Double_t x0,Double_t y0,
367                                             Double_t z0,Int_t ix0,Int_t iz0,
368                                             Double_t el,Double_t sig,Int_t t,
369                                             Int_t ti,Int_t hi,Int_t mod){
370     // Spreads the charge over neighboring cells. Assume charge is distributed
371     // as charge(x,z) = (el/2*pi*sig*sig)*exp(-arg)
372     // arg=((x-x0)*(x-x0)/2*sig*sig)+((z-z0*z-z0)/2*sig*sig)
373     // Defined this way, the integral over all x and z is el.
374     const Int_t knx = 3,knz = 2;
375     const Double_t kRoot2 = 1.414213562; // Sqrt(2).
376     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
377     Int_t ix,iz,ixs,ixe,izs,ize;
378     Float_t x,z;
379     Double_t x1,x2,z1,z2,s,sp;
380
381     if(sig<=0.0) {
382         fpList->AddSignal(iz0,ix0,t,hi,mod,el);
383         return;
384     } // end if
385     sp = 1.0/(sig*kRoot2);
386 #ifdef DEBUG
387     cout << "sig=" << sig << " sp=" << sp << endl;
388 #endif
389     ixs = TMath::Max(-knx+ix0,0);
390     ixe = TMath::Min(knx+ix0,GetSeg()->Npx()-1);
391     izs = TMath::Max(-knz+iz0,0);
392     ize = TMath::Min(knz+iz0,GetSeg()->Npz()-1);
393     for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
394         GetSeg()->DetToLocal(ix,iz,x,z); // pixel center
395         x1 = x;
396         z1 = z;
397         x2  = x1 + 0.5*kmictocm*GetSeg()->Dpx(ix); // Upper
398         x1 -= 0.5*kmictocm*GetSeg()->Dpx(ix);  // Lower
399         z2  = z1 + 0.5*kmictocm*GetSeg()->Dpz(iz); // Upper
400         z1 -= 0.5*kmictocm*GetSeg()->Dpz(iz);  // Lower
401         x1 -= x0; // Distance from where track traveled
402         x2 -= x0; // Distance from where track traveled
403         z1 -= z0; // Distance from where track traveled
404         z2 -= z0; // Distance from where track traveled
405         s = 0.25; // Correction based on definision of Erfc
406         s *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2);
407 #ifdef DEBUG
408         cout << "el=" << el << " ix0=" << ix0 << " ix=" << ix << " x0="<< x <<
409             " iz0=" << iz0 << " iz=" << iz << " z0=" << z  << 
410             " sp*x1=" << sp*x1 <<" sp*x2=" << sp*x2 << " s=" << s;
411 #endif
412         s *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2);
413 #ifdef DEBUG
414         cout << " sp*z1=" << sp*z1 <<" sp*z2=" << sp*z2 << " s=" << s << endl;
415 #endif
416         fpList->AddSignal(iz,ix,t,hi,mod,s*el);
417     } // end for ix, iz
418 }
419 //______________________________________________________________________
420 void AliITSsimulationSPDdubna::HitToSDigitOld(AliITSmodule *mod, Int_t module,
421                                            Int_t dummy, AliITSpList *pList){
422     // digitize module 
423     const Float_t kEnToEl = 2.778e+8; // GeV->charge in electrons 
424                                       // for 3.6 eV/pair 
425     const Float_t kconv = 10000.;     // cm -> microns
426
427     Float_t spdLength = fSegmentation->Dz();
428     Float_t spdWidth = fSegmentation->Dx();
429     Float_t spdThickness = fSegmentation->Dy();
430     Float_t difCoef, dum;       
431     fResponse->DiffCoeff(difCoef,dum); 
432     if(spdThickness > 290) difCoef = 0.00613;  
433
434     Float_t zPix0 = 1e+6;
435     Float_t xPix0 = 1e+6;
436     Float_t yPrev = 1e+6;   
437
438     Float_t zPitch = fSegmentation->Dpz(0);
439     Float_t xPitch = fSegmentation->Dpx(0);
440   
441     TObjArray *fHits = mod->GetHits();
442     module = mod->GetIndex();
443     Int_t nhits = fHits->GetEntriesFast();
444     if (!nhits) return;
445 #ifdef DEBUG
446     cout<<"len,wid,thickness,nx,nz,pitchx,pitchz,difcoef ="<<spdLength<<","
447         <<spdWidth<<","<<spdThickness<<","<<fNPixelsX<<","<<fNPixelsZ<<","
448         <<xPitch<<","<<zPitch<<","<<difCoef<<endl;
449 #endif
450     //  Array of pointers to the label-signal list
451     Int_t indexRange[4] = {0,0,0,0};
452
453     // Fill detector maps with GEANT hits
454     // loop over hits in the module
455     static Bool_t first;
456     Int_t lasttrack=-2;
457     Int_t hit, iZi, jz, jx;
458     Int_t idhit=-1; //!
459 #ifdef DEBUG
460     cout<<"SPDdubna: module,nhits ="<<module<<","<<nhits<<endl;
461 #endif
462     for (hit=0;hit<nhits;hit++) {
463         AliITShit *iHit = (AliITShit*) fHits->At(hit);
464 #ifdef DEBUG
465         cout << "Hits=" << hit << "," << *iHit << endl;
466 #endif
467         //Int_t layer = iHit->GetLayer();
468         Float_t yPix0 = -spdThickness/2; 
469
470         // work with the idtrack=entry number in the TreeH
471         //Int_t idhit,idtrack; //!
472         //mod->GetHitTrackAndHitIndex(hit,idtrack,idhit);  //!    
473         //Int_t idtrack=mod->GetHitTrackIndex(hit);  
474         // or store straight away the particle position in the array
475         // of particles : 
476         if(iHit->StatusEntering()) idhit=hit;
477         Int_t itrack = iHit->GetTrack();
478         Int_t dray = 0;
479    
480         if (lasttrack != itrack || hit==(nhits-1)) first = kTRUE; 
481
482         //Int_t parent = iHit->GetParticle()->GetFirstMother();
483         Int_t partcode = iHit->GetParticle()->GetPdgCode();
484
485         //  partcode (pdgCode): 11 - e-, 13 - mu-, 22 - gamma, 111 - pi0,
486         // 211 - pi+,  310 - K0s, 321 - K+, 2112 - n, 2212 - p, 3122 - lambda
487
488         Float_t pmod = iHit->GetParticle()->P(); // total momentum at the
489                                                    // vertex
490         pmod *= 1000;
491
492         if(partcode == 11 && pmod < 6) dray = 1; // delta ray is e-
493                                                  // at p < 6 MeV/c
494
495         //  Get hit z and x(r*phi) cordinates for each module (detector)
496         //  in local system.
497
498         Float_t zPix = kconv*iHit->GetZL();
499         Float_t xPix = kconv*iHit->GetXL();
500         Float_t yPix = kconv*iHit->GetYL();
501
502         // Get track status
503         Int_t status = iHit->GetTrackStatus();      
504
505         // Check boundaries
506         if(zPix  > spdLength/2) {
507 #ifdef DEBUG
508             cout<<"!!! SPD: z outside ="<<zPix<<endl;
509 #endif
510             zPix = spdLength/2 - 10;
511         }
512         if(zPix  < 0 && zPix < -spdLength/2) {
513 #ifdef DEBUG
514             cout<<"!!! SPD: z outside ="<<zPix<<endl;
515 #endif
516             zPix = -spdLength/2 + 10;
517         }
518         if(xPix  > spdWidth/2) {
519 #ifdef DEBUG
520             cout<<"!!! SPD: x outside ="<<xPix<<endl;
521 #endif
522             xPix = spdWidth/2 - 10;
523         }
524         if(xPix  < 0 && xPix < -spdWidth/2) {
525 #ifdef DEBUG
526             cout<<"!!! SPD: x outside ="<<xPix<<endl;
527 #endif
528             xPix = -spdWidth/2 + 10;
529         }
530         Int_t trdown = 0;
531
532         // enter Si or after event in Si
533         if (status == 66 ) {  
534             zPix0 = zPix;
535             xPix0 = xPix;
536             yPrev = yPix; 
537         } // end if status == 66
538
539         Float_t depEnergy = iHit->GetIonization();
540         // skip if the input point to Si       
541
542         if(depEnergy <= 0.) continue;        
543
544         // if track returns to the opposite direction:
545         if (yPix < yPrev) {
546             trdown = 1;
547         } // end if yPix < yPrev
548
549         // take into account the holes diffusion inside the Silicon
550         // the straight line between the entrance and exit points in Si is
551         // divided into the several steps; the diffusion is considered 
552         // for each end point of step and charge
553         // is distributed between the pixels through the diffusion.
554
555         //  ---------- the diffusion in Z (beam) direction -------
556         Float_t charge = depEnergy*kEnToEl;         // charge in e-
557         Float_t drPath = 0.;   
558         Float_t tang = 0.;
559         Float_t sigmaDif = 0.; 
560         Float_t zdif = zPix - zPix0;
561         Float_t xdif = xPix - xPix0;
562         Float_t ydif = TMath::Abs(yPix - yPrev);
563         Float_t ydif0 = TMath::Abs(yPrev - yPix0);
564
565         if(ydif < 1) continue; // ydif is not zero
566
567         Float_t projDif = sqrt(xdif*xdif + zdif*zdif);
568
569         Int_t ndZ = (Int_t)TMath::Abs(zdif/zPitch) + 1;
570         Int_t ndX = (Int_t)TMath::Abs(xdif/xPitch) + 1; 
571
572         // number of the steps along the track:
573         Int_t nsteps = ndZ;
574         if(ndX > ndZ) nsteps = ndX;
575         if(nsteps < 20) nsteps = 20;  // minimum number of the steps 
576
577         if (projDif < 5 ) {
578             drPath = (yPix-yPix0)*1.e-4;  
579             drPath = TMath::Abs(drPath);        // drift path in cm
580             sigmaDif = difCoef*sqrt(drPath);    // sigma diffusion in cm
581             sigmaDif = sigmaDif*kconv;         // sigma diffusion in microns
582             nsteps = 1;
583         }  // end if projDif < 5
584
585         if(projDif > 5) tang = ydif/projDif;
586         Float_t dCharge = charge/nsteps;       // charge in e- for one step
587         Float_t dZ = zdif/nsteps;
588         Float_t dX = xdif/nsteps;
589
590         for (iZi = 1; iZi <= nsteps;iZi++) {
591             Float_t dZn = iZi*dZ;
592             Float_t dXn = iZi*dX;
593             Float_t zPixn = zPix0 + dZn;
594             Float_t xPixn = xPix0 + dXn;
595
596             if(projDif >= 5) {
597                 Float_t dProjn = sqrt(dZn*dZn+dXn*dXn);
598                 drPath = dProjn*tang*1.e-4; // drift path for iZi+1 step in cm 
599                 if(trdown == 0) {
600                     drPath = TMath::Abs(drPath) + ydif0*1.e-4;
601                 }// end if trdow ==0
602                 if(trdown == 1) {
603                     drPath = ydif0*1.e-4 - TMath::Abs(drPath);
604                     drPath = TMath::Abs(drPath);
605                 } // end if trdown == 1
606                 sigmaDif = difCoef*sqrt(drPath);    
607                 sigmaDif = sigmaDif*kconv;       // sigma diffusion in microns
608             } // end if projdif >= 5
609
610             zPixn = (zPixn + spdLength/2.);  
611             xPixn = (xPixn + spdWidth/2.);  
612             Int_t nZpix, nXpix;
613             fSegmentation->GetPadIxz(xPixn,zPixn,nXpix,nZpix);
614             zPitch = fSegmentation->Dpz(nZpix);
615             fSegmentation->GetPadTxz(xPixn,zPixn);
616             // set the window for the integration
617             Int_t jzmin = 1;  
618             Int_t jzmax = 3; 
619             if(nZpix == 1) jzmin =2;
620             if(nZpix == fNPixelsZ) jzmax = 2; 
621
622             Int_t jxmin = 1;  
623             Int_t jxmax = 3; 
624             if(nXpix == 1) jxmin =2;
625             if(nXpix == fNPixelsX) jxmax = 2; 
626
627             Float_t zpix = nZpix; 
628             Float_t dZright = zPitch*(zpix - zPixn);
629             Float_t dZleft = zPitch - dZright;
630
631             Float_t xpix = nXpix; 
632             Float_t dXright = xPitch*(xpix - xPixn);
633             Float_t dXleft = xPitch - dXright;
634
635             Float_t dZprev = 0.;
636             Float_t dZnext = 0.;
637             Float_t dXprev = 0.;
638             Float_t dXnext = 0.;
639
640             for(jz=jzmin; jz <=jzmax; jz++) {
641                 if(jz == 1) {
642                     dZprev = -zPitch - dZleft;
643                     dZnext = -dZleft;
644                 } else if(jz == 2) {
645                     dZprev = -dZleft;
646                     dZnext = dZright;
647                 } else if(jz == 3) {
648                     dZprev = dZright;
649                     dZnext = dZright + zPitch;
650                 } // end if jz
651                 // kz changes from 1 to the fNofPixels(270)  
652                 Int_t kz = nZpix + jz -2; 
653
654                 Float_t zArg1 = dZprev/sigmaDif;
655                 Float_t zArg2 = dZnext/sigmaDif;
656                 Float_t zProb1 = TMath::Erfc(zArg1);
657                 Float_t zProb2 = TMath::Erfc(zArg2);
658                 Float_t dZCharge =0.5*(zProb1-zProb2)*dCharge; 
659
660
661                 // ----------- holes diffusion in X(r*phi) direction  --------
662
663                 if(dZCharge > 1.) { 
664                     for(jx=jxmin; jx <=jxmax; jx++) {
665                         if(jx == 1) {
666                             dXprev = -xPitch - dXleft;
667                             dXnext = -dXleft;
668                         } else if(jx == 2) {
669                             dXprev = -dXleft;
670                             dXnext = dXright;
671                         } else if(jx == 3) {
672                             dXprev = dXright;
673                             dXnext = dXright + xPitch;
674                         }  // end if jx
675                         Int_t kx = nXpix + jx -2;
676                         Float_t xArg1 = dXprev/sigmaDif;
677                         Float_t xArg2 = dXnext/sigmaDif;
678                         Float_t xProb1 = TMath::Erfc(xArg1);
679                         Float_t xProb2 = TMath::Erfc(xArg2);
680                         Float_t dXCharge =0.5*(xProb1-xProb2)*dZCharge;
681
682                         if(dXCharge > 1.) {
683                             if (first) {
684                                 indexRange[0]=indexRange[1]=kz-1;
685                                 indexRange[2]=indexRange[3]=kx-1;
686                                 first=kFALSE;
687                             } // end if first
688                             indexRange[0]=TMath::Min(indexRange[0],kz-1);
689                             indexRange[1]=TMath::Max(indexRange[1],kz-1);
690                             indexRange[2]=TMath::Min(indexRange[2],kx-1);
691                             indexRange[3]=TMath::Max(indexRange[3],kx-1);
692 /*
693                             // build the list of digits for this module 
694                             Double_t signal = fMapA2->GetSignal(kz-1,kx-1);
695                             signal+=dXCharge;
696                             fMapA2->SetHit(kz-1,kx-1,(double)signal);
697 */
698                             // The calling sequence for UpdateMapSignal was 
699                             // moved into the (dx > 1 e-) loop because it 
700                             // needs to call signal which is defined inside 
701                             // this loop
702                             fModule   = module;//Defined because functions 
703                                                // called by UpdateMapSignal 
704                                                // expect module to be an 
705                                                // integer
706                             UpdateMapSignal(kz-1,kx-1,
707 //                                          mod->GetHitTrackIndex(hit),
708                              ((AliITShit*)(mod->GetHit(hit)))->GetTrack(),
709                                             hit,fModule,dXCharge,pList);
710                         }      // dXCharge > 1 e-
711                     }       // jx loop
712                 }       // dZCharge > 1 e-
713             }        // jz loop
714         }         // iZi loop
715         if (status == 65) {   // the step is inside of Si
716             zPix0 = zPix;
717             xPix0 = xPix;
718         } // end if status == 65
719         yPrev = yPix;
720     }   // hit loop inside the module
721 }
722 //______________________________________________________________________
723 void AliITSsimulationSPDdubna::ChargeToSignal(AliITSpList *pList){
724     // add noise and electronics, perform the zero suppression and add the
725     // digit to the list
726     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
727     Float_t threshold = (float)fResponse->MinVal();
728     Int_t j;
729 //    Int_t    digits[3], tracks[3], hits[3];
730 //    Float_t  charges[3];
731     Float_t  electronics;
732 //    Float_t  phys; 
733     Double_t sig;
734     const Int_t    nmaxtrk=AliITSdigitSPD::GetNTracks();
735     static AliITSdigitSPD dig;
736
737     for(Int_t iz=0; iz<fNPixelsZ; iz++){
738         for(Int_t ix=0; ix<fNPixelsX; ix++){
739             electronics = fBaseline + fNoise*gRandom->Gaus();
740             sig = pList->GetSignalOnly(iz,ix);
741             UpdateMapNoise(iz,ix,fModule,sig,electronics,pList);
742 #ifdef DEBUG
743 //          cout << sig << "+" << electronics <<">threshold=" << threshold 
744 //               << endl;
745 #endif
746             if (sig+electronics > threshold) {
747                 dig.fCoord1 = iz;
748                 dig.fCoord2 = ix;
749                 dig.fSignal = 1;
750                 dig.fSignalSPD = (Int_t) pList->GetSignal(iz,ix);
751                 /*
752                 digits[0] = iz;
753                 digits[1] = ix;
754                 digits[2] = 1; */
755                 for(j=0;j<nmaxtrk;j++){
756 //                  charges[j] = 0.0;
757                     if (j<pList->GetNEnteries()) {
758                         dig.fTracks[j] = pList->GetTrack(iz,ix,j);
759                         dig.fHits[j]   = pList->GetHit(iz,ix,j);
760                         /*
761                         tracks[j] = pList->GetTrack(iz,ix,j);
762                         hits[j]   = pList->GetHit(iz,ix,j);
763                         */
764                     }else { // Default values
765                         dig.fTracks[j] = -3;
766                         dig.fHits[j]   = -1;
767 /*                      tracks[j] = -2; //noise
768                         hits[j]   = -1;  */
769                     } // end if pList
770                 } // end for j
771 //              charges[0] = (Float_t) pList->GetSumSignal(iz,ix);
772 /*
773                 if(tracks[0] == tracks[1] && tracks[0] == tracks[2]) {
774                     tracks[1] = -3;
775                     hits[1] = -1;
776                     tracks[2] = -3;
777                     hits[2] = -1;
778                 } else if(tracks[0] == tracks[1] && tracks[0] != tracks[2]) {
779                     tracks[1] = -3;
780                     hits[1] = -1;   
781                 } else if(tracks[0] == tracks[2] && tracks[0] != tracks[1]) {
782                     tracks[2] = -3;
783                     hits[2] = -1;   
784                 } else if(tracks[1] == tracks[2] && tracks[0] != tracks[1]) {
785                     tracks[2] = -3;
786                     hits[2] = -1;   
787                 } // end if
788 */
789 //              phys = 0.0;
790 #ifdef DEBUG
791                 cout << iz << "," << ix << "," << 
792                     *(pList->GetpListItem(iz,ix)) << endl;
793 #endif
794 //              aliITS->AddSimDigit(0, phys, digits, tracks, hits, charges);
795                 aliITS->AddSimDigit(0,&dig);
796             } // 
797         } // 
798     } //
799 }
800 //______________________________________________________________________
801 void AliITSsimulationSPDdubna::CreateHistograms(){
802     // create 1D histograms for tests
803
804     printf("SPD - create histograms\n");
805
806     fHis=new TObjArray(fNPixelsZ);
807     TString spdName("spd_");
808     for (Int_t i=0;i<fNPixelsZ;i++) {
809         Char_t pixelz[4];
810         sprintf(pixelz,"%d",i+1);
811         spdName.Append(pixelz);
812         //PH       (*fHis)[i] = new TH1F(spdName.Data(),"SPD maps",
813         //PH                       fNPixelsX,0.,(Float_t) fNPixelsX);
814         fHis->AddAt(new TH1F(spdName.Data(),"SPD maps",
815                              fNPixelsX,0.,(Float_t) fNPixelsX), i);
816     } // end for i
817 }
818 //______________________________________________________________________
819 void AliITSsimulationSPDdubna::ResetHistograms(){
820     //
821     // Reset histograms for this detector
822     //
823
824     for ( int i=0;i<fNPixelsZ;i++ ) {
825         //PH    if ((*fHis)[i])    ((TH1F*)(*fHis)[i])->Reset();
826         if (fHis->At(i))    ((TH1F*)fHis->At(i))->Reset();
827     } // end for i
828 }