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