]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSsimulationSPDdubna.cxx
SSD calibration
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationSPDdubna.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 <TMath.h>
23 #include <TParticle.h>
24 #include <TRandom.h>
25 #include <TString.h>
26 #include "AliITS.h"
27 #include "AliITSMapA2.h" 
28 #include "AliITSdigitSPD.h"
29 #include "AliITSgeom.h"
30 #include "AliITShit.h"
31 #include "AliITSmodule.h"
32 #include "AliITSpList.h"
33 #include "AliITSCalibrationSPD.h"
34 #include "AliITSsegmentationSPD.h"
35 #include "AliITSsimulationSPDdubna.h"
36 #include "AliLog.h"
37 #include "AliRun.h"
38
39 //#define DEBUG
40
41 ClassImp(AliITSsimulationSPDdubna)
42 ////////////////////////////////////////////////////////////////////////
43 // Version: 1
44 // Modified by Bjorn S. Nilsen
45 // Version: 0
46 // Written by Boris Batyunya
47 // December 20 1999
48 //
49 // AliITSsimulationSPDdubna is to do the simulation of SPDs.
50 //______________________________________________________________________
51 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna():
52 AliITSsimulation(),
53 fHis(0),
54 fSPDname(),
55 fCoupling(0){
56     // Default constructor.
57     // Inputs:
58     //    none.
59     // Outputs:
60     //    none.
61     // Return:
62     //    A default constructed AliITSsimulationSPDdubna class.
63
64     AliDebug(1,Form("Calling degault constructor"));
65 }
66 //______________________________________________________________________
67 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(AliITSDetTypeSim *dettyp, Int_t cup):
68 AliITSsimulation(dettyp),
69 fHis(0),
70 fSPDname(),
71 fCoupling(cup){
72     // standard constructor
73     // Inputs:
74     //    AliITSsegmentation *seg  A pointer to the segmentation class
75     //                             to be used for this simulation
76     //    AliITSCalibration     *resp A pointer to the responce class to
77     //                             be used for this simulation
78     //    Int_t              cup   The type of coupling to be used
79     //                             =1 uses SetCoupling, =2 uses SetCouplingOld
80     //                             With diffusion tured off
81     //                             =3 uses SetCoupling, =4 uses SetCouplingOld
82     //                             with diffusion on other, no coupling.
83     // Outputs:
84     //    none.
85     // Return:
86     //    A default constructed AliITSsimulationSPDdubna class.
87
88     AliDebug(1,
89              Form("Calling degault constructor cup=%d",cup));
90     if(cup==1||cup==2){ // For the moment, remove defusion if Coupling is
91         // set.
92       AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
93       res->SetTemperature(0.0);
94       res->SetDistanceOverVoltage(0.0);
95     } // end if
96     Init();
97 }
98 //______________________________________________________________________
99 void AliITSsimulationSPDdubna::Init(){
100     // Initilization
101     // Inputs:
102     //    none.
103     // Outputs:
104     //    none.
105     // Return:
106     //    none.
107     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
108
109     SetModuleNumber(0);
110     SetEventNumber(0);
111     SetMap(new AliITSpList(GetNPixelsZ(),GetNPixelsX()));
112     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
113     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
114     res->SetDistanceOverVoltage(kmictocm*seg->Dy(),50.0);
115 }
116 //______________________________________________________________________
117 AliITSsimulationSPDdubna::~AliITSsimulationSPDdubna(){
118     // destructor
119     // Inputs:
120     //    none.
121     // Outputs:
122     //    none.
123     // Return:
124     //     none.
125
126     if (fHis) {
127         fHis->Delete(); 
128         delete fHis;     
129     } // end if fHis
130 }
131 //______________________________________________________________________
132 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(const 
133                                                    AliITSsimulationSPDdubna 
134                                                    &s) : AliITSsimulation(s){
135     //     Copy Constructor
136     // Inputs:
137     //    AliITSsimulationSPDdubna &s The original class for which
138     //                                this class is a copy of
139     // Outputs:
140     //    none.
141     // Return:
142
143     *this = s;
144     return;
145 }
146 //______________________________________________________________________
147 AliITSsimulationSPDdubna&  AliITSsimulationSPDdubna::operator=(const 
148                                            AliITSsimulationSPDdubna &s){
149     //    Assignment operator
150     // Inputs:
151     //    AliITSsimulationSPDdubna &s The original class for which
152     //                                this class is a copy of
153     // Outputs:
154     //    none.
155     // Return:
156
157     if(&s == this) return *this;
158     this->fHis = s.fHis;
159     fCoupling  = s.fCoupling;
160     fSPDname   = s.fSPDname;
161     return *this;
162 }
163 //______________________________________________________________________
164 AliITSsimulation&  AliITSsimulationSPDdubna::operator=(const 
165                                            AliITSsimulation &s){
166     //    Assignment operator
167     // Inputs:
168     //    AliITSsimulationSPDdubna &s The original class for which
169     //                                this class is a copy of
170     // Outputs:
171     //    none.
172     // Return:
173
174     if(&s == this) return *this;
175     Error("AliITSsimulationSPDdubna","Not allowed to make a = with "
176           "AliITSsimulationSPDdubna","Using default creater instead");
177
178     return *this;
179 }
180 //______________________________________________________________________
181 void AliITSsimulationSPDdubna::InitSimulationModule(Int_t module, Int_t event){
182     //  This function creates maps to build the list of tracks for each
183     //  summable digit. Inputs defined by base class.
184     //  Inputs:
185     //    Int_t module   // Module number to be simulated
186     //    Int_t event    // Event number to be simulated
187     //  Outputs:
188     //    none
189     //  Returns:
190     //    none
191
192     AliDebug(1,Form("(module=%d,event=%d)",module,event));
193     SetModuleNumber(module);
194     SetEventNumber(event);
195     ClearMap();
196 }
197 //_____________________________________________________________________
198 void AliITSsimulationSPDdubna::SDigitiseModule(AliITSmodule *mod,Int_t,
199                                                Int_t event){
200     //  This function begins the work of creating S-Digits.  Inputs defined
201     //  by base class.
202     //  Inputs:
203     //    AliITSmodule *mod  //  module
204     //    Int_t              //  not used
205     //    Int_t event        //  Event number
206     //  Outputs:
207     //    none
208     //  Return:
209     //    test              //  test returns kTRUE if the module contained hits
210     //                      //  test returns kFALSE if it did not contain hits
211
212     AliDebug(1,Form("(mod=%p, ,event=%d)",mod,event));
213     if(!(mod->GetNhits())){
214         AliDebug(1,Form("In event %d module %d there are %d hits returning.",
215                         event, mod->GetIndex(),mod->GetNhits()));
216         return;// if module has no hits don't create Sdigits
217     } // end if
218     SetModuleNumber(mod->GetIndex());
219     SetEventNumber(event);
220     HitToSDigit(mod);
221     WriteSDigits();
222     ClearMap();
223 }
224 //______________________________________________________________________
225 void AliITSsimulationSPDdubna::WriteSDigits(){
226     //  This function adds each S-Digit to pList
227     //  Inputs:
228     //    none.
229     //  Outputs:
230     //    none.
231     //  Return:
232     //    none
233     Int_t ix, nix, iz, niz;
234     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
235
236     AliDebug(1,Form("Writing SDigits for module %d",GetModuleNumber()));
237     GetMap()->GetMaxMapIndex(niz, nix);
238     for(iz=0; iz<niz; iz++)for(ix=0; ix<nix; ix++){
239         if(GetMap()->GetSignalOnly(iz,ix)>0.0){
240             aliITS->AddSumDigit(*(GetMap()->GetpListItem(iz,ix)));
241             if(AliDebugLevel()>0) {
242               AliDebug(1,Form("%d, %d",iz,ix));
243               cout << *(GetMap()->GetpListItem(iz,ix)) << endl;
244             } // end if GetDebug
245         } // end if GetMap()->GetSignalOnly(iz,ix)>0.0
246     } // end for iz,ix
247     return; 
248 }
249 //______________________________________________________________________
250 void AliITSsimulationSPDdubna::FinishSDigitiseModule(){
251     //  This function calls SDigitsToDigits which creates Digits from SDigits
252     //  Inputs:
253     //    none
254     //  Outputs:
255     //    none
256     //  Return
257     //    none
258
259     AliDebug(1,"()");
260     pListToDigits(); // Charge To Signal both adds noise and
261     ClearMap();
262     return;
263 }
264 //______________________________________________________________________
265 void AliITSsimulationSPDdubna::DigitiseModule(AliITSmodule *mod,Int_t,
266                                               Int_t){
267     //  This function creates Digits straight from the hits and then adds
268     //  electronic noise to the digits before adding them to pList
269     //  Each of the input variables is passed along to HitToSDigit
270     //  Inputs:
271     //    AliITSmodule *mod     module
272     //    Int_t                 Dummy.
273     //    Int_t                 Dummy
274     //  Outputs:
275     //     none.
276     //  Return:
277     //    none.
278
279     AliDebug(1,Form("(mod=%p,,)",mod));
280     HitToSDigit(mod);
281     pListToDigits();
282     ClearMap();
283 }
284 //______________________________________________________________________
285 void AliITSsimulationSPDdubna::HitToSDigit(AliITSmodule *mod){
286     // Does the charge distributions using Gaussian diffusion charge charing.
287     // Inputs:
288     //    AliITSmodule *mod  Pointer to this module
289     // Output:
290     //    none.
291     // Return:
292     //    none.
293     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
294     TObjArray *hits = mod->GetHits();
295     Int_t nhits = hits->GetEntriesFast();
296     Int_t h,ix,iz,i;
297     Int_t idtrack;
298     Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
299     Double_t x,y,z,t,tp,st,dt=0.2,el,sig;
300     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
301     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
302     Double_t thick = kmictocm*seg->Dy();
303
304     AliDebug(1,Form("(mod=%p) fCoupling=%d",mod,fCoupling));
305     if(nhits<=0) return;
306     for(h=0;h<nhits;h++){
307       if(AliDebugLevel()>0) {
308         AliDebug(1,Form("Hits, %d", h));
309         cout << *(mod->GetHit(h)) << endl;
310       } // end if GetDebug
311         if(!mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)) continue;
312         st = TMath::Sqrt(x1*x1+y1*y1+z1*z1);
313         if(st>0.0){
314             st = (Double_t)((Int_t)(st/kmictocm)); // number of microns
315             if(st<=1.0) st = 1.0;
316             dt = 1.0/st;
317             for(t=0.0;t<1.0;t+=dt){ // Integrate over t
318                 tp  = t+0.5*dt;
319                 x   = x0+x1*tp;
320                 y   = y0+y1*tp;
321                 z   = z0+z1*tp;
322                 if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
323                 el  = res->GeVToCharge((Double_t)(dt*de));
324                 if(GetDebug(1)){
325                     if(el<=0.0) cout<<"el="<<el<<" dt="<<dt
326                                     <<" de="<<de<<endl;
327                 } // end if GetDebug
328                 sig = res->SigmaDiffusion1D(thick + y);
329                 SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
330             } // end for t
331         } else { // st == 0.0 deposit it at this point
332             x   = x0;
333             y   = y0;
334             z   = z0;
335             if(!(seg->LocalToDet(x,z,ix,iz))) continue; // outside
336             el  = res->GeVToCharge((Double_t)de);
337             sig = res->SigmaDiffusion1D(thick + y);
338             SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
339         } // end if st>0.0
340         // Coupling
341         switch (fCoupling) {
342         default:
343             break;
344         case 1: case 3:
345             // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
346             for(i=0;i<GetMap()->GetEntries();i++) 
347                 if(GetMap()->GetpListItem(i)==0) continue;
348                 else{
349                     GetMap()->GetMapIndex(
350                               GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
351                     SetCoupling(iz,ix,idtrack,h);
352                 } // end for i
353             break;
354         case 2: case 4:
355             // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
356             for(i=0;i<GetMap()->GetEntries();i++) 
357                 if(GetMap()->GetpListItem(i)==0) continue;
358                 else{
359                     GetMap()->GetMapIndex(
360                                 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
361                     SetCouplingOld(iz,ix,idtrack,h);
362                 } // end for i
363             break;
364         } // end switch
365     } // Loop over all hits h
366     if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
367 }
368 //______________________________________________________________________
369 void AliITSsimulationSPDdubna::SpreadCharge(Double_t x0,Double_t z0,
370                                             Int_t ix0,Int_t iz0,
371                                             Double_t el,Double_t sig,Int_t t,
372                                             Int_t hi){
373     // Spreads the charge over neighboring cells. Assume charge is distributed
374     // as charge(x,z) = (el/2*pi*sig*sig)*exp(-arg)
375     // arg=((x-x0)*(x-x0)/2*sig*sig)+((z-z0*z-z0)/2*sig*sig)
376     // Defined this way, the integral over all x and z is el.
377     // Inputs:
378     //    Double_t x0   x position of point where charge is liberated
379     //    Double_t y0   y position of point where charge is liberated
380     //    Double_t z0   z position of point where charge is liberated
381     //    Int_t    ix0  row of cell corresponding to point x0
382     //    Int_t    iz0  columb of cell corresponding to point z0
383     //    Double_t el   number of electrons liberated in this step
384     //    Double_t sig  Sigma difusion for this step (y0 dependent)
385     //    Int_t    t    track number
386     //    Int_t    ti   hit track index number
387     //    Int_t    hi   hit "hit" index number
388     // Outputs:
389     //     none.
390     // Return:
391     //     none.
392     const Int_t knx = 3,knz = 2;
393     const Double_t kRoot2 = 1.414213562; // Sqrt(2).
394     const Double_t kmictocm = 1.0e-4; // convert microns to cm.
395     Int_t ix,iz,ixs,ixe,izs,ize;
396     Float_t x,z;
397     Double_t x1,x2,z1,z2,s,sp;
398     AliITSsegmentationSPD* seg = (AliITSsegmentationSPD*)GetSegmentationModel(0);
399
400
401     if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
402                          "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sig,t,hi);
403     if(sig<=0.0) { // if sig<=0 No diffusion to simulate.
404         GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
405         if(GetDebug(2)){
406             cout << "sig<=0.0=" << sig << endl;
407         } // end if GetDebug
408         return;
409     } // end if
410     sp = 1.0/(sig*kRoot2);
411     if(GetDebug(2)){
412         cout << "sig=" << sig << " sp=" << sp << endl;
413     } // end if GetDebug
414     ixs = TMath::Max(-knx+ix0,0);
415     ixe = TMath::Min(knx+ix0,seg->Npx()-1);
416     izs = TMath::Max(-knz+iz0,0);
417     ize = TMath::Min(knz+iz0,seg->Npz()-1);
418     for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
419         seg->DetToLocal(ix,iz,x,z); // pixel center
420         x1  = x;
421         z1  = z;
422         x2  = x1 + 0.5*kmictocm*seg->Dpx(ix); // Upper
423         x1 -= 0.5*kmictocm*seg->Dpx(ix);  // Lower
424         z2  = z1 + 0.5*kmictocm*seg->Dpz(iz); // Upper
425         z1 -= 0.5*kmictocm*seg->Dpz(iz);  // Lower
426         x1 -= x0; // Distance from where track traveled
427         x2 -= x0; // Distance from where track traveled
428         z1 -= z0; // Distance from where track traveled
429         z2 -= z0; // Distance from where track traveled
430         s   = 0.25; // Correction based on definision of Erfc
431         s  *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2);
432         if(GetDebug(3)){
433             cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
434                 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<< 
435                 " sp*x1="<<sp*x1<<" sp*x2="<<sp*x2<<" s="<<s;
436         } // end if GetDebug
437         s  *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2);
438         if(GetDebug(3)){
439             cout<<" sp*z1="<<sp*z1<<" sp*z2="<<sp*z2<<" s="<<s<< endl;
440         } // end if GetDebug
441         GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
442     } // end for ix, iz
443 }
444 //______________________________________________________________________
445 void AliITSsimulationSPDdubna::pListToDigits(){
446     // add noise and electronics, perform the zero suppression and add the
447     // digit to the list
448     // Inputs:
449     //    none.
450     // Outputs:
451     //    none.
452     // Return:
453     //    none.
454     static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
455     Int_t j,ix,iz;
456     Double_t  electronics;
457     Double_t sig;
458     const Int_t    nmaxtrk=AliITSdigitSPD::GetNTracks();
459     static AliITSdigitSPD dig;
460     AliITSCalibrationSPD* res = (AliITSCalibrationSPD*)GetCalibrationModel(fDetType->GetITSgeom()->GetStartSPD());
461     if(GetDebug(1)) Info("pListToDigits","()");
462     for(iz=0; iz<GetNPixelsZ(); iz++) for(ix=0; ix<GetNPixelsX(); ix++){
463         // Apply Noise/Dead channals and the like
464         if(res->IsPixelDead(GetModuleNumber(),ix,iz)) continue;
465         electronics = res->ApplyBaselineAndNoise();
466         UpdateMapNoise(ix,iz,electronics);
467         //
468         // Apply Threshold and write Digits.
469         sig = GetMap()->GetSignalOnly(iz,ix);
470         FillHistograms(ix,iz,sig+electronics);
471         if(GetDebug(3)){
472             cout<<sig<<"+"<<electronics<<">threshold("<<ix<<","<<iz
473                 <<")="<<GetThreshold() <<endl;
474         } // end if GetDebug
475         if (sig+electronics <= GetThreshold()) continue;
476         dig.SetCoord1(iz);
477         dig.SetCoord2(ix);
478         dig.SetSignal(1);
479         dig.SetSignalSPD((Int_t) GetMap()->GetSignal(iz,ix));
480         for(j=0;j<nmaxtrk;j++){
481             if (j<GetMap()->GetNEntries()) {
482                 dig.SetTrack(j,GetMap()->GetTrack(iz,ix,j));
483                 dig.SetHit(j,GetMap()->GetHit(iz,ix,j));
484             }else { // Default values
485                 dig.SetTrack(j,-3);
486                 dig.SetHit(j,-1);
487             } // end if GetMap()
488         } // end for j
489         if(GetDebug(3)){
490             cout<<iz<<","<<ix<<","<<*(GetMap()->GetpListItem(iz,ix))<<endl;
491         } // end if GetDebug
492         aliITS->AddSimDigit(0,&dig);
493     } //  for ix/iz
494 }
495 //______________________________________________________________________
496 void AliITSsimulationSPDdubna::CreateHistograms(){
497     // create 1D histograms for tests
498     // Inputs:
499     //    none.
500     // Outputs:
501     //    none.
502     // Return:
503     //     none.
504
505     if(GetDebug(1)) Info("CreateHistograms","create histograms");
506
507     fHis = new TObjArray(GetNPixelsZ());
508     TString fSPDname("spd_");
509     for(Int_t i=0;i<GetNPixelsZ();i++) {
510         Char_t pixelz[4];
511         sprintf(pixelz,"%d",i);
512         fSPDname.Append(pixelz);
513         fHis->AddAt(new TH1F(fSPDname.Data(),"SPD maps",
514                              GetNPixelsX(),0.,(Double_t)GetNPixelsX()),i);
515     } // end for i
516 }
517 //______________________________________________________________________
518 void AliITSsimulationSPDdubna::FillHistograms(Int_t ix,Int_t iz,Double_t v){
519     // Fill the histogram
520     // Inputs:
521     //    none.
522     // Outputs:
523     //    none.
524     // Return:
525     //     none.
526
527     if(!GetHistArray()) return; // Only fill if setup.
528     if(GetDebug(2)) Info("FillHistograms","fill histograms");
529     GetHistogram(iz)->Fill(ix,v);
530 }
531 //______________________________________________________________________
532 void AliITSsimulationSPDdubna::ResetHistograms(){
533     // Reset histograms for this detector
534     // Inputs:
535     //    none.
536     // Outputs:
537     //    none.
538     // Return:
539     //     none.
540
541     if(!GetHistArray()) return; // Only fill if setup.
542     if(GetDebug(2)) Info("FillHistograms","fill histograms");
543     for ( int i=0;i<GetNPixelsZ();i++ ) {
544         if (fHis->At(i))    ((TH1F*)fHis->At(i))->Reset();
545     } // end for i
546 }
547
548 //______________________________________________________________________
549 void AliITSsimulationSPDdubna::SetCoupling(Int_t row, Int_t col, Int_t ntrack,
550                                       Int_t idhit) {
551     //  Take into account the coupling between adiacent pixels.
552     //  The parameters probcol and probrow are the probability of the
553     //  signal in one pixel shared in the two adjacent pixels along
554     //  the column and row direction, respectively.
555     //  Note pList is goten via GetMap() and module is not need any more.
556     //  Otherwise it is identical to that coded by Tiziano Virgili (BSN).
557     //Begin_Html
558     /*
559       <img src="picts/ITS/barimodel_3.gif">
560       </pre>
561       <br clear=left>
562       <font size=+2 color=red>
563       <a href="mailto:tiziano.virgili@cern.ch"></a>.
564       </font>
565       <pre>
566     */
567     //End_Html
568     // Inputs:
569     //    Int_t row            z cell index
570     //    Int_t col            x cell index
571     //    Int_t ntrack         track incex number
572     //    Int_t idhit          hit index number
573     // Outputs:
574     //    none.
575     // Return:
576     //     none.
577     Int_t j1,j2,flag=0;
578     Double_t pulse1,pulse2;
579     Double_t couplR=0.0,couplC=0.0;
580     Double_t xr=0.;
581
582     GetCouplings(couplR,couplC);
583     if(GetDebug(3)) Info("SetCoupling","(row=%d,col=%d,ntrack=%d,idhit=%d) "
584                          "Calling SetCoupling couplR=%e couplC=%e",
585                          row,col,ntrack,idhit,couplR,couplC);
586     j1 = row;
587     j2 = col;
588     pulse1 = GetMap()->GetSignalOnly(row,col);
589     pulse2 = pulse1;
590     for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
591         do{
592             j1 += isign;
593             //   pulse1 *= couplR; 
594             xr = gRandom->Rndm();
595             //if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold())){
596             if ((j1<0) || (j1>GetNPixelsZ()-1) || (xr>couplR)){
597                 j1 = row;
598                 flag = 1;
599             }else{
600                 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
601                 //  flag = 0;
602                 flag = 1; // only first next!!
603             } // end if
604         } while(flag == 0);
605         // loop in column direction
606         do{
607             j2 += isign;
608             // pulse2 *= couplC; 
609             xr = gRandom->Rndm();
610             //if((j2<0)||j2>(GetNPixelsX()-1)||pulse2<GetThreshold()){
611             if ((j2<0) || (j2>GetNPixelsX()-1) || (xr>couplC)){
612                 j2 = col;
613                 flag = 1;
614             }else{
615                 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
616                 //  flag = 0;
617                 flag = 1; // only first next!!
618             } // end if
619         } while(flag == 0);
620     } // for isign
621 }
622 //______________________________________________________________________
623 void AliITSsimulationSPDdubna::SetCouplingOld(Int_t row, Int_t col,
624                 Int_t ntrack,Int_t idhit) {
625     //  Take into account the coupling between adiacent pixels.
626     //  The parameters probcol and probrow are the fractions of the
627     //  signal in one pixel shared in the two adjacent pixels along
628     //  the column and row direction, respectively.
629     //Begin_Html
630     /*
631       <img src="picts/ITS/barimodel_3.gif">
632       </pre>
633       <br clear=left>
634       <font size=+2 color=red>
635       <a href="mailto:Rocco.Caliandro@ba.infn.it"></a>.
636       </font>
637       <pre>
638     */
639     //End_Html
640     // Inputs:
641     //    Int_t row            z cell index
642     //    Int_t col            x cell index
643     //    Int_t ntrack         track incex number
644     //    Int_t idhit          hit index number
645     //    Int_t module         module number
646     // Outputs:
647     //    none.
648     // Return:
649     //     none.
650     Int_t j1,j2,flag=0;
651     Double_t pulse1,pulse2;
652     Double_t couplR=0.0,couplC=0.0;
653
654     GetCouplings(couplR,couplC);
655     if(GetDebug(3)) Info("SetCouplingOld","(row=%d,col=%d,ntrack=%d,idhit=%d) "
656                          "Calling SetCoupling couplR=%e couplC=%e",
657                          row,col,ntrack,idhit,couplR,couplC);
658     j1 = row;
659     j2 = col;
660     pulse1 = GetMap()->GetSignalOnly(row,col);
661     pulse2 = pulse1;
662     for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
663         do{
664             j1 += isign;
665             pulse1 *= couplR;
666             if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold())){
667                 pulse1 = GetMap()->GetSignalOnly(row,col);
668                 j1 = row;
669                 flag = 1;
670             }else{
671                 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
672                 flag = 0;
673             } // end if
674         } while(flag == 0);
675         // loop in column direction
676         do{
677             j2 += isign;
678             pulse2 *= couplC;
679             if((j2<0)||(j2>(GetNPixelsX()-1))||(pulse2<GetThreshold())){
680                 pulse2 = GetMap()->GetSignalOnly(row,col);
681                 j2 = col;
682                 flag = 1;
683             }else{
684                 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
685                 flag = 0;
686             } // end if
687         } while(flag == 0);
688     } // for isign
689 }