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