]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCParamSR.cxx
Corrections needed when AliH2F was moved to CONTAINERS
[u/mrichter/AliRoot.git] / TPC / AliTPCParamSR.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 $Log$
18 Revision 1.10  2002/10/14 14:57:43  hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
20
21 Revision 1.7.4.1  2002/07/24 10:09:01  alibrary
22 Updating VirtualMC
23
24 Revision 1.9  2002/06/07 06:13:07  kowal2
25 Corrected calculation of the cross talk in the outer sectors
26
27 Revision 1.8  2002/06/05 15:37:31  kowal2
28 Added cross-talk from the wires beyond the first and the last rows
29
30 Revision 1.7  2002/03/18 17:59:13  kowal2
31 Chnges in the pad geometry - 3 pad lengths introduced.
32
33 Revision 1.6  2002/02/25 11:02:56  kowal2
34 Changes towards speeding up the code. Thanks to Marian Ivanov.
35
36 Revision 1.5  2001/12/06 07:49:30  kowal2
37 corrected number of pads calculation
38
39 Revision 1.4  2000/11/02 07:33:15  kowal2
40 Improvements of the code.
41
42 Revision 1.3  2000/06/30 12:07:50  kowal2
43 Updated from the TPC-PreRelease branch
44
45 Revision 1.2.4.2  2000/06/14 16:48:24  kowal2
46 Parameter setting improved. Removed compiler warnings
47
48 Revision 1.2.4.1  2000/06/09 07:55:39  kowal2
49
50 Updated defaults
51
52 Revision 1.2  2000/04/17 09:37:33  kowal2
53 removed obsolete AliTPCDigitsDisplay.C
54
55 Revision 1.1.4.2  2000/04/10 11:36:13  kowal2
56
57 New Detector parameters handling class
58
59 */
60
61 ///////////////////////////////////////////////////////////////////////
62 //  Manager and of geomety  classes for set: TPC                     //
63 //                                                                   //
64 //  !sectors are numbered from  0                                     //
65 //  !pad rows are numbered from 0                                     //
66 //  
67 //  27.7.   - AliTPCPaaramSr object for TPC 
68 //            TPC with straight pad rows 
69 //  Origin:  Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk // 
70 //                                                                   //  
71 ///////////////////////////////////////////////////////////////////////
72
73
74 #include <Riostream.h>
75 #include <TMath.h>
76 #include <TObject.h>
77 #include <AliTPCParamSR.h>
78 #include "AliTPCPRF2D.h"
79 #include "AliTPCRF1D.h"
80 #include "TH1.h"
81
82
83 ClassImp(AliTPCParamSR)
84 const static  Int_t kMaxRows=600;
85 const static  Float_t  kEdgeSectorSpace = 2.5;
86 const static Float_t kFacSigmaPadRow=3.;
87 const static Float_t kFacSigmaPad=3.;
88 const static Float_t kFacSigmaTime=3.;
89
90
91 AliTPCParamSR::AliTPCParamSR()
92 {   
93   //
94   //constructor set the default parameters
95   fInnerPRF=0;
96   fOuter1PRF=0;
97   fOuter2PRF=0;
98   fTimeRF = 0;
99   fFacSigmaPadRow = Float_t(kFacSigmaPadRow);
100   fFacSigmaPad = Float_t(kFacSigmaPad);
101   fFacSigmaTime = Float_t(kFacSigmaTime);
102   SetDefault();
103   Update();
104 }
105
106 AliTPCParamSR::~AliTPCParamSR()
107 {
108   //
109   //destructor destroy some dynmicaly alocated variables
110   if (fInnerPRF != 0) delete fInnerPRF;
111   if (fOuter1PRF != 0) delete fOuter1PRF;
112   if (fOuter2PRF != 0) delete fOuter2PRF;
113   if (fTimeRF != 0) delete fTimeRF;
114 }
115
116 void AliTPCParamSR::SetDefault()
117 {
118   //set default TPC param   
119   fbStatus = kFALSE;
120   AliTPCParam::SetDefault();  
121 }  
122
123 Int_t  AliTPCParamSR::CalcResponse(Float_t* xyz, Int_t * index, Int_t row)
124 {
125   //
126   //calculate bin response as function of the input position -x 
127   //return number of valid response bin
128   //
129   //we suppose that coordinate is expressed in float digits 
130   // it's mean coordinate system 8
131   //xyz[0] - float padrow xyz[1] is float pad  (center pad is number 0) and xyz[2] is float time bin
132   if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){ 
133     Error("AliTPCParamSR", "response function was not adjusted");
134     return -1;
135   }
136   
137   Float_t sfpadrow;   // sigma of response function
138   Float_t sfpad;      // sigma  of 
139   Float_t sftime= fFacSigmaTime*fTimeRF->GetSigma()/fZWidth;     //3 sigma of time response
140   if (index[1]<fNInnerSector){
141     sfpadrow =fFacSigmaPadRow*fInnerPRF->GetSigmaY()/fInnerPadPitchLength;
142     sfpad    =fFacSigmaPad*fInnerPRF->GetSigmaX()/fInnerPadPitchWidth;
143   }
144   else{
145   if(row<fNRowUp1){
146     sfpadrow =fFacSigmaPadRow*fOuter1PRF->GetSigmaY()/fOuter1PadPitchLength;
147     sfpad    =fFacSigmaPad*fOuter1PRF->GetSigmaX()/fOuterPadPitchWidth;}
148     else{
149       sfpadrow =fFacSigmaPadRow*fOuter2PRF->GetSigmaY()/fOuter2PadPitchLength;
150       sfpad    =fFacSigmaPad*fOuter2PRF->GetSigmaX()/fOuterPadPitchWidth;
151     }   
152   }
153
154   Int_t fpadrow = TMath::Max(TMath::Nint(index[2]+xyz[0]-sfpadrow),0);  //"first" padrow
155   Int_t fpad    = TMath::Nint(xyz[1]-sfpad);     //first pad
156   Int_t ftime   = TMath::Max(TMath::Nint(xyz[2]+GetZOffset()/GetZWidth()-sftime),0);  // first time
157   Int_t lpadrow = TMath::Min(TMath::Nint(index[2]+xyz[0]+sfpadrow),fpadrow+19);  //"last" padrow
158   lpadrow       = TMath::Min(GetNRow(index[1])-1,lpadrow);
159   Int_t lpad    = TMath::Min(TMath::Nint(xyz[1]+sfpad),fpad+19);     //last pad
160   Int_t ltime   = TMath::Min(TMath::Nint(xyz[2]+GetZOffset()/GetZWidth()+sftime),ftime+19);    // last time
161   ltime         = TMath::Min(ltime,GetMaxTBin()-1); 
162   // 
163   Int_t npads = GetNPads(index[1],row);
164   if (fpad<-npads/2) 
165     fpad = -npads/2;
166   if (lpad>npads/2) 
167     lpad= npads/2;
168   if (ftime<0) ftime=0;
169   // 
170   if (row>=0) { //if we are interesting about given pad row
171     if (fpadrow<=row) fpadrow =row;
172     else 
173       return 0;
174     if (lpadrow>=row) lpadrow = row;
175     else 
176       return 0;
177   }
178
179  
180   Float_t  padres[20][20];  //I don't expect bigger number of bins
181   Float_t  timeres[20];     
182   Int_t cindex3=0;
183   Int_t cindex=0;
184   Float_t cweight = 0;
185   if (fpadrow>=0) {
186   //calculate padresponse function    
187   Int_t padrow, pad;
188   for (padrow = fpadrow;padrow<=lpadrow;padrow++)
189     for (pad = fpad;pad<=lpad;pad++){
190       Float_t dy = (xyz[0]+Float_t(index[2]-padrow));
191       Float_t dx = (xyz[1]+Float_t(pad));
192       if (index[1]<fNInnerSector)
193         padres[padrow-fpadrow][pad-fpad]=fInnerPRF->GetPRF(dx*fInnerPadPitchWidth,dy*fInnerPadPitchLength);
194       else{
195         if(row<fNRowUp1){
196         padres[padrow-fpadrow][pad-fpad]=fOuter1PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter1PadPitchLength);}
197         else{
198           padres[padrow-fpadrow][pad-fpad]=fOuter2PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter2PadPitchLength);}}}
199   //calculate time response function
200   Int_t time;
201   for (time = ftime;time<=ltime;time++) 
202     timeres[time-ftime]= fTimeRF->GetRF((-xyz[2]+Float_t(time))*fZWidth);     
203   //write over threshold values to stack
204   for (padrow = fpadrow;padrow<=lpadrow;padrow++)
205     for (pad = fpad;pad<=lpad;pad++)
206       for (time = ftime;time<=ltime;time++){
207         cweight = timeres[time-ftime]*padres[padrow-fpadrow][pad-fpad];
208         if (cweight>fResponseThreshold) {
209           fResponseBin[cindex3]=padrow;
210           fResponseBin[cindex3+1]=pad;
211           fResponseBin[cindex3+2]=time;
212           cindex3+=3;  
213           fResponseWeight[cindex]=cweight;
214           cindex++;
215         }
216       }
217   }
218   fCurrentMax=cindex;   
219   return fCurrentMax;
220 }
221
222 void AliTPCParamSR::TransformTo8(Float_t *xyz, Int_t *index) const
223 {
224   //
225   // transformate point to digit coordinate
226   //
227   if (index[0]==0) Transform0to1(xyz,index);
228   if (index[0]==1) Transform1to2(xyz,index);
229   if (index[0]==2) Transform2to3(xyz,index);
230   if (index[0]==3) Transform3to4(xyz,index);
231   if (index[0]==4) Transform4to8(xyz,index);
232 }
233
234 void AliTPCParamSR::TransformTo2(Float_t *xyz, Int_t *index) const
235 {
236   //
237   //transformate point to rotated coordinate
238   //
239   //we suppose that   
240   if (index[0]==0) Transform0to1(xyz,index);
241   if (index[0]==1) Transform1to2(xyz,index);
242   if (index[0]==4) Transform4to3(xyz,index);
243   if (index[0]==8) {  //if we are in digit coordinate system transform to global
244     Transform8to4(xyz,index);
245     Transform4to3(xyz,index);  
246   }
247 }
248
249 void AliTPCParamSR::CRXYZtoXYZ(Float_t *xyz,
250                const Int_t &sector, const Int_t & padrow, Int_t option) const  
251 {  
252   //transform relative coordinates to absolute
253   Bool_t rel = ( (option&2)!=0);
254   Int_t index[2]={sector,padrow};
255   if (rel==kTRUE)      Transform4to3(xyz,index);//if the position is relative to pad row  
256   Transform2to1(xyz,index);
257 }
258
259 void AliTPCParamSR::XYZtoCRXYZ(Float_t *xyz,
260                              Int_t &sector, Int_t & padrow, Int_t option) const
261 {
262    //transform global position to the position relative to the sector padrow
263   //if option=0  X calculate absolute            calculate sector
264   //if option=1  X           absolute            use input sector
265   //if option=2  X           relative to pad row calculate sector
266   //if option=3  X           relative            use input sector
267   //!!!!!!!!! WE start to calculate rows from row = 0
268   Int_t index[2];
269   Bool_t rel = ( (option&2)!=0);  
270
271   //option 0 and 2  means that we don't have information about sector
272   if ((option&1)==0)   Transform0to1(xyz,index);  //we calculate sector number 
273   else
274     index[0]=sector;
275   Transform1to2(xyz,index);
276   Transform2to3(xyz,index);
277   //if we store relative position calculate position relative to pad row
278   if (rel==kTRUE) Transform3to4(xyz,index);
279   sector = index[0];
280   padrow = index[1];
281 }
282
283 Float_t AliTPCParamSR::GetPrimaryLoss(Float_t *x, Int_t *index, Float_t *angle)
284 {
285   //
286   //
287   Float_t padlength=GetPadPitchLength(index[1]);
288   Float_t a1=TMath::Sin(angle[0]);
289   a1*=a1;
290   Float_t a2=TMath::Sin(angle[1]);
291   a2*=a2;
292   Float_t length =padlength*TMath::Sqrt(1+a1+a2);
293   return length*fNPrimLoss;
294 }
295
296 Float_t AliTPCParamSR::GetTotalLoss(Float_t *x, Int_t *index, Float_t *angle)
297 {
298   //
299   //
300   Float_t padlength=GetPadPitchLength(index[1]);
301   Float_t a1=TMath::Sin(angle[0]);
302   a1*=a1;
303   Float_t a2=TMath::Sin(angle[1]);
304   a2*=a2;
305   Float_t length =padlength*TMath::Sqrt(1+a1+a2);
306   return length*fNTotalLoss;
307   
308 }
309
310
311 void AliTPCParamSR::GetClusterSize(Float_t *x, Int_t *index, Float_t *angle, Int_t mode, Float_t *sigma)
312 {
313   //
314   //return cluster sigma2 (x,y) for particle at position x
315   // in this case x coordinata is in drift direction
316   //and y in pad row direction
317   //we suppose that input coordinate system is digit system
318    
319   Float_t  xx;
320   Float_t lx[3] = {x[0],x[1],x[2]};
321   Int_t   li[3] = {index[0],index[1],index[2]};
322   TransformTo2(lx,li);
323   //  Float_t  sigmadiff;
324   sigma[0]=0;
325   sigma[1]=0;
326   
327   xx = lx[2];  //calculate drift length in cm
328   if (xx>0) {
329     sigma[0]+= xx*GetDiffL()*GetDiffL();
330     sigma[1]+= xx*GetDiffT()*GetDiffT(); 
331   }
332
333
334   //sigma[0]=sigma[1]=0;
335   if (GetTimeRF()!=0) sigma[0]+=GetTimeRF()->GetSigma()*GetTimeRF()->GetSigma();
336   if ( (index[1]<fNInnerSector) &&(GetInnerPRF()!=0))   
337     sigma[1]+=GetInnerPRF()->GetSigmaX()*GetInnerPRF()->GetSigmaX();
338   if ( (index[1]>=fNInnerSector) &&(index[2]<fNRowUp1) && (GetOuter1PRF()!=0))
339     sigma[1]+=GetOuter1PRF()->GetSigmaX()*GetOuter1PRF()->GetSigmaX();
340   if( (index[1]>=fNInnerSector) &&(index[2]>=fNRowUp1) && (GetOuter2PRF()!=0))
341     sigma[1]+=GetOuter2PRF()->GetSigmaX()*GetOuter2PRF()->GetSigmaX();
342
343
344   sigma[0]/= GetZWidth()*GetZWidth();
345   sigma[1]/=GetPadPitchWidth(index[0])*GetPadPitchWidth(index[0]);
346 }
347
348
349
350
351 void AliTPCParamSR::GetSpaceResolution(Float_t *x, Int_t *index, Float_t *angle, 
352                                        Float_t amplitude, Int_t mode, Float_t *sigma)
353 {
354   //
355   //
356   //
357   
358 }
359 Float_t  AliTPCParamSR::GetAmp(Float_t *x, Int_t *index, Float_t *angle)
360 {
361   //
362   //
363   //
364   return 0;
365 }
366
367 Float_t * AliTPCParamSR::GetAnglesAccMomentum(Float_t *x, Int_t * index, Float_t* momentum, Float_t *angle)
368 {
369   //
370   //calculate angle of track to padrow at given position
371   // for given magnetic field and momentum of the particle
372   //
373
374   TransformTo2(x,index);
375   AliDetectorParam::GetAnglesAccMomentum(x,index,momentum,angle);    
376   Float_t addangle = TMath::ASin(x[1]/GetPadRowRadii(index[1],index[2]));
377   angle[1] +=addangle;
378   return angle;                          
379 }
380
381          
382 Bool_t AliTPCParamSR::Update()
383 {
384   Int_t i;
385   if (AliTPCParam::Update()==kFALSE) return kFALSE;
386   fbStatus = kFALSE;
387
388  Float_t firstrow = fInnerRadiusLow + 2.225 ;   
389  for( i= 0;i<fNRowLow;i++)
390    {
391      Float_t x = firstrow + fInnerPadPitchLength*(Float_t)i;  
392      fPadRowLow[i]=x;
393      // number of pads per row
394      Float_t y = (x-0.5*fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount-
395        fInnerPadPitchWidth/2.;
396      // 0 and fNRowLow+1 reserved for cross talk rows
397      fYInner[i+1]  = x*tan(fInnerAngle/2.)-fInnerWireMount;
398      fNPadsLow[i] = 1+2*(Int_t)(y/fInnerPadPitchWidth) ;
399    }
400  // cross talk rows
401  fYInner[0]=(fPadRowLow[0]-fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
402  fYInner[fNRowLow+1]=(fPadRowLow[fNRowLow-1]+fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount; 
403  firstrow = fOuterRadiusLow + 1.6;
404  for(i=0;i<fNRowUp;i++)
405    {
406      if(i<fNRowUp1){
407        Float_t x = firstrow + fOuter1PadPitchLength*(Float_t)i; 
408        fPadRowUp[i]=x;
409     Float_t y =(x-0.5*fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
410           fOuterPadPitchWidth/2.;
411      fYOuter[i+1]= x*tan(fOuterAngle/2.)-fOuterWireMount;
412      fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
413      if(i==fNRowUp1-1) {
414        fLastWireUp1=fPadRowUp[i] +0.375;
415        firstrow = fPadRowUp[i] + 0.5*(fOuter1PadPitchLength+fOuter2PadPitchLength);
416      }
417      }
418      else
419        {
420          Float_t x = firstrow + fOuter2PadPitchLength*(Float_t)(i-64);
421          fPadRowUp[i]=x;
422 Float_t y =(x-0.5*fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
423            fOuterPadPitchWidth/2.;
424          fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ; 
425        }
426      fYOuter[i+1]  = fPadRowUp[i]*tan(fOuterAngle/2.)-fOuterWireMount;
427    }
428  // cross talk rows
429  fYOuter[0]=(fPadRowUp[0]-fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
430  fYOuter[fNRowUp+1]=(fPadRowUp[fNRowUp-1]+fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
431  fNtRows = fNInnerSector*fNRowLow+fNOuterSector*fNRowUp;
432  fbStatus = kTRUE;
433  return kTRUE;
434 }
435 Float_t AliTPCParamSR::GetYInner(Int_t irow) const
436 {
437   return fYInner[irow];
438 }
439 Float_t AliTPCParamSR::GetYOuter(Int_t irow) const
440 {
441   return fYOuter[irow];
442 }
443
444 void AliTPCParamSR::Streamer(TBuffer &R__b)
445 {
446    // Stream an object of class AliTPC.
447
448    if (R__b.IsReading()) {
449       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
450       //      TObject::Streamer(R__b);
451       AliTPCParam::Streamer(R__b);
452       //      if (R__v < 2) return;
453        Update();
454    } else {
455       R__b.WriteVersion(AliTPCParamSR::IsA());
456       //TObject::Streamer(R__b);  
457       AliTPCParam::Streamer(R__b);    
458    }
459 }
460 Int_t  AliTPCParamSR::CalcResponseFast(Float_t* xyz, Int_t * index, Int_t row)
461 {
462   //
463   //calculate bin response as function of the input position -x 
464   //return number of valid response bin
465   //
466   //we suppose that coordinate is expressed in float digits 
467   // it's mean coordinate system 8
468   //xyz[0] - electron position w.r.t. pad center, normalized to pad length,
469   //xyz[1] is float pad  (center pad is number 0) and xyz[2] is float time bin
470   if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){ 
471     Error("AliTPCParamSR", "response function was not adjusted");
472     return -1;
473   }
474   
475   const Int_t padn =  500;
476   const Float_t fpadn =  500.;
477   const Int_t timen = 500;
478   const Float_t ftimen = 500.;
479   const Int_t padrn = 500;
480   const Float_t fpadrn = 500.;
481
482  
483
484   static Float_t prfinner[2*padrn][5*padn];  //pad divided by 50
485   static Float_t prfouter1[2*padrn][5*padn];  //prfouter division
486   static Float_t prfouter2[2*padrn][5*padn];
487
488   static Float_t rftime[5*timen];         //time division
489   static Int_t blabla=0;
490   static Float_t zoffset=0;
491   static Float_t zwidth=0;
492   static Float_t zoffset2=0;
493   static TH1F * hdiff=0;
494   static TH1F * hdiff1=0;
495   static TH1F * hdiff2=0;
496   
497   if (blabla==0) {  //calculate Response function - only at the begginning
498     hdiff =new TH1F("prf_diff","prf_diff",10000,-1,1);
499     hdiff1 =new TH1F("no_repsonse1","no_response1",10000,-1,1);
500     hdiff2 =new TH1F("no_response2","no_response2",10000,-1,1);
501     
502     blabla=1;
503     zoffset = GetZOffset();
504     zwidth  = fZWidth;
505     zoffset2 = zoffset/zwidth;
506     for (Int_t i=0;i<5*timen;i++){
507       rftime[i] = fTimeRF->GetRF(((i-2.5*ftimen)/ftimen)*zwidth+zoffset);
508     }
509     for (Int_t i=0;i<5*padn;i++){    
510       for (Int_t j=0;j<2*padrn;j++){
511         prfinner[j][i] =
512           fInnerPRF->GetPRF((i-2.5*fpadn)/fpadn
513                             *fInnerPadPitchWidth,(j-fpadrn)/fpadrn*fInnerPadPitchLength);
514         prfouter1[j][i] =
515           fOuter1PRF->GetPRF((i-2.5*fpadn)/fpadn
516                             *fOuterPadPitchWidth,(j-fpadrn)/fpadrn*fOuter1PadPitchLength);
517
518         //
519         prfouter2[j][i] =
520           fOuter2PRF->GetPRF((i-2.5*fpadn)/fpadn
521                             *fOuterPadPitchWidth,(j-fpadrn)/fpadrn*fOuter2PadPitchLength);
522       }
523     }      
524   } // the above is calculated only once
525
526   // calculate central padrow, pad, time
527   Int_t npads = GetNPads(index[1],index[3]-1);
528   Int_t cpadrow = index[2]; // electrons are here
529   Int_t cpad    = TMath::Nint(xyz[1]);
530   Int_t ctime   = TMath::Nint(xyz[2]+zoffset2);
531   //calulate deviation
532   Float_t dpadrow = xyz[0];
533   Float_t dpad    = xyz[1]-cpad;
534   Float_t dtime   = xyz[2]+zoffset2-ctime;
535   Int_t cindex =0;
536   Int_t cindex3 =0;
537   Int_t maxt =GetMaxTBin();
538
539   Int_t fpadrow;
540   Int_t lpadrow;
541
542   if (row>=0) { //if we are interesting about given pad row
543     fpadrow = row-cpadrow;
544     lpadrow = row-cpadrow;
545   }else{
546     fpadrow = (index[2]>1) ? -1 :0;
547     lpadrow = (index[2]<GetNRow(index[1])-1) ? 1:0;
548   }
549   Int_t fpad =  (cpad > -npads/2+1) ? -2: -npads/2-cpad;
550   Int_t lpad =  (cpad < npads/2-1)  ?  2: npads/2-cpad;
551   Int_t ftime =  (ctime>1) ? -2: -ctime;
552   Int_t ltime =  (ctime<maxt-2) ? 2: maxt-ctime-1;
553
554   // cross talk from long pad to short one
555   if(row==fNRowUp1 && fpadrow==-1) {
556     dpadrow *= fOuter2PadPitchLength;
557     dpadrow += fOuterWWPitch;
558     dpadrow /= fOuter1PadPitchLength;
559   }    
560   // cross talk from short pad to long one
561   if(row==fNRowUp1+1 && fpadrow==1){ 
562     dpadrow *= fOuter1PadPitchLength;
563     if(dpadrow < 0.) dpadrow = -1.; //protection against 3rd wire
564     dpadrow += fOuterWWPitch;
565     dpadrow /= fOuter2PadPitchLength;
566     
567   }
568   // "normal"
569   Int_t apadrow = TMath::Nint((dpadrow-fpadrow)*fpadrn+fpadrn);
570   for (Int_t ipadrow = fpadrow; ipadrow<=lpadrow;ipadrow++){
571     if ( (apadrow<0) || (apadrow>=2*padrn)) 
572       continue;
573     Int_t apad= TMath::Nint((dpad-fpad)*fpadn+2.5*fpadn);
574     for (Int_t ipad = fpad; ipad<=lpad;ipad++){
575         Float_t cweight;
576         if (index[1]<fNInnerSector)
577           cweight=prfinner[apadrow][apad];
578         else{
579           if(row < fNRowUp1+1)
580             cweight=prfouter1[apadrow][apad];
581           else cweight=prfouter2[apadrow][apad];
582         }
583
584         //      if (cweight<fResponseThreshold) continue;
585         Int_t atime = TMath::Nint((dtime-ftime)*ftimen+2.5*ftimen);
586         for (Int_t itime = ftime;itime<=ltime;itime++){ 
587           Float_t cweight2 = cweight*rftime[atime];
588           if (cweight2>fResponseThreshold) {
589             fResponseBin[cindex3++]=cpadrow+ipadrow;
590             fResponseBin[cindex3++]=cpad+ipad;
591             fResponseBin[cindex3++]=ctime+itime;
592             fResponseWeight[cindex++]=cweight2;
593             
594             if (cweight2>100) 
595               {
596                 printf("Pici pici %d %f %d\n",ipad,dpad,apad);
597               }
598             
599           }
600           atime-=timen;
601         }
602         apad-= padn;    
603     }
604     apadrow-=padrn;
605   }
606   fCurrentMax=cindex;   
607   return fCurrentMax;    
608   
609 }
610
611
612
613
614
615
616