]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSsimulationSPD.cxx
Using references in the arguments of the methods, removing warnings
[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         // Coupling
402         switch (fCoupling) {
403         default:
404             break;
405         case 1: //case 3:
406             for(i=0;i<GetMap()->GetEntries();i++) 
407                 if(GetMap()->GetpListItem(i)==0) continue;
408                 else{
409                     GetMap()->GetMapIndex(GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
410                     SetCoupling(iz,ix,idtrack,h);
411                 } // end for i
412             break;
413         case 2: // case 4:
414             for(i=0;i<GetMap()->GetEntries();i++) 
415                 if(GetMap()->GetpListItem(i)==0) continue;
416                 else{
417                     GetMap()->GetMapIndex(GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
418                     SetCouplingOld(iz,ix,idtrack,h);
419                 } // end for i
420             break;
421         } // end switch
422     } // Loop over all hits h
423     if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
424 }
425 //______________________________________________________________________
426 void AliITSsimulationSPD::HitToSDigitFast(AliITSmodule *mod){
427     // Does the charge distributions using Gaussian diffusion charge charing.    // Inputs:
428     //    AliITSmodule *mod  Pointer to this module
429     // Output:
430     //    none.
431     // Return:
432     //    none.
433     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
434     const Int_t kn10=10;
435     const Double_t kti[kn10]={7.443716945e-3,2.166976971e-1,3.397047841e-1,
436                             4.325316833e-1,4.869532643e-1,5.130467358e-1,
437                             5.674683167e-1,6.602952159e-1,7.833023029e-1,
438                             9.255628306e-1};
439     const Double_t kwi[kn10]={1.477621124e-1,1.346333597e-1,1.095431813e-1,
440                             7.472567455e-2,3.333567215e-2,3.333567215e-2,
441                             7.472567455e-2,1.095431813e-1,1.346333597e-1,
442                             1.477621124e-1};
443     TObjArray *hits = mod->GetHits();
444     Int_t nhits = hits->GetEntriesFast();
445     Int_t h,ix,iz,i;
446     Int_t idtrack;
447     Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
448     Double_t x,y,z,t,st,el,sig,sigx,sigz,fda;
449     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
450     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
451     Double_t thick = 0.5*kmictocm*seg->Dy();  // Half thickness
452     res->GetSigmaDiffusionAsymmetry(fda);
453 //    cout << "Half Thickness " << thick << endl;  // dom
454 //    cout << "Diffusion asymm " << fda << endl;  // dom
455
456     AliDebug(1,Form("(mod=%p) fCoupling=%d",mod,fCoupling));
457     if(nhits<=0) return;
458     for(h=0;h<nhits;h++){
459       if(AliDebugLevel()>0) {
460         AliDebug(1,Form("Hits, %d", h));
461         cout << *(mod->GetHit(h)) << endl;
462       } // end if GetDebug
463         if(!mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)) continue;
464         st = TMath::Sqrt(x1*x1+y1*y1+z1*z1);
465         if(st>0.0) for(i=0;i<kn10;i++){ // Integrate over t
466             t   = kti[i];
467             x   = x0+x1*t;
468             y   = y0+y1*t;
469             z   = z0+z1*t;
470                 if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
471                 // el  = res->GeVToCharge((Double_t)(dt*de));
472                 // el  = 1./kn10*res->GeVToCharge((Double_t)de);
473                 el  = kwi[i]*res->GeVToCharge((Double_t)de); 
474                 if(GetDebug(1)){
475                     if(el<=0.0) cout<<"el="<<el<<" kwi["<<i<<"]="<<kwi[i]
476                                     <<" de="<<de<<endl;
477                 } // end if GetDebug
478                 sig = res->SigmaDiffusion1D(TMath::Abs(thick + y));
479                 sigx=sig;
480                 sigz=sig*fda;
481                 //SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
482                 SpreadChargeAsym(x,z,ix,iz,el,sigx,sigz,idtrack,h);
483                 // cout << "sigx sigz " << sigx << " " << sigz << endl; // dom
484             } // end for i // End Integrate over t
485         else { // st == 0.0 deposit it at this point
486             x   = x0;
487             y   = y0;
488             z   = z0;
489             if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
490             el  = res->GeVToCharge((Double_t)de);
491             sig = res->SigmaDiffusion1D(TMath::Abs(thick + y));
492             //SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
493             sigx=sig;
494             sigz=sig*fda;
495             SpreadChargeAsym(x,z,ix,iz,el,sigx,sigz,idtrack,h);
496         } // end if st>0.0
497         // Coupling
498         switch (fCoupling) {
499         default:
500             break;
501         case 1: // case 3:
502             for(i=0;i<GetMap()->GetEntries();i++)
503                 if(GetMap()->GetpListItem(i)==0) continue;
504                 else{
505                     GetMap()->GetMapIndex(GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
506                     SetCoupling(iz,ix,idtrack,h);
507                 } // end for i
508             break;
509         case 2: // case 4:
510             for(i=0;i<GetMap()->GetEntries();i++)
511                 if(GetMap()->GetpListItem(i)==0) continue;
512                 else{
513                     GetMap()->GetMapIndex(GetMap()->GetpListItem(i)->GetIndex(),iz,ix);  
514                     SetCouplingOld(iz,ix,idtrack,h);
515                 } // end for i
516             break;
517         } // end switch
518     } // Loop over all hits h
519     if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
520 }
521 //______________________________________________________________________
522 void AliITSsimulationSPD::SpreadCharge(Double_t x0,Double_t z0,
523                                             Int_t ix0,Int_t iz0,
524                                             Double_t el,Double_t sig,Int_t t,
525                                             Int_t hi){
526     // Spreads the charge over neighboring cells. Assume charge is distributed
527     // as charge(x,z) = (el/2*pi*sig*sig)*exp(-arg)
528     // arg=((x-x0)*(x-x0)/2*sig*sig)+((z-z0*z-z0)/2*sig*sig)
529     // Defined this way, the integral over all x and z is el.
530     // Inputs:
531     //    Double_t x0   x position of point where charge is liberated
532     //    Double_t y0   y position of point where charge is liberated
533     //    Double_t z0   z position of point where charge is liberated
534     //    Int_t    ix0  row of cell corresponding to point x0
535     //    Int_t    iz0  columb of cell corresponding to point z0
536     //    Double_t el   number of electrons liberated in this step
537     //    Double_t sig  Sigma difusion for this step (y0 dependent)
538     //    Int_t    t    track number
539     //    Int_t    ti   hit track index number
540     //    Int_t    hi   hit "hit" index number
541     // Outputs:
542     //     none.
543     // Return:
544     //     none.
545     const Int_t knx = 3,knz = 2;
546     const Double_t kRoot2 = 1.414213562; // Sqrt(2).
547     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
548     Int_t ix,iz,ixs,ixe,izs,ize;
549     Float_t x,z;
550     Double_t x1,x2,z1,z2,s,sp;
551     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
552
553
554     if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
555                          "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sig,t,hi);
556     if(sig<=0.0) { // if sig<=0 No diffusion to simulate.
557         GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
558         if(GetDebug(2)){
559             cout << "sig<=0.0=" << sig << endl;
560         } // end if GetDebug
561         return;
562     } // end if
563     sp = 1.0/(sig*kRoot2);
564     if(GetDebug(2)){
565         cout << "sig=" << sig << " sp=" << sp << endl;
566     } // end if GetDebug
567     ixs = TMath::Max(-knx+ix0,0);
568     ixe = TMath::Min(knx+ix0,seg->Npx()-1);
569     izs = TMath::Max(-knz+iz0,0);
570     ize = TMath::Min(knz+iz0,seg->Npz()-1);
571     for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
572         seg->DetToLocal(ix,iz,x,z); // pixel center
573         x1  = x;
574         z1  = z;
575         x2  = x1 + 0.5*kmictocm*seg->Dpx(ix); // Upper
576         x1 -= 0.5*kmictocm*seg->Dpx(ix);  // Lower
577         z2  = z1 + 0.5*kmictocm*seg->Dpz(iz); // Upper
578         z1 -= 0.5*kmictocm*seg->Dpz(iz);  // Lower
579         x1 -= x0; // Distance from where track traveled
580         x2 -= x0; // Distance from where track traveled
581         z1 -= z0; // Distance from where track traveled
582         z2 -= z0; // Distance from where track traveled
583         s   = 0.25; // Correction based on definision of Erfc
584         s  *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2);
585         if(GetDebug(3)){
586             cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
587                 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<< 
588                 " sp*x1="<<sp*x1<<" sp*x2="<<sp*x2<<" s="<<s;
589         } // end if GetDebug
590         s  *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2);
591         if(GetDebug(3)){
592             cout<<" sp*z1="<<sp*z1<<" sp*z2="<<sp*z2<<" s="<<s<< endl;
593         } // end if GetDebug
594         GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
595     } // end for ix, iz
596 }
597 //______________________________________________________________________
598 void AliITSsimulationSPD::SpreadChargeAsym(Double_t x0,Double_t z0,
599                                             Int_t ix0,Int_t iz0,
600                                             Double_t el,Double_t sigx,Double_t sigz,
601                                             Int_t t,Int_t hi){
602     // Spreads the charge over neighboring cells. Assume charge is distributed
603     // as charge(x,z) = (el/2*pi*sigx*sigz)*exp(-arg)
604     // arg=((x-x0)*(x-x0)/2*sigx*sigx)+((z-z0*z-z0)/2*sigz*sigz)
605     // Defined this way, the integral over all x and z is el.
606     // Inputs:
607     //    Double_t x0   x position of point where charge is liberated
608     //    Double_t y0   y position of point where charge is liberated
609     //    Double_t z0   z position of point where charge is liberated
610     //    Int_t    ix0  row of cell corresponding to point x0
611     //    Int_t    iz0  columb of cell corresponding to point z0
612     //    Double_t el   number of electrons liberated in this step
613     //    Double_t sigx Sigma difusion along x for this step (y0 dependent)
614     //    Double_t sigz Sigma difusion along z for this step (y0 dependent)
615     //    Int_t    t    track number
616     //    Int_t    ti   hit track index number
617     //    Int_t    hi   hit "hit" index number
618     // Outputs:
619     //     none.
620     // Return:
621     //     none.
622     const Int_t knx = 3,knz = 2;
623     const Double_t kRoot2 = 1.414213562; // Sqrt(2).
624     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
625     Int_t ix,iz,ixs,ixe,izs,ize;
626     Float_t x,z;
627     Double_t x1,x2,z1,z2,s,spx,spz;
628     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
629
630
631     if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
632                          "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sigx,sigz,t,hi);
633     if(sigx<=0.0 || sigz<=0.0) { // if sig<=0 No diffusion to simulate.
634         GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
635         if(GetDebug(2)){
636             cout << "sigx<=0.0=" << sigx << endl;
637             cout << "sigz<=0.0=" << sigz << endl;
638         } // end if GetDebug
639         return;
640     } // end if
641     spx = 1.0/(sigx*kRoot2);     spz = 1.0/(sigz*kRoot2);
642     if(GetDebug(2)){
643         cout << "sigx=" << sigx << " spx=" << spx << endl;
644         cout << "sigz=" << sigz << " spz=" << spz << endl;
645     } // end if GetDebug
646     ixs = TMath::Max(-knx+ix0,0);
647     ixe = TMath::Min(knx+ix0,seg->Npx()-1);
648     izs = TMath::Max(-knz+iz0,0);
649     ize = TMath::Min(knz+iz0,seg->Npz()-1);
650     for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
651         seg->DetToLocal(ix,iz,x,z); // pixel center
652         x1  = x;
653         z1  = z;
654         x2  = x1 + 0.5*kmictocm*seg->Dpx(ix); // Upper
655         x1 -= 0.5*kmictocm*seg->Dpx(ix);  // Lower
656         z2  = z1 + 0.5*kmictocm*seg->Dpz(iz); // Upper
657         z1 -= 0.5*kmictocm*seg->Dpz(iz);  // Lower
658         x1 -= x0; // Distance from where track traveled
659         x2 -= x0; // Distance from where track traveled
660         z1 -= z0; // Distance from where track traveled
661         z2 -= z0; // Distance from where track traveled
662         s   = 0.25; // Correction based on definision of Erfc
663         s  *= TMath::Erfc(spx*x1) - TMath::Erfc(spx*x2);
664         if(GetDebug(3)){
665             cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
666                 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<< 
667                 " spx*x1="<<spx*x1<<" spx*x2="<<spx*x2<<" s="<<s;
668         } // end if GetDebug
669         s  *= TMath::Erfc(spz*z1) - TMath::Erfc(spz*z2);
670         if(GetDebug(3)){
671             cout<<" spz*z1="<<spz*z1<<" spz*z2="<<spz*z2<<" s="<<s<< endl;
672         } // end if GetDebug
673         GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
674     } // end for ix, iz
675 }
676 //______________________________________________________________________
677 void AliITSsimulationSPD::RemoveDeadPixels(AliITSmodule *mod){
678     //    Removes dead pixels on each module (ladder)
679     // Inputs:
680     //    Module Index (0,239)
681     // Outputs:
682     //    none.
683     // Return:
684     //    none.
685
686   Int_t moduleNr = mod->GetIndex();
687   Int_t nrDead = fCalObj[moduleNr]->GetNrDead();
688   for (Int_t i=0; i<nrDead; i++) {
689     GetMap()->DeleteHit(fCalObj[moduleNr]->GetDeadColAt(i),fCalObj[moduleNr]->GetDeadRowAt(i));
690   }
691 }
692 //______________________________________________________________________
693 void AliITSsimulationSPD::FrompListToDigits(){
694     // add noise and electronics, perform the zero suppression and add the
695     // digit to the list
696     // Inputs:
697     //    none.
698     // Outputs:
699     //    none.
700     // Return:
701     //    none.
702     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
703     Int_t j,ix,iz;
704     Double_t  electronics;
705     Double_t sig;
706     const Int_t    knmaxtrk=AliITSdigit::GetNTracks();
707     static AliITSdigitSPD dig;
708     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
709     if(GetDebug(1)) Info("FrompListToDigits","()");
710     for(iz=0; iz<GetNPixelsZ(); iz++) for(ix=0; ix<GetNPixelsX(); ix++){
711 // NEW (for the moment plugged by hand, in the future possibly read from Data Base)
712 // here parametrize the efficiency of the pixel along the row for the test columns (1,9,17,25)
713 //        if(iz==1 || iz == 9 || iz == 17 || iz == 25) {
714 //        Double_t eff,p1=0.,p2=0.; 
715 //        Double_t x=ix;
716 //        switch (iz) {
717 //          case 1:   p1=0.63460;p2=0.42438E-01;break;  
718 //          case 9:   p1=0.41090;p2=0.75914E-01;break;
719 //        case 17:  p1=0.31883;p2=0.91502E-01;break;
720 //        case 25:  p1=0.48828;p2=0.57975E-01;break;
721 //         } // end switch
722 //          eff=1.-p1*exp(-p2*x);
723 //          if (gRandom->Rndm() >= eff) continue;
724 //        } // end  if 
725  // END parametrize the efficiency
726  // 
727         electronics = res->ApplyBaselineAndNoise();
728         UpdateMapNoise(ix,iz,electronics);
729         //
730         // Apply Threshold and write Digits.
731         sig = GetMap()->GetSignalOnly(iz,ix);
732         FillHistograms(ix,iz,sig+electronics);
733         if(GetDebug(3)){
734             cout<<sig<<"+"<<electronics<<">threshold("<<ix<<","<<iz
735                 <<")="<<GetThreshold() <<endl;
736         } // end if GetDebug
737         if (sig+electronics <= GetThreshold()) continue;
738         dig.SetCoord1(iz);
739         dig.SetCoord2(ix);
740         dig.SetSignal(1);
741
742 //        dig.SetSignalSPD((Int_t) GetMap()->GetSignal(iz,ix));
743         Double_t aSignal =  GetMap()->GetSignal(iz,ix);
744         if (TMath::Abs(aSignal)>2147483647.0) {
745           //PH 2147483647 is the max. integer
746           //PH This apparently is a problem which needs investigation
747           AliWarning(Form("Too big or too small signal value %f",aSignal));
748           aSignal = TMath::Sign((Double_t)2147483647,aSignal);
749         }
750         dig.SetSignalSPD((Int_t)aSignal);
751
752         for(j=0;j<knmaxtrk;j++){
753             if (j<GetMap()->GetNEntries()) {
754                 dig.SetTrack(j,GetMap()->GetTrack(iz,ix,j));
755                 dig.SetHit(j,GetMap()->GetHit(iz,ix,j));
756             }else { // Default values
757                 dig.SetTrack(j,-3);
758                 dig.SetHit(j,-1);
759             } // end if GetMap()
760         } // end for j
761         if(GetDebug(3)){
762             cout<<iz<<","<<ix<<","<<*(GetMap()->GetpListItem(iz,ix))<<endl;
763         } // end if GetDebug
764         aliITS->AddSimDigit(0,&dig);
765     } //  for ix/iz
766 }
767 //______________________________________________________________________
768 void AliITSsimulationSPD::CreateHistograms(){
769     // create 1D histograms for tests
770     // Inputs:
771     //    none.
772     // Outputs:
773     //    none.
774     // Return:
775     //     none.
776
777     if(GetDebug(1)) Info("CreateHistograms","create histograms");
778
779     fHis = new TObjArray(GetNPixelsZ());
780     TString fSPDname("spd_");
781     for(Int_t i=0;i<GetNPixelsZ();i++) {
782         Char_t pixelz[4];
783         sprintf(pixelz,"%d",i);
784         fSPDname.Append(pixelz);
785         fHis->AddAt(new TH1F(fSPDname.Data(),"SPD maps",
786                              GetNPixelsX(),0.,(Double_t)GetNPixelsX()),i);
787     } // end for i
788 }
789 //______________________________________________________________________
790 void AliITSsimulationSPD::FillHistograms(Int_t ix,Int_t iz,Double_t v){
791     // Fill the histogram
792     // Inputs:
793     //    none.
794     // Outputs:
795     //    none.
796     // Return:
797     //     none.
798
799     if(!GetHistArray()) return; // Only fill if setup.
800     if(GetDebug(2)) Info("FillHistograms","fill histograms");
801     GetHistogram(iz)->Fill(ix,v);
802 }
803 //______________________________________________________________________
804 void AliITSsimulationSPD::ResetHistograms(){
805     // Reset histograms for this detector
806     // Inputs:
807     //    none.
808     // Outputs:
809     //    none.
810     // Return:
811     //     none.
812
813     if(!GetHistArray()) return; // Only fill if setup.
814     if(GetDebug(2)) Info("FillHistograms","fill histograms");
815     for ( int i=0;i<GetNPixelsZ();i++ ) {
816         if (fHis->At(i))    ((TH1F*)fHis->At(i))->Reset();
817     } // end for i
818 }
819
820 //______________________________________________________________________
821 void AliITSsimulationSPD::SetCoupling(Int_t col, Int_t row, Int_t ntrack,
822                                       Int_t idhit) {
823     //  Take into account the coupling between adiacent pixels.
824     //  The parameters probcol and probrow are the probability of the
825     //  signal in one pixel shared in the two adjacent pixels along
826     //  the column and row direction, respectively.
827     //  Note pList is goten via GetMap() and module is not need any more.
828     //  Otherwise it is identical to that coded by Tiziano Virgili (BSN).
829     //Begin_Html
830     /*
831       <img src="picts/ITS/barimodel_3.gif">
832       </pre>
833       <br clear=left>
834       <font size=+2 color=red>
835       <a href="mailto:tiziano.virgili@cern.ch"></a>.
836       </font>
837       <pre>
838     */
839     //End_Html
840     // Inputs:
841     //    Int_t col            z cell index
842     //    Int_t row            x cell index
843     //    Int_t ntrack         track incex number
844     //    Int_t idhit          hit index number
845     // Outputs:
846     //    none.
847     // Return:
848     //     none.
849     Int_t j1,j2,flag=0;
850     Double_t pulse1,pulse2;
851     Double_t couplR=0.0,couplC=0.0;
852     Double_t xr=0.;
853
854     GetCouplings(couplC,couplR);
855     if(GetDebug(3)) Info("SetCoupling","(col=%d,row=%d,ntrack=%d,idhit=%d) "
856                          "Calling SetCoupling couplC=%e couplR=%e",
857                          col,row,ntrack,idhit,couplC,couplR);
858     j1 = col;
859     j2 = row;
860     pulse1 = GetMap()->GetSignalOnly(col,row);
861     pulse2 = pulse1;
862     for (Int_t isign=-1;isign<=1;isign+=2){// loop in col direction
863         do{
864             j1 += isign;
865             xr = gRandom->Rndm();
866             if ((j1<0) || (j1>GetNPixelsZ()-1) || (xr>couplC)){
867                 j1 = col;
868                 flag = 1;
869             }else{
870                 UpdateMapSignal(row,j1,ntrack,idhit,pulse1);
871                 //  flag = 0;
872                 flag = 1; // only first next!!
873             } // end if
874         } while(flag == 0);
875         // loop in row direction
876         do{
877             j2 += isign;
878             xr = gRandom->Rndm();
879             if ((j2<0) || (j2>GetNPixelsX()-1) || (xr>couplR)){
880                 j2 = row;
881                 flag = 1;
882             }else{
883                 UpdateMapSignal(j2,col,ntrack,idhit,pulse2);
884                 //  flag = 0;
885                 flag = 1; // only first next!!
886             } // end if
887         } while(flag == 0);
888     } // for isign
889 }
890 //______________________________________________________________________
891 void AliITSsimulationSPD::SetCouplingOld(Int_t col, Int_t row,
892                 Int_t ntrack,Int_t idhit) {
893     //  Take into account the coupling between adiacent pixels.
894     //  The parameters probcol and probrow are the fractions of the
895     //  signal in one pixel shared in the two adjacent pixels along
896     //  the column and row direction, respectively.
897     //Begin_Html
898     /*
899       <img src="picts/ITS/barimodel_3.gif">
900       </pre>
901       <br clear=left>
902       <font size=+2 color=red>
903       <a href="mailto:Rocco.Caliandro@ba.infn.it"></a>.
904       </font>
905       <pre>
906     */
907     //End_Html
908     // Inputs:
909     //    Int_t col            z cell index
910     //    Int_t row            x cell index
911     //    Int_t ntrack         track incex number
912     //    Int_t idhit          hit index number
913     //    Int_t module         module number
914     // Outputs:
915     //    none.
916     // Return:
917     //     none.
918     Int_t j1,j2,flag=0;
919     Double_t pulse1,pulse2;
920     Double_t couplR=0.0,couplC=0.0;
921
922     GetCouplings(couplC,couplR);
923
924     //  Debugging ...
925 //    cout << "Threshold --> " << GetThreshold() << endl;  // dom
926 //    cout << "Couplings --> " << couplC << " " << couplR << endl;  //dom
927
928
929     if(GetDebug(3)) Info("SetCouplingOld","(col=%d,row=%d,ntrack=%d,idhit=%d) "
930                          "Calling SetCoupling couplC=%e couplR=%e",
931                          col,row,ntrack,idhit,couplC,couplR);
932     for (Int_t isign=-1;isign<=1;isign+=2){// loop in col direction
933     pulse1 = GetMap()->GetSignalOnly(col,row);
934     pulse2 = pulse1;
935     j1 = col;
936     j2 = row;
937         do{
938             j1 += isign;
939             pulse1 *= couplC;
940             if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold())){
941                 pulse1 = GetMap()->GetSignalOnly(col,row);
942                 j1 = col;
943                 flag = 1;
944             }else{
945                 UpdateMapSignal(row,j1,ntrack,idhit,pulse1);
946                 // flag = 0;
947                 flag = 1;  // only first next !!
948             } // end if
949         } while(flag == 0);
950         // loop in row direction
951         do{
952             j2 += isign;
953             pulse2 *= couplR;
954             if((j2<0)||(j2>(GetNPixelsX()-1))||(pulse2<GetThreshold())){
955                 pulse2 = GetMap()->GetSignalOnly(col,row);
956                 j2 = row;
957                 flag = 1;
958             }else{
959                 UpdateMapSignal(j2,col,ntrack,idhit,pulse2);
960                 // flag = 0;
961                 flag = 1; // only first next!!
962             } // end if
963         } while(flag == 0);
964     } // for isign
965 }