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