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