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