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