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