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