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