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