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