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