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