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