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