]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSsimulationSPD.cxx
New SPD simulation (Massimo):
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationSPD.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 $Id$
18 */
19
20 #include <Riostream.h>
21 #include <TH1.h>
22 #include <TMath.h>
23 #include <TParticle.h>
24 #include <TRandom.h>
25 #include <TString.h>
26 #include "AliITS.h"
27 #include "AliITSMapA2.h" 
28 #include "AliITSdigitSPD.h"
29 #include "AliITSgeom.h"
30 #include "AliITShit.h"
31 #include "AliITSmodule.h"
32 #include "AliITSpList.h"
33 #include "AliITSCalibrationSPD.h"
34 #include "AliITSsegmentationSPD.h"
35 #include "AliITSsimulationSPD.h"
36 #include "AliLog.h"
37 #include "AliRun.h"
38
39 //#define DEBUG
40
41 ClassImp(AliITSsimulationSPD)
42 ////////////////////////////////////////////////////////////////////////
43 //  Version: 1
44 //  Modified by Bjorn S. Nilsen, G.E. Bruno, H. Tydesjo
45 //  March-April 2006
46 //
47 //  Version: 0
48 //  Written by Boris Batyunya
49 //  December 20 1999
50 //
51 //
52 // AliITSsimulationSPD is to do the simulation of SPDs.
53 //
54 ////////////////////////////////////////////////////////////////////////
55
56 //______________________________________________________________________
57 AliITSsimulationSPD::AliITSsimulationSPD():
58 AliITSsimulation(),
59 fHis(0),
60 fSPDname(),
61 fCoupling(){
62     // Default constructor.
63     // Inputs:
64     //    none.
65     // Outputs:
66     //    none.
67     // Return:
68     //    A default constructed AliITSsimulationSPD class.
69
70     AliDebug(1,Form("Calling default constructor"));
71 //    Init();
72 }
73 //______________________________________________________________________
74 AliITSsimulationSPD::AliITSsimulationSPD(AliITSDetTypeSim *dettyp):
75 AliITSsimulation(dettyp),
76 fHis(0),
77 fSPDname(),
78 fCoupling(){
79     // standard constructor
80     // Inputs:
81     //    AliITSsegmentation *seg  A pointer to the segmentation class
82     //                             to be used for this simulation
83     //    AliITSCalibration     *resp A pointer to the responce class to
84     //                             be used for this simulation
85     // Outputs:
86     //    none.
87     // Return:
88     //    A default constructed AliITSsimulationSPD class.
89
90     AliDebug(1,Form("Calling standard constructor "));
91 //      AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
92 //      res->SetTemperature(0.0);
93 //      res->SetDistanceOverVoltage(0.0);
94     Init();
95 }
96 //______________________________________________________________________
97 void AliITSsimulationSPD::Init(){
98     // Initilization
99     // Inputs:
100     //    none.
101     // Outputs:
102     //    none.
103     // Return:
104     //    none.
105     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
106
107     SetModuleNumber(0);
108     SetEventNumber(0);
109     SetMap(new AliITSpList(GetNPixelsZ(),GetNPixelsX()));
110     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
111     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
112     Double_t bias = res->GetBiasVoltage();
113 //    cout << "Bias Voltage --> " << bias << endl; // dom    
114     res->SetDistanceOverVoltage(kmictocm*seg->Dy(),bias);
115 // set kind of coupling ("old" or "new")
116     char opt[20];
117     res->GetCouplingOption(opt);
118     char *old = strstr(opt,"old");
119     if (old) {
120         fCoupling=2;
121     } else {
122         fCoupling=1;
123     } // end if
124
125     // Get the calibration objects for each module(ladder)
126     GetCalibrationObjects(0); //RunNr 0 hard coded for now
127
128 }
129 //______________________________________________________________________
130 AliITSsimulationSPD::~AliITSsimulationSPD(){
131     // destructor
132     // Inputs:
133     //    none.
134     // Outputs:
135     //    none.
136     // Return:
137     //     none.
138
139     if (fHis) {
140         fHis->Delete(); 
141         delete fHis;     
142     } // end if fHis
143 }
144 //______________________________________________________________________
145 AliITSsimulationSPD::AliITSsimulationSPD(const 
146                                                    AliITSsimulationSPD 
147                                                    &s) : AliITSsimulation(s){
148     //     Copy Constructor
149     // Inputs:
150     //    AliITSsimulationSPD &s The original class for which
151     //                                this class is a copy of
152     // Outputs:
153     //    none.
154     // Return:
155
156     *this = s;
157     return;
158 }
159 //______________________________________________________________________
160 AliITSsimulationSPD&  AliITSsimulationSPD::operator=(const 
161                                            AliITSsimulationSPD &s){
162     //    Assignment operator
163     // Inputs:
164     //    AliITSsimulationSPD &s The original class for which
165     //                                this class is a copy of
166     // Outputs:
167     //    none.
168     // Return:
169
170     if(&s == this) return *this;
171     this->fHis = s.fHis;
172     fCoupling  = s.fCoupling;
173     fSPDname   = s.fSPDname;
174     return *this;
175 }
176 //______________________________________________________________________
177 AliITSsimulation&  AliITSsimulationSPD::operator=(const 
178                                            AliITSsimulation &s){
179     //    Assignment operator
180     // Inputs:
181     //    AliITSsimulationSPD &s The original class for which
182     //                                this class is a copy of
183     // Outputs:
184     //    none.
185     // Return:
186
187     if(&s == this) return *this;
188     Error("AliITSsimulationSPD","Not allowed to make a = with "
189           "AliITSsimulationSPD","Using default creater instead");
190
191     return *this;
192 }
193
194 //______________________________________________________________________
195 void AliITSsimulationSPD::GetCalibrationObjects(Int_t RunNr) {
196   AliCDBManager* man = AliCDBManager::Instance();
197   if(!man->IsDefaultStorageSet()) {
198     man->SetDefaultStorage("local://$ALICE_ROOT");
199   }
200   AliCDBEntry *entrySPD = man->Get("ITS/Calib/CalibSPD", RunNr);
201   if(!entrySPD){
202     AliWarning("Cannot find SPD calibration entry!");
203     return;
204   }
205   TObjArray *respSPD = (TObjArray *)entrySPD->GetObject();
206   if ((! respSPD)) {
207     AliWarning("Cannot get data from SPD database entry!");
208     return;
209   }
210   for (Int_t mod=0; mod<240; mod++) {
211     calObj[mod] = (AliITSCalibrationSPD*) respSPD->At(mod);
212   }
213 }
214
215 //______________________________________________________________________
216 void AliITSsimulationSPD::InitSimulationModule(Int_t module, Int_t event){
217     //  This function creates maps to build the list of tracks for each
218     //  summable digit. Inputs defined by base class.
219     //  Inputs:
220     //    Int_t module   // Module number to be simulated
221     //    Int_t event    // Event number to be simulated
222     //  Outputs:
223     //    none
224     //  Returns:
225     //    none
226
227     AliDebug(1,Form("(module=%d,event=%d)",module,event));
228     SetModuleNumber(module);
229     SetEventNumber(event);
230     ClearMap();
231 }
232 //_____________________________________________________________________
233 void AliITSsimulationSPD::SDigitiseModule(AliITSmodule *mod,Int_t,
234                                                Int_t event){
235     //  This function begins the work of creating S-Digits.  Inputs defined
236     //  by base class.
237     //  Inputs:
238     //    AliITSmodule *mod  //  module
239     //    Int_t              //  not used
240     //    Int_t event        //  Event number
241     //  Outputs:
242     //    none
243     //  Return:
244     //    test              //  test returns kTRUE if the module contained hits
245     //                      //  test returns kFALSE if it did not contain hits
246
247     AliDebug(1,Form("(mod=%p, ,event=%d)",mod,event));
248     if(!(mod->GetNhits())){
249         AliDebug(1,Form("In event %d module %d there are %d hits returning.",
250                         event, mod->GetIndex(),mod->GetNhits()));
251         return;// if module has no hits don't create Sdigits
252     } // end if
253     SetModuleNumber(mod->GetIndex());
254     SetEventNumber(event);
255     // HitToSDigit(mod);
256     HitToSDigitFast(mod);
257     WriteSDigits();
258     ClearMap();
259 }
260 //______________________________________________________________________
261 void AliITSsimulationSPD::WriteSDigits(){
262     //  This function adds each S-Digit to pList
263     //  Inputs:
264     //    none.
265     //  Outputs:
266     //    none.
267     //  Return:
268     //    none
269     Int_t ix, nix, iz, niz;
270     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
271
272     AliDebug(1,Form("Writing SDigits for module %d",GetModuleNumber()));
273     GetMap()->GetMaxMapIndex(niz, nix);
274     for(iz=0; iz<niz; iz++)for(ix=0; ix<nix; ix++){
275         if(GetMap()->GetSignalOnly(iz,ix)>0.0){
276             aliITS->AddSumDigit(*(GetMap()->GetpListItem(iz,ix)));
277             if(AliDebugLevel()>0) {
278               AliDebug(1,Form("%d, %d",iz,ix));
279               cout << *(GetMap()->GetpListItem(iz,ix)) << endl;
280             } // end if GetDebug
281         } // end if GetMap()->GetSignalOnly(iz,ix)>0.0
282     } // end for iz,ix
283     return; 
284 }
285 //______________________________________________________________________
286 void AliITSsimulationSPD::FinishSDigitiseModule(){
287     //  This function calls SDigitsToDigits which creates Digits from SDigits
288     //  Inputs:
289     //    none
290     //  Outputs:
291     //    none
292     //  Return
293     //    none
294
295     AliDebug(1,"()");
296     pListToDigits(); // Charge To Signal both adds noise and
297     ClearMap();
298     return;
299 }
300 //______________________________________________________________________
301 void AliITSsimulationSPD::DigitiseModule(AliITSmodule *mod,Int_t,
302                                               Int_t){
303     //  This function creates Digits straight from the hits and then adds
304     //  electronic noise to the digits before adding them to pList
305     //  Each of the input variables is passed along to HitToSDigit
306     //  Inputs:
307     //    AliITSmodule *mod     module
308     //    Int_t                 Dummy.
309     //    Int_t                 Dummy
310     //  Outputs:
311     //     none.
312     //  Return:
313     //    none.
314
315     AliDebug(1,Form("(mod=%p,,)",mod));
316     // HitToSDigit(mod);
317     HitToSDigitFast(mod);
318     RemoveDeadPixels(mod);
319     pListToDigits();
320     ClearMap();
321 }
322 //______________________________________________________________________
323 void AliITSsimulationSPD::HitToSDigit(AliITSmodule *mod){
324     // Does the charge distributions using Gaussian diffusion charge charing.
325     // Inputs:
326     //    AliITSmodule *mod  Pointer to this module
327     // Output:
328     //    none.
329     // Return:
330     //    none.
331     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
332     TObjArray *hits = mod->GetHits();
333     Int_t nhits = hits->GetEntriesFast();
334     Int_t h,ix,iz,i;
335     Int_t idtrack;
336     Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
337     Double_t x,y,z,t,tp,st,dt=0.2,el,sig,sigx,sigz,fda;
338     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
339     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
340     Double_t thick = 0.5*kmictocm*seg->Dy();  // Half Thickness
341     res->GetSigmaDiffusionAsymmetry(fda);
342
343     AliDebug(1,Form("(mod=%p) fCoupling=%d",mod,fCoupling));
344     if(nhits<=0) return;
345     for(h=0;h<nhits;h++){
346       if(AliDebugLevel()>0) {
347         AliDebug(1,Form("Hits, %d", h));
348         cout << *(mod->GetHit(h)) << endl;
349       } // end if GetDebug
350         if(!mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)) continue;
351         st = TMath::Sqrt(x1*x1+y1*y1+z1*z1);
352         if(st>0.0){
353             st = (Double_t)((Int_t)(st/kmictocm)); // number of microns
354             if(st<=1.0) st = 1.0;
355             dt = 1.0/st;
356             for(t=0.0;t<1.0;t+=dt){ // Integrate over t
357                 tp  = t+0.5*dt;
358                 x   = x0+x1*tp;
359                 y   = y0+y1*tp;
360                 z   = z0+z1*tp;
361                 if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
362                 el  = res->GeVToCharge((Double_t)(dt*de));
363                 if(GetDebug(1)){
364                     if(el<=0.0) cout<<"el="<<el<<" dt="<<dt
365                                     <<" de="<<de<<endl;
366                 } // end if GetDebug
367                 sig = res->SigmaDiffusion1D(TMath::Abs(thick + y)); 
368                 //  SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
369                 sigx=sig;
370                 sigz=sig*fda;
371                 SpreadChargeAsym(x,z,ix,iz,el,sigx,sigz,idtrack,h);
372                 cout << "sigx, sigz, y "<< sigx << " " << sigz<< " " << TMath::Abs(thick + y) << endl;// ciccio
373             } // end for t
374         } else { // st == 0.0 deposit it at this point
375             x   = x0;
376             y   = y0;
377             z   = z0;
378             if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
379             el  = res->GeVToCharge((Double_t)de);
380             sig = res->SigmaDiffusion1D(TMath::Abs(thick + y));
381             // SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
382             sigx=sig;
383             sigz=sig*fda;
384             SpreadChargeAsym(x,z,ix,iz,el,sigx,sigz,idtrack,h);
385         } // end if st>0.0
386         // Coupling
387         switch (fCoupling) {
388         default:
389             break;
390         case 1: //case 3:
391             // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
392             for(i=0;i<GetMap()->GetEntries();i++) 
393                 if(GetMap()->GetpListItem(i)==0) continue;
394                 else{
395                     GetMap()->GetMapIndex(
396                               GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
397                     SetCoupling(iz,ix,idtrack,h);
398                 } // end for i
399             break;
400         case 2: // case 4:
401             // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
402             for(i=0;i<GetMap()->GetEntries();i++) 
403                 if(GetMap()->GetpListItem(i)==0) continue;
404                 else{
405                     GetMap()->GetMapIndex(
406                                 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
407                     SetCouplingOld(iz,ix,idtrack,h);
408                 } // end for i
409             break;
410         } // end switch
411     } // Loop over all hits h
412     if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
413 }
414 //______________________________________________________________________
415 void AliITSsimulationSPD::HitToSDigitFast(AliITSmodule *mod){
416     // Does the charge distributions using Gaussian diffusion charge charing.    // Inputs:
417     //    AliITSmodule *mod  Pointer to this module
418     // Output:
419     //    none.
420     // Return:
421     //    none.
422     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
423     const Int_t n10=10;
424     const Double_t ti[n10]={7.443716945e-3,2.166976971e-1,3.397047841e-1,
425                             4.325316833e-1,4.869532643e-1,5.130467358e-1,
426                             5.674683167e-1,6.602952159e-1,7.833023029e-1,
427                             9.255628306e-1};
428     const Double_t wi[n10]={1.477621124e-1,1.346333597e-1,1.095431813e-1,
429                             7.472567455e-2,3.333567215e-2,3.333567215e-2,
430                             7.472567455e-2,1.095431813e-1,1.346333597e-1,
431                             1.477621124e-1};
432     TObjArray *hits = mod->GetHits();
433     Int_t nhits = hits->GetEntriesFast();
434     Int_t h,ix,iz,i;
435     Int_t idtrack;
436     Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
437     Double_t x,y,z,t,st,el,sig,sigx,sigz,fda;
438     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
439     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
440     Double_t thick = 0.5*kmictocm*seg->Dy();  // Half thickness
441     res->GetSigmaDiffusionAsymmetry(fda);
442 //    cout << "Half Thickness " << thick << endl;  // dom
443 //    cout << "Diffusion asymm " << fda << endl;  // dom
444
445     AliDebug(1,Form("(mod=%p) fCoupling=%d",mod,fCoupling));
446     if(nhits<=0) return;
447     for(h=0;h<nhits;h++){
448       if(AliDebugLevel()>0) {
449         AliDebug(1,Form("Hits, %d", h));
450         cout << *(mod->GetHit(h)) << endl;
451       } // end if GetDebug
452         if(!mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)) continue;
453         st = TMath::Sqrt(x1*x1+y1*y1+z1*z1);
454         if(st>0.0) for(i=0;i<n10;i++){ // Integrate over t
455             t   = ti[i];
456             x   = x0+x1*t;
457             y   = y0+y1*t;
458             z   = z0+z1*t;
459                 if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
460                 // el  = res->GeVToCharge((Double_t)(dt*de));
461                 // el  = 1./n10*res->GeVToCharge((Double_t)de);
462                 el  = wi[i]*res->GeVToCharge((Double_t)de); 
463                 if(GetDebug(1)){
464                     if(el<=0.0) cout<<"el="<<el<<" wi["<<i<<"]="<<wi[i]
465                                     <<" de="<<de<<endl;
466                 } // end if GetDebug
467                 sig = res->SigmaDiffusion1D(TMath::Abs(thick + y));
468                 sigx=sig;
469                 sigz=sig*fda;
470                 //SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
471                 SpreadChargeAsym(x,z,ix,iz,el,sigx,sigz,idtrack,h);
472                 // cout << "sigx sigz " << sigx << " " << sigz << endl; // dom
473             } // end for i // End Integrate over t
474         else { // st == 0.0 deposit it at this point
475             x   = x0;
476             y   = y0;
477             z   = z0;
478             if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
479             el  = res->GeVToCharge((Double_t)de);
480             sig = res->SigmaDiffusion1D(TMath::Abs(thick + y));
481             //SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
482             sigx=sig;
483             sigz=sig*fda;
484             SpreadChargeAsym(x,z,ix,iz,el,sigx,sigz,idtrack,h);
485         } // end if st>0.0
486         // Coupling
487         switch (fCoupling) {
488         default:
489             break;
490         case 1: // case 3:
491             // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
492             for(i=0;i<GetMap()->GetEntries();i++)
493                 if(GetMap()->GetpListItem(i)==0) continue;
494                 else{
495                     GetMap()->GetMapIndex(
496                               GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
497                     SetCoupling(iz,ix,idtrack,h);
498                 } // end for i
499             break;
500         case 2: // case 4:
501             // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
502             for(i=0;i<GetMap()->GetEntries();i++)
503                 if(GetMap()->GetpListItem(i)==0) continue;
504                 else{
505                     GetMap()->GetMapIndex(
506                                 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);                    SetCouplingOld(iz,ix,idtrack,h);
507                 } // end for i
508             break;
509         } // end switch
510     } // Loop over all hits h
511     if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
512 }
513 //______________________________________________________________________
514 void AliITSsimulationSPD::SpreadCharge(Double_t x0,Double_t z0,
515                                             Int_t ix0,Int_t iz0,
516                                             Double_t el,Double_t sig,Int_t t,
517                                             Int_t hi){
518     // Spreads the charge over neighboring cells. Assume charge is distributed
519     // as charge(x,z) = (el/2*pi*sig*sig)*exp(-arg)
520     // arg=((x-x0)*(x-x0)/2*sig*sig)+((z-z0*z-z0)/2*sig*sig)
521     // Defined this way, the integral over all x and z is el.
522     // Inputs:
523     //    Double_t x0   x position of point where charge is liberated
524     //    Double_t y0   y position of point where charge is liberated
525     //    Double_t z0   z position of point where charge is liberated
526     //    Int_t    ix0  row of cell corresponding to point x0
527     //    Int_t    iz0  columb of cell corresponding to point z0
528     //    Double_t el   number of electrons liberated in this step
529     //    Double_t sig  Sigma difusion for this step (y0 dependent)
530     //    Int_t    t    track number
531     //    Int_t    ti   hit track index number
532     //    Int_t    hi   hit "hit" index number
533     // Outputs:
534     //     none.
535     // Return:
536     //     none.
537     const Int_t knx = 3,knz = 2;
538     const Double_t kRoot2 = 1.414213562; // Sqrt(2).
539     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
540     Int_t ix,iz,ixs,ixe,izs,ize;
541     Float_t x,z;
542     Double_t x1,x2,z1,z2,s,sp;
543     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
544
545
546     if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
547                          "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sig,t,hi);
548     if(sig<=0.0) { // if sig<=0 No diffusion to simulate.
549         GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
550         if(GetDebug(2)){
551             cout << "sig<=0.0=" << sig << endl;
552         } // end if GetDebug
553         return;
554     } // end if
555     sp = 1.0/(sig*kRoot2);
556     if(GetDebug(2)){
557         cout << "sig=" << sig << " sp=" << sp << endl;
558     } // end if GetDebug
559     ixs = TMath::Max(-knx+ix0,0);
560     ixe = TMath::Min(knx+ix0,seg->Npx()-1);
561     izs = TMath::Max(-knz+iz0,0);
562     ize = TMath::Min(knz+iz0,seg->Npz()-1);
563     for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
564         seg->DetToLocal(ix,iz,x,z); // pixel center
565         x1  = x;
566         z1  = z;
567         x2  = x1 + 0.5*kmictocm*seg->Dpx(ix); // Upper
568         x1 -= 0.5*kmictocm*seg->Dpx(ix);  // Lower
569         z2  = z1 + 0.5*kmictocm*seg->Dpz(iz); // Upper
570         z1 -= 0.5*kmictocm*seg->Dpz(iz);  // Lower
571         x1 -= x0; // Distance from where track traveled
572         x2 -= x0; // Distance from where track traveled
573         z1 -= z0; // Distance from where track traveled
574         z2 -= z0; // Distance from where track traveled
575         s   = 0.25; // Correction based on definision of Erfc
576         s  *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2);
577         if(GetDebug(3)){
578             cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
579                 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<< 
580                 " sp*x1="<<sp*x1<<" sp*x2="<<sp*x2<<" s="<<s;
581         } // end if GetDebug
582         s  *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2);
583         if(GetDebug(3)){
584             cout<<" sp*z1="<<sp*z1<<" sp*z2="<<sp*z2<<" s="<<s<< endl;
585         } // end if GetDebug
586         GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
587     } // end for ix, iz
588 }
589 //______________________________________________________________________
590 void AliITSsimulationSPD::SpreadChargeAsym(Double_t x0,Double_t z0,
591                                             Int_t ix0,Int_t iz0,
592                                             Double_t el,Double_t sigx,Double_t sigz,
593                                             Int_t t,Int_t hi){
594     // Spreads the charge over neighboring cells. Assume charge is distributed
595     // as charge(x,z) = (el/2*pi*sigx*sigz)*exp(-arg)
596     // arg=((x-x0)*(x-x0)/2*sigx*sigx)+((z-z0*z-z0)/2*sigz*sigz)
597     // Defined this way, the integral over all x and z is el.
598     // Inputs:
599     //    Double_t x0   x position of point where charge is liberated
600     //    Double_t y0   y position of point where charge is liberated
601     //    Double_t z0   z position of point where charge is liberated
602     //    Int_t    ix0  row of cell corresponding to point x0
603     //    Int_t    iz0  columb of cell corresponding to point z0
604     //    Double_t el   number of electrons liberated in this step
605     //    Double_t sigx Sigma difusion along x for this step (y0 dependent)
606     //    Double_t sigz Sigma difusion along z for this step (y0 dependent)
607     //    Int_t    t    track number
608     //    Int_t    ti   hit track index number
609     //    Int_t    hi   hit "hit" index number
610     // Outputs:
611     //     none.
612     // Return:
613     //     none.
614     const Int_t knx = 3,knz = 2;
615     const Double_t kRoot2 = 1.414213562; // Sqrt(2).
616     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
617     Int_t ix,iz,ixs,ixe,izs,ize;
618     Float_t x,z;
619     Double_t x1,x2,z1,z2,s,spx,spz;
620     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
621
622
623     if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
624                          "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sigx,sigz,t,hi);
625     if(sigx<=0.0 || sigz<=0.0) { // if sig<=0 No diffusion to simulate.
626         GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
627         if(GetDebug(2)){
628             cout << "sigx<=0.0=" << sigx << endl;
629             cout << "sigz<=0.0=" << sigz << endl;
630         } // end if GetDebug
631         return;
632     } // end if
633     spx = 1.0/(sigx*kRoot2);     spz = 1.0/(sigz*kRoot2);
634     if(GetDebug(2)){
635         cout << "sigx=" << sigx << " spx=" << spx << endl;
636         cout << "sigz=" << sigz << " spz=" << spz << endl;
637     } // end if GetDebug
638     ixs = TMath::Max(-knx+ix0,0);
639     ixe = TMath::Min(knx+ix0,seg->Npx()-1);
640     izs = TMath::Max(-knz+iz0,0);
641     ize = TMath::Min(knz+iz0,seg->Npz()-1);
642     for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
643         seg->DetToLocal(ix,iz,x,z); // pixel center
644         x1  = x;
645         z1  = z;
646         x2  = x1 + 0.5*kmictocm*seg->Dpx(ix); // Upper
647         x1 -= 0.5*kmictocm*seg->Dpx(ix);  // Lower
648         z2  = z1 + 0.5*kmictocm*seg->Dpz(iz); // Upper
649         z1 -= 0.5*kmictocm*seg->Dpz(iz);  // Lower
650         x1 -= x0; // Distance from where track traveled
651         x2 -= x0; // Distance from where track traveled
652         z1 -= z0; // Distance from where track traveled
653         z2 -= z0; // Distance from where track traveled
654         s   = 0.25; // Correction based on definision of Erfc
655         s  *= TMath::Erfc(spx*x1) - TMath::Erfc(spx*x2);
656         if(GetDebug(3)){
657             cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
658                 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<< 
659                 " spx*x1="<<spx*x1<<" spx*x2="<<spx*x2<<" s="<<s;
660         } // end if GetDebug
661         s  *= TMath::Erfc(spz*z1) - TMath::Erfc(spz*z2);
662         if(GetDebug(3)){
663             cout<<" spz*z1="<<spz*z1<<" spz*z2="<<spz*z2<<" s="<<s<< endl;
664         } // end if GetDebug
665         GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
666     } // end for ix, iz
667 }
668 //______________________________________________________________________
669 void AliITSsimulationSPD::RemoveDeadPixels(AliITSmodule *mod){
670   Int_t module_nr = mod->GetIndex();
671   Int_t nr_dead = calObj[module_nr]->GetNrDead();
672   for (Int_t i=0; i<nr_dead; i++) {
673     GetMap()->DeleteHit(calObj[module_nr]->GetDeadColAt(i),calObj[module_nr]->GetDeadRowAt(i));
674   }
675 }
676 //______________________________________________________________________
677 void AliITSsimulationSPD::pListToDigits(){
678     // add noise and electronics, perform the zero suppression and add the
679     // digit to the list
680     // Inputs:
681     //    none.
682     // Outputs:
683     //    none.
684     // Return:
685     //    none.
686     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
687     Int_t j,ix,iz;
688     Double_t  electronics;
689     Double_t sig;
690     const Int_t    nmaxtrk=AliITSdigitSPD::GetNTracks();
691     static AliITSdigitSPD dig;
692     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
693     if(GetDebug(1)) Info("pListToDigits","()");
694     for(iz=0; iz<GetNPixelsZ(); iz++) for(ix=0; ix<GetNPixelsX(); ix++){
695 // NEW (for the moment plugged by hand, in the future possibly read from Data Base)
696 // here parametrize the efficiency of the pixel along the row for the test columns (1,9,17,25)
697 //        if(iz==1 || iz == 9 || iz == 17 || iz == 25) {
698 //        Double_t eff,p1=0.,p2=0.; 
699 //        Double_t x=ix;
700 //        switch (iz) {
701 //          case 1:   p1=0.63460;p2=0.42438E-01;break;  
702 //          case 9:   p1=0.41090;p2=0.75914E-01;break;
703 //        case 17:  p1=0.31883;p2=0.91502E-01;break;
704 //        case 25:  p1=0.48828;p2=0.57975E-01;break;
705 //         } // end switch
706 //          eff=1.-p1*exp(-p2*x);
707 //          if (gRandom->Rndm() >= eff) continue;
708 //        } // end  if 
709  // END parametrize the efficiency
710  // 
711         electronics = res->ApplyBaselineAndNoise();
712         UpdateMapNoise(ix,iz,electronics);
713         //
714         // Apply Threshold and write Digits.
715         sig = GetMap()->GetSignalOnly(iz,ix);
716         FillHistograms(ix,iz,sig+electronics);
717         if(GetDebug(3)){
718             cout<<sig<<"+"<<electronics<<">threshold("<<ix<<","<<iz
719                 <<")="<<GetThreshold() <<endl;
720         } // end if GetDebug
721         if (sig+electronics <= GetThreshold()) continue;
722         dig.SetCoord1(iz);
723         dig.SetCoord2(ix);
724         dig.SetSignal(1);
725         dig.SetSignalSPD((Int_t) GetMap()->GetSignal(iz,ix));
726         for(j=0;j<nmaxtrk;j++){
727             if (j<GetMap()->GetNEntries()) {
728                 dig.SetTrack(j,GetMap()->GetTrack(iz,ix,j));
729                 dig.SetHit(j,GetMap()->GetHit(iz,ix,j));
730             }else { // Default values
731                 dig.SetTrack(j,-3);
732                 dig.SetHit(j,-1);
733             } // end if GetMap()
734         } // end for j
735         if(GetDebug(3)){
736             cout<<iz<<","<<ix<<","<<*(GetMap()->GetpListItem(iz,ix))<<endl;
737         } // end if GetDebug
738         aliITS->AddSimDigit(0,&dig);
739     } //  for ix/iz
740 }
741 //______________________________________________________________________
742 void AliITSsimulationSPD::CreateHistograms(){
743     // create 1D histograms for tests
744     // Inputs:
745     //    none.
746     // Outputs:
747     //    none.
748     // Return:
749     //     none.
750
751     if(GetDebug(1)) Info("CreateHistograms","create histograms");
752
753     fHis = new TObjArray(GetNPixelsZ());
754     TString fSPDname("spd_");
755     for(Int_t i=0;i<GetNPixelsZ();i++) {
756         Char_t pixelz[4];
757         sprintf(pixelz,"%d",i);
758         fSPDname.Append(pixelz);
759         fHis->AddAt(new TH1F(fSPDname.Data(),"SPD maps",
760                              GetNPixelsX(),0.,(Double_t)GetNPixelsX()),i);
761     } // end for i
762 }
763 //______________________________________________________________________
764 void AliITSsimulationSPD::FillHistograms(Int_t ix,Int_t iz,Double_t v){
765     // Fill the histogram
766     // Inputs:
767     //    none.
768     // Outputs:
769     //    none.
770     // Return:
771     //     none.
772
773     if(!GetHistArray()) return; // Only fill if setup.
774     if(GetDebug(2)) Info("FillHistograms","fill histograms");
775     GetHistogram(iz)->Fill(ix,v);
776 }
777 //______________________________________________________________________
778 void AliITSsimulationSPD::ResetHistograms(){
779     // Reset histograms for this detector
780     // Inputs:
781     //    none.
782     // Outputs:
783     //    none.
784     // Return:
785     //     none.
786
787     if(!GetHistArray()) return; // Only fill if setup.
788     if(GetDebug(2)) Info("FillHistograms","fill histograms");
789     for ( int i=0;i<GetNPixelsZ();i++ ) {
790         if (fHis->At(i))    ((TH1F*)fHis->At(i))->Reset();
791     } // end for i
792 }
793
794 //______________________________________________________________________
795 void AliITSsimulationSPD::SetCoupling(Int_t row, Int_t col, Int_t ntrack,
796                                       Int_t idhit) {
797     //  Take into account the coupling between adiacent pixels.
798     //  The parameters probcol and probrow are the probability of the
799     //  signal in one pixel shared in the two adjacent pixels along
800     //  the column and row direction, respectively.
801     //  Note pList is goten via GetMap() and module is not need any more.
802     //  Otherwise it is identical to that coded by Tiziano Virgili (BSN).
803     //Begin_Html
804     /*
805       <img src="picts/ITS/barimodel_3.gif">
806       </pre>
807       <br clear=left>
808       <font size=+2 color=red>
809       <a href="mailto:tiziano.virgili@cern.ch"></a>.
810       </font>
811       <pre>
812     */
813     //End_Html
814     // Inputs:
815     //    Int_t row            z cell index
816     //    Int_t col            x cell index
817     //    Int_t ntrack         track incex number
818     //    Int_t idhit          hit index number
819     // Outputs:
820     //    none.
821     // Return:
822     //     none.
823     Int_t j1,j2,flag=0;
824     Double_t pulse1,pulse2;
825     Double_t couplR=0.0,couplC=0.0;
826     Double_t xr=0.;
827
828     GetCouplings(couplC,couplR);
829     if(GetDebug(3)) Info("SetCoupling","(row=%d,col=%d,ntrack=%d,idhit=%d) "
830                          "Calling SetCoupling couplR=%e couplC=%e",
831                          row,col,ntrack,idhit,couplR,couplC);
832     j1 = row;
833     j2 = col;
834     pulse1 = GetMap()->GetSignalOnly(row,col);
835     pulse2 = pulse1;
836     for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
837         do{
838             j1 += isign;
839             //   pulse1 *= couplR; 
840             xr = gRandom->Rndm();
841             //if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold())){
842             if ((j1<0) || (j1>GetNPixelsZ()-1) || (xr>couplR)){
843                 j1 = row;
844                 flag = 1;
845             }else{
846                 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
847                 //  flag = 0;
848                 flag = 1; // only first next!!
849             } // end if
850         } while(flag == 0);
851         // loop in column direction
852         do{
853             j2 += isign;
854             // pulse2 *= couplC; 
855             xr = gRandom->Rndm();
856             //if((j2<0)||j2>(GetNPixelsX()-1)||pulse2<GetThreshold()){
857             if ((j2<0) || (j2>GetNPixelsX()-1) || (xr>couplC)){
858                 j2 = col;
859                 flag = 1;
860             }else{
861                 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
862                 //  flag = 0;
863                 flag = 1; // only first next!!
864             } // end if
865         } while(flag == 0);
866     } // for isign
867 }
868 //______________________________________________________________________
869 void AliITSsimulationSPD::SetCouplingOld(Int_t row, Int_t col,
870                 Int_t ntrack,Int_t idhit) {
871     //  Take into account the coupling between adiacent pixels.
872     //  The parameters probcol and probrow are the fractions of the
873     //  signal in one pixel shared in the two adjacent pixels along
874     //  the column and row direction, respectively.
875     //Begin_Html
876     /*
877       <img src="picts/ITS/barimodel_3.gif">
878       </pre>
879       <br clear=left>
880       <font size=+2 color=red>
881       <a href="mailto:Rocco.Caliandro@ba.infn.it"></a>.
882       </font>
883       <pre>
884     */
885     //End_Html
886     // Inputs:
887     //    Int_t row            z cell index
888     //    Int_t col            x cell index
889     //    Int_t ntrack         track incex number
890     //    Int_t idhit          hit index number
891     //    Int_t module         module number
892     // Outputs:
893     //    none.
894     // Return:
895     //     none.
896     Int_t j1,j2,flag=0;
897     Double_t pulse1,pulse2;
898     Double_t couplR=0.0,couplC=0.0;
899
900     GetCouplings(couplC,couplR);
901
902     //  Debugging ...
903 //    cout << "Threshold --> " << GetThreshold() << endl;  // dom
904 //    cout << "Couplings --> " << couplC << " " << couplR << endl;  //dom
905
906
907     if(GetDebug(3)) Info("SetCouplingOld","(row=%d,col=%d,ntrack=%d,idhit=%d) "
908                          "Calling SetCoupling couplR=%e couplC=%e",
909                          row,col,ntrack,idhit,couplR,couplC);
910     for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
911     pulse1 = GetMap()->GetSignalOnly(row,col);
912     pulse2 = pulse1;
913     j1 = row;
914     j2 = col;
915         do{
916             j1 += isign;
917             pulse1 *= couplR;
918             if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold())){
919                 pulse1 = GetMap()->GetSignalOnly(row,col);
920                 j1 = row;
921                 flag = 1;
922             }else{
923                 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
924                 // flag = 0;
925                 flag = 1;  // only first next !!
926             } // end if
927         } while(flag == 0);
928         // loop in column direction
929         do{
930             j2 += isign;
931             pulse2 *= couplC;
932             if((j2<0)||(j2>(GetNPixelsX()-1))||(pulse2<GetThreshold())){
933                 pulse2 = GetMap()->GetSignalOnly(row,col);
934                 j2 = col;
935                 flag = 1;
936             }else{
937                 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
938                 // flag = 0;
939                 flag = 1; // only first next!!
940             } // end if
941         } while(flag == 0);
942     } // for isign
943 }