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