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