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